Go to the documentation of this file.
53 Info<<
"ddtScheme<Type>::New(const fvMesh&, Istream&) : "
54 "constructing ddtScheme<Type>"
63 ) <<
"Ddt scheme not specified" <<
endl <<
endl
64 <<
"Valid ddt schemes are :" <<
endl
65 << IstreamConstructorTablePtr_->sortedToc()
69 const word schemeName(schemeData);
71 typename IstreamConstructorTable::iterator cstrIter =
72 IstreamConstructorTablePtr_->find(schemeName);
74 if (cstrIter == IstreamConstructorTablePtr_->end())
79 ) <<
"Unknown ddt scheme " << schemeName <<
nl <<
nl
80 <<
"Valid ddt schemes are :" <<
endl
81 << IstreamConstructorTablePtr_->sortedToc()
85 return cstrIter()(
mesh, schemeData);
159 U.boundaryField()[
patchi].fixesValue()
169 Info<<
"ddtCouplingCoeff mean max min = "
176 return tddtCouplingCoeff;
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Type gAverage(const FieldField< Field, Type > &f)
bool eof() const
Return true if end of input seen.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
faceListList boundary(nPatches)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual ~ddtScheme()
Destructor.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
tmp< surfaceScalarField > fvcDdtPhiCoeff(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi, const fluxFieldType &phiCorr)
InternalField & internalField()
Return internal field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)=0
const dimensionSet dimVol(dimVolume)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Type gMin(const FieldField< Field, Type > &f)
const dimensionSet & dimensions() const
Return const reference to dimensions.
Generic GeometricField class.
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)=0