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 > &) |
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 41 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and d2dt2Scheme< Type >::New().
tmp< fvMatrix< Type > > d2dt2 | ( | const dimensionedScalar & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 56 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), d2dt2Scheme< Type >::New(), and rho.
tmp< fvMatrix< Type > > d2dt2 | ( | const volScalarField & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 72 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), d2dt2Scheme< Type >::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 41 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and ddtScheme< Type >::New().
Referenced by multiphaseMangrovesSource::addSup(), relaxation::correct(), IATE::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), LienLeschziner::correct(), kOmega< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), thixotropicViscosity::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), SSG< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), EBRSM< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), kkLOmega::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), ddt(), thermo::evolveRegion(), scalarTransport::execute(), energyTransport::execute(), twoPhaseSystem::solve(), populationBalanceModel::solve(), reactingOneDim::solveContinuity(), kinematicSingleLayer::solveContinuity(), thermalBaffle::solveEnergy(), reactingOneDim::solveEnergy(), thermoSingleLayer::solveEnergy(), kinematicSingleLayer::solveMomentum(), reactingOneDim::solveSpeciesMass(), and kinematicSingleLayer::solveThickness().
tmp< fvMatrix< Type > > ddt | ( | const Foam::one | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > ddt | ( | const dimensionedScalar & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 68 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), ddtScheme< Type >::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 84 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), ddtScheme< Type >::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 100 of file fvmDdt.C.
References Foam::constant::atomic::alpha, DimensionedField< Type, GeoMesh >::mesh(), dimensioned< Type >::name(), IOobject::name(), ddtScheme< Type >::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const Foam::one | , |
const Foam::one | , | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > ddt | ( | const Foam::one | , |
const volScalarField & | rho, | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | alpha, |
const Foam::one | , | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 149 of file fvmDdt.C.
References Foam::constant::atomic::alpha, and ddt().
tmp< fvMatrix< Type > > div | ( | const surfaceScalarField & | flux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 41 of file fvmDiv.C.
References Foam::fvc::flux(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and convectionScheme< Type >::New().
Referenced by ATCUaGradU::addATC(), ATCstandard::addATC(), relaxation::correct(), IATE::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), LienLeschziner::correct(), kOmega< BasicTurbulenceModel >::correct(), radiativeIntensityRay::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), thixotropicViscosity::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), SSG< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), EBRSM< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), kineticTheoryModel::correct(), kkLOmega::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), div(), age::execute(), scalarTransport::execute(), energyTransport::execute(), adjointSimple::mainIter(), simple::mainIter(), adjointEikonalSolver::solve(), 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 58 of file fvmDiv.C.
References tmp< T >::clear(), div(), and Foam::name().
tmp< fvMatrix< Type > > div | ( | const surfaceScalarField & | flux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 73 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 84 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 41 of file fvmLaplacian.C.
References Foam::dimless, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), IOobjectOption::NO_READ, and IOobject::time().
Referenced by isotropic::addRegularisationTerm(), jouleHeatingSource::addSup(), P1::calculate(), hydrostaticPressure::calculateAndWrite(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), Poisson::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), LienLeschziner::correct(), kOmega< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), SSG< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), EBRSM< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), kineticTheoryModel::correct(), kkLOmega::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), adjointLaminar::divDevReff(), adjointSpalartAllmaras::divDevReff(), adjointkOmegaSST::divDevReff(), Maxwell< BasicTurbulenceModel >::divDevRhoReff(), kineticTheoryModel::divDevRhoReff(), thermo::evolveRegion(), age::execute(), scalarTransport::execute(), energyTransport::execute(), electricPotential::execute(), laplacian(), adjointSimple::mainIter(), simple::mainIter(), velocityComponentLaplacianFvMotionSolver::solve(), velocityLaplacianFvMotionSolver::solve(), laplacianMotionSolver::solve(), displacementComponentLaplacianFvMotionSolver::solve(), displacementSBRStressFvMotionSolver::solve(), adjointMeshMovementSolver::solve(), solidBodyDisplacementLaplacianFvMotionSolver::solve(), displacementLaplacianFvMotionSolver::solve(), pLaplacianMotionSolver::solve(), elasticityMotionSolver::solve(), twoPhaseSystem::solve(), surfaceAlignedSBRStressFvMotionSolver::solve(), adjointEikonalSolver::solve(), thermalBaffle::solveEnergy(), reactingOneDim::solveEnergy(), shapeDesignVariables::solveMeshMovementEqn(), and kinematicSingleLayer::solveThickness().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 66 of file fvmLaplacian.C.
References Foam::dimless, laplacian(), DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), IOobjectOption::NO_READ, and IOobject::time().
tmp< fvMatrix< Type > > laplacian | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 95 of file fvmLaplacian.C.
tmp< fvMatrix< Type > > laplacian | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 111 of file fvmLaplacian.C.
tmp< fvMatrix< Type > > laplacian | ( | const Foam::one | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 126 of file fvmLaplacian.C.
References laplacian(), and Foam::name().
tmp< fvMatrix< Type > > laplacian | ( | const Foam::one | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 139 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 151 of file fvmLaplacian.C.
References gamma, IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and IOobjectOption::NO_READ.
tmp< fvMatrix< Type > > laplacian | ( | const dimensioned< GType > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 177 of file fvmLaplacian.C.
References gamma, IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), and IOobjectOption::NO_READ.
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvPatchField, volMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 204 of file fvmLaplacian.C.
References gamma, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), 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 221 of file fvmLaplacian.C.
References laplacian(), and Foam::name().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvPatchField, volMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 236 of file fvmLaplacian.C.
References gamma, laplacian(), and IOobject::name().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvPatchField, volMesh >> & | tgamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 253 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 269 of file fvmLaplacian.C.
References gamma, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), 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 286 of file fvmLaplacian.C.
References laplacian(), and Foam::name().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< GType, fvsPatchField, surfaceMesh > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 301 of file fvmLaplacian.C.
References gamma, laplacian(), and IOobject::name().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvsPatchField, surfaceMesh >> & | tGamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 318 of file fvmLaplacian.C.
References laplacian().
zeroField Foam::fvm::Su | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A no-op source.
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 dimensioned< Type > & | su, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A uniform source (no-op for small values)
tmp<fvMatrix<Type> > Foam::fvm::Su | ( | const DimensionedField< Type, volMesh > & | su, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
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 > & | |||
) |
zeroField Foam::fvm::Sp | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A no-op source.
Referenced by topOSource::addSup(), multiphaseMangrovesSource::addSup(), atmPlantCanopyUSource::addSup(), acousticDampingSource::addSup(), interRegionHeatTransferModel::addSup(), P1::calculate(), IATE::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), LienLeschziner::correct(), radiativeIntensityRay::correct(), kOmega< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), thixotropicViscosity::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), SSG< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), SpalartAllmarasBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), EBRSM< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), kkLOmega::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), adjointkOmegaSST::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), continuousGasKEpsilon< BasicTurbulenceModel >::epsilonSource(), LaheyKEpsilon< BasicTurbulenceModel >::epsilonSource(), energyTransport::execute(), boundedDdtScheme< Type >::fvmDdt(), boundedConvectionScheme< Type >::fvmDiv(), MassTransferPhaseSystem< BasePhaseSystem >::heatTransfer(), continuousGasKEqn< BasicTurbulenceModel >::kSource(), continuousGasKEpsilon< BasicTurbulenceModel >::kSource(), NicenoKEqn< BasicTurbulenceModel >::kSource(), LaheyKEpsilon< BasicTurbulenceModel >::kSource(), thermoSingleLayer::q(), singleStepCombustion< ReactionThermo, ThermoType >::R(), radiationModel::Sh(), populationBalanceModel::solve(), multiphaseSystem::solveAlphas(), radiationModel::ST(), laminar::Su(), and MassTransferPhaseSystem< BasePhaseSystem >::volTransfer().
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const dimensionedScalar & | sp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const DimensionedField< scalar, volMesh > & | sp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const tmp< DimensionedField< scalar, volMesh >> & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const tmp< volScalarField > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
zeroField Foam::fvm::SuSp | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A no-op source.
Referenced by relaxation::correct(), IATE::correct(), kEqn< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), kineticTheoryModel::correct(), kkLOmega::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), buoyantKEpsilon< BasicTurbulenceModel >::epsilonSource(), buoyantKEpsilon< BasicTurbulenceModel >::kSource(), ThermoCloud< Foam::DSMCCloud >::Sh(), adjointEikonalSolver::solve(), populationBalanceModel::solve(), and adjointkOmegaSST::waEqnSourceFromCDkOmega().
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const dimensionedScalar & | susp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A uniform source (no-op for small values)
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const DimensionedField< scalar, volMesh > & | susp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const tmp< DimensionedField< scalar, volMesh >> & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::SuSp | ( | const tmp< volScalarField > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |