28 #include "surfaceInterpolate.H"
60 ) <<
"Ddt scheme not specified" <<
endl <<
endl
61 <<
"Valid ddt schemes are :" <<
endl
62 << IstreamConstructorTablePtr_->sortedToc()
66 const word schemeName(schemeData);
68 typename IstreamConstructorTable::iterator cstrIter =
69 IstreamConstructorTablePtr_->find(schemeName);
71 if (cstrIter == IstreamConstructorTablePtr_->end())
76 ) <<
"Unknown ddt scheme " << schemeName <<
nl <<
nl
77 <<
"Valid ddt schemes are :" <<
endl
78 << IstreamConstructorTablePtr_->sortedToc()
82 return cstrIter()(mesh, schemeData);
184 if (!
U.boundaryField()[
patchi].coupled())
193 <<
"ddtCouplingCoeff mean max min = "
200 return tddtCouplingCoeff;
213 return fvcDdtPhiCoeff(
U, phi, phiCorr);
236 return fvcDdtPhiCoeff
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
bool eof() const
Return true if end of input seen.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)=0
virtual ~ddtScheme()
Destructor.
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)=0
virtual tmp< surfaceScalarField > fvcDdtPhiCoeff(const VolField< Type > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
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.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
volScalarField sf(fieldObject, mesh)
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define InfoInFunction
Report an information message using Foam::Info.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
errorManipArg< error, int > exit(error &err, const int errNo=1)
static const scalar SMALL
Ostream & endl(Ostream &os)
Add newline and flush stream.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gAverage(const FieldField< Field, Type > &f)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.