00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include <finiteVolume/volFields.H>
00027 #include <finiteVolume/surfaceFields.H>
00028 #include <finiteVolume/calculatedFvPatchFields.H>
00029 #include <finiteVolume/zeroGradientFvPatchFields.H>
00030 #include <finiteVolume/coupledFvPatchFields.H>
00031
00032
00033
00034 template<class Type>
00035 template<class Type2>
00036 void Foam::fvMatrix<Type>::addToInternalField
00037 (
00038 const unallocLabelList& addr,
00039 const Field<Type2>& pf,
00040 Field<Type2>& intf
00041 ) const
00042 {
00043 if (addr.size() != pf.size())
00044 {
00045 FatalErrorIn
00046 (
00047 "fvMatrix<Type>::addToInternalField(const unallocLabelList&, "
00048 "const Field&, Field&)"
00049 ) << "sizes of addressing and field are different"
00050 << abort(FatalError);
00051 }
00052
00053 forAll(addr, faceI)
00054 {
00055 intf[addr[faceI]] += pf[faceI];
00056 }
00057 }
00058
00059
00060 template<class Type>
00061 template<class Type2>
00062 void Foam::fvMatrix<Type>::addToInternalField
00063 (
00064 const unallocLabelList& addr,
00065 const tmp<Field<Type2> >& tpf,
00066 Field<Type2>& intf
00067 ) const
00068 {
00069 addToInternalField(addr, tpf(), intf);
00070 tpf.clear();
00071 }
00072
00073
00074 template<class Type>
00075 template<class Type2>
00076 void Foam::fvMatrix<Type>::subtractFromInternalField
00077 (
00078 const unallocLabelList& addr,
00079 const Field<Type2>& pf,
00080 Field<Type2>& intf
00081 ) const
00082 {
00083 if (addr.size() != pf.size())
00084 {
00085 FatalErrorIn
00086 (
00087 "fvMatrix<Type>::addToInternalField(const unallocLabelList&, "
00088 "const Field&, Field&)"
00089 ) << "sizes of addressing and field are different"
00090 << abort(FatalError);
00091 }
00092
00093 forAll(addr, faceI)
00094 {
00095 intf[addr[faceI]] -= pf[faceI];
00096 }
00097 }
00098
00099
00100 template<class Type>
00101 template<class Type2>
00102 void Foam::fvMatrix<Type>::subtractFromInternalField
00103 (
00104 const unallocLabelList& addr,
00105 const tmp<Field<Type2> >& tpf,
00106 Field<Type2>& intf
00107 ) const
00108 {
00109 subtractFromInternalField(addr, tpf(), intf);
00110 tpf.clear();
00111 }
00112
00113
00114 template<class Type>
00115 void Foam::fvMatrix<Type>::addBoundaryDiag
00116 (
00117 scalarField& diag,
00118 const direction solveCmpt
00119 ) const
00120 {
00121 forAll(internalCoeffs_, patchI)
00122 {
00123 addToInternalField
00124 (
00125 lduAddr().patchAddr(patchI),
00126 internalCoeffs_[patchI].component(solveCmpt),
00127 diag
00128 );
00129 }
00130 }
00131
00132
00133 template<class Type>
00134 void Foam::fvMatrix<Type>::addCmptAvBoundaryDiag(scalarField& diag) const
00135 {
00136 forAll(internalCoeffs_, patchI)
00137 {
00138 addToInternalField
00139 (
00140 lduAddr().patchAddr(patchI),
00141 cmptAv(internalCoeffs_[patchI]),
00142 diag
00143 );
00144 }
00145 }
00146
00147
00148 template<class Type>
00149 void Foam::fvMatrix<Type>::addBoundarySource
00150 (
00151 Field<Type>& source,
00152 const bool couples
00153 ) const
00154 {
00155 forAll(psi_.boundaryField(), patchI)
00156 {
00157 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
00158 const Field<Type>& pbc = boundaryCoeffs_[patchI];
00159
00160 if (!ptf.coupled())
00161 {
00162 addToInternalField(lduAddr().patchAddr(patchI), pbc, source);
00163 }
00164 else if (couples)
00165 {
00166 tmp<Field<Type> > tpnf = ptf.patchNeighbourField();
00167 const Field<Type>& pnf = tpnf();
00168
00169 const unallocLabelList& addr = lduAddr().patchAddr(patchI);
00170
00171 forAll(addr, facei)
00172 {
00173 source[addr[facei]] += cmptMultiply(pbc[facei], pnf[facei]);
00174 }
00175 }
00176 }
00177 }
00178
00179
00180
00181
00182 template<class Type>
00183 Foam::fvMatrix<Type>::fvMatrix
00184 (
00185 GeometricField<Type, fvPatchField, volMesh>& psi,
00186 const dimensionSet& ds
00187 )
00188 :
00189 lduMatrix(psi.mesh()),
00190 psi_(psi),
00191 dimensions_(ds),
00192 source_(psi.size(), pTraits<Type>::zero),
00193 internalCoeffs_(psi.mesh().boundary().size()),
00194 boundaryCoeffs_(psi.mesh().boundary().size()),
00195 faceFluxCorrectionPtr_(NULL)
00196 {
00197 if (debug)
00198 {
00199 Info<< "fvMatrix<Type>(GeometricField<Type, fvPatchField, volMesh>&,"
00200 " const dimensionSet&) : "
00201 "constructing fvMatrix<Type> for field " << psi_.name()
00202 << endl;
00203 }
00204
00205
00206 forAll (psi.mesh().boundary(), patchI)
00207 {
00208 internalCoeffs_.set
00209 (
00210 patchI,
00211 new Field<Type>
00212 (
00213 psi.mesh().boundary()[patchI].size(),
00214 pTraits<Type>::zero
00215 )
00216 );
00217
00218 boundaryCoeffs_.set
00219 (
00220 patchI,
00221 new Field<Type>
00222 (
00223 psi.mesh().boundary()[patchI].size(),
00224 pTraits<Type>::zero
00225 )
00226 );
00227 }
00228
00229 psi_.boundaryField().updateCoeffs();
00230 }
00231
00232
00233 template<class Type>
00234 Foam::fvMatrix<Type>::fvMatrix(const fvMatrix<Type>& fvm)
00235 :
00236 refCount(),
00237 lduMatrix(fvm),
00238 psi_(fvm.psi_),
00239 dimensions_(fvm.dimensions_),
00240 source_(fvm.source_),
00241 internalCoeffs_(fvm.internalCoeffs_),
00242 boundaryCoeffs_(fvm.boundaryCoeffs_),
00243 faceFluxCorrectionPtr_(NULL)
00244 {
00245 if (debug)
00246 {
00247 Info<< "fvMatrix<Type>::fvMatrix(const fvMatrix<Type>&) : "
00248 << "copying fvMatrix<Type> for field " << psi_.name()
00249 << endl;
00250 }
00251
00252 if (fvm.faceFluxCorrectionPtr_)
00253 {
00254 faceFluxCorrectionPtr_ = new
00255 GeometricField<Type, fvsPatchField, surfaceMesh>
00256 (
00257 *(fvm.faceFluxCorrectionPtr_)
00258 );
00259 }
00260 }
00261
00262
00263 #ifdef ConstructFromTmp
00264 template<class Type>
00265 Foam::fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type> >& tfvm)
00266 :
00267 refCount(),
00268 lduMatrix
00269 (
00270 const_cast<fvMatrix<Type>&>(tfvm()),
00271 tfvm.isTmp()
00272 ),
00273 psi_(tfvm().psi_),
00274 dimensions_(tfvm().dimensions_),
00275 source_
00276 (
00277 const_cast<fvMatrix<Type>&>(tfvm()).source_,
00278 tfvm.isTmp()
00279 ),
00280 internalCoeffs_
00281 (
00282 const_cast<fvMatrix<Type>&>(tfvm()).internalCoeffs_,
00283 tfvm.isTmp()
00284 ),
00285 boundaryCoeffs_
00286 (
00287 const_cast<fvMatrix<Type>&>(tfvm()).boundaryCoeffs_,
00288 tfvm.isTmp()
00289 ),
00290 faceFluxCorrectionPtr_(NULL)
00291 {
00292 if (debug)
00293 {
00294 Info<< "fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type> >&) : "
00295 << "copying fvMatrix<Type> for field " << psi_.name()
00296 << endl;
00297 }
00298
00299 if (tfvm().faceFluxCorrectionPtr_)
00300 {
00301 if (tfvm.isTmp())
00302 {
00303 faceFluxCorrectionPtr_ = tfvm().faceFluxCorrectionPtr_;
00304 tfvm().faceFluxCorrectionPtr_ = NULL;
00305 }
00306 else
00307 {
00308 faceFluxCorrectionPtr_ = new
00309 GeometricField<Type, fvsPatchField, surfaceMesh>
00310 (
00311 *(tfvm().faceFluxCorrectionPtr_)
00312 );
00313 }
00314 }
00315
00316 tfvm.clear();
00317 }
00318 #endif
00319
00320
00321 template<class Type>
00322 Foam::fvMatrix<Type>::fvMatrix
00323 (
00324 GeometricField<Type, fvPatchField, volMesh>& psi,
00325 Istream& is
00326 )
00327 :
00328 lduMatrix(psi.mesh()),
00329 psi_(psi),
00330 dimensions_(is),
00331 source_(is),
00332 internalCoeffs_(psi.mesh().boundary().size()),
00333 boundaryCoeffs_(psi.mesh().boundary().size()),
00334 faceFluxCorrectionPtr_(NULL)
00335 {
00336 if (debug)
00337 {
00338 Info<< "fvMatrix<Type>"
00339 "(GeometricField<Type, fvPatchField, volMesh>&, Istream&) : "
00340 "constructing fvMatrix<Type> for field " << psi_.name()
00341 << endl;
00342 }
00343
00344
00345 forAll (psi.mesh().boundary(), patchI)
00346 {
00347 internalCoeffs_.set
00348 (
00349 patchI,
00350 new Field<Type>
00351 (
00352 psi.mesh().boundary()[patchI].size(),
00353 pTraits<Type>::zero
00354 )
00355 );
00356
00357 boundaryCoeffs_.set
00358 (
00359 patchI,
00360 new Field<Type>
00361 (
00362 psi.mesh().boundary()[patchI].size(),
00363 pTraits<Type>::zero
00364 )
00365 );
00366 }
00367
00368 }
00369
00370
00371 template<class Type>
00372 Foam::fvMatrix<Type>::~fvMatrix()
00373 {
00374 if (debug)
00375 {
00376 Info<< "fvMatrix<Type>::~fvMatrix<Type>() : "
00377 << "destroying fvMatrix<Type> for field " << psi_.name()
00378 << endl;
00379 }
00380
00381 if (faceFluxCorrectionPtr_)
00382 {
00383 delete faceFluxCorrectionPtr_;
00384 }
00385 }
00386
00387
00388
00389
00390
00391
00392 template<class Type>
00393 void Foam::fvMatrix<Type>::setValues
00394 (
00395 const labelList& cellLabels,
00396 const Field<Type>& values
00397 )
00398 {
00399 const fvMesh& mesh = psi_.mesh();
00400
00401 const cellList& cells = mesh.cells();
00402 const unallocLabelList& own = mesh.owner();
00403 const unallocLabelList& nei = mesh.neighbour();
00404
00405 scalarField& Diag = diag();
00406
00407 forAll(cellLabels, i)
00408 {
00409 label celli = cellLabels[i];
00410
00411 psi_[celli] = values[i];
00412 source_[celli] = values[i]*Diag[celli];
00413
00414 if (symmetric() || asymmetric())
00415 {
00416 const cell& c = cells[celli];
00417
00418 forAll(c, j)
00419 {
00420 label facei = c[j];
00421
00422 if (mesh.isInternalFace(facei))
00423 {
00424 if (symmetric())
00425 {
00426 if (celli == own[facei])
00427 {
00428 source_[nei[facei]] -= upper()[facei]*values[i];
00429 }
00430 else
00431 {
00432 source_[own[facei]] -= upper()[facei]*values[i];
00433 }
00434
00435 upper()[facei] = 0.0;
00436 }
00437 else
00438 {
00439 if (celli == own[facei])
00440 {
00441 source_[nei[facei]] -= lower()[facei]*values[i];
00442 }
00443 else
00444 {
00445 source_[own[facei]] -= upper()[facei]*values[i];
00446 }
00447
00448 upper()[facei] = 0.0;
00449 lower()[facei] = 0.0;
00450 }
00451 }
00452 else
00453 {
00454 label patchi = mesh.boundaryMesh().whichPatch(facei);
00455
00456 if (internalCoeffs_[patchi].size())
00457 {
00458 label patchFacei =
00459 mesh.boundaryMesh()[patchi].whichFace(facei);
00460
00461 internalCoeffs_[patchi][patchFacei] =
00462 pTraits<Type>::zero;
00463
00464 boundaryCoeffs_[patchi][patchFacei] =
00465 pTraits<Type>::zero;
00466 }
00467 }
00468 }
00469 }
00470 }
00471 }
00472
00473
00474 template<class Type>
00475 void Foam::fvMatrix<Type>::setReference
00476 (
00477 const label celli,
00478 const Type& value,
00479 const bool forceReference
00480 )
00481 {
00482 if ((forceReference || psi_.needReference()) && celli >= 0)
00483 {
00484 source()[celli] += diag()[celli]*value;
00485 diag()[celli] += diag()[celli];
00486 }
00487 }
00488
00489
00490 template<class Type>
00491 void Foam::fvMatrix<Type>::relax(const scalar alpha)
00492 {
00493 if (alpha <= 0)
00494 {
00495 return;
00496 }
00497
00498 Field<Type>& S = source();
00499 scalarField& D = diag();
00500
00501
00502 scalarField D0(D);
00503
00504
00505 scalarField sumOff(D.size(), 0.0);
00506 sumMagOffDiag(sumOff);
00507
00508
00509 forAll(psi_.boundaryField(), patchI)
00510 {
00511 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
00512
00513 if (ptf.size())
00514 {
00515 const unallocLabelList& pa = lduAddr().patchAddr(patchI);
00516 Field<Type>& iCoeffs = internalCoeffs_[patchI];
00517
00518 if (ptf.coupled())
00519 {
00520 const Field<Type>& pCoeffs = boundaryCoeffs_[patchI];
00521
00522
00523
00524 forAll(pa, face)
00525 {
00526 D[pa[face]] += component(iCoeffs[face], 0);
00527 sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
00528 }
00529 }
00530 else
00531 {
00532
00533
00534
00535
00536 forAll(pa, face)
00537 {
00538 Type iCoeff0 = iCoeffs[face];
00539 iCoeffs[face] = cmptMag(iCoeffs[face]);
00540 sumOff[pa[face]] -= cmptMin(iCoeffs[face]);
00541 iCoeffs[face] /= alpha;
00542 S[pa[face]] +=
00543 cmptMultiply(iCoeffs[face] - iCoeff0, psi_[pa[face]]);
00544 }
00545 }
00546 }
00547 }
00548
00549
00550 max(D, D, sumOff);
00551
00552
00553 D /= alpha;
00554
00555
00556 forAll(psi_.boundaryField(), patchI)
00557 {
00558 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
00559
00560 if (ptf.size())
00561 {
00562 const unallocLabelList& pa = lduAddr().patchAddr(patchI);
00563 Field<Type>& iCoeffs = internalCoeffs_[patchI];
00564
00565 if (ptf.coupled())
00566 {
00567 forAll(pa, face)
00568 {
00569 D[pa[face]] -= component(iCoeffs[face], 0);
00570 }
00571 }
00572 }
00573 }
00574
00575
00576 S += (D - D0)*psi_.internalField();
00577 }
00578
00579
00580 template<class Type>
00581 void Foam::fvMatrix<Type>::relax()
00582 {
00583 if (psi_.mesh().relax(psi_.name()))
00584 {
00585 relax(psi_.mesh().relaxationFactor(psi_.name()));
00586 }
00587 }
00588
00589
00590 template<class Type>
00591 void Foam::fvMatrix<Type>::boundaryManipulate
00592 (
00593 typename GeometricField<Type, fvPatchField, volMesh>::
00594 GeometricBoundaryField& bFields
00595 )
00596 {
00597 forAll(bFields, patchI)
00598 {
00599 bFields[patchI].manipulateMatrix(*this);
00600 }
00601 }
00602
00603
00604 template<class Type>
00605 Foam::tmp<Foam::scalarField> Foam::fvMatrix<Type>::D() const
00606 {
00607 tmp<scalarField> tdiag(new scalarField(diag()));
00608 addCmptAvBoundaryDiag(tdiag());
00609 return tdiag;
00610 }
00611
00612
00613 template<class Type>
00614 Foam::tmp<Foam::Field<Type> > Foam::fvMatrix<Type>::DD() const
00615 {
00616 tmp<Field<Type> > tdiag(pTraits<Type>::one*diag());
00617
00618 forAll(psi_.boundaryField(), patchI)
00619 {
00620 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
00621
00622 if (!ptf.coupled() && ptf.size())
00623 {
00624 addToInternalField
00625 (
00626 lduAddr().patchAddr(patchI),
00627 internalCoeffs_[patchI],
00628 tdiag()
00629 );
00630 }
00631 }
00632
00633 return tdiag;
00634 }
00635
00636
00637 template<class Type>
00638 Foam::tmp<Foam::volScalarField> Foam::fvMatrix<Type>::A() const
00639 {
00640 tmp<volScalarField> tAphi
00641 (
00642 new volScalarField
00643 (
00644 IOobject
00645 (
00646 "A("+psi_.name()+')',
00647 psi_.instance(),
00648 psi_.mesh(),
00649 IOobject::NO_READ,
00650 IOobject::NO_WRITE
00651 ),
00652 psi_.mesh(),
00653 dimensions_/psi_.dimensions()/dimVol,
00654 zeroGradientFvPatchScalarField::typeName
00655 )
00656 );
00657
00658 tAphi().internalField() = D()/psi_.mesh().V();
00659 tAphi().correctBoundaryConditions();
00660
00661 return tAphi;
00662 }
00663
00664
00665 template<class Type>
00666 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
00667 Foam::fvMatrix<Type>::H() const
00668 {
00669 tmp<GeometricField<Type, fvPatchField, volMesh> > tHphi
00670 (
00671 new GeometricField<Type, fvPatchField, volMesh>
00672 (
00673 IOobject
00674 (
00675 "H("+psi_.name()+')',
00676 psi_.instance(),
00677 psi_.mesh(),
00678 IOobject::NO_READ,
00679 IOobject::NO_WRITE
00680 ),
00681 psi_.mesh(),
00682 dimensions_/dimVol,
00683 zeroGradientFvPatchScalarField::typeName
00684 )
00685 );
00686 GeometricField<Type, fvPatchField, volMesh>& Hphi = tHphi();
00687
00688
00689 for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
00690 {
00691 scalarField psiCmpt = psi_.internalField().component(cmpt);
00692
00693 scalarField boundaryDiagCmpt(psi_.size(), 0.0);
00694 addBoundaryDiag(boundaryDiagCmpt, cmpt);
00695 boundaryDiagCmpt.negate();
00696 addCmptAvBoundaryDiag(boundaryDiagCmpt);
00697
00698 Hphi.internalField().replace(cmpt, boundaryDiagCmpt*psiCmpt);
00699 }
00700
00701 Hphi.internalField() += lduMatrix::H(psi_.internalField()) + source_;
00702 addBoundarySource(Hphi.internalField());
00703
00704 Hphi.internalField() /= psi_.mesh().V();
00705 Hphi.correctBoundaryConditions();
00706
00707 typename Type::labelType validComponents
00708 (
00709 pow
00710 (
00711 psi_.mesh().solutionD(),
00712 pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero
00713 )
00714 );
00715
00716 for(direction cmpt=0; cmpt<Type::nComponents; cmpt++)
00717 {
00718 if (validComponents[cmpt] == -1)
00719 {
00720 Hphi.replace
00721 (
00722 cmpt,
00723 dimensionedScalar("0", Hphi.dimensions(), 0.0)
00724 );
00725 }
00726 }
00727
00728 return tHphi;
00729 }
00730
00731
00732 template<class Type>
00733 Foam::tmp<Foam::volScalarField> Foam::fvMatrix<Type>::H1() const
00734 {
00735 tmp<volScalarField> tH1
00736 (
00737 new volScalarField
00738 (
00739 IOobject
00740 (
00741 "H(1)",
00742 psi_.instance(),
00743 psi_.mesh(),
00744 IOobject::NO_READ,
00745 IOobject::NO_WRITE
00746 ),
00747 psi_.mesh(),
00748 dimensions_/(dimVol*psi_.dimensions()),
00749 zeroGradientFvPatchScalarField::typeName
00750 )
00751 );
00752 volScalarField& H1_ = tH1();
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771 H1_.internalField() = lduMatrix::H1();
00772
00773 H1_.internalField() /= psi_.mesh().V();
00774 H1_.correctBoundaryConditions();
00775
00776 return tH1;
00777 }
00778
00779
00780 template<class Type>
00781 Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh> >
00782 Foam::fvMatrix<Type>::
00783 flux() const
00784 {
00785 if (!psi_.mesh().fluxRequired(psi_.name()))
00786 {
00787 FatalErrorIn("fvMatrix<Type>::flux()")
00788 << "flux requested but " << psi_.name()
00789 << " not specified in the fluxRequired sub-dictionary"
00790 " of fvSchemes."
00791 << abort(FatalError);
00792 }
00793
00794
00795 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tfieldFlux
00796 (
00797 new GeometricField<Type, fvsPatchField, surfaceMesh>
00798 (
00799 IOobject
00800 (
00801 "flux("+psi_.name()+')',
00802 psi_.instance(),
00803 psi_.mesh(),
00804 IOobject::NO_READ,
00805 IOobject::NO_WRITE
00806 ),
00807 psi_.mesh(),
00808 dimensions()
00809 )
00810 );
00811 GeometricField<Type, fvsPatchField, surfaceMesh>& fieldFlux = tfieldFlux();
00812
00813 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
00814 {
00815 fieldFlux.internalField().replace
00816 (
00817 cmpt,
00818 lduMatrix::faceH(psi_.internalField().component(cmpt))
00819 );
00820 }
00821
00822 FieldField<Field, Type> InternalContrib = internalCoeffs_;
00823
00824 forAll(InternalContrib, patchI)
00825 {
00826 InternalContrib[patchI] =
00827 cmptMultiply
00828 (
00829 InternalContrib[patchI],
00830 psi_.boundaryField()[patchI].patchInternalField()
00831 );
00832 }
00833
00834 FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
00835
00836 forAll(NeighbourContrib, patchI)
00837 {
00838 if (psi_.boundaryField()[patchI].coupled())
00839 {
00840 NeighbourContrib[patchI] =
00841 cmptMultiply
00842 (
00843 NeighbourContrib[patchI],
00844 psi_.boundaryField()[patchI].patchNeighbourField()
00845 );
00846 }
00847 }
00848
00849 forAll(fieldFlux.boundaryField(), patchI)
00850 {
00851 fieldFlux.boundaryField()[patchI] =
00852 InternalContrib[patchI] - NeighbourContrib[patchI];
00853 }
00854
00855 if (faceFluxCorrectionPtr_)
00856 {
00857 fieldFlux += *faceFluxCorrectionPtr_;
00858 }
00859
00860 return tfieldFlux;
00861 }
00862
00863
00864
00865
00866 template<class Type>
00867 void Foam::fvMatrix<Type>::operator=(const fvMatrix<Type>& fvmv)
00868 {
00869 if (this == &fvmv)
00870 {
00871 FatalErrorIn("fvMatrix<Type>::operator=(const fvMatrix<Type>&)")
00872 << "attempted assignment to self"
00873 << abort(FatalError);
00874 }
00875
00876 if (&psi_ != &(fvmv.psi_))
00877 {
00878 FatalErrorIn("fvMatrix<Type>::operator=(const fvMatrix<Type>&)")
00879 << "different fields"
00880 << abort(FatalError);
00881 }
00882
00883 lduMatrix::operator=(fvmv);
00884 source_ = fvmv.source_;
00885 internalCoeffs_ = fvmv.internalCoeffs_;
00886 boundaryCoeffs_ = fvmv.boundaryCoeffs_;
00887
00888 if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
00889 {
00890 *faceFluxCorrectionPtr_ = *fvmv.faceFluxCorrectionPtr_;
00891 }
00892 else if (fvmv.faceFluxCorrectionPtr_)
00893 {
00894 faceFluxCorrectionPtr_ =
00895 new GeometricField<Type, fvsPatchField, surfaceMesh>
00896 (*fvmv.faceFluxCorrectionPtr_);
00897 }
00898 }
00899
00900
00901 template<class Type>
00902 void Foam::fvMatrix<Type>::operator=(const tmp<fvMatrix<Type> >& tfvmv)
00903 {
00904 operator=(tfvmv());
00905 tfvmv.clear();
00906 }
00907
00908
00909 template<class Type>
00910 void Foam::fvMatrix<Type>::negate()
00911 {
00912 lduMatrix::negate();
00913 source_.negate();
00914 internalCoeffs_.negate();
00915 boundaryCoeffs_.negate();
00916
00917 if (faceFluxCorrectionPtr_)
00918 {
00919 faceFluxCorrectionPtr_->negate();
00920 }
00921 }
00922
00923
00924 template<class Type>
00925 void Foam::fvMatrix<Type>::operator+=(const fvMatrix<Type>& fvmv)
00926 {
00927 checkMethod(*this, fvmv, "+=");
00928
00929 dimensions_ += fvmv.dimensions_;
00930 lduMatrix::operator+=(fvmv);
00931 source_ += fvmv.source_;
00932 internalCoeffs_ += fvmv.internalCoeffs_;
00933 boundaryCoeffs_ += fvmv.boundaryCoeffs_;
00934
00935 if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
00936 {
00937 *faceFluxCorrectionPtr_ += *fvmv.faceFluxCorrectionPtr_;
00938 }
00939 else if (fvmv.faceFluxCorrectionPtr_)
00940 {
00941 faceFluxCorrectionPtr_ = new
00942 GeometricField<Type, fvsPatchField, surfaceMesh>
00943 (
00944 *fvmv.faceFluxCorrectionPtr_
00945 );
00946 }
00947 }
00948
00949
00950 template<class Type>
00951 void Foam::fvMatrix<Type>::operator+=(const tmp<fvMatrix<Type> >& tfvmv)
00952 {
00953 operator+=(tfvmv());
00954 tfvmv.clear();
00955 }
00956
00957
00958 template<class Type>
00959 void Foam::fvMatrix<Type>::operator-=(const fvMatrix<Type>& fvmv)
00960 {
00961 checkMethod(*this, fvmv, "+=");
00962
00963 dimensions_ -= fvmv.dimensions_;
00964 lduMatrix::operator-=(fvmv);
00965 source_ -= fvmv.source_;
00966 internalCoeffs_ -= fvmv.internalCoeffs_;
00967 boundaryCoeffs_ -= fvmv.boundaryCoeffs_;
00968
00969 if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
00970 {
00971 *faceFluxCorrectionPtr_ -= *fvmv.faceFluxCorrectionPtr_;
00972 }
00973 else if (fvmv.faceFluxCorrectionPtr_)
00974 {
00975 faceFluxCorrectionPtr_ =
00976 new GeometricField<Type, fvsPatchField, surfaceMesh>
00977 (-*fvmv.faceFluxCorrectionPtr_);
00978 }
00979 }
00980
00981
00982 template<class Type>
00983 void Foam::fvMatrix<Type>::operator-=(const tmp<fvMatrix<Type> >& tfvmv)
00984 {
00985 operator-=(tfvmv());
00986 tfvmv.clear();
00987 }
00988
00989
00990 template<class Type>
00991 void Foam::fvMatrix<Type>::operator+=
00992 (
00993 const DimensionedField<Type, volMesh>& su
00994 )
00995 {
00996 checkMethod(*this, su, "+=");
00997 source() -= su.mesh().V()*su.field();
00998 }
00999
01000
01001 template<class Type>
01002 void Foam::fvMatrix<Type>::operator+=
01003 (
01004 const tmp<DimensionedField<Type, volMesh> >& tsu
01005 )
01006 {
01007 operator+=(tsu());
01008 tsu.clear();
01009 }
01010
01011
01012 template<class Type>
01013 void Foam::fvMatrix<Type>::operator+=
01014 (
01015 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01016 )
01017 {
01018 operator+=(tsu());
01019 tsu.clear();
01020 }
01021
01022
01023 template<class Type>
01024 void Foam::fvMatrix<Type>::operator-=
01025 (
01026 const DimensionedField<Type, volMesh>& su
01027 )
01028 {
01029 checkMethod(*this, su, "-=");
01030 source() += su.mesh().V()*su.field();
01031 }
01032
01033
01034 template<class Type>
01035 void Foam::fvMatrix<Type>::operator-=
01036 (
01037 const tmp<DimensionedField<Type, volMesh> >& tsu
01038 )
01039 {
01040 operator-=(tsu());
01041 tsu.clear();
01042 }
01043
01044
01045 template<class Type>
01046 void Foam::fvMatrix<Type>::operator-=
01047 (
01048 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01049 )
01050 {
01051 operator-=(tsu());
01052 tsu.clear();
01053 }
01054
01055
01056 template<class Type>
01057 void Foam::fvMatrix<Type>::operator+=
01058 (
01059 const dimensioned<Type>& su
01060 )
01061 {
01062 source() -= psi().mesh().V()*su;
01063 }
01064
01065
01066 template<class Type>
01067 void Foam::fvMatrix<Type>::operator-=
01068 (
01069 const dimensioned<Type>& su
01070 )
01071 {
01072 source() += psi().mesh().V()*su;
01073 }
01074
01075
01076 template<class Type>
01077 void Foam::fvMatrix<Type>::operator+=
01078 (
01079 const zeroField&
01080 )
01081 {}
01082
01083
01084 template<class Type>
01085 void Foam::fvMatrix<Type>::operator-=
01086 (
01087 const zeroField&
01088 )
01089 {}
01090
01091
01092 template<class Type>
01093 void Foam::fvMatrix<Type>::operator*=
01094 (
01095 const DimensionedField<scalar, volMesh>& dsf
01096 )
01097 {
01098 dimensions_ *= dsf.dimensions();
01099 lduMatrix::operator*=(dsf.field());
01100 source_ *= dsf.field();
01101
01102 forAll(boundaryCoeffs_, patchI)
01103 {
01104 scalarField psf
01105 (
01106 dsf.mesh().boundary()[patchI].patchInternalField(dsf.field())
01107 );
01108
01109 internalCoeffs_[patchI] *= psf;
01110 boundaryCoeffs_[patchI] *= psf;
01111 }
01112
01113 if (faceFluxCorrectionPtr_)
01114 {
01115 FatalErrorIn
01116 (
01117 "fvMatrix<Type>::operator*="
01118 "(const DimensionedField<scalar, volMesh>&)"
01119 ) << "cannot scale a matrix containing a faceFluxCorrection"
01120 << abort(FatalError);
01121 }
01122 }
01123
01124
01125 template<class Type>
01126 void Foam::fvMatrix<Type>::operator*=
01127 (
01128 const tmp<DimensionedField<scalar, volMesh> >& tdsf
01129 )
01130 {
01131 operator*=(tdsf());
01132 tdsf.clear();
01133 }
01134
01135
01136 template<class Type>
01137 void Foam::fvMatrix<Type>::operator*=
01138 (
01139 const tmp<volScalarField>& tvsf
01140 )
01141 {
01142 operator*=(tvsf());
01143 tvsf.clear();
01144 }
01145
01146
01147 template<class Type>
01148 void Foam::fvMatrix<Type>::operator*=
01149 (
01150 const dimensioned<scalar>& ds
01151 )
01152 {
01153 dimensions_ *= ds.dimensions();
01154 lduMatrix::operator*=(ds.value());
01155 source_ *= ds.value();
01156 internalCoeffs_ *= ds.value();
01157 boundaryCoeffs_ *= ds.value();
01158
01159 if (faceFluxCorrectionPtr_)
01160 {
01161 *faceFluxCorrectionPtr_ *= ds.value();
01162 }
01163 }
01164
01165
01166
01167
01168 template<class Type>
01169 void Foam::checkMethod
01170 (
01171 const fvMatrix<Type>& fvm1,
01172 const fvMatrix<Type>& fvm2,
01173 const char* op
01174 )
01175 {
01176 if (&fvm1.psi() != &fvm2.psi())
01177 {
01178 FatalErrorIn
01179 (
01180 "checkMethod(const fvMatrix<Type>&, const fvMatrix<Type>&)"
01181 ) << "incompatible fields for operation "
01182 << endl << " "
01183 << "[" << fvm1.psi().name() << "] "
01184 << op
01185 << " [" << fvm2.psi().name() << "]"
01186 << abort(FatalError);
01187 }
01188
01189 if (dimensionSet::debug && fvm1.dimensions() != fvm2.dimensions())
01190 {
01191 FatalErrorIn
01192 (
01193 "checkMethod(const fvMatrix<Type>&, const fvMatrix<Type>&)"
01194 ) << "incompatible dimensions for operation "
01195 << endl << " "
01196 << "[" << fvm1.psi().name() << fvm1.dimensions()/dimVolume << " ] "
01197 << op
01198 << " [" << fvm2.psi().name() << fvm2.dimensions()/dimVolume << " ]"
01199 << abort(FatalError);
01200 }
01201 }
01202
01203
01204 template<class Type>
01205 void Foam::checkMethod
01206 (
01207 const fvMatrix<Type>& fvm,
01208 const DimensionedField<Type, volMesh>& df,
01209 const char* op
01210 )
01211 {
01212 if (dimensionSet::debug && fvm.dimensions()/dimVolume != df.dimensions())
01213 {
01214 FatalErrorIn
01215 (
01216 "checkMethod(const fvMatrix<Type>&, const GeometricField<Type, "
01217 "fvPatchField, volMesh>&)"
01218 ) << "incompatible dimensions for operation "
01219 << endl << " "
01220 << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
01221 << op
01222 << " [" << df.name() << df.dimensions() << " ]"
01223 << abort(FatalError);
01224 }
01225 }
01226
01227
01228 template<class Type>
01229 void Foam::checkMethod
01230 (
01231 const fvMatrix<Type>& fvm,
01232 const dimensioned<Type>& dt,
01233 const char* op
01234 )
01235 {
01236 if (dimensionSet::debug && fvm.dimensions()/dimVolume != dt.dimensions())
01237 {
01238 FatalErrorIn
01239 (
01240 "checkMethod(const fvMatrix<Type>&, const dimensioned<Type>&)"
01241 ) << "incompatible dimensions for operation "
01242 << endl << " "
01243 << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
01244 << op
01245 << " [" << dt.name() << dt.dimensions() << " ]"
01246 << abort(FatalError);
01247 }
01248 }
01249
01250
01251 template<class Type>
01252 Foam::lduMatrix::solverPerformance Foam::solve
01253 (
01254 fvMatrix<Type>& fvm,
01255 const dictionary& solverControls
01256 )
01257 {
01258 return fvm.solve(solverControls);
01259 }
01260
01261 template<class Type>
01262 Foam::lduMatrix::solverPerformance Foam::solve
01263 (
01264 const tmp<fvMatrix<Type> >& tfvm,
01265 const dictionary& solverControls
01266 )
01267 {
01268 lduMatrix::solverPerformance solverPerf =
01269 const_cast<fvMatrix<Type>&>(tfvm()).solve(solverControls);
01270
01271 tfvm.clear();
01272
01273 return solverPerf;
01274 }
01275
01276
01277 template<class Type>
01278 Foam::lduMatrix::solverPerformance Foam::solve(fvMatrix<Type>& fvm)
01279 {
01280 return fvm.solve();
01281 }
01282
01283 template<class Type>
01284 Foam::lduMatrix::solverPerformance Foam::solve(const tmp<fvMatrix<Type> >& tfvm)
01285 {
01286 lduMatrix::solverPerformance solverPerf =
01287 const_cast<fvMatrix<Type>&>(tfvm()).solve();
01288
01289 tfvm.clear();
01290
01291 return solverPerf;
01292 }
01293
01294
01295 template<class Type>
01296 Foam::tmp<Foam::fvMatrix<Type> > Foam::correction
01297 (
01298 const fvMatrix<Type>& A
01299 )
01300 {
01301 tmp<Foam::fvMatrix<Type> > tAcorr = A - (A & A.psi());
01302
01303 if
01304 (
01305 (A.hasUpper() || A.hasLower())
01306 && A.psi().mesh().fluxRequired(A.psi().name())
01307 )
01308 {
01309 tAcorr().faceFluxCorrectionPtr() = (-A.flux()).ptr();
01310 }
01311
01312 return tAcorr;
01313 }
01314
01315
01316 template<class Type>
01317 Foam::tmp<Foam::fvMatrix<Type> > Foam::correction
01318 (
01319 const tmp<fvMatrix<Type> >& tA
01320 )
01321 {
01322 tmp<Foam::fvMatrix<Type> > tAcorr = tA - (tA() & tA().psi());
01323
01324
01325 const fvMatrix<Type>& A = tAcorr();
01326
01327 if
01328 (
01329 (A.hasUpper() || A.hasLower())
01330 && A.psi().mesh().fluxRequired(A.psi().name())
01331 )
01332 {
01333 tAcorr().faceFluxCorrectionPtr() = (-A.flux()).ptr();
01334 }
01335
01336 return tAcorr;
01337 }
01338
01339
01340
01341
01342 template<class Type>
01343 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01344 (
01345 const fvMatrix<Type>& A,
01346 const fvMatrix<Type>& B
01347 )
01348 {
01349 checkMethod(A, B, "==");
01350 return (A - B);
01351 }
01352
01353 template<class Type>
01354 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01355 (
01356 const tmp<fvMatrix<Type> >& tA,
01357 const fvMatrix<Type>& B
01358 )
01359 {
01360 checkMethod(tA(), B, "==");
01361 return (tA - B);
01362 }
01363
01364 template<class Type>
01365 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01366 (
01367 const fvMatrix<Type>& A,
01368 const tmp<fvMatrix<Type> >& tB
01369 )
01370 {
01371 checkMethod(A, tB(), "==");
01372 return (A - tB);
01373 }
01374
01375 template<class Type>
01376 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01377 (
01378 const tmp<fvMatrix<Type> >& tA,
01379 const tmp<fvMatrix<Type> >& tB
01380 )
01381 {
01382 checkMethod(tA(), tB(), "==");
01383 return (tA - tB);
01384 }
01385
01386 template<class Type>
01387 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01388 (
01389 const fvMatrix<Type>& A,
01390 const DimensionedField<Type, volMesh>& su
01391 )
01392 {
01393 checkMethod(A, su, "==");
01394 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01395 tC().source() += su.mesh().V()*su.field();
01396 return tC;
01397 }
01398
01399 template<class Type>
01400 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01401 (
01402 const fvMatrix<Type>& A,
01403 const tmp<DimensionedField<Type, volMesh> >& tsu
01404 )
01405 {
01406 checkMethod(A, tsu(), "==");
01407 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01408 tC().source() += tsu().mesh().V()*tsu().field();
01409 tsu.clear();
01410 return tC;
01411 }
01412
01413 template<class Type>
01414 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01415 (
01416 const fvMatrix<Type>& A,
01417 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01418 )
01419 {
01420 checkMethod(A, tsu(), "==");
01421 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01422 tC().source() += tsu().mesh().V()*tsu().internalField();
01423 tsu.clear();
01424 return tC;
01425 }
01426
01427 template<class Type>
01428 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01429 (
01430 const tmp<fvMatrix<Type> >& tA,
01431 const DimensionedField<Type, volMesh>& su
01432 )
01433 {
01434 checkMethod(tA(), su, "==");
01435 tmp<fvMatrix<Type> > tC(tA.ptr());
01436 tC().source() += su.mesh().V()*su.field();
01437 return tC;
01438 }
01439
01440 template<class Type>
01441 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01442 (
01443 const tmp<fvMatrix<Type> >& tA,
01444 const tmp<DimensionedField<Type, volMesh> >& tsu
01445 )
01446 {
01447 checkMethod(tA(), tsu(), "==");
01448 tmp<fvMatrix<Type> > tC(tA.ptr());
01449 tC().source() += tsu().mesh().V()*tsu().field();
01450 tsu.clear();
01451 return tC;
01452 }
01453
01454 template<class Type>
01455 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01456 (
01457 const tmp<fvMatrix<Type> >& tA,
01458 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01459 )
01460 {
01461 checkMethod(tA(), tsu(), "==");
01462 tmp<fvMatrix<Type> > tC(tA.ptr());
01463 tC().source() += tsu().mesh().V()*tsu().internalField();
01464 tsu.clear();
01465 return tC;
01466 }
01467
01468 template<class Type>
01469 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01470 (
01471 const fvMatrix<Type>& A,
01472 const dimensioned<Type>& su
01473 )
01474 {
01475 checkMethod(A, su, "==");
01476 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01477 tC().source() += A.psi().mesh().V()*su.value();
01478 return tC;
01479 }
01480
01481 template<class Type>
01482 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01483 (
01484 const tmp<fvMatrix<Type> >& tA,
01485 const dimensioned<Type>& su
01486 )
01487 {
01488 checkMethod(tA(), su, "==");
01489 tmp<fvMatrix<Type> > tC(tA.ptr());
01490 tC().source() += tC().psi().mesh().V()*su.value();
01491 return tC;
01492 }
01493
01494 template<class Type>
01495 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01496 (
01497 const fvMatrix<Type>& A,
01498 const zeroField&
01499 )
01500 {
01501 return A;
01502 }
01503
01504
01505 template<class Type>
01506 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
01507 (
01508 const tmp<fvMatrix<Type> >& tA,
01509 const zeroField&
01510 )
01511 {
01512 return tA;
01513 }
01514
01515
01516 template<class Type>
01517 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01518 (
01519 const fvMatrix<Type>& A
01520 )
01521 {
01522 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01523 tC().negate();
01524 return tC;
01525 }
01526
01527 template<class Type>
01528 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01529 (
01530 const tmp<fvMatrix<Type> >& tA
01531 )
01532 {
01533 tmp<fvMatrix<Type> > tC(tA.ptr());
01534 tC().negate();
01535 return tC;
01536 }
01537
01538
01539 template<class Type>
01540 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01541 (
01542 const fvMatrix<Type>& A,
01543 const fvMatrix<Type>& B
01544 )
01545 {
01546 checkMethod(A, B, "+");
01547 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01548 tC() += B;
01549 return tC;
01550 }
01551
01552 template<class Type>
01553 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01554 (
01555 const tmp<fvMatrix<Type> >& tA,
01556 const fvMatrix<Type>& B
01557 )
01558 {
01559 checkMethod(tA(), B, "+");
01560 tmp<fvMatrix<Type> > tC(tA.ptr());
01561 tC() += B;
01562 return tC;
01563 }
01564
01565 template<class Type>
01566 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01567 (
01568 const fvMatrix<Type>& A,
01569 const tmp<fvMatrix<Type> >& tB
01570 )
01571 {
01572 checkMethod(A, tB(), "+");
01573 tmp<fvMatrix<Type> > tC(tB.ptr());
01574 tC() += A;
01575 return tC;
01576 }
01577
01578 template<class Type>
01579 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01580 (
01581 const tmp<fvMatrix<Type> >& tA,
01582 const tmp<fvMatrix<Type> >& tB
01583 )
01584 {
01585 checkMethod(tA(), tB(), "+");
01586 tmp<fvMatrix<Type> > tC(tA.ptr());
01587 tC() += tB();
01588 tB.clear();
01589 return tC;
01590 }
01591
01592 template<class Type>
01593 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01594 (
01595 const fvMatrix<Type>& A,
01596 const DimensionedField<Type, volMesh>& su
01597 )
01598 {
01599 checkMethod(A, su, "+");
01600 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01601 tC().source() -= su.mesh().V()*su.field();
01602 return tC;
01603 }
01604
01605 template<class Type>
01606 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01607 (
01608 const fvMatrix<Type>& A,
01609 const tmp<DimensionedField<Type, volMesh> >& tsu
01610 )
01611 {
01612 checkMethod(A, tsu(), "+");
01613 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01614 tC().source() -= tsu().mesh().V()*tsu().field();
01615 tsu.clear();
01616 return tC;
01617 }
01618
01619 template<class Type>
01620 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01621 (
01622 const fvMatrix<Type>& A,
01623 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01624 )
01625 {
01626 checkMethod(A, tsu(), "+");
01627 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01628 tC().source() -= tsu().mesh().V()*tsu().internalField();
01629 tsu.clear();
01630 return tC;
01631 }
01632
01633 template<class Type>
01634 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01635 (
01636 const tmp<fvMatrix<Type> >& tA,
01637 const DimensionedField<Type, volMesh>& su
01638 )
01639 {
01640 checkMethod(tA(), su, "+");
01641 tmp<fvMatrix<Type> > tC(tA.ptr());
01642 tC().source() -= su.mesh().V()*su.field();
01643 return tC;
01644 }
01645
01646 template<class Type>
01647 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01648 (
01649 const tmp<fvMatrix<Type> >& tA,
01650 const tmp<DimensionedField<Type, volMesh> >& tsu
01651 )
01652 {
01653 checkMethod(tA(), tsu(), "+");
01654 tmp<fvMatrix<Type> > tC(tA.ptr());
01655 tC().source() -= tsu().mesh().V()*tsu().field();
01656 tsu.clear();
01657 return tC;
01658 }
01659
01660 template<class Type>
01661 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01662 (
01663 const tmp<fvMatrix<Type> >& tA,
01664 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01665 )
01666 {
01667 checkMethod(tA(), tsu(), "+");
01668 tmp<fvMatrix<Type> > tC(tA.ptr());
01669 tC().source() -= tsu().mesh().V()*tsu().internalField();
01670 tsu.clear();
01671 return tC;
01672 }
01673
01674 template<class Type>
01675 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01676 (
01677 const DimensionedField<Type, volMesh>& su,
01678 const fvMatrix<Type>& A
01679 )
01680 {
01681 checkMethod(A, su, "+");
01682 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01683 tC().source() -= su.mesh().V()*su.field();
01684 return tC;
01685 }
01686
01687 template<class Type>
01688 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01689 (
01690 const tmp<DimensionedField<Type, volMesh> >& tsu,
01691 const fvMatrix<Type>& A
01692 )
01693 {
01694 checkMethod(A, tsu(), "+");
01695 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01696 tC().source() -= tsu().mesh().V()*tsu().field();
01697 tsu.clear();
01698 return tC;
01699 }
01700
01701 template<class Type>
01702 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01703 (
01704 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
01705 const fvMatrix<Type>& A
01706 )
01707 {
01708 checkMethod(A, tsu(), "+");
01709 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01710 tC().source() -= tsu().mesh().V()*tsu().internalField();
01711 tsu.clear();
01712 return tC;
01713 }
01714
01715 template<class Type>
01716 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01717 (
01718 const DimensionedField<Type, volMesh>& su,
01719 const tmp<fvMatrix<Type> >& tA
01720 )
01721 {
01722 checkMethod(tA(), su, "+");
01723 tmp<fvMatrix<Type> > tC(tA.ptr());
01724 tC().source() -= su.mesh().V()*su.field();
01725 return tC;
01726 }
01727
01728 template<class Type>
01729 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01730 (
01731 const tmp<DimensionedField<Type, volMesh> >& tsu,
01732 const tmp<fvMatrix<Type> >& tA
01733 )
01734 {
01735 checkMethod(tA(), tsu(), "+");
01736 tmp<fvMatrix<Type> > tC(tA.ptr());
01737 tC().source() -= tsu().mesh().V()*tsu().field();
01738 tsu.clear();
01739 return tC;
01740 }
01741
01742 template<class Type>
01743 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01744 (
01745 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
01746 const tmp<fvMatrix<Type> >& tA
01747 )
01748 {
01749 checkMethod(tA(), tsu(), "+");
01750 tmp<fvMatrix<Type> > tC(tA.ptr());
01751 tC().source() -= tsu().mesh().V()*tsu().internalField();
01752 tsu.clear();
01753 return tC;
01754 }
01755
01756
01757 template<class Type>
01758 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01759 (
01760 const fvMatrix<Type>& A,
01761 const fvMatrix<Type>& B
01762 )
01763 {
01764 checkMethod(A, B, "-");
01765 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01766 tC() -= B;
01767 return tC;
01768 }
01769
01770 template<class Type>
01771 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01772 (
01773 const tmp<fvMatrix<Type> >& tA,
01774 const fvMatrix<Type>& B
01775 )
01776 {
01777 checkMethod(tA(), B, "-");
01778 tmp<fvMatrix<Type> > tC(tA.ptr());
01779 tC() -= B;
01780 return tC;
01781 }
01782
01783 template<class Type>
01784 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01785 (
01786 const fvMatrix<Type>& A,
01787 const tmp<fvMatrix<Type> >& tB
01788 )
01789 {
01790 checkMethod(A, tB(), "-");
01791 tmp<fvMatrix<Type> > tC(tB.ptr());
01792 tC() -= A;
01793 tC().negate();
01794 return tC;
01795 }
01796
01797 template<class Type>
01798 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01799 (
01800 const tmp<fvMatrix<Type> >& tA,
01801 const tmp<fvMatrix<Type> >& tB
01802 )
01803 {
01804 checkMethod(tA(), tB(), "-");
01805 tmp<fvMatrix<Type> > tC(tA.ptr());
01806 tC() -= tB();
01807 tB.clear();
01808 return tC;
01809 }
01810
01811 template<class Type>
01812 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01813 (
01814 const fvMatrix<Type>& A,
01815 const DimensionedField<Type, volMesh>& su
01816 )
01817 {
01818 checkMethod(A, su, "-");
01819 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01820 tC().source() += su.mesh().V()*su.field();
01821 return tC;
01822 }
01823
01824 template<class Type>
01825 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01826 (
01827 const fvMatrix<Type>& A,
01828 const tmp<DimensionedField<Type, volMesh> >& tsu
01829 )
01830 {
01831 checkMethod(A, tsu(), "-");
01832 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01833 tC().source() += tsu().mesh().V()*tsu().field();
01834 tsu.clear();
01835 return tC;
01836 }
01837
01838 template<class Type>
01839 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01840 (
01841 const fvMatrix<Type>& A,
01842 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01843 )
01844 {
01845 checkMethod(A, tsu(), "-");
01846 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01847 tC().source() += tsu().mesh().V()*tsu().internalField();
01848 tsu.clear();
01849 return tC;
01850 }
01851
01852 template<class Type>
01853 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01854 (
01855 const tmp<fvMatrix<Type> >& tA,
01856 const DimensionedField<Type, volMesh>& su
01857 )
01858 {
01859 checkMethod(tA(), su, "-");
01860 tmp<fvMatrix<Type> > tC(tA.ptr());
01861 tC().source() += su.mesh().V()*su.field();
01862 return tC;
01863 }
01864
01865 template<class Type>
01866 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01867 (
01868 const tmp<fvMatrix<Type> >& tA,
01869 const tmp<DimensionedField<Type, volMesh> >& tsu
01870 )
01871 {
01872 checkMethod(tA(), tsu(), "-");
01873 tmp<fvMatrix<Type> > tC(tA.ptr());
01874 tC().source() += tsu().mesh().V()*tsu().field();
01875 tsu.clear();
01876 return tC;
01877 }
01878
01879 template<class Type>
01880 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01881 (
01882 const tmp<fvMatrix<Type> >& tA,
01883 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
01884 )
01885 {
01886 checkMethod(tA(), tsu(), "-");
01887 tmp<fvMatrix<Type> > tC(tA.ptr());
01888 tC().source() += tsu().mesh().V()*tsu().internalField();
01889 tsu.clear();
01890 return tC;
01891 }
01892
01893 template<class Type>
01894 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01895 (
01896 const DimensionedField<Type, volMesh>& su,
01897 const fvMatrix<Type>& A
01898 )
01899 {
01900 checkMethod(A, su, "-");
01901 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01902 tC().negate();
01903 tC().source() -= su.mesh().V()*su.field();
01904 return tC;
01905 }
01906
01907 template<class Type>
01908 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01909 (
01910 const tmp<DimensionedField<Type, volMesh> >& tsu,
01911 const fvMatrix<Type>& A
01912 )
01913 {
01914 checkMethod(A, tsu(), "-");
01915 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01916 tC().negate();
01917 tC().source() -= tsu().mesh().V()*tsu().field();
01918 tsu.clear();
01919 return tC;
01920 }
01921
01922 template<class Type>
01923 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01924 (
01925 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
01926 const fvMatrix<Type>& A
01927 )
01928 {
01929 checkMethod(A, tsu(), "-");
01930 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01931 tC().negate();
01932 tC().source() -= tsu().mesh().V()*tsu().internalField();
01933 tsu.clear();
01934 return tC;
01935 }
01936
01937 template<class Type>
01938 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01939 (
01940 const DimensionedField<Type, volMesh>& su,
01941 const tmp<fvMatrix<Type> >& tA
01942 )
01943 {
01944 checkMethod(tA(), su, "-");
01945 tmp<fvMatrix<Type> > tC(tA.ptr());
01946 tC().negate();
01947 tC().source() -= su.mesh().V()*su.field();
01948 return tC;
01949 }
01950
01951 template<class Type>
01952 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01953 (
01954 const tmp<DimensionedField<Type, volMesh> >& tsu,
01955 const tmp<fvMatrix<Type> >& tA
01956 )
01957 {
01958 checkMethod(tA(), tsu(), "-");
01959 tmp<fvMatrix<Type> > tC(tA.ptr());
01960 tC().negate();
01961 tC().source() -= tsu().mesh().V()*tsu().field();
01962 tsu.clear();
01963 return tC;
01964 }
01965
01966 template<class Type>
01967 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
01968 (
01969 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
01970 const tmp<fvMatrix<Type> >& tA
01971 )
01972 {
01973 checkMethod(tA(), tsu(), "-");
01974 tmp<fvMatrix<Type> > tC(tA.ptr());
01975 tC().negate();
01976 tC().source() -= tsu().mesh().V()*tsu().internalField();
01977 tsu.clear();
01978 return tC;
01979 }
01980
01981 template<class Type>
01982 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01983 (
01984 const fvMatrix<Type>& A,
01985 const dimensioned<Type>& su
01986 )
01987 {
01988 checkMethod(A, su, "+");
01989 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
01990 tC().source() -= su.value()*A.psi().mesh().V();
01991 return tC;
01992 }
01993
01994 template<class Type>
01995 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
01996 (
01997 const tmp<fvMatrix<Type> >& tA,
01998 const dimensioned<Type>& su
01999 )
02000 {
02001 checkMethod(tA(), su, "+");
02002 tmp<fvMatrix<Type> > tC(tA.ptr());
02003 tC().source() -= su.value()*tC().psi().mesh().V();
02004 return tC;
02005 }
02006
02007 template<class Type>
02008 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
02009 (
02010 const dimensioned<Type>& su,
02011 const fvMatrix<Type>& A
02012 )
02013 {
02014 checkMethod(A, su, "+");
02015 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02016 tC().source() -= su.value()*A.psi().mesh().V();
02017 return tC;
02018 }
02019
02020 template<class Type>
02021 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
02022 (
02023 const dimensioned<Type>& su,
02024 const tmp<fvMatrix<Type> >& tA
02025 )
02026 {
02027 checkMethod(tA(), su, "+");
02028 tmp<fvMatrix<Type> > tC(tA.ptr());
02029 tC().source() -= su.value()*tC().psi().mesh().V();
02030 return tC;
02031 }
02032
02033 template<class Type>
02034 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
02035 (
02036 const fvMatrix<Type>& A,
02037 const dimensioned<Type>& su
02038 )
02039 {
02040 checkMethod(A, su, "-");
02041 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02042 tC().source() += su.value()*tC().psi().mesh().V();
02043 return tC;
02044 }
02045
02046 template<class Type>
02047 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
02048 (
02049 const tmp<fvMatrix<Type> >& tA,
02050 const dimensioned<Type>& su
02051 )
02052 {
02053 checkMethod(tA(), su, "-");
02054 tmp<fvMatrix<Type> > tC(tA.ptr());
02055 tC().source() += su.value()*tC().psi().mesh().V();
02056 return tC;
02057 }
02058
02059 template<class Type>
02060 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
02061 (
02062 const dimensioned<Type>& su,
02063 const fvMatrix<Type>& A
02064 )
02065 {
02066 checkMethod(A, su, "-");
02067 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02068 tC().negate();
02069 tC().source() -= su.value()*A.psi().mesh().V();
02070 return tC;
02071 }
02072
02073 template<class Type>
02074 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
02075 (
02076 const dimensioned<Type>& su,
02077 const tmp<fvMatrix<Type> >& tA
02078 )
02079 {
02080 checkMethod(tA(), su, "-");
02081 tmp<fvMatrix<Type> > tC(tA.ptr());
02082 tC().negate();
02083 tC().source() -= su.value()*tC().psi().mesh().V();
02084 return tC;
02085 }
02086
02087
02088 template<class Type>
02089 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02090 (
02091 const DimensionedField<scalar, volMesh>& dsf,
02092 const fvMatrix<Type>& A
02093 )
02094 {
02095 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02096 tC() *= dsf;
02097 return tC;
02098 }
02099
02100 template<class Type>
02101 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02102 (
02103 const tmp< DimensionedField<scalar, volMesh> >& tdsf,
02104 const fvMatrix<Type>& A
02105 )
02106 {
02107 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02108 tC() *= tdsf;
02109 return tC;
02110 }
02111
02112 template<class Type>
02113 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02114 (
02115 const tmp<volScalarField>& tvsf,
02116 const fvMatrix<Type>& A
02117 )
02118 {
02119 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02120 tC() *= tvsf;
02121 return tC;
02122 }
02123
02124 template<class Type>
02125 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02126 (
02127 const DimensionedField<scalar, volMesh>& dsf,
02128 const tmp<fvMatrix<Type> >& tA
02129 )
02130 {
02131 tmp<fvMatrix<Type> > tC(tA.ptr());
02132 tC() *= dsf;
02133 return tC;
02134 }
02135
02136 template<class Type>
02137 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02138 (
02139 const tmp<DimensionedField<scalar, volMesh> >& tdsf,
02140 const tmp<fvMatrix<Type> >& tA
02141 )
02142 {
02143 tmp<fvMatrix<Type> > tC(tA.ptr());
02144 tC() *= tdsf;
02145 return tC;
02146 }
02147
02148 template<class Type>
02149 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02150 (
02151 const tmp<volScalarField>& tvsf,
02152 const tmp<fvMatrix<Type> >& tA
02153 )
02154 {
02155 tmp<fvMatrix<Type> > tC(tA.ptr());
02156 tC() *= tvsf;
02157 return tC;
02158 }
02159
02160 template<class Type>
02161 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02162 (
02163 const dimensioned<scalar>& ds,
02164 const fvMatrix<Type>& A
02165 )
02166 {
02167 tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
02168 tC() *= ds;
02169 return tC;
02170 }
02171
02172 template<class Type>
02173 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
02174 (
02175 const dimensioned<scalar>& ds,
02176 const tmp<fvMatrix<Type> >& tA
02177 )
02178 {
02179 tmp<fvMatrix<Type> > tC(tA.ptr());
02180 tC() *= ds;
02181 return tC;
02182 }
02183
02184
02185 template<class Type>
02186 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
02187 Foam::operator&
02188 (
02189 const fvMatrix<Type>& M,
02190 const DimensionedField<Type, volMesh>& psi
02191 )
02192 {
02193 tmp<GeometricField<Type, fvPatchField, volMesh> > tMphi
02194 (
02195 new GeometricField<Type, fvPatchField, volMesh>
02196 (
02197 IOobject
02198 (
02199 "M&" + psi.name(),
02200 psi.instance(),
02201 psi.mesh(),
02202 IOobject::NO_READ,
02203 IOobject::NO_WRITE
02204 ),
02205 psi.mesh(),
02206 M.dimensions()/dimVol,
02207 zeroGradientFvPatchScalarField::typeName
02208 )
02209 );
02210 GeometricField<Type, fvPatchField, volMesh>& Mphi = tMphi();
02211
02212
02213 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
02214 {
02215 scalarField psiCmpt = psi.field().component(cmpt);
02216 scalarField boundaryDiagCmpt(M.diag());
02217 M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
02218 Mphi.internalField().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
02219 }
02220
02221 Mphi.internalField() += M.lduMatrix::H(psi.field()) + M.source();
02222 M.addBoundarySource(Mphi.internalField());
02223
02224 Mphi.internalField() /= -psi.mesh().V();
02225 Mphi.correctBoundaryConditions();
02226
02227 return tMphi;
02228 }
02229
02230 template<class Type>
02231 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
02232 Foam::operator&
02233 (
02234 const fvMatrix<Type>& M,
02235 const tmp<DimensionedField<Type, volMesh> >& tpsi
02236 )
02237 {
02238 tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = M & tpsi();
02239 tpsi.clear();
02240 return tMpsi;
02241 }
02242
02243 template<class Type>
02244 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
02245 Foam::operator&
02246 (
02247 const fvMatrix<Type>& M,
02248 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tpsi
02249 )
02250 {
02251 tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = M & tpsi();
02252 tpsi.clear();
02253 return tMpsi;
02254 }
02255
02256 template<class Type>
02257 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
02258 Foam::operator&
02259 (
02260 const tmp<fvMatrix<Type> >& tM,
02261 const DimensionedField<Type, volMesh>& psi
02262 )
02263 {
02264 tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & psi;
02265 tM.clear();
02266 return tMpsi;
02267 }
02268
02269 template<class Type>
02270 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
02271 Foam::operator&
02272 (
02273 const tmp<fvMatrix<Type> >& tM,
02274 const tmp<DimensionedField<Type, volMesh> >& tpsi
02275 )
02276 {
02277 tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & tpsi();
02278 tM.clear();
02279 tpsi.clear();
02280 return tMpsi;
02281 }
02282
02283 template<class Type>
02284 Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
02285 Foam::operator&
02286 (
02287 const tmp<fvMatrix<Type> >& tM,
02288 const tmp<GeometricField<Type, fvPatchField, volMesh> >& tpsi
02289 )
02290 {
02291 tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & tpsi();
02292 tM.clear();
02293 tpsi.clear();
02294 return tMpsi;
02295 }
02296
02297
02298
02299
02300 template<class Type>
02301 Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
02302 {
02303 os << static_cast<const lduMatrix&>(fvm) << nl
02304 << fvm.dimensions_ << nl
02305 << fvm.source_ << nl
02306 << fvm.internalCoeffs_ << nl
02307 << fvm.boundaryCoeffs_ << endl;
02308
02309 os.check("Ostream& operator<<(Ostream&, fvMatrix<Type>&");
02310
02311 return os;
02312 }
02313
02314
02315
02316
02317 #include <finiteVolume/fvMatrixSolve.C>
02318
02319