Foam::fvm Namespace Reference

Namespace of functions to calculate implicit derivatives returning a matrix. More...

Functions

template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const dimensionedScalar &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const Foam::one, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const dimensionedScalar &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const volScalarField &alpha, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const Foam::one, const Foam::one, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const Foam::one, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > ddt (const volScalarField &alpha, const Foam::one, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > div (const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type >
tmp< fvMatrix< Type > > div (const tmp< surfaceScalarField > &tflux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type >
tmp< fvMatrix< Type > > div (const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > div (const tmp< surfaceScalarField > &tflux, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > laplacian (const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type >
tmp< fvMatrix< Type > > laplacian (const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > laplacian (const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type >
tmp< fvMatrix< Type > > laplacian (const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
tmp< fvMatrix< Type > > laplacian (const Foam::one, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type >
tmp< fvMatrix< Type > > laplacian (const Foam::one, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const dimensioned< GType > &gamma, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const dimensioned< GType > &gamma, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const GeometricField< GType, fvPatchField, volMesh > &gamma, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const tmp< GeometricField< GType, fvPatchField, volMesh >> &tgamma, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const GeometricField< GType, fvPatchField, volMesh > &gamma, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const tmp< GeometricField< GType, fvPatchField, volMesh >> &tgamma, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const GeometricField< GType, fvsPatchField, surfaceMesh > &gamma, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> &tgamma, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const GeometricField< GType, fvsPatchField, surfaceMesh > &gamma, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type , class GType >
tmp< fvMatrix< Type > > laplacian (const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> &tGamma, const GeometricField< Type, fvPatchField, volMesh > &vf)
 
template<class Type >
zeroField Su (const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
 A no-op source. More...
 
template<class Type >
tmp< fvMatrix< Type > > Su (const dimensioned< Type > &su, const GeometricField< Type, fvPatchField, volMesh > &)
 A uniform source (no-op for small values) More...
 
template<class Type >
tmp< fvMatrix< Type > > Su (const DimensionedField< Type, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
tmp< fvMatrix< Type > > Su (const tmp< DimensionedField< Type, volMesh >> &, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
tmp< fvMatrix< Type > > Su (const tmp< GeometricField< Type, fvPatchField, volMesh >> &, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
zeroField Sp (const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
 A no-op source. More...
 
template<class Type >
tmp< fvMatrix< Type > > Sp (const dimensionedScalar &sp, const GeometricField< Type, fvPatchField, volMesh > &)
 A uniform source (no-op for small values) More...
 
template<class Type >
tmp< fvMatrix< Type > > Sp (const DimensionedField< scalar, volMesh > &sp, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
tmp< fvMatrix< Type > > Sp (const tmp< DimensionedField< scalar, volMesh >> &, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
tmp< fvMatrix< Type > > Sp (const tmp< volScalarField > &, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
zeroField SuSp (const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
 A no-op source. More...
 
template<class Type >
tmp< fvMatrix< Type > > SuSp (const dimensionedScalar &susp, const GeometricField< Type, fvPatchField, volMesh > &)
 A uniform source (no-op for small values) More...
 
template<class Type >
tmp< fvMatrix< Type > > SuSp (const DimensionedField< scalar, volMesh > &susp, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
tmp< fvMatrix< Type > > SuSp (const tmp< DimensionedField< scalar, volMesh >> &, const GeometricField< Type, fvPatchField, volMesh > &)
 
template<class Type >
tmp< fvMatrix< Type > > SuSp (const tmp< volScalarField > &, const GeometricField< Type, fvPatchField, volMesh > &)
 

Detailed Description

Namespace of functions to calculate implicit derivatives returning a matrix.

Temporal derivatives are calculated using Euler-implicit, backward differencing or Crank-Nicolson. Spatial derivatives are calculated using Gauss' Theorem.

Function Documentation

◆ d2dt2() [1/3]

tmp< fvMatrix< Type > > d2dt2 ( const GeometricField< Type, fvPatchField, volMesh > &  vf)

Definition at line 40 of file fvmD2dt2.C.

References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and d2dt2Scheme< Type >::New().

Here is the call graph for this function:

◆ d2dt2() [2/3]

tmp< fvMatrix< Type > > d2dt2 ( const dimensionedScalar rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 55 of file fvmD2dt2.C.

References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), d2dt2Scheme< Type >::New(), and rho.

Here is the call graph for this function:

◆ d2dt2() [3/3]

tmp< fvMatrix< Type > > d2dt2 ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 71 of file fvmD2dt2.C.

References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), d2dt2Scheme< Type >::New(), and rho.

Here is the call graph for this function:

◆ ddt() [1/8]

tmp< fvMatrix< Type > > ddt ( const GeometricField< Type, fvPatchField, volMesh > &  vf)

Definition at line 40 of file fvmDdt.C.

References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and ddtScheme< Type >::New().

Referenced by multiphaseMangrovesSource::addSup(), Implicit< CloudType >::cacheFields(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kL< BasicTurbulenceModel >::correct(), adjointSpalartAllmaras::correct(), kineticTheoryModel::correct(), IATE::correct(), kkLOmega::correct(), LamBremhorstKE::correct(), LienCubicKE::correct(), LienLeschziner::correct(), qZeta::correct(), ShihQuadraticKE::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kOmegaSSTBase< BasicEddyViscosityModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), EBRSM< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), relaxation::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), ddt(), thermo::evolveRegion(), energyTransport::execute(), scalarTransport::execute(), AnisothermalPhaseModel< BasePhaseModel >::heEqn(), MomentumTransferPhaseSystem< BasePhaseSystem >::momentumTransfer(), populationBalanceModel::solve(), twoPhaseSystem::solve(), reactingOneDim::solveContinuity(), kinematicSingleLayer::solveContinuity(), reactingOneDim::solveEnergy(), thermoSingleLayer::solveEnergy(), thermalBaffle::solveEnergy(), kinematicSingleLayer::solveMomentum(), reactingOneDim::solveSpeciesMass(), kinematicSingleLayer::solveThickness(), MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi(), MovingPhaseModel< BasePhaseModel >::UEqn(), and MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::YiEqn().

Here is the call graph for this function:

◆ ddt() [2/8]

tmp< fvMatrix< Type > > ddt ( const Foam::one  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 55 of file fvmDdt.C.

References ddt().

Here is the call graph for this function:

◆ ddt() [3/8]

tmp< fvMatrix< Type > > ddt ( const dimensionedScalar rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 67 of file fvmDdt.C.

References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), ddtScheme< Type >::New(), and rho.

Here is the call graph for this function:

◆ ddt() [4/8]

tmp< fvMatrix< Type > > ddt ( const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 83 of file fvmDdt.C.

References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), ddtScheme< Type >::New(), and rho.

Here is the call graph for this function:

◆ ddt() [5/8]

tmp< fvMatrix< Type > > ddt ( const volScalarField alpha,
const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 99 of file fvmDdt.C.

References Foam::constant::atomic::alpha, DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), dimensioned< Type >::name(), ddtScheme< Type >::New(), and rho.

Here is the call graph for this function:

◆ ddt() [6/8]

tmp< fvMatrix< Type > > ddt ( const Foam::one  ,
const Foam::one  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 122 of file fvmDdt.C.

References ddt().

Here is the call graph for this function:

◆ ddt() [7/8]

tmp< fvMatrix< Type > > ddt ( const Foam::one  ,
const volScalarField rho,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 135 of file fvmDdt.C.

References ddt(), and rho.

Here is the call graph for this function:

◆ ddt() [8/8]

tmp< fvMatrix< Type > > ddt ( const volScalarField alpha,
const Foam::one  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 148 of file fvmDdt.C.

References Foam::constant::atomic::alpha, and ddt().

Here is the call graph for this function:

◆ div() [1/4]

tmp<fvMatrix<Type> > Foam::fvm::div ( const surfaceScalarField flux,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 40 of file fvmDiv.C.

References DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and convectionScheme< Type >::New().

Referenced by ATCstandard::addATC(), ATCUaGradU::addATC(), Implicit< CloudType >::cacheFields(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kL< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), adjointSpalartAllmaras::correct(), kineticTheoryModel::correct(), IATE::correct(), radiativeIntensityRay::correct(), kkLOmega::correct(), LamBremhorstKE::correct(), LienCubicKE::correct(), LienLeschziner::correct(), qZeta::correct(), ShihQuadraticKE::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kOmegaSSTBase< BasicEddyViscosityModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), EBRSM< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), relaxation::correct(), advectionDiffusion::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), div(), age::execute(), energyTransport::execute(), scalarTransport::execute(), AnisothermalPhaseModel< BasePhaseModel >::heEqn(), adjointSimple::mainIter(), simple::mainIter(), MomentumTransferPhaseSystem< BasePhaseSystem >::momentumTransfer(), MomentumTransferPhaseSystem< BasePhaseSystem >::momentumTransferf(), adjointEikonalSolver::solve(), thermoSingleLayer::solveEnergy(), kinematicSingleLayer::solveMomentum(), kinematicSingleLayer::solveThickness(), MovingPhaseModel< BasePhaseModel >::UEqn(), MovingPhaseModel< BasePhaseModel >::UfEqn(), and MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::YiEqn().

Here is the call graph for this function:

◆ div() [2/4]

tmp<fvMatrix<Type> > Foam::fvm::div ( const tmp< surfaceScalarField > &  tflux,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 57 of file fvmDiv.C.

References tmp< T >::clear(), div(), and Foam::name().

Here is the call graph for this function:

◆ div() [3/4]

tmp<fvMatrix<Type> > Foam::fvm::div ( const surfaceScalarField flux,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 72 of file fvmDiv.C.

References div(), and IOobject::name().

Here is the call graph for this function:

◆ div() [4/4]

tmp<fvMatrix<Type> > Foam::fvm::div ( const tmp< surfaceScalarField > &  tflux,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 83 of file fvmDiv.C.

References tmp< T >::clear(), and div().

Here is the call graph for this function:

◆ laplacian() [1/16]

tmp< fvMatrix< Type > > laplacian ( const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 40 of file fvmLaplacian.C.

References TimePaths::constant(), Foam::dimless, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), IOobjectOption::NO_READ, and IOobject::time().

Referenced by isotropic::addRegularisationTerm(), jouleHeatingSource::addSup(), Implicit< CloudType >::cacheFields(), P1::calculate(), hydrostaticPressure::calculateAndWrite(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kL< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), adjointSpalartAllmaras::correct(), kineticTheoryModel::correct(), kkLOmega::correct(), LamBremhorstKE::correct(), LienCubicKE::correct(), LienLeschziner::correct(), qZeta::correct(), ShihQuadraticKE::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kOmegaSSTBase< BasicEddyViscosityModel >::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), EBRSM< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), Poisson::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), adjointkOmegaSST::divDevReff(), adjointLaminar::divDevReff(), adjointSpalartAllmaras::divDevReff(), kineticTheoryModel::divDevRhoReff(), Maxwell< BasicTurbulenceModel >::divDevRhoReff(), thermo::evolveRegion(), age::execute(), electricPotential::execute(), energyTransport::execute(), scalarTransport::execute(), AnisothermalPhaseModel< BasePhaseModel >::heEqn(), laplacian(), adjointSimple::mainIter(), simple::mainIter(), displacementComponentLaplacianFvMotionSolver::solve(), velocityComponentLaplacianFvMotionSolver::solve(), displacementLaplacianFvMotionSolver::solve(), displacementSBRStressFvMotionSolver::solve(), solidBodyDisplacementLaplacianFvMotionSolver::solve(), surfaceAlignedSBRStressFvMotionSolver::solve(), velocityLaplacianFvMotionSolver::solve(), elasticityMotionSolver::solve(), laplacianMotionSolver::solve(), pLaplacianMotionSolver::solve(), adjointEikonalSolver::solve(), adjointMeshMovementSolver::solve(), twoPhaseSystem::solve(), reactingOneDim::solveEnergy(), thermalBaffle::solveEnergy(), shapeDesignVariables::solveMeshMovementEqn(), kinematicSingleLayer::solveThickness(), MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi(), and MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::YiEqn().

Here is the call graph for this function:

◆ laplacian() [2/16]

tmp< fvMatrix< Type > > laplacian ( const GeometricField< Type, fvPatchField, volMesh > &  vf)

◆ laplacian() [3/16]

tmp< fvMatrix< Type > > laplacian ( const Foam::zero  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 94 of file fvmLaplacian.C.

◆ laplacian() [4/16]

tmp< fvMatrix< Type > > laplacian ( const Foam::zero  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 110 of file fvmLaplacian.C.

◆ laplacian() [5/16]

tmp< fvMatrix< Type > > laplacian ( const Foam::one  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 125 of file fvmLaplacian.C.

References laplacian(), and Foam::name().

Here is the call graph for this function:

◆ laplacian() [6/16]

tmp< fvMatrix< Type > > laplacian ( const Foam::one  ,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 138 of file fvmLaplacian.C.

References laplacian().

Here is the call graph for this function:

◆ laplacian() [7/16]

tmp< fvMatrix< Type > > laplacian ( const dimensioned< GType > &  gamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 150 of file fvmLaplacian.C.

References gamma, IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and IOobjectOption::NO_READ.

Here is the call graph for this function:

◆ laplacian() [8/16]

tmp< fvMatrix< Type > > laplacian ( const dimensioned< GType > &  gamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 176 of file fvmLaplacian.C.

References gamma, IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), and IOobjectOption::NO_READ.

Here is the call graph for this function:

◆ laplacian() [9/16]

tmp< fvMatrix< Type > > laplacian ( const GeometricField< GType, fvPatchField, volMesh > &  gamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 203 of file fvmLaplacian.C.

References gamma, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and laplacianScheme< Type, GType >::New().

Here is the call graph for this function:

◆ laplacian() [10/16]

tmp< fvMatrix< Type > > laplacian ( const tmp< GeometricField< GType, fvPatchField, volMesh >> &  tgamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 220 of file fvmLaplacian.C.

References laplacian(), and Foam::name().

Here is the call graph for this function:

◆ laplacian() [11/16]

tmp< fvMatrix< Type > > laplacian ( const GeometricField< GType, fvPatchField, volMesh > &  gamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 235 of file fvmLaplacian.C.

References gamma, laplacian(), and IOobject::name().

Here is the call graph for this function:

◆ laplacian() [12/16]

tmp< fvMatrix< Type > > laplacian ( const tmp< GeometricField< GType, fvPatchField, volMesh >> &  tgamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 252 of file fvmLaplacian.C.

References laplacian().

Here is the call graph for this function:

◆ laplacian() [13/16]

tmp< fvMatrix< Type > > laplacian ( const GeometricField< GType, fvsPatchField, surfaceMesh > &  gamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 268 of file fvmLaplacian.C.

References gamma, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and laplacianScheme< Type, GType >::New().

Here is the call graph for this function:

◆ laplacian() [14/16]

tmp< fvMatrix< Type > > laplacian ( const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> &  tgamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf,
const word name 
)

Definition at line 285 of file fvmLaplacian.C.

References laplacian(), and Foam::name().

Here is the call graph for this function:

◆ laplacian() [15/16]

tmp< fvMatrix< Type > > laplacian ( const GeometricField< GType, fvsPatchField, surfaceMesh > &  gamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 300 of file fvmLaplacian.C.

References gamma, laplacian(), and IOobject::name().

Here is the call graph for this function:

◆ laplacian() [16/16]

tmp< fvMatrix< Type > > laplacian ( const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> &  tGamma,
const GeometricField< Type, fvPatchField, volMesh > &  vf 
)

Definition at line 317 of file fvmLaplacian.C.

References laplacian().

Here is the call graph for this function:

◆ Su() [1/5]

◆ Su() [2/5]

tmp<fvMatrix<Type> > Foam::fvm::Su ( const dimensioned< Type > &  su,
const GeometricField< Type, fvPatchField, volMesh > &   
)

A uniform source (no-op for small values)

◆ Su() [3/5]

tmp<fvMatrix<Type> > Foam::fvm::Su ( const DimensionedField< Type, volMesh > &  su,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ Su() [4/5]

tmp<fvMatrix<Type> > Foam::fvm::Su ( const tmp< DimensionedField< Type, volMesh >> &  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ Su() [5/5]

tmp<fvMatrix<Type> > Foam::fvm::Su ( const tmp< GeometricField< Type, fvPatchField, volMesh >> &  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ Sp() [1/5]

zeroField Foam::fvm::Sp ( const Foam::zero  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)

A no-op source.

Referenced by PhaseLimitStabilization< Type >::addSup(), interRegionHeatTransferModel::addSup(), multiphaseMangrovesTurbulenceModel::addSup(), atmPlantCanopyUSource::addSup(), acousticDampingSource::addSup(), topOSource::addSup(), multiphaseMangrovesSource::addSup(), P1::calculate(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), adjointSpalartAllmaras::correct(), kineticTheoryModel::correct(), IATE::correct(), radiativeIntensityRay::correct(), kkLOmega::correct(), LamBremhorstKE::correct(), LienCubicKE::correct(), LienLeschziner::correct(), qZeta::correct(), ShihQuadraticKE::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kOmegaSSTBase< BasicEddyViscosityModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), EBRSM< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), continuousGasKEpsilon< BasicTurbulenceModel >::epsilonSource(), LaheyKEpsilon< BasicTurbulenceModel >::epsilonSource(), energyTransport::execute(), boundedDdtScheme< Type >::fvmDdt(), boundedConvectionScheme< Type >::fvmDiv(), OneResistanceHeatTransferPhaseSystem< BasePhaseSystem >::heatTransfer(), TwoResistanceHeatTransferPhaseSystem< BasePhaseSystem >::heatTransfer(), MassTransferPhaseSystem< BasePhaseSystem >::heatTransfer(), AnisothermalPhaseModel< BasePhaseModel >::heEqn(), continuousGasKEqn< BasicTurbulenceModel >::kSource(), NicenoKEqn< BasicTurbulenceModel >::kSource(), continuousGasKEpsilon< BasicTurbulenceModel >::kSource(), LaheyKEpsilon< BasicTurbulenceModel >::kSource(), InterfaceCompositionPhaseChangePhaseSystem< BasePhaseSystem >::massTransfer(), PhaseTransferPhaseSystem< BasePhaseSystem >::massTransfer(), PopulationBalancePhaseSystem< BasePhaseSystem >::massTransfer(), MomentumTransferPhaseSystem< BasePhaseSystem >::momentumTransfer(), MomentumTransferPhaseSystem< BasePhaseSystem >::momentumTransferf(), thermoSingleLayer::q(), randomCoalescence::R(), wallBoiling::R(), singleStepCombustion< ReactionThermo, ThermoType >::R(), radiationModel::Sh(), populationBalanceModel::solve(), multiphaseSystem::solveAlphas(), Helmholtz::solveEqn(), MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi(), radiationModel::ST(), laminar::Su(), MovingPhaseModel< BasePhaseModel >::UfEqn(), and MassTransferPhaseSystem< BasePhaseSystem >::volTransfer().

◆ Sp() [2/5]

tmp<fvMatrix<Type> > Foam::fvm::Sp ( const dimensionedScalar sp,
const GeometricField< Type, fvPatchField, volMesh > &   
)

A uniform source (no-op for small values)

◆ Sp() [3/5]

tmp<fvMatrix<Type> > Foam::fvm::Sp ( const DimensionedField< scalar, volMesh > &  sp,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ Sp() [4/5]

tmp<fvMatrix<Type> > Foam::fvm::Sp ( const tmp< DimensionedField< scalar, volMesh >> &  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ Sp() [5/5]

tmp<fvMatrix<Type> > Foam::fvm::Sp ( const tmp< volScalarField > &  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ SuSp() [1/5]

◆ SuSp() [2/5]

tmp<fvMatrix<Type> > Foam::fvm::SuSp ( const dimensionedScalar susp,
const GeometricField< Type, fvPatchField, volMesh > &   
)

A uniform source (no-op for small values)

◆ SuSp() [3/5]

tmp<fvMatrix<Type> > Foam::fvm::SuSp ( const DimensionedField< scalar, volMesh > &  susp,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ SuSp() [4/5]

tmp<fvMatrix<Type> > Foam::fvm::SuSp ( const tmp< DimensionedField< scalar, volMesh >> &  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)

◆ SuSp() [5/5]

tmp<fvMatrix<Type> > Foam::fvm::SuSp ( const tmp< volScalarField > &  ,
const GeometricField< Type, fvPatchField, volMesh > &   
)