51 class polyTopoChangeMap;
53 class polyDistributionMap;
67 const word modelType_;
111 template<
class Type,
class ... AlphaRhoFieldTypes>
116 const AlphaRhoFieldTypes& ... alphaRhoFields
135 const word& modelType,
146 template<
class AlphaRhoFieldType,
class ... AlphaRhoFieldTypes>
150 const AlphaRhoFieldType& alphaRhoField,
151 const AlphaRhoFieldTypes& ... alphaRhoFields
159 template<
class AlphaRhoFieldType,
class ... AlphaRhoFieldTypes>
162 const AlphaRhoFieldType& alphaRhoField,
163 const AlphaRhoFieldTypes& ... alphaRhoFields
168 template<
class AlphaRhoFieldType>
169 static const word&
fieldName(
const AlphaRhoFieldType& alphaRhoField);
182 const word& modelType,
407 virtual
bool write(const
bool write = true) const;
417 #include "fvModelI.H"
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Dimension set for the base types.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
Return pointer to new fvModel object created.
autoPtr< fvModel > operator()(Istream &is) const
iNew(const fvMesh &mesh, const word &name)
Finite volume model abstract base class.
tmp< fvMatrix< Type > > sourceTerm(const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
Return a source for an equation.
virtual void mapMesh(const polyMeshMap &)=0
Update from another mesh using the given map.
static autoPtr< fvModel > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected fvModel.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual bool movePoints()=0
Update for mesh motion.
TypeName("fvModel")
Runtime type information.
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &field) const
Return source for an equation with a second time derivative.
virtual void distribute(const polyDistributionMap &)=0
Redistribute or update using the given distribution map.
void addSupType(const VolField< Type > &field, fvMatrix< Type > &eqn) const
Add a source term to an equation.
declareRunTimeSelectionTable(autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
virtual ~fvModel()
Destructor.
virtual void correct()
Correct the fvModel.
const dictionary & coeffs() const
Return dictionary.
static dimensionSet sourceDims(const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
Return the dimensions of the matrix of a source term.
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void preUpdateMesh()
Prepare for mesh update.
tmp< fvMatrix< Type > > sourceProxy(const VolField< Type > &eqnField) const
Add a source term to an equation.
virtual bool write(const bool write=true) const
Write fvModel data.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
const fvMesh & mesh() const
Return const access to the mesh database.
fvModel(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
autoPtr< fvModel > clone() const
Return clone.
const word & name() const
Return const access to the source name.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
static const word & fieldName()
Return the name of the field associated with a source term. Special.
virtual void topoChange(const polyTopoChangeMap &)=0
Update topology using the given map.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
A class for handling words, derived from string.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Forward declarations of fvMatrix specialisations.
#define DEFINE_FV_MODEL_ADD_ALPHA_RHO_FIELD_SUP(Type, nullArg)
#define DEFINE_FV_MODEL_ADD_RHO_FIELD_SUP(Type, nullArg)
#define DEFINE_FV_MODEL_ADD_FIELD_SUP(Type, nullArg)
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)