65 void Foam::fv::codedFvModel::readCoeffs(
const dictionary&
dict)
67 fieldName_ =
dict.lookup<word>(
"field");
71 Foam::word Foam::fv::codedFvModel::fieldPrimitiveTypeName()
const
73 #define fieldPrimitiveTypeNameTernary(Type, nullArg) \
74 mesh().foundObject<VolField<Type>>(fieldName_) \
75 ? pTraits<Type>::typeName \
82 void Foam::fv::codedFvModel::prepare
85 const dynamicCodeContext& context
88 const word primitiveTypeName = fieldPrimitiveTypeName();
91 dynCode.setFilterVariable(
"typeName",
name());
92 dynCode.setFilterVariable(
"TemplateType", primitiveTypeName);
93 dynCode.setFilterVariable(
"SourceType", primitiveTypeName +
"Source");
96 dynCode.addCompileFile(
"codedFvModelTemplate.C");
99 dynCode.addCopyFile(
"codedFvModelTemplate.H");
102 dynCode.setFilterVariable(
"verbose",
Foam::name(
bool(debug)));
105 dynCode.setMakeOptions
108 "-I$(LIB_SRC)/finiteVolume/lnInclude \\\n"
109 "-I$(LIB_SRC)/meshTools/lnInclude \\\n"
110 "-I$(LIB_SRC)/sampling/lnInclude \\\n"
111 "-I$(LIB_SRC)/fvModels/general/lnInclude \\\n"
113 +
"\n\nLIB_LIBS = \\\n"
114 +
" -lmeshTools \\\n"
117 +
" -lfiniteVolume \\\n"
123 Foam::fvModel& Foam::fv::codedFvModel::redirectFvModel()
const
125 if (!redirectFvModelPtr_.valid())
127 updateLibrary(coeffsDict_);
129 dictionary constructDict(coeffsDict_);
130 constructDict.set(
"type",
name());
139 return redirectFvModelPtr_();
144 void Foam::fv::codedFvModel::addSupType
146 const VolField<Type>& field,
154 Info<<
"codedFvModel::addSup for source " <<
name() <<
endl;
157 redirectFvModel().addSup(field, eqn);
163 void Foam::fv::codedFvModel::addSupType
166 const VolField<Type>& field,
174 Info<<
"codedFvModel::addSup for source " <<
name() <<
endl;
177 redirectFvModel().addSup(
rho, field, eqn);
183 void Foam::fv::codedFvModel::addSupType
187 const VolField<Type>& field,
195 Info<<
"codedFvModel::addSup for source " <<
name() <<
endl;
198 redirectFvModel().addSup(
alpha,
rho, field, eqn);
208 const word& modelType,
215 fieldName_(
word::null),
216 coeffsDict_(coeffs(
dict))
218 readCoeffs(coeffsDict_);
243 bool Foam::
fv::codedFvModel::movePoints()
245 return redirectFvModel().movePoints();
251 redirectFvModel().topoChange(map);
257 redirectFvModel().mapMesh(map);
263 redirectFvModel().distribute(map);
271 redirectFvModelPtr_.clear();
272 coeffsDict_ = coeffs(
dict);
273 readCoeffs(coeffsDict_);
275 updateLibrary(coeffsDict_);
Macros for easy insertion into run-time selection tables.
Base class for function objects and boundary conditions using dynamic code.
void read(const dictionary &dict)
Read the dictionary and update the code.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Mesh data needed to do the Finite Volume discretisation.
Finite volume model abstract base class.
static autoPtr< fvModel > New(const word &name, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected fvModel.
virtual bool read(const dictionary &dict)
Read source dictionary.
Constructs on-the-fly fvModel source from user-supplied code.
codedFvModel(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
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 handling words, derived from string.
static const word null
An empty word.
#define fieldPrimitiveTypeNameTernary(Type, nullArg)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
#define IMPLEMENT_FV_MODEL_ADD_FIELD_SUP(Type, modelType)
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_FIELD_SUP(Type, modelType)
#define IMPLEMENT_FV_MODEL_ADD_RHO_FIELD_SUP(Type, modelType)
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
List< word > wordList
A List of words.
Ostream & endl(Ostream &os)
Add newline and flush stream.
VolField< scalar > volScalarField
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)