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 > | |
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 47 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and Time::New().
tmp< fvMatrix< Type > > d2dt2 | ( | const dimensionedScalar & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 62 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Time::New(), and rho.
tmp< fvMatrix< Type > > d2dt2 | ( | const volScalarField & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 78 of file fvmD2dt2.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Time::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 47 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), and Time::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(), SpalartAllmarasDES< BasicTurbulenceModel >::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(), SpalartAllmaras< 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().
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 74 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Time::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | rho, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 90 of file fvmDdt.C.
References DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Time::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const volScalarField & | alpha, |
const volScalarField & | rho, | ||
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 106 of file fvmDdt.C.
References alpha, DimensionedField< Type, GeoMesh >::mesh(), IOobject::name(), Time::New(), and rho.
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 47 of file fvmDiv.C.
References DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and Time::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(), SpalartAllmarasDES< BasicTurbulenceModel >::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(), SpalartAllmaras< 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().
tmp< fvMatrix< Type > > div | ( | const tmp< surfaceScalarField > & | tflux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 64 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 79 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 90 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 47 of file fvmLaplacian.C.
References TimePaths::constant(), Foam::dimless, laplacian(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), IOobject::NO_READ, and IOobject::time().
Referenced by 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(), SpalartAllmarasDES< BasicTurbulenceModel >::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(), SpalartAllmaras< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), linearViscousStress< LESModel< BasicTurbulenceModel > >::correct(), advectionDiffusion::correct(), Poisson::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), adjointkOmegaSST::divDevReff(), adjointLaminar::divDevReff(), adjointSpalartAllmaras::divDevReff(), linearViscousStress< BasicTurbulenceModel >::divDevRhoReff(), 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(), adjointEikonalSolver::solve(), adjointMeshMovementSolver::solve(), twoPhaseSystem::solve(), reactingOneDim::solveEnergy(), thermalBaffle::solveEnergy(), sensitivityBezierFI::solveMeshMovementEqn(), kinematicSingleLayer::solveThickness(), MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi(), and MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::YiEqn().
tmp< fvMatrix< Type > > laplacian | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 72 of file fvmLaplacian.C.
References TimePaths::constant(), 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 101 of file fvmLaplacian.C.
tmp< fvMatrix< Type > > laplacian | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 117 of file fvmLaplacian.C.
tmp< fvMatrix< Type > > laplacian | ( | const one & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 132 of file fvmLaplacian.C.
References laplacian(), and Foam::name().
tmp< fvMatrix< Type > > laplacian | ( | const one & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 145 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 157 of file fvmLaplacian.C.
References gamma, IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and IOobject::NO_READ.
tmp< fvMatrix< Type > > laplacian | ( | const dimensioned< GType > & | gamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 183 of file fvmLaplacian.C.
References gamma, IOobject::instance(), laplacian(), DimensionedField< Type, GeoMesh >::mesh(), 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 210 of file fvmLaplacian.C.
References gamma, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and Time::New().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvPatchField, volMesh > > & | tgamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 227 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 242 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 259 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 275 of file fvmLaplacian.C.
References gamma, DimensionedField< Type, GeoMesh >::mesh(), Foam::name(), and Time::New().
tmp< fvMatrix< Type > > laplacian | ( | const tmp< GeometricField< GType, fvsPatchField, surfaceMesh > > & | tgamma, |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 292 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 307 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 324 of file fvmLaplacian.C.
References laplacian().
zeroField Su | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A no-op source.
Referenced by mixtureKEpsilon< BasicTurbulenceModel >::epsilonSource(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::epsilonSource(), mixtureKEpsilon< BasicTurbulenceModel >::kSource(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::kSource(), kOmegaSSTSAS< BasicTurbulenceModel >::Qsas(), dummy::R(), and turbulentBreakUp::R().
tmp< fvMatrix< Type > > Su | ( | const dimensioned< Type > & | su, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A uniform source (no-op for small values)
tmp< fvMatrix< Type > > Su | ( | const DimensionedField< Type, volMesh > & | su, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp< fvMatrix< Type > > Su | ( | const tmp< DimensionedField< Type, volMesh > > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp< fvMatrix< Type > > Su | ( | const tmp< GeometricField< Type, fvPatchField, volMesh > > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
zeroField Sp | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A no-op source.
Referenced by solver::addOptimisationTypeSource(), PhaseLimitStabilization< Type >::addSup(), interRegionHeatTransferModel::addSup(), multiphaseMangrovesTurbulenceModel::addSup(), atmPlantCanopyUSource::addSup(), acousticDampingSource::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(), SpalartAllmarasDES< BasicTurbulenceModel >::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(), SpalartAllmaras< BasicTurbulenceModel >::correct(), SSG< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), advectionDiffusion::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(), MultiComponentPhaseModel< BasePhaseModel, phaseThermo >::solveYi(), radiationModel::ST(), laminar::Su(), MovingPhaseModel< BasePhaseModel >::UfEqn(), and MassTransferPhaseSystem< BasePhaseSystem >::volTransfer().
tmp< fvMatrix< Type > > Sp | ( | const dimensionedScalar & | sp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A uniform source (no-op for small values)
tmp< fvMatrix< Type > > Sp | ( | const DimensionedField< scalar, volMesh > & | sp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp< fvMatrix< Type > > Sp | ( | const tmp< DimensionedField< scalar, volMesh > > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp< fvMatrix< Type > > Sp | ( | const tmp< volScalarField > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
zeroField SuSp | ( | const Foam::zero | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A no-op source.
Referenced by SemiImplicitSource< Type >::addSup(), multiphaseStabilizedTurbulence::addSup(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kL< BasicTurbulenceModel >::correct(), adjointkOmegaSST::correct(), adjointSpalartAllmaras::correct(), kineticTheoryModel::correct(), IATE::correct(), kkLOmega::correct(), qZeta::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kOmegaSSTBase< BasicEddyViscosityModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), relaxation::correct(), buoyantKEpsilon< BasicTurbulenceModel >::epsilonSource(), buoyantKEpsilon< BasicTurbulenceModel >::kSource(), phaseChange::R(), wakeEntrainmentCoalescence::R(), adjointEikonalSolver::solve(), populationBalanceModel::solve(), MovingPhaseModel< BasePhaseModel >::UEqn(), MovingPhaseModel< BasePhaseModel >::UfEqn(), and adjointkOmegaSST::waEqnSourceFromCDkOmega().
tmp< fvMatrix< Type > > SuSp | ( | const dimensionedScalar & | susp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
A uniform source (no-op for small values)
tmp< fvMatrix< Type > > SuSp | ( | const DimensionedField< scalar, volMesh > & | susp, |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp< fvMatrix< Type > > SuSp | ( | const tmp< DimensionedField< scalar, volMesh > > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
tmp< fvMatrix< Type > > SuSp | ( | const tmp< volScalarField > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |