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 41 of file fvmD2dt2.C.
References 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 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 d2dt2Scheme< Type >::New(), and rho.
tmp< fvMatrix< Type > > ddt | ( | const GeometricField< Type, fvPatchField, volMesh > & | vf | ) |
Definition at line 41 of file fvmDdt.C.
References 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(), kOmega< BasicTurbulenceModel >::correct(), LienLeschziner::correct(), thixotropicViscosity::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kkLOmega::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 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 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 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, dimensioned::name(), ddtScheme< Type >::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 | ||
) |
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::name(), and convectionScheme< Type >::New().
Referenced by ATCstandard::addATC(), ATCUaGradU::addATC(), relaxation::correct(), IATE::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), radiativeIntensityRay::correct(), kEqn< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), kOmega< BasicTurbulenceModel >::correct(), LienLeschziner::correct(), thixotropicViscosity::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kkLOmega::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), div(), age::execute(), scalarTransport::execute(), energyTransport::execute(), simple::mainIter(), adjointSimple::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::clear(), div(), and Foam::name().
tmp< fvMatrix< Type > > div | ( | const surfaceScalarField & | flux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
tmp< fvMatrix< Type > > div | ( | const tmp< surfaceScalarField > & | tflux, |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 84 of file fvmDiv.C.
References tmp::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, Foam::name(), and IOobject::NO_READ.
Referenced by jouleHeatingSource::addSup(), P1::calculate(), hydrostaticPressure::calculateAndWrite(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), kEqn< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), Poisson::correct(), LamBremhorstKE::correct(), LienLeschziner::correct(), kOmega< BasicTurbulenceModel >::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kkLOmega::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), incompressiblePrimalSolver::correctBoundaryConditions(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), adjointLaminar::divDevReff(), adjointSpalartAllmaras::divDevReff(), Maxwell< BasicTurbulenceModel >::divDevRhoReff(), kineticTheoryModel::divDevRhoReff(), thermo::evolveRegion(), age::execute(), scalarTransport::execute(), electricPotential::execute(), energyTransport::execute(), laplacian(), simple::mainIter(), adjointSimple::mainIter(), velocityComponentLaplacianFvMotionSolver::solve(), velocityLaplacianFvMotionSolver::solve(), adjointMeshMovementSolver::solve(), displacementComponentLaplacianFvMotionSolver::solve(), displacementSBRStressFvMotionSolver::solve(), solidBodyDisplacementLaplacianFvMotionSolver::solve(), laplacianMotionSolver::solve(), displacementLaplacianFvMotionSolver::solve(), twoPhaseSystem::solve(), elasticityMotionSolver::solve(), surfaceAlignedSBRStressFvMotionSolver::solve(), adjointEikonalSolver::solve(), thermalBaffle::solveEnergy(), reactingOneDim::solveEnergy(), sensitivityBezierFI::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(), and IOobject::NO_READ.
tmp< fvMatrix< Type > > laplacian | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf, | ||
const word & | name | ||
) |
Definition at line 95 of file fvmLaplacian.C.
tmp< fvMatrix< Type > > laplacian | ( | const zero & | , |
const GeometricField< Type, fvPatchField, volMesh > & | vf | ||
) |
Definition at line 111 of file fvmLaplacian.C.
tmp< fvMatrix< Type > > laplacian | ( | const 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 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, laplacian(), Foam::name(), and IOobject::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, laplacian(), 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 204 of file fvmLaplacian.C.
References gamma, Foam::name(), and laplacianScheme< Type, GType >::New().
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, and laplacian().
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, Foam::name(), and laplacianScheme< Type, GType >::New().
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, and laplacian().
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().
tmp<fvMatrix<Type> > Foam::fvm::Su | ( | const DimensionedField< Type, volMesh > & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
Referenced by ATCstandard::addATC(), ATCUaGradU::addATC(), 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 > & | |||
) |
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 multiphaseMangrovesSource::addSup(), acousticDampingSource::addSup(), interRegionHeatTransferModel::addSup(), P1::calculate(), IATE::correct(), DeardorffDiffStress< BasicTurbulenceModel >::correct(), radiativeIntensityRay::correct(), kEqn< BasicTurbulenceModel >::correct(), Maxwell< BasicTurbulenceModel >::correct(), dynamicLagrangian< BasicTurbulenceModel >::correct(), ShihQuadraticKE::correct(), LamBremhorstKE::correct(), kOmega< BasicTurbulenceModel >::correct(), LienLeschziner::correct(), thixotropicViscosity::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), LienCubicKE::correct(), SSG< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), LRR< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), SpalartAllmarasDES< BasicTurbulenceModel >::correct(), advectionDiffusion::correct(), SpalartAllmaras< BasicTurbulenceModel >::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kkLOmega::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), waxSolventEvaporation::correctModel(), kOmegaSSTLM< BasicTurbulenceModel >::correctReThetatGammaInt(), continuousGasKEpsilon< BasicTurbulenceModel >::epsilonSource(), LaheyKEpsilon< BasicTurbulenceModel >::epsilonSource(), energyTransport::execute(), 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(), populationBalanceModel::solve(), 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 > & | |||
) |
tmp<fvMatrix<Type> > Foam::fvm::Sp | ( | const dimensionedScalar & | , |
const GeometricField< Type, fvPatchField, volMesh > & | |||
) |
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(), IATE::correct(), kEqn< BasicTurbulenceModel >::correct(), kOmega< BasicTurbulenceModel >::correct(), thixotropicViscosity::correct(), dynamicKEqn< BasicTurbulenceModel >::correct(), realizableKE< BasicTurbulenceModel >::correct(), kEpsilon< EddyDiffusivity< compressible::turbulenceModel > >::correct(), LaunderSharmaKE< BasicTurbulenceModel >::correct(), RNGkEpsilon< BasicTurbulenceModel >::correct(), qZeta::correct(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), kineticTheoryModel::correct(), mixtureKEpsilon< BasicTurbulenceModel >::correct(), kEpsilonPhitF< BasicTurbulenceModel >::correct(), kkLOmega::correct(), adjointSpalartAllmaras::correct(), kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::correct(), kL< BasicTurbulenceModel >::correct(), buoyantKEpsilon< BasicTurbulenceModel >::epsilonSource(), buoyantKEpsilon< BasicTurbulenceModel >::kSource(), ThermoCloud< CloudType >::Sh(), adjointEikonalSolver::solve(), and populationBalanceModel::solve().
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 > & | |||
) |
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.