27 #include "surfaceInterpolate.H"
43 void SLTSDdtScheme<Type>::relaxedDiag
57 diag[owner[facei]] += phi[facei];
58 rD[neighbour[facei]] += phi[facei];
62 diag[neighbour[facei]] -= phi[facei];
63 rD[owner[facei]] -= phi[facei];
70 const labelUList& faceCells = pphi.patch().patch().faceCells();
74 if (pphi[patchFacei] > 0.0)
76 diag[faceCells[patchFacei]] += pphi[patchFacei];
80 rD[faceCells[patchFacei]] -= pphi[patchFacei];
85 rD += (1.0/alpha_ - 2.0)*
diag;
90 tmp<volScalarField> SLTSDdtScheme<Type>::SLrDeltaT()
const
93 mesh().objectRegistry::template
94 lookupObject<surfaceScalarField>(phiName_);
98 tmp<volScalarField> trDeltaT
105 extrapolatedCalculatedFvPatchScalarField::typeName
111 relaxedDiag(rDeltaT, phi);
113 if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
115 rDeltaT.primitiveFieldRef() =
max
117 rDeltaT.primitiveField()/
mesh().V(),
118 scalar(1)/deltaT.value()
121 else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
124 mesh().objectRegistry::template lookupObject<volScalarField>
129 rDeltaT.primitiveFieldRef() =
max
131 rDeltaT.primitiveField()/(
rho.primitiveField()*
mesh().V()),
132 scalar(1)/deltaT.value()
138 <<
"Incorrect dimensions of phi: " << phi.dimensions()
142 rDeltaT.correctBoundaryConditions();
157 const word ddtName(
"ddt("+dt.
name()+
')');
176 tdtdt.
ref().primitiveFieldRef() =
211 const word ddtName(
"ddt("+vf.
name()+
')');
252 const word ddtName(
"ddt("+
rho.name()+
','+vf.
name()+
')');
293 const word ddtName(
"ddt("+
rho.name()+
','+vf.
name()+
')');
310 -
rho.oldTime().boundaryField()
358 alpha.boundaryField()
362 -
alpha.oldTime().boundaryField()
363 *
rho.oldTime().boundaryField()
404 const scalarField rDeltaT(SLrDeltaT()().primitiveField());
439 const scalarField rDeltaT(SLrDeltaT()().primitiveField());
476 const scalarField rDeltaT(SLrDeltaT()().primitiveField());
483 *
rho.oldTime().primitiveField()
489 *
rho.oldTime().primitiveField()
520 const scalarField rDeltaT(SLrDeltaT()().primitiveField());
528 *
alpha.oldTime().primitiveField()
529 *
rho.oldTime().primitiveField()
535 *
alpha.oldTime().primitiveField()
536 *
rho.oldTime().primitiveField()
562 "ddtCorr(" +
U.name() +
',' + Uf.
name() +
')',
563 this->fvcDdtPhiCoeff(
U.oldTime(), phiUf0, phiCorr)
586 "ddtCorr(" +
U.name() +
',' + phi.
name() +
')',
587 this->fvcDdtPhiCoeff(
U.oldTime(), phi.
oldTime(), phiCorr)
612 rho.oldTime()*
U.oldTime()
620 "ddtCorr(" +
rho.name() +
',' +
U.name() +
',' + rhoUf.
name() +
')',
621 this->fvcDdtPhiCoeff(rhoU0, phiUf0, phiCorr,
rho.oldTime())
639 "ddtCorr(" +
rho.name() +
',' +
U.name() +
',' + rhoUf.
name() +
')',
652 <<
"dimensions of rhoUf are not correct"
655 return fluxFieldType::null();
679 rho.oldTime()*
U.oldTime()
689 "ddtCorr(" +
rho.name() +
',' +
U.name() +
',' + phi.
name() +
')',
712 "ddtCorr(" +
rho.name() +
',' +
U.name() +
',' + phi.
name() +
')',
725 <<
"dimensions of phi are not correct"
728 return fluxFieldType::null();
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const word & name() const
Return name.
const Field0Type & oldTime() const
Return the old-time field.
dimensionedScalar deltaT() const
Return time step.
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const Time & time() const
Return the top-level database.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const labelUList & owner() const
Internal face owner.
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
const labelUList & neighbour() const
Internal face neighbour.
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionSet dimVolumetricFlux
errorManip< error > abort(error &err)
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void diag(LagrangianPatchField< vector > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimVelocity
UList< label > labelUList
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
faceListList boundary(nPatches)