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 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 one &, const one &, const GeometricField< Type, fvPatchField, volMesh > &vf) |
template<class Type > | |
tmp< fvMatrix< Type > > | ddt (const one &, const volScalarField &rho, const GeometricField< Type, fvPatchField, volMesh > &vf) |
template<class Type > | |
tmp< fvMatrix< Type > > | ddt (const volScalarField &alpha, const 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 zero &, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name) |
template<class Type > | |
tmp< fvMatrix< Type > > | laplacian (const zero &, const GeometricField< Type, fvPatchField, volMesh > &vf) |
template<class Type > | |
tmp< fvMatrix< Type > > | laplacian (const one &, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name) |
template<class Type > | |
tmp< fvMatrix< Type > > | laplacian (const 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 > | |
tmp< fvMatrix< Type > > | Su (const DimensionedField< Type, volMesh > &, 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 | Su (const zero &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | Sp (const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | Sp (const tmp< volScalarField::Internal > &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | Sp (const tmp< volScalarField > &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | Sp (const dimensionedScalar &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
zeroField | Sp (const zero &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | SuSp (const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | SuSp (const tmp< volScalarField::Internal > &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
tmp< fvMatrix< Type > > | SuSp (const tmp< volScalarField > &, const GeometricField< Type, fvPatchField, volMesh > &) |
template<class Type > | |
zeroField | SuSp (const zero &, const GeometricField< Type, fvPatchField, volMesh > &) |
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.
tmp< fvMatrix< Type > > d2dt2 | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 46 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and d2dt2Scheme< Type >::New().
Referenced by d2dt2().
tmp< fvMatrix< Type > > d2dt2 | ( | const dimensionedScalar & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 61 of file fvmD2dt2.C.
References d2dt2(), DimensionedField< Type, GeoMesh >::mesh(), dimensioned< Type >::name(), IOobject::name(), and d2dt2Scheme< Type >::New().
tmp< fvMatrix< Type > > d2dt2 | ( | const volScalarField & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 77 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and d2dt2Scheme< Type >::New().
tmp< fvMatrix< Type > > ddt | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 46 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and ddtScheme< Type >::New().
Referenced by relaxation::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), LamBremhorstKE::correct(), ShihQuadraticKE::correct(), LienLeschziner::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), v2f< BasicTurbulenceModel >::correct(), kkLOmega::correct(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), ddt(), scalarTransport::execute(), reactingOneDim::solveContinuity(), kinematicSingleLayer::solveContinuity(), thermalBaffle::solveEnergy(), reactingOneDim::solveEnergy(), thermoSingleLayer::solveEnergy(), kinematicSingleLayer::solveMomentum(), reactingOneDim::solveSpeciesMass(), and kinematicSingleLayer::solveThickness().
tmp< fvMatrix< Type > > ddt | ( | const one & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > ddt | ( | const dimensionedScalar & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 73 of file fvmDdt.C.
References ddt(), DimensionedField< Type, GeoMesh >::mesh(), dimensioned< Type >::name(), IOobject::name(), and ddtScheme< Type >::New().
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 89 of file fvmDdt.C.
References ddt(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and ddtScheme< Type >::New().
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 105 of file fvmDdt.C.
References ddt(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and ddtScheme< Type >::New().
tmp< fvMatrix< Type > > ddt | ( | const one & | , |
const one & | , | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > ddt | ( | const one & | , |
const volScalarField & | rho, | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | alpha, |
const one & | , | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > div | ( | const surfaceScalarField & | flux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 46 of file fvmDiv.C.
References Foam::fvc::flux(), DimensionedField< Type, GeoMesh >::mesh(), and convectionScheme< Type >::New().
Referenced by relaxation::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), radiativeIntensityRay::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), LienLeschziner::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), v2f< BasicTurbulenceModel >::correct(), kkLOmega::correct(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), div(), scalarTransport::execute(), thermoSingleLayer::solveEnergy(), kinematicSingleLayer::solveMomentum(), and kinematicSingleLayer::solveThickness().
tmp< fvMatrix< Type > > div | ( | const tmp< surfaceScalarField > & | tflux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 63 of file fvmDiv.C.
References tmp< T >::clear(), and div().
tmp< fvMatrix< Type > > div | ( | const surfaceScalarField & | flux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 78 of file fvmDiv.C.
References div(), and IOobject::name().
tmp< fvMatrix< Type > > div | ( | const tmp< surfaceScalarField > & | tflux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 89 of file fvmDiv.C.
References tmp< T >::clear(), and div().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf, |
const word & | name | ||
) |
Definition at line 46 of file fvmLaplacian.C.
References Foam::dimless, DimensionedField< Type, GeoMesh >::mesh(), IOobject::NO_READ, and IOobject::time().
Referenced by P1::calculate(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), LamBremhorstKE::correct(), ShihQuadraticKE::correct(), Poisson::correct(), LienLeschziner::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), kOmega< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), v2f< BasicTurbulenceModel >::correct(), kkLOmega::correct(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::correct(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), Maxwell< BasicTurbulenceModel >::divDevRhoReff(), scalarTransport::execute(), laplacian(), velocityComponentLaplacianFvMotionSolver::solve(), velocityLaplacianFvMotionSolver::solve(), displacementComponentLaplacianFvMotionSolver::solve(), displacementSBRStressFvMotionSolver::solve(), displacementLaplacianFvMotionSolver::solve(), thermalBaffle::solveEnergy(), reactingOneDim::solveEnergy(), and kinematicSingleLayer::solveThickness().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 71 of file fvmLaplacian.C.
References Foam::dimless, laplacian(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), IOobject::NO_READ, and IOobject::time().
tmp< fvMatrix< Type > > laplacian | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 100 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 116 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const one & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 131 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const one & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 144 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const dimensioned< GType > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 156 of file fvmLaplacian.C.
References IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), dimensioned< Type >::name(), and IOobject::NO_READ.
tmp< fvMatrix< Type > > laplacian | ( | const dimensioned< GType > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 182 of file fvmLaplacian.C.
References IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), dimensioned< Type >::name(), and IOobject::NO_READ.
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvPatchField, volMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 209 of file fvmLaplacian.C.
References laplacian(), DimensionedField< Type, GeoMesh >::mesh(), laplacianScheme< Type, GType >::New(), and ref().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvPatchField, volMesh >> & | tgamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 226 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvPatchField, volMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 241 of file fvmLaplacian.C.
References laplacian(), and IOobject::name().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvPatchField, volMesh >> & | tgamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 258 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvsPatchField, surfaceMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 274 of file fvmLaplacian.C.
References laplacian(), DimensionedField< Type, GeoMesh >::mesh(), laplacianScheme< Type, GType >::New(), and ref().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> & | tgamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 291 of file fvmLaplacian.C.
References laplacian().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvsPatchField, surfaceMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 306 of file fvmLaplacian.C.
References laplacian(), and IOobject::name().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> & | tGamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 323 of file fvmLaplacian.C.
References laplacian().
tmp<fvMatrix<Type> > Foam::fvm::Su | ( | const DimensionedField< Type, volMesh > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Referenced by kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilonSource(), mixtureKEpsilon< BasicTurbulenceModel >::epsilonSource(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::kSource(), mixtureKEpsilon< BasicTurbulenceModel >::kSource(), and kOmegaSSTSAS< BasicTurbulenceModel >::Qsas().
tmp<fvMatrix<Type> > Foam::fvm::Su | ( | const tmp< DimensionedField< Type, volMesh >> & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Su | ( | const tmp< GeometricField< Type, fvPatchField, volMesh >> & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Definition at line 55 of file fvcSup.C.
References Foam::fvc::Sp().
zeroField Foam::fvm::Su | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const volScalarField::Internal & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Referenced by interRegionHeatTransferModel::addSup(), P1::calculate(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), radiativeIntensityRay::correct(), kEqn< BasicTurbulenceModel >::correct(), LamBremhorstKE::correct(), ShihQuadraticKE::correct(), LienLeschziner::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), kOmega< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), v2f< BasicTurbulenceModel >::correct(), kkLOmega::correct(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), continuousGasKEpsilon< BasicTurbulenceModel >::epsilonSource(), LaheyKEpsilon< BasicTurbulenceModel >::epsilonSource(), boundedDdtScheme< Type >::fvmDdt(), boundedConvectionScheme< Type >::fvmDiv(), continuousGasKEqn< BasicTurbulenceModel >::kSource(), continuousGasKEpsilon< BasicTurbulenceModel >::kSource(), NicenoKEqn< BasicTurbulenceModel >::kSource(), LaheyKEpsilon< BasicTurbulenceModel >::kSource(), thermoSingleLayer::q(), singleStepCombustion< ReactionThermo, ThermoType >::R(), radiationModel::Sh(), reactingOneDim::solveEnergy(), radiationModel::ST(), and laminar::Su().
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const tmp< volScalarField::Internal > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const tmp< volScalarField > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Definition at line 78 of file fvcSup.C.
References Foam::fvc::Sp().
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const dimensionedScalar & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Definition at line 90 of file fvcSup.C.
References Foam::fvc::SuSp().
zeroField Foam::fvm::Sp | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const volScalarField::Internal & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Referenced by relaxation::correct(), kEqn< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), v2f< BasicTurbulenceModel >::correct(), kkLOmega::correct(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::correct(), buoyantKEpsilon< BasicTurbulenceModel >::epsilonSource(), buoyantKEpsilon< BasicTurbulenceModel >::kSource(), and ThermoCloud< Foam::DSMCCloud >::Sh().
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const tmp< volScalarField::Internal > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const tmp< volScalarField > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
zeroField Foam::fvm::SuSp | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |