45 void Foam::fv::massSource::readCoeffs()
47 phaseName_ = coeffs().lookupOrDefault<word>(
"phase",
word::null);
50 coeffs().lookupOrDefault<word>
58 mesh().foundObject<basicThermo>
64 const basicThermo& thermo =
65 mesh().lookupObject<basicThermo>
69 heName_ = thermo.he().name();
70 TName_ = thermo.T().name();
74 const dictionary& fieldCoeffs = coeffs().subDict(
"fieldValues");
80 new unknownTypeFunction1(iter().keyword(), fieldCoeffs)
89 void Foam::fv::massSource::addGeneralSupType
95 const scalar t =
mesh().time().userTimeValue();
96 const scalar massFlowRate = massFlowRate_->value(t);
97 const Type value = fieldValues_[fieldName]->value<Type>(t);
103 eqn.source()[cells[i]] -=
104 mesh().V()[cells[i]]/set_.V()*massFlowRate*value;
110 void Foam::fv::massSource::addSupType
113 const word& fieldName
116 addGeneralSupType(eqn, fieldName);
120 void Foam::fv::massSource::addSupType
122 fvMatrix<scalar>& eqn,
123 const word& fieldName
128 if (fieldName == rhoName_)
130 const scalar t =
mesh().time().userTimeValue();
131 const scalar massFlowRate = massFlowRate_->value(t);
135 eqn.source()[cells[i]] -=
136 mesh().V()[cells[i]]/set_.V()*massFlowRate;
139 else if (fieldName == heName_ && fieldValues_.found(TName_))
141 if (fieldValues_.found(heName_))
144 <<
"Source " <<
name() <<
" defined for both field " << heName_
145 <<
" and " << TName_ <<
". Only one of these should be present." 149 const scalar t =
mesh().time().userTimeValue();
150 const scalar massFlowRate = massFlowRate_->value(t);
151 const scalar T = fieldValues_[TName_]->value<scalar>(t);
152 const basicThermo& thermo =
153 mesh().lookupObject<basicThermo>
164 eqn.source()[cells[i]] -=
165 mesh().V()[cells[i]]/set_.V()*massFlowRate*hs[i];
170 addGeneralSupType(eqn, fieldName);
176 void Foam::fv::massSource::addSupType
180 const word& fieldName
183 addSupType(eqn, fieldName);
188 void Foam::fv::massSource::addSupType
193 const word& fieldName
196 addSupType(eqn, fieldName);
205 const word& modelType,
210 fvModel(name, modelType, dict, mesh),
211 set_(coeffs(), mesh),
232 && !(fieldName == rhoName_)
233 && !(fieldName == heName_ && fieldValues_.found(TName_))
234 && !fieldValues_.found(fieldName)
238 <<
"No value supplied for field " << fieldName <<
" in " 252 if (fieldValues_.found(TName_))
254 fieldNames[
findIndex(fieldNames, TName_)] = heName_;
279 set_.topoChange(map);
291 set_.distribute(map);
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
#define forAll(list, i)
Loop across all elements in list.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
A list of keyword definitions, which are a keyword followed by any number of values (e...
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
virtual bool read(const dictionary &dict)
Read source dictionary.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Finite volume model abstract base class.
word group() const
Return group (extension part of name)
Macros for easy insertion into run-time selection tables.
static List< word > fieldNames
GeometricField< scalar, fvPatchField, volMesh > volScalarField
A class for handling words, derived from string.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
massSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
static const word null
An empty word.
List< label > labelList
A List of labels.
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
#define WarningInFunction
Report a warning using Foam::Warning.
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Mesh data needed to do the Finite Volume discretisation.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
virtual bool movePoints()
Update for mesh motion.
Class containing mesh-to-mesh mapping information.
virtual bool read(const dictionary &dict)
Read source dictionary.
This fvModel applies a mass source to the continuity equation and to all field equations.
static autoPtr< Function1< scalar > > New(const word &name, const dictionary &dict)
Selector.