43 tmp<GeometricField<Type, fvPatchField, volMesh>>
51 4.0/
sqr(
mesh().time().deltaT() +
mesh().time().deltaT0())
54 const word d2dt2name(
"d2dt2("+vf.
name()+
')');
56 const scalar deltaT =
mesh().time().deltaTValue();
57 const scalar deltaT0 =
mesh().time().deltaT0Value();
59 const scalar coefft = (deltaT + deltaT0)/(2*deltaT);
60 const scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
61 const scalar coefft0 = coefft + coefft00;
77 - (coefft*VV0 + coefft00*V0V00)
80 + (coefft00*V0V00)*vf.
oldTime().oldTime()()
85 - coefft0*vf.
oldTime().boundaryField()
86 + coefft00*vf.
oldTime().oldTime().boundaryField()
99 + coefft00*vf.
oldTime().oldTime()
116 4.0/
sqr(
mesh().time().deltaT() +
mesh().time().deltaT0())
119 const word d2dt2name(
"d2dt2("+rho.
name()+
','+vf.
name()+
')');
121 const scalar deltaT =
mesh().time().deltaTValue();
122 const scalar deltaT0 =
mesh().time().deltaT0Value();
124 const scalar coefft = (deltaT + deltaT0)/(2*deltaT);
125 const scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
151 coefft*VV0rhoRho0*vf()
153 - (coefft*VV0rhoRho0 + coefft00*V0V00rho0Rho00)
156 + (coefft00*V0V00rho0Rho00)
159 halfRdeltaT2.
value()*
182 )*vf.
oldTime().oldTime().boundaryField()
199 - (coefft*rhoRho0 + coefft00*rho0Rho00)*vf.
oldTime()
200 + coefft00*rho0Rho00*vf.
oldTime().oldTime()
225 const scalar deltaT =
mesh().time().deltaTValue();
226 const scalar deltaT0 =
mesh().time().deltaT0Value();
228 const scalar coefft = (deltaT + deltaT0)/(2*deltaT);
229 const scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
230 const scalar coefft0 = coefft + coefft00;
232 const scalar rDeltaT2 = 4.0/
sqr(deltaT + deltaT0);
236 const scalar halfRdeltaT2 = rDeltaT2/2.0;
241 fvm.
diag() = (coefft*halfRdeltaT2)*VV0;
243 fvm.
source() = halfRdeltaT2*
245 (coefft*VV0 + coefft00*V0V00)
248 - (coefft00*V0V00)*vf.
oldTime().oldTime().primitiveField()
253 fvm.
diag() = (coefft*rDeltaT2)*
mesh().V();
257 coefft0*vf.
oldTime().primitiveField()
258 - coefft00*vf.
oldTime().oldTime().primitiveField()
286 const scalar deltaT =
mesh().time().deltaTValue();
287 const scalar deltaT0 =
mesh().time().deltaT0Value();
289 const scalar coefft = (deltaT + deltaT0)/(2*deltaT);
290 const scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
292 const scalar rDeltaT2 = 4.0/
sqr(deltaT + deltaT0);
296 const scalar halfRdeltaT2 = 0.5*rDeltaT2;
301 fvm.
diag() = rho.
value()*(coefft*halfRdeltaT2)*VV0;
305 (coefft*VV0 + coefft00*V0V00)
308 - (coefft00*V0V00)*vf.
oldTime().oldTime().primitiveField()
317 (coefft + coefft00)*vf.
oldTime().primitiveField()
318 - coefft00*vf.
oldTime().oldTime().primitiveField()
346 const scalar deltaT =
mesh().time().deltaTValue();
347 const scalar deltaT0 =
mesh().time().deltaT0Value();
349 const scalar coefft = (deltaT + deltaT0)/(2*deltaT);
350 const scalar coefft00 = (deltaT + deltaT0)/(2*deltaT0);
352 const scalar rDeltaT2 = 4.0/
sqr(deltaT + deltaT0);
356 const scalar quarterRdeltaT2 = 0.25*rDeltaT2;
373 fvm.
diag() = (coefft*quarterRdeltaT2)*VV0rhoRho0;
375 fvm.
source() = quarterRdeltaT2*
377 (coefft*VV0rhoRho0 + coefft00*V0V00rho0Rho00)
380 - (coefft00*V0V00rho0Rho00)
381 *vf.
oldTime().oldTime().primitiveField()
386 const scalar halfRdeltaT2 = 0.5*rDeltaT2;
400 fvm.
diag() = (coefft*halfRdeltaT2)*
mesh().V()*rhoRho0;
404 (coefft*rhoRho0 + coefft00*rho0Rho00)
407 - (coefft00*rho0Rho00)
408 *vf.
oldTime().oldTime().primitiveField()
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const word & name() const
Return name.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
tmp< fvMatrix< Type > > fvmD2dt2(const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
A class for handling words, derived from string.
const dimensionSet dimVol
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Calculate the divergence of the given field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet & dimensions() const
Return const reference to dimensions.
A class for managing temporary objects.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcD2dt2(const GeometricField< Type, fvPatchField, volMesh > &)