Files
openfoam/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.C
2009-02-24 19:20:55 +00:00

2280 lines
52 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "volFields.H"
#include "surfaceFields.H"
#include "calculatedFvPatchFields.H"
#include "zeroGradientFvPatchFields.H"
#include "coupledFvPatchFields.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class Type>
template<class Type2>
void Foam::fvMatrix<Type>::addToInternalField
(
const unallocLabelList& addr,
const Field<Type2>& pf,
Field<Type2>& intf
) const
{
if (addr.size() != pf.size())
{
FatalErrorIn
(
"fvMatrix<Type>::addToInternalField(const unallocLabelList&, "
"const Field&, Field&)"
) << "sizes of addressing and field are different"
<< abort(FatalError);
}
forAll(addr, faceI)
{
intf[addr[faceI]] += pf[faceI];
}
}
template<class Type>
template<class Type2>
void Foam::fvMatrix<Type>::addToInternalField
(
const unallocLabelList& addr,
const tmp<Field<Type2> >& tpf,
Field<Type2>& intf
) const
{
addToInternalField(addr, tpf(), intf);
tpf.clear();
}
template<class Type>
template<class Type2>
void Foam::fvMatrix<Type>::subtractFromInternalField
(
const unallocLabelList& addr,
const Field<Type2>& pf,
Field<Type2>& intf
) const
{
if (addr.size() != pf.size())
{
FatalErrorIn
(
"fvMatrix<Type>::addToInternalField(const unallocLabelList&, "
"const Field&, Field&)"
) << "sizes of addressing and field are different"
<< abort(FatalError);
}
forAll(addr, faceI)
{
intf[addr[faceI]] -= pf[faceI];
}
}
template<class Type>
template<class Type2>
void Foam::fvMatrix<Type>::subtractFromInternalField
(
const unallocLabelList& addr,
const tmp<Field<Type2> >& tpf,
Field<Type2>& intf
) const
{
subtractFromInternalField(addr, tpf(), intf);
tpf.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::addBoundaryDiag
(
scalarField& diag,
const direction solveCmpt
) const
{
forAll(internalCoeffs_, patchI)
{
addToInternalField
(
lduAddr().patchAddr(patchI),
internalCoeffs_[patchI].component(solveCmpt),
diag
);
}
}
template<class Type>
void Foam::fvMatrix<Type>::addCmptAvBoundaryDiag(scalarField& diag) const
{
forAll(internalCoeffs_, patchI)
{
addToInternalField
(
lduAddr().patchAddr(patchI),
cmptAv(internalCoeffs_[patchI]),
diag
);
}
}
template<class Type>
void Foam::fvMatrix<Type>::addBoundarySource
(
Field<Type>& source,
const bool couples
) const
{
forAll(psi_.boundaryField(), patchI)
{
const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
const Field<Type>& pbc = boundaryCoeffs_[patchI];
if (!ptf.coupled())
{
addToInternalField(lduAddr().patchAddr(patchI), pbc, source);
}
else if (couples)
{
tmp<Field<Type> > tpnf = ptf.patchNeighbourField();
const Field<Type>& pnf = tpnf();
const unallocLabelList& addr = lduAddr().patchAddr(patchI);
forAll(addr, facei)
{
source[addr[facei]] += cmptMultiply(pbc[facei], pnf[facei]);
}
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class Type>
Foam::fvMatrix<Type>::fvMatrix
(
GeometricField<Type, fvPatchField, volMesh>& psi,
const dimensionSet& ds
)
:
lduMatrix(psi.mesh()),
psi_(psi),
dimensions_(ds),
source_(psi.size(), pTraits<Type>::zero),
internalCoeffs_(psi.mesh().boundary().size()),
boundaryCoeffs_(psi.mesh().boundary().size()),
faceFluxCorrectionPtr_(NULL)
{
if (debug)
{
Info<< "fvMatrix<Type>(GeometricField<Type, fvPatchField, volMesh>&,"
" const dimensionSet&) : "
"constructing fvMatrix<Type> for field " << psi_.name()
<< endl;
}
// Initialise coupling coefficients
forAll (psi.mesh().boundary(), patchI)
{
internalCoeffs_.set
(
patchI,
new Field<Type>
(
psi.mesh().boundary()[patchI].size(),
pTraits<Type>::zero
)
);
boundaryCoeffs_.set
(
patchI,
new Field<Type>
(
psi.mesh().boundary()[patchI].size(),
pTraits<Type>::zero
)
);
}
psi_.boundaryField().updateCoeffs();
}
template<class Type>
Foam::fvMatrix<Type>::fvMatrix(const fvMatrix<Type>& fvm)
:
refCount(),
lduMatrix(fvm),
psi_(fvm.psi_),
dimensions_(fvm.dimensions_),
source_(fvm.source_),
internalCoeffs_(fvm.internalCoeffs_),
boundaryCoeffs_(fvm.boundaryCoeffs_),
faceFluxCorrectionPtr_(NULL)
{
if (debug)
{
Info<< "fvMatrix<Type>::fvMatrix(const fvMatrix<Type>&) : "
<< "copying fvMatrix<Type> for field " << psi_.name()
<< endl;
}
if (fvm.faceFluxCorrectionPtr_)
{
faceFluxCorrectionPtr_ = new
GeometricField<Type, fvsPatchField, surfaceMesh>
(
*(fvm.faceFluxCorrectionPtr_)
);
}
}
#ifdef ConstructFromTmp
template<class Type>
Foam::fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type> >& tfvm)
:
refCount(),
lduMatrix
(
const_cast<fvMatrix<Type>&>(tfvm()),
tfvm.isTmp()
),
psi_(tfvm().psi_),
dimensions_(tfvm().dimensions_),
source_
(
const_cast<fvMatrix<Type>&>(tfvm()).source_,
tfvm.isTmp()
),
internalCoeffs_
(
const_cast<fvMatrix<Type>&>(tfvm()).internalCoeffs_,
tfvm.isTmp()
),
boundaryCoeffs_
(
const_cast<fvMatrix<Type>&>(tfvm()).boundaryCoeffs_,
tfvm.isTmp()
),
faceFluxCorrectionPtr_(NULL)
{
if (debug)
{
Info<< "fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type> >&) : "
<< "copying fvMatrix<Type> for field " << psi_.name()
<< endl;
}
if (tfvm().faceFluxCorrectionPtr_)
{
if (tfvm.isTmp())
{
faceFluxCorrectionPtr_ = tfvm().faceFluxCorrectionPtr_;
tfvm().faceFluxCorrectionPtr_ = NULL;
}
else
{
faceFluxCorrectionPtr_ = new
GeometricField<Type, fvsPatchField, surfaceMesh>
(
*(tfvm().faceFluxCorrectionPtr_)
);
}
}
tfvm.clear();
}
#endif
template<class Type>
Foam::fvMatrix<Type>::fvMatrix
(
GeometricField<Type, fvPatchField, volMesh>& psi,
Istream& is
)
:
lduMatrix(psi.mesh()),
psi_(psi),
dimensions_(is),
source_(is),
internalCoeffs_(psi.mesh().boundary().size()),
boundaryCoeffs_(psi.mesh().boundary().size()),
faceFluxCorrectionPtr_(NULL)
{
if (debug)
{
Info<< "fvMatrix<Type>"
"(GeometricField<Type, fvPatchField, volMesh>&, Istream&) : "
"constructing fvMatrix<Type> for field " << psi_.name()
<< endl;
}
// Initialise coupling coefficients
forAll (psi.mesh().boundary(), patchI)
{
internalCoeffs_.set
(
patchI,
new Field<Type>
(
psi.mesh().boundary()[patchI].size(),
pTraits<Type>::zero
)
);
boundaryCoeffs_.set
(
patchI,
new Field<Type>
(
psi.mesh().boundary()[patchI].size(),
pTraits<Type>::zero
)
);
}
}
template<class Type>
Foam::fvMatrix<Type>::~fvMatrix()
{
if (debug)
{
Info<< "fvMatrix<Type>::~fvMatrix<Type>() : "
<< "destroying fvMatrix<Type> for field " << psi_.name()
<< endl;
}
if (faceFluxCorrectionPtr_)
{
delete faceFluxCorrectionPtr_;
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
// Set solution in given cells and eliminate corresponding
// equations from the matrix
template<class Type>
void Foam::fvMatrix<Type>::setValues
(
const labelList& cellLabels,
const Field<Type>& values
)
{
const fvMesh& mesh = psi_.mesh();
const cellList& cells = mesh.cells();
const unallocLabelList& own = mesh.owner();
const unallocLabelList& nei = mesh.neighbour();
scalarField& Diag = diag();
forAll(cellLabels, i)
{
label celli = cellLabels[i];
psi_[celli] = values[i];
source_[celli] = values[i]*Diag[celli];
if (symmetric() || asymmetric())
{
const cell& c = cells[celli];
forAll(c, j)
{
label facei = c[j];
if (mesh.isInternalFace(facei))
{
if (symmetric())
{
if (celli == own[facei])
{
source_[nei[facei]] -= upper()[facei]*values[i];
}
else
{
source_[own[facei]] -= upper()[facei]*values[i];
}
upper()[facei] = 0.0;
}
else
{
if (celli == own[facei])
{
source_[nei[facei]] -= lower()[facei]*values[i];
}
else
{
source_[own[facei]] -= upper()[facei]*values[i];
}
upper()[facei] = 0.0;
lower()[facei] = 0.0;
}
}
else
{
label patchi = mesh.boundaryMesh().whichPatch(facei);
if (internalCoeffs_[patchi].size())
{
label patchFacei =
mesh.boundaryMesh()[patchi].whichFace(facei);
internalCoeffs_[patchi][patchFacei] =
pTraits<Type>::zero;
boundaryCoeffs_[patchi][patchFacei] =
pTraits<Type>::zero;
}
}
}
}
}
}
// Set reference level for solution
template<class Type>
void Foam::fvMatrix<Type>::setReference
(
const label cell,
const Type& value,
const bool forceReference
)
{
if (psi_.needReference() || forceReference)
{
if (cell >= 0)
{
source()[cell] += diag()[cell]*value;
diag()[cell] += diag()[cell];
}
}
}
template<class Type>
void Foam::fvMatrix<Type>::relax(const scalar alpha)
{
if (alpha <= 0)
{
return;
}
Field<Type>& S = source();
scalarField& D = diag();
// Store the current unrelaxed diagonal for use in updating the source
scalarField D0(D);
// Calculate the sum-mag off-diagonal from the interior faces
scalarField sumOff(D.size(), 0.0);
sumMagOffDiag(sumOff);
// Handle the boundary contributions to the diagonal
forAll(psi_.boundaryField(), patchI)
{
const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
if (ptf.size())
{
const unallocLabelList& pa = lduAddr().patchAddr(patchI);
Field<Type>& iCoeffs = internalCoeffs_[patchI];
if (ptf.coupled())
{
const Field<Type>& pCoeffs = boundaryCoeffs_[patchI];
// For coupled boundaries add the diagonal and
// off-diagonal contributions
forAll(pa, face)
{
D[pa[face]] += component(iCoeffs[face], 0);
sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
}
}
else
{
// For non-coupled boundaries subtract the diagonal
// contribution off-diagonal sum which avoids having to remove
// it from the diagonal later.
// Also add the source contribution from the relaxation
forAll(pa, face)
{
Type iCoeff0 = iCoeffs[face];
iCoeffs[face] = cmptMag(iCoeffs[face]);
sumOff[pa[face]] -= cmptMin(iCoeffs[face]);
iCoeffs[face] /= alpha;
S[pa[face]] +=
cmptMultiply(iCoeffs[face] - iCoeff0, psi_[pa[face]]);
}
}
}
}
// Ensure the matrix is diagonally dominant...
max(D, D, sumOff);
// ... then relax
D /= alpha;
// Now remove the diagonal contribution from coupled boundaries
forAll(psi_.boundaryField(), patchI)
{
const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
if (ptf.size())
{
const unallocLabelList& pa = lduAddr().patchAddr(patchI);
Field<Type>& iCoeffs = internalCoeffs_[patchI];
if (ptf.coupled())
{
forAll(pa, face)
{
D[pa[face]] -= component(iCoeffs[face], 0);
}
}
}
}
// Finally add the relaxation contribution to the source.
S += (D - D0)*psi_.internalField();
}
template<class Type>
void Foam::fvMatrix<Type>::relax()
{
if (psi_.mesh().relax(psi_.name()))
{
relax(psi_.mesh().relaxationFactor(psi_.name()));
}
}
template<class Type>
void Foam::fvMatrix<Type>::boundaryManipulate
(
typename GeometricField<Type, fvPatchField, volMesh>::
GeometricBoundaryField& bFields
)
{
forAll(bFields, patchI)
{
bFields[patchI].manipulateMatrix(*this);
}
}
template<class Type>
Foam::tmp<Foam::scalarField> Foam::fvMatrix<Type>::D() const
{
tmp<scalarField> tdiag(new scalarField(diag()));
addCmptAvBoundaryDiag(tdiag());
return tdiag;
}
template<class Type>
Foam::tmp<Foam::Field<Type> > Foam::fvMatrix<Type>::DD() const
{
tmp<Field<Type> > tdiag(pTraits<Type>::one*diag());
forAll(psi_.boundaryField(), patchI)
{
const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
if (!ptf.coupled() && ptf.size())
{
addToInternalField
(
lduAddr().patchAddr(patchI),
internalCoeffs_[patchI],
tdiag()
);
}
}
return tdiag;
}
template<class Type>
Foam::tmp<Foam::volScalarField> Foam::fvMatrix<Type>::A() const
{
tmp<volScalarField> tAphi
(
new volScalarField
(
IOobject
(
"A("+psi_.name()+')',
psi_.instance(),
psi_.mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
psi_.mesh(),
dimensions_/psi_.dimensions()/dimVol,
zeroGradientFvPatchScalarField::typeName
)
);
tAphi().internalField() = D()/psi_.mesh().V();
tAphi().correctBoundaryConditions();
return tAphi;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::fvMatrix<Type>::H() const
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tHphi
(
new GeometricField<Type, fvPatchField, volMesh>
(
IOobject
(
"H("+psi_.name()+')',
psi_.instance(),
psi_.mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
psi_.mesh(),
dimensions_/dimVol,
zeroGradientFvPatchScalarField::typeName
)
);
GeometricField<Type, fvPatchField, volMesh>& Hphi = tHphi();
// Loop over field components
for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
{
scalarField psiCmpt = psi_.internalField().component(cmpt);
scalarField boundaryDiagCmpt(psi_.size(), 0.0);
addBoundaryDiag(boundaryDiagCmpt, cmpt);
boundaryDiagCmpt.negate();
addCmptAvBoundaryDiag(boundaryDiagCmpt);
Hphi.internalField().replace(cmpt, boundaryDiagCmpt*psiCmpt);
}
Hphi.internalField() += lduMatrix::H(psi_.internalField()) + source_;
addBoundarySource(Hphi.internalField());
Hphi.internalField() /= psi_.mesh().V();
Hphi.correctBoundaryConditions();
typename Type::labelType validComponents
(
pow
(
psi_.mesh().solutionD(),
pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero
)
);
for(direction cmpt=0; cmpt<Type::nComponents; cmpt++)
{
if (validComponents[cmpt] == -1)
{
Hphi.replace
(
cmpt,
dimensionedScalar("0", Hphi.dimensions(), 0.0)
);
}
}
return tHphi;
}
template<class Type>
Foam::tmp<Foam::volScalarField> Foam::fvMatrix<Type>::H1() const
{
tmp<volScalarField> tH1
(
new volScalarField
(
IOobject
(
"H(1)",
psi_.instance(),
psi_.mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
psi_.mesh(),
dimensions_/(dimVol*psi_.dimensions()),
zeroGradientFvPatchScalarField::typeName
)
);
volScalarField& H1_ = tH1();
// Loop over field components
/*
for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
{
scalarField psiCmpt = psi_.internalField().component(cmpt);
scalarField boundaryDiagCmpt(psi_.size(), 0.0);
addBoundaryDiag(boundaryDiagCmpt, cmpt);
boundaryDiagCmpt.negate();
addCmptAvBoundaryDiag(boundaryDiagCmpt);
H1_.internalField().replace(cmpt, boundaryDiagCmpt);
}
H1_.internalField() += lduMatrix::H1();
*/
H1_.internalField() = lduMatrix::H1();
H1_.internalField() /= psi_.mesh().V();
H1_.correctBoundaryConditions();
return tH1;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh> >
Foam::fvMatrix<Type>::
flux() const
{
if (!psi_.mesh().fluxRequired(psi_.name()))
{
FatalErrorIn("fvMatrix<Type>::flux()")
<< "flux requested but " << psi_.name()
<< " not specified in the fluxRequired sub-dictionary"
" of fvSchemes."
<< abort(FatalError);
}
// construct GeometricField<Type, fvsPatchField, surfaceMesh>
tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tfieldFlux
(
new GeometricField<Type, fvsPatchField, surfaceMesh>
(
IOobject
(
"flux("+psi_.name()+')',
psi_.instance(),
psi_.mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
psi_.mesh(),
dimensions()
)
);
GeometricField<Type, fvsPatchField, surfaceMesh>& fieldFlux = tfieldFlux();
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
{
fieldFlux.internalField().replace
(
cmpt,
lduMatrix::faceH(psi_.internalField().component(cmpt))
);
}
FieldField<Field, Type> InternalContrib = internalCoeffs_;
forAll(InternalContrib, patchI)
{
InternalContrib[patchI] =
cmptMultiply
(
InternalContrib[patchI],
psi_.boundaryField()[patchI].patchInternalField()
);
}
FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
forAll(NeighbourContrib, patchI)
{
if (psi_.boundaryField()[patchI].coupled())
{
NeighbourContrib[patchI] =
cmptMultiply
(
NeighbourContrib[patchI],
psi_.boundaryField()[patchI].patchNeighbourField()
);
}
}
forAll(fieldFlux.boundaryField(), patchI)
{
fieldFlux.boundaryField()[patchI] =
InternalContrib[patchI] - NeighbourContrib[patchI];
}
if (faceFluxCorrectionPtr_)
{
fieldFlux += *faceFluxCorrectionPtr_;
}
return tfieldFlux;
}
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
template<class Type>
void Foam::fvMatrix<Type>::operator=(const fvMatrix<Type>& fvmv)
{
if (this == &fvmv)
{
FatalErrorIn("fvMatrix<Type>::operator=(const fvMatrix<Type>&)")
<< "attempted assignment to self"
<< abort(FatalError);
}
if (&psi_ != &(fvmv.psi_))
{
FatalErrorIn("fvMatrix<Type>::operator=(const fvMatrix<Type>&)")
<< "different fields"
<< abort(FatalError);
}
lduMatrix::operator=(fvmv);
source_ = fvmv.source_;
internalCoeffs_ = fvmv.internalCoeffs_;
boundaryCoeffs_ = fvmv.boundaryCoeffs_;
if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
{
*faceFluxCorrectionPtr_ = *fvmv.faceFluxCorrectionPtr_;
}
else if (fvmv.faceFluxCorrectionPtr_)
{
faceFluxCorrectionPtr_ =
new GeometricField<Type, fvsPatchField, surfaceMesh>
(*fvmv.faceFluxCorrectionPtr_);
}
}
template<class Type>
void Foam::fvMatrix<Type>::operator=(const tmp<fvMatrix<Type> >& tfvmv)
{
operator=(tfvmv());
tfvmv.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::negate()
{
lduMatrix::negate();
source_.negate();
internalCoeffs_.negate();
boundaryCoeffs_.negate();
if (faceFluxCorrectionPtr_)
{
faceFluxCorrectionPtr_->negate();
}
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=(const fvMatrix<Type>& fvmv)
{
checkMethod(*this, fvmv, "+=");
dimensions_ += fvmv.dimensions_;
lduMatrix::operator+=(fvmv);
source_ += fvmv.source_;
internalCoeffs_ += fvmv.internalCoeffs_;
boundaryCoeffs_ += fvmv.boundaryCoeffs_;
if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
{
*faceFluxCorrectionPtr_ += *fvmv.faceFluxCorrectionPtr_;
}
else if (fvmv.faceFluxCorrectionPtr_)
{
faceFluxCorrectionPtr_ = new
GeometricField<Type, fvsPatchField, surfaceMesh>
(
*fvmv.faceFluxCorrectionPtr_
);
}
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=(const tmp<fvMatrix<Type> >& tfvmv)
{
operator+=(tfvmv());
tfvmv.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator-=(const fvMatrix<Type>& fvmv)
{
checkMethod(*this, fvmv, "+=");
dimensions_ -= fvmv.dimensions_;
lduMatrix::operator-=(fvmv);
source_ -= fvmv.source_;
internalCoeffs_ -= fvmv.internalCoeffs_;
boundaryCoeffs_ -= fvmv.boundaryCoeffs_;
if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
{
*faceFluxCorrectionPtr_ -= *fvmv.faceFluxCorrectionPtr_;
}
else if (fvmv.faceFluxCorrectionPtr_)
{
faceFluxCorrectionPtr_ =
new GeometricField<Type, fvsPatchField, surfaceMesh>
(-*fvmv.faceFluxCorrectionPtr_);
}
}
template<class Type>
void Foam::fvMatrix<Type>::operator-=(const tmp<fvMatrix<Type> >& tfvmv)
{
operator-=(tfvmv());
tfvmv.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=
(
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(*this, su, "+=");
source() -= su.mesh().V()*su.field();
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=
(
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
operator+=(tsu());
tsu.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=
(
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
operator+=(tsu());
tsu.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator-=
(
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(*this, su, "-=");
source() += su.mesh().V()*su.field();
}
template<class Type>
void Foam::fvMatrix<Type>::operator-=
(
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
operator-=(tsu());
tsu.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator-=
(
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
operator-=(tsu());
tsu.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=
(
const dimensioned<Type>& su
)
{
source() -= psi().mesh().V()*su;
}
template<class Type>
void Foam::fvMatrix<Type>::operator-=
(
const dimensioned<Type>& su
)
{
source() += psi().mesh().V()*su;
}
template<class Type>
void Foam::fvMatrix<Type>::operator+=
(
const zeroField&
)
{}
template<class Type>
void Foam::fvMatrix<Type>::operator-=
(
const zeroField&
)
{}
template<class Type>
void Foam::fvMatrix<Type>::operator*=
(
const DimensionedField<scalar, volMesh>& dsf
)
{
dimensions_ *= dsf.dimensions();
lduMatrix::operator*=(dsf.field());
source_ *= dsf.field();
forAll(boundaryCoeffs_, patchI)
{
scalarField psf
(
dsf.mesh().boundary()[patchI].patchInternalField(dsf.field())
);
internalCoeffs_[patchI] *= psf;
boundaryCoeffs_[patchI] *= psf;
}
if (faceFluxCorrectionPtr_)
{
FatalErrorIn
(
"fvMatrix<Type>::operator*="
"(const DimensionedField<scalar, volMesh>&)"
) << "cannot scale a matrix containing a faceFluxCorrection"
<< abort(FatalError);
}
}
template<class Type>
void Foam::fvMatrix<Type>::operator*=
(
const tmp<DimensionedField<scalar, volMesh> >& tdsf
)
{
operator*=(tdsf());
tdsf.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator*=
(
const tmp<volScalarField>& tvsf
)
{
operator*=(tvsf());
tvsf.clear();
}
template<class Type>
void Foam::fvMatrix<Type>::operator*=
(
const dimensioned<scalar>& ds
)
{
dimensions_ *= ds.dimensions();
lduMatrix::operator*=(ds.value());
source_ *= ds.value();
internalCoeffs_ *= ds.value();
boundaryCoeffs_ *= ds.value();
if (faceFluxCorrectionPtr_)
{
*faceFluxCorrectionPtr_ *= ds.value();
}
}
// * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::checkMethod
(
const fvMatrix<Type>& fvm1,
const fvMatrix<Type>& fvm2,
const char* op
)
{
if (&fvm1.psi() != &fvm2.psi())
{
FatalErrorIn
(
"checkMethod(const fvMatrix<Type>&, const fvMatrix<Type>&)"
) << "incompatible fields for operation "
<< endl << " "
<< "[" << fvm1.psi().name() << "] "
<< op
<< " [" << fvm2.psi().name() << "]"
<< abort(FatalError);
}
if (dimensionSet::debug && fvm1.dimensions() != fvm2.dimensions())
{
FatalErrorIn
(
"checkMethod(const fvMatrix<Type>&, const fvMatrix<Type>&)"
) << "incompatible dimensions for operation "
<< endl << " "
<< "[" << fvm1.psi().name() << fvm1.dimensions()/dimVolume << " ] "
<< op
<< " [" << fvm2.psi().name() << fvm2.dimensions()/dimVolume << " ]"
<< abort(FatalError);
}
}
template<class Type>
void Foam::checkMethod
(
const fvMatrix<Type>& fvm,
const DimensionedField<Type, volMesh>& df,
const char* op
)
{
if (dimensionSet::debug && fvm.dimensions()/dimVolume != df.dimensions())
{
FatalErrorIn
(
"checkMethod(const fvMatrix<Type>&, const GeometricField<Type, "
"fvPatchField, volMesh>&)"
) << "incompatible dimensions for operation "
<< endl << " "
<< "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
<< op
<< " [" << df.name() << df.dimensions() << " ]"
<< abort(FatalError);
}
}
template<class Type>
void Foam::checkMethod
(
const fvMatrix<Type>& fvm,
const dimensioned<Type>& dt,
const char* op
)
{
if (dimensionSet::debug && fvm.dimensions()/dimVolume != dt.dimensions())
{
FatalErrorIn
(
"checkMethod(const fvMatrix<Type>&, const dimensioned<Type>&)"
) << "incompatible dimensions for operation "
<< endl << " "
<< "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
<< op
<< " [" << dt.name() << dt.dimensions() << " ]"
<< abort(FatalError);
}
}
template<class Type>
Foam::lduMatrix::solverPerformance Foam::solve
(
fvMatrix<Type>& fvm,
const dictionary& solverControls
)
{
return fvm.solve(solverControls);
}
template<class Type>
Foam::lduMatrix::solverPerformance Foam::solve
(
const tmp<fvMatrix<Type> >& tfvm,
const dictionary& solverControls
)
{
lduMatrix::solverPerformance solverPerf =
const_cast<fvMatrix<Type>&>(tfvm()).solve(solverControls);
tfvm.clear();
return solverPerf;
}
template<class Type>
Foam::lduMatrix::solverPerformance Foam::solve(fvMatrix<Type>& fvm)
{
return fvm.solve();
}
template<class Type>
Foam::lduMatrix::solverPerformance Foam::solve(const tmp<fvMatrix<Type> >& tfvm)
{
lduMatrix::solverPerformance solverPerf =
const_cast<fvMatrix<Type>&>(tfvm()).solve();
tfvm.clear();
return solverPerf;
}
// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const fvMatrix<Type>& B
)
{
checkMethod(A, B, "==");
return (A - B);
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const fvMatrix<Type>& B
)
{
checkMethod(tA(), B, "==");
return (tA - B);
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const tmp<fvMatrix<Type> >& tB
)
{
checkMethod(A, tB(), "==");
return (A - tB);
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const tmp<fvMatrix<Type> >& tB
)
{
checkMethod(tA(), tB(), "==");
return (tA - tB);
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(A, su, "==");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
checkMethod(A, tsu(), "==");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
checkMethod(A, tsu(), "==");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(tA(), su, "==");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
checkMethod(tA(), tsu(), "==");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
checkMethod(tA(), tsu(), "==");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const dimensioned<Type>& su
)
{
checkMethod(A, su, "==");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += A.psi().mesh().V()*su.value();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const dimensioned<Type>& su
)
{
checkMethod(tA(), su, "==");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += tC().psi().mesh().V()*su.value();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const fvMatrix<Type>& A,
const zeroField&
)
{
return A;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
(
const tmp<fvMatrix<Type> >& tA,
const zeroField&
)
{
return tA;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A
)
{
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().negate();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA
)
{
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().negate();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const fvMatrix<Type>& A,
const fvMatrix<Type>& B
)
{
checkMethod(A, B, "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC() += B;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<fvMatrix<Type> >& tA,
const fvMatrix<Type>& B
)
{
checkMethod(tA(), B, "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() += B;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const fvMatrix<Type>& A,
const tmp<fvMatrix<Type> >& tB
)
{
checkMethod(A, tB(), "+");
tmp<fvMatrix<Type> > tC(tB.ptr());
tC() += A;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<fvMatrix<Type> >& tA,
const tmp<fvMatrix<Type> >& tB
)
{
checkMethod(tA(), tB(), "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() += tB();
tB.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const fvMatrix<Type>& A,
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(A, su, "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const fvMatrix<Type>& A,
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
checkMethod(A, tsu(), "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const fvMatrix<Type>& A,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
checkMethod(A, tsu(), "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<fvMatrix<Type> >& tA,
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(tA(), su, "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<fvMatrix<Type> >& tA,
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
checkMethod(tA(), tsu(), "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<fvMatrix<Type> >& tA,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
checkMethod(tA(), tsu(), "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const DimensionedField<Type, volMesh>& su,
const fvMatrix<Type>& A
)
{
checkMethod(A, su, "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<DimensionedField<Type, volMesh> >& tsu,
const fvMatrix<Type>& A
)
{
checkMethod(A, tsu(), "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
const fvMatrix<Type>& A
)
{
checkMethod(A, tsu(), "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const DimensionedField<Type, volMesh>& su,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), su, "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<DimensionedField<Type, volMesh> >& tsu,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), tsu(), "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), tsu(), "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A,
const fvMatrix<Type>& B
)
{
checkMethod(A, B, "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC() -= B;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA,
const fvMatrix<Type>& B
)
{
checkMethod(tA(), B, "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() -= B;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A,
const tmp<fvMatrix<Type> >& tB
)
{
checkMethod(A, tB(), "-");
tmp<fvMatrix<Type> > tC(tB.ptr());
tC() -= A;
tC().negate();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA,
const tmp<fvMatrix<Type> >& tB
)
{
checkMethod(tA(), tB(), "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() -= tB();
tB.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A,
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(A, su, "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A,
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
checkMethod(A, tsu(), "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
checkMethod(A, tsu(), "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA,
const DimensionedField<Type, volMesh>& su
)
{
checkMethod(tA(), su, "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA,
const tmp<DimensionedField<Type, volMesh> >& tsu
)
{
checkMethod(tA(), tsu(), "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
)
{
checkMethod(tA(), tsu(), "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const DimensionedField<Type, volMesh>& su,
const fvMatrix<Type>& A
)
{
checkMethod(A, su, "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().negate();
tC().source() -= su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<DimensionedField<Type, volMesh> >& tsu,
const fvMatrix<Type>& A
)
{
checkMethod(A, tsu(), "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().negate();
tC().source() -= tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
const fvMatrix<Type>& A
)
{
checkMethod(A, tsu(), "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().negate();
tC().source() -= tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const DimensionedField<Type, volMesh>& su,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), su, "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().negate();
tC().source() -= su.mesh().V()*su.field();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<DimensionedField<Type, volMesh> >& tsu,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), tsu(), "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().negate();
tC().source() -= tsu().mesh().V()*tsu().field();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), tsu(), "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().negate();
tC().source() -= tsu().mesh().V()*tsu().internalField();
tsu.clear();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const fvMatrix<Type>& A,
const dimensioned<Type>& su
)
{
checkMethod(A, su, "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= su.value()*A.psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const tmp<fvMatrix<Type> >& tA,
const dimensioned<Type>& su
)
{
checkMethod(tA(), su, "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= su.value()*tC().psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const dimensioned<Type>& su,
const fvMatrix<Type>& A
)
{
checkMethod(A, su, "+");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() -= su.value()*A.psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
(
const dimensioned<Type>& su,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), su, "+");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() -= su.value()*tC().psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const fvMatrix<Type>& A,
const dimensioned<Type>& su
)
{
checkMethod(A, su, "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().source() += su.value()*tC().psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const tmp<fvMatrix<Type> >& tA,
const dimensioned<Type>& su
)
{
checkMethod(tA(), su, "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().source() += su.value()*tC().psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const dimensioned<Type>& su,
const fvMatrix<Type>& A
)
{
checkMethod(A, su, "-");
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC().negate();
tC().source() -= su.value()*A.psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
(
const dimensioned<Type>& su,
const tmp<fvMatrix<Type> >& tA
)
{
checkMethod(tA(), su, "-");
tmp<fvMatrix<Type> > tC(tA.ptr());
tC().negate();
tC().source() -= su.value()*tC().psi().mesh().V();
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const DimensionedField<scalar, volMesh>& dsf,
const fvMatrix<Type>& A
)
{
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC() *= dsf;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const tmp< DimensionedField<scalar, volMesh> >& tdsf,
const fvMatrix<Type>& A
)
{
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC() *= tdsf;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const tmp<volScalarField>& tvsf,
const fvMatrix<Type>& A
)
{
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC() *= tvsf;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const DimensionedField<scalar, volMesh>& dsf,
const tmp<fvMatrix<Type> >& tA
)
{
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() *= dsf;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const tmp<DimensionedField<scalar, volMesh> >& tdsf,
const tmp<fvMatrix<Type> >& tA
)
{
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() *= tdsf;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const tmp<volScalarField>& tvsf,
const tmp<fvMatrix<Type> >& tA
)
{
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() *= tvsf;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const dimensioned<scalar>& ds,
const fvMatrix<Type>& A
)
{
tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
tC() *= ds;
return tC;
}
template<class Type>
Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
(
const dimensioned<scalar>& ds,
const tmp<fvMatrix<Type> >& tA
)
{
tmp<fvMatrix<Type> > tC(tA.ptr());
tC() *= ds;
return tC;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::operator&
(
const fvMatrix<Type>& M,
const DimensionedField<Type, volMesh>& psi
)
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tMphi
(
new GeometricField<Type, fvPatchField, volMesh>
(
IOobject
(
"M&" + psi.name(),
psi.instance(),
psi.mesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
psi.mesh(),
M.dimensions()/dimVol,
zeroGradientFvPatchScalarField::typeName
)
);
GeometricField<Type, fvPatchField, volMesh>& Mphi = tMphi();
// Loop over field components
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
{
scalarField psiCmpt = psi.field().component(cmpt);
scalarField boundaryDiagCmpt(M.diag());
M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
Mphi.internalField().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
}
Mphi.internalField() += M.lduMatrix::H(psi.field()) + M.source();
M.addBoundarySource(Mphi.internalField());
Mphi.internalField() /= -psi.mesh().V();
Mphi.correctBoundaryConditions();
return tMphi;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::operator&
(
const fvMatrix<Type>& M,
const tmp<DimensionedField<Type, volMesh> >& tpsi
)
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = M & tpsi();
tpsi.clear();
return tMpsi;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::operator&
(
const fvMatrix<Type>& M,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tpsi
)
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = M & tpsi();
tpsi.clear();
return tMpsi;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::operator&
(
const tmp<fvMatrix<Type> >& tM,
const DimensionedField<Type, volMesh>& psi
)
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & psi;
tM.clear();
return tMpsi;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::operator&
(
const tmp<fvMatrix<Type> >& tM,
const tmp<DimensionedField<Type, volMesh> >& tpsi
)
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & tpsi();
tM.clear();
tpsi.clear();
return tMpsi;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvPatchField, Foam::volMesh> >
Foam::operator&
(
const tmp<fvMatrix<Type> >& tM,
const tmp<GeometricField<Type, fvPatchField, volMesh> >& tpsi
)
{
tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & tpsi();
tM.clear();
tpsi.clear();
return tMpsi;
}
// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
template<class Type>
Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
{
os << static_cast<const lduMatrix&>(fvm) << nl
<< fvm.dimensions_ << nl
<< fvm.source_ << nl
<< fvm.internalCoeffs_ << nl
<< fvm.boundaryCoeffs_ << endl;
os.check("Ostream& operator<<(Ostream&, fvMatrix<Type>&");
return os;
}
// * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
#include "fvMatrixSolve.C"
// ************************************************************************* //