27 #include "surfaceInterpolate.H" 44 scalar backwardDdtScheme<Type>::deltaT_()
const 46 return mesh().time().deltaTValue();
51 scalar backwardDdtScheme<Type>::deltaT0_()
const 53 return mesh().time().deltaT0Value();
58 template<
class GeoField>
59 scalar backwardDdtScheme<Type>::deltaT0_(
const GeoField& vf)
const 61 if (vf.nOldTimes() < 2)
75 tmp<GeometricField<Type, fvPatchField, volMesh>>
86 mesh().time().timeName(),
90 scalar deltaT = deltaT_();
91 scalar deltaT0 = deltaT0_();
93 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
94 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
95 scalar coefft0 = coefft + coefft00;
114 tdtdt.
ref().primitiveFieldRef() = rDeltaT.
value()*dt.
value()*
116 coefft - (coefft0*
mesh().V0() - coefft00*
mesh().V00())/
mesh().V()
153 "ddt("+vf.
name()+
')',
154 mesh().time().timeName(),
158 scalar deltaT = deltaT_();
159 scalar deltaT0 = deltaT0_(vf);
161 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
162 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
163 scalar coefft0 = coefft + coefft00;
179 - coefft00*vf.
oldTime().oldTime().primitiveField()
187 coefft0*vf.
oldTime().boundaryField()
188 - coefft00*vf.
oldTime().oldTime().boundaryField()
205 + coefft00*vf.
oldTime().oldTime()
225 "ddt("+rho.
name()+
','+vf.
name()+
')',
226 mesh().time().timeName(),
230 scalar deltaT = deltaT_();
231 scalar deltaT0 = deltaT0_(vf);
233 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
234 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
235 scalar coefft0 = coefft + coefft00;
251 - coefft00*vf.
oldTime().oldTime().primitiveField()
259 coefft0*vf.
oldTime().boundaryField()
260 - coefft00*vf.
oldTime().oldTime().boundaryField()
277 + coefft00*vf.
oldTime().oldTime()
297 "ddt("+rho.
name()+
','+vf.
name()+
')',
298 mesh().time().timeName(),
302 scalar deltaT = deltaT_();
303 scalar deltaT0 = deltaT0_(vf);
305 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
306 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
307 scalar coefft0 = coefft + coefft00;
325 *vf.
oldTime().oldTime().primitiveField()*
mesh().V00()
335 *vf.
oldTime().oldTime().boundaryField()
374 mesh().time().timeName(),
378 scalar deltaT = deltaT_();
379 scalar deltaT0 = deltaT0_(vf);
381 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
382 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
383 scalar coefft0 = coefft + coefft00;
410 *vf.
oldTime().oldTime().primitiveField()*
mesh().V00()
428 *vf.
oldTime().oldTime().boundaryField()
472 scalar rDeltaT = 1.0/deltaT_();
474 scalar deltaT = deltaT_();
475 scalar deltaT0 = deltaT0_(vf);
477 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
478 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
479 scalar coefft0 = coefft + coefft00;
481 fvm.
diag() = (coefft*rDeltaT)*
mesh().V();
488 - coefft00*vf.
oldTime().oldTime().primitiveField()
496 coefft0*vf.
oldTime().primitiveField()
497 - coefft00*vf.
oldTime().oldTime().primitiveField()
523 scalar rDeltaT = 1.0/deltaT_();
525 scalar deltaT = deltaT_();
526 scalar deltaT0 = deltaT0_(vf);
528 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
529 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
530 scalar coefft0 = coefft + coefft00;
539 - coefft00*vf.
oldTime().oldTime().primitiveField()
547 coefft0*vf.
oldTime().primitiveField()
548 - coefft00*vf.
oldTime().oldTime().primitiveField()
574 scalar rDeltaT = 1.0/deltaT_();
576 scalar deltaT = deltaT_();
577 scalar deltaT0 = deltaT0_(vf);
579 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
580 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
581 scalar coefft0 = coefft + coefft00;
592 *vf.
oldTime().oldTime().primitiveField()*
mesh().V00()
602 *vf.
oldTime().oldTime().primitiveField()
629 scalar rDeltaT = 1.0/deltaT_();
631 scalar deltaT = deltaT_();
632 scalar deltaT0 = deltaT0_(vf);
634 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
635 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
636 scalar coefft0 = coefft + coefft00;
653 *vf.
oldTime().oldTime().primitiveField()*
mesh().V00()
668 *vf.
oldTime().oldTime().primitiveField()
686 scalar deltaT = deltaT_();
687 scalar deltaT0 = deltaT0_(U);
689 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
690 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
691 scalar coefft0 = coefft + coefft00;
699 "ddtCorr(" + U.
name() +
',' + Uf.
name() +
')',
700 mesh().time().timeName(),
730 scalar deltaT = deltaT_();
731 scalar deltaT0 = deltaT0_(U);
733 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
734 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
735 scalar coefft0 = coefft + coefft00;
743 "ddtCorr(" + U.
name() +
',' + phi.
name() +
')',
744 mesh().time().timeName(),
779 scalar deltaT = deltaT_();
780 scalar deltaT0 = deltaT0_(U);
782 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
783 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
784 scalar coefft0 = coefft + coefft00;
804 mesh().time().timeName(),
807 this->fvcDdtPhiCoeff(rhoU0,
mesh().Sf() & Uf.
oldTime())
825 return fvcDdtUfCorr(U, Uf);
830 <<
"dimensions of phi are not correct" 833 return fluxFieldType::null();
855 scalar deltaT = deltaT_();
856 scalar deltaT0 = deltaT0_(U);
858 scalar coefft = 1 + deltaT/(deltaT + deltaT0);
859 scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
860 scalar coefft0 = coefft + coefft00;
880 mesh().time().timeName(),
883 this->fvcDdtPhiCoeff(rhoU0, phi.
oldTime())
890 coefft0*rhoU0 - coefft00*rhoU00
902 return fvcDdtPhiCorr(U, phi);
907 <<
"dimensions of phi are not correct" 910 return fluxFieldType::null();
921 scalar deltaT = deltaT_();
922 scalar deltaT0 = deltaT0_(vf);
925 scalar coefft0_00 = deltaT/(deltaT + deltaT0);
928 scalar coefftn_0 = 1 + coefft0_00;
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
const dimensionSet & dimensions() const
Return const reference to dimensions.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
const dimensionSet dimVol(dimVolume)
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Type & value() const
Return const reference to value.
Generic GeometricField class.
Generic dimensioned Type class.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
const word & name() const
Return const reference to name.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
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. Face addressing is used to make all matrix assembly and solution loops vectorise.
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
errorManip< error > abort(error &err)
const dimensionSet & dimensions() const
Return dimensions.
Calculate the divergence of the given field.
word name(const complex &)
Return a string representation of a complex.
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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
A class for managing temporary objects.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
T & ref() const
Return non-const reference or generate a fatal error.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const word & name() const
Return name.
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity