Go to the documentation of this file.
29 #include "surfaceInterpolate.H"
45 const volScalarField& localEulerDdtScheme<Type>::localRDeltaT()
const
59 tmp<GeometricField<Type, fvPatchField, volMesh>>
67 "ddt(" + dt.
name() +
')',
68 mesh().time().timeName(),
93 "ddt(" + vf.name() +
')',
94 mesh().time().timeName(),
121 "ddt(" +
rho.name() +
',' + vf.name() +
')',
122 mesh().time().timeName(),
149 "ddt(" +
rho.name() +
',' + vf.name() +
')',
150 mesh().time().timeName(),
179 mesh().time().timeName(),
209 "ddt("+sf.name()+
')',
210 mesh().time().timeName(),
245 fvm.diag() = rDeltaT*
mesh().Vsc();
272 fvm.diag() = rDeltaT*
rho.value()*
mesh().Vsc();
301 fvm.diag() = rDeltaT*
rho.primitiveField()*
mesh().Vsc();
304 *
rho.oldTime().primitiveField()
333 rDeltaT*
alpha.primitiveField()*
rho.primitiveField()*
mesh().Vsc();
336 *
alpha.oldTime().primitiveField()
337 *
rho.oldTime().primitiveField()
366 "ddtCorr(" +
U.name() +
',' +
Uf.name() +
')',
367 mesh().time().timeName(),
370 this->fvcDdtPhiCoeff(
U.oldTime(), phiUf0, phiCorr)
398 "ddtCorr(" +
U.name() +
',' +
phi.name() +
')',
399 mesh().time().timeName(),
402 this->fvcDdtPhiCoeff(
U.oldTime(),
phi.oldTime(), phiCorr)
428 rho.oldTime()*
U.oldTime()
441 +
rho.name() +
',' +
U.name() +
',' +
Uf.name() +
')',
442 mesh().time().timeName(),
445 this->fvcDdtPhiCoeff(rhoU0, phiUf0, phiCorr,
rho.oldTime())
469 +
rho.name() +
',' +
U.name() +
',' +
Uf.name() +
')',
470 mesh().time().timeName(),
486 <<
"dimensions of Uf are not correct"
489 return fluxFieldType::null();
513 rho.oldTime()*
U.oldTime()
528 +
rho.name() +
',' +
U.name() +
',' +
phi.name() +
')',
529 mesh().time().timeName(),
560 +
rho.name() +
',' +
U.name() +
',' +
phi.name() +
')',
561 mesh().time().timeName(),
577 <<
"dimensions of phi are not correct"
580 return fluxFieldType::null();
609 tmeshPhi.
ref().setOriented();
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimVelocity
const dimensionSet dimDensity
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
autoPtr< surfaceVectorField > Uf
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
const word & name() const
Return const reference to name.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Generic dimensioned Type class.
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
errorManip< error > abort(error &err)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const dimensionSet dimVolume(pow3(dimLength))