massSource.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "massSource.H"
27 #include "fvMatrices.H"
28 #include "basicThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37  defineTypeNameAndDebug(massSource, 0);
38  addToRunTimeSelectionTable(fvModel, massSource, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::fv::massSource::readCoeffs()
46 {
47  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
48 
49  rhoName_ =
50  coeffs().lookupOrDefault<word>
51  (
52  "rho",
53  IOobject::groupName("rho", phaseName_)
54  );
55 
56  if
57  (
58  mesh().foundObject<basicThermo>
59  (
60  IOobject::groupName(physicalProperties::typeName, phaseName_)
61  )
62  )
63  {
64  const basicThermo& thermo =
65  mesh().lookupObject<basicThermo>
66  (
67  IOobject::groupName(physicalProperties::typeName, phaseName_)
68  );
69  heName_ = thermo.he().name();
70  TName_ = thermo.T().name();
71  }
72 
73  fieldValues_.clear();
74  const dictionary& fieldCoeffs = coeffs().subDict("fieldValues");
75  forAllConstIter(dictionary, fieldCoeffs, iter)
76  {
77  fieldValues_.set
78  (
79  iter().keyword(),
80  new unknownTypeFunction1(iter().keyword(), fieldCoeffs)
81  );
82  }
83 
84  massFlowRate_.reset(Function1<scalar>::New("massFlowRate", coeffs()).ptr());
85 }
86 
87 
88 template<class Type>
89 void Foam::fv::massSource::addGeneralSupType
90 (
91  fvMatrix<Type>& eqn,
92  const word& fieldName
93 ) const
94 {
95  const scalar t = mesh().time().userTimeValue();
96  const scalar massFlowRate = massFlowRate_->value(t);
97  const Type value = fieldValues_[fieldName]->value<Type>(t);
98 
99  const labelList& cells = set_.cells();
100 
101  forAll(cells, i)
102  {
103  eqn.source()[cells[i]] -=
104  mesh().V()[cells[i]]/set_.V()*massFlowRate*value;
105  }
106 }
107 
108 
109 template<class Type>
110 void Foam::fv::massSource::addSupType
111 (
112  fvMatrix<Type>& eqn,
113  const word& fieldName
114 ) const
115 {
116  addGeneralSupType(eqn, fieldName);
117 }
118 
119 
120 void Foam::fv::massSource::addSupType
121 (
122  fvMatrix<scalar>& eqn,
123  const word& fieldName
124 ) const
125 {
126  const labelList& cells = set_.cells();
127 
128  if (fieldName == rhoName_)
129  {
130  const scalar t = mesh().time().userTimeValue();
131  const scalar massFlowRate = massFlowRate_->value(t);
132 
133  forAll(cells, i)
134  {
135  eqn.source()[cells[i]] -=
136  mesh().V()[cells[i]]/set_.V()*massFlowRate;
137  }
138  }
139  else if (fieldName == heName_ && fieldValues_.found(TName_))
140  {
141  if (fieldValues_.found(heName_))
142  {
144  << "Source " << name() << " defined for both field " << heName_
145  << " and " << TName_ << ". Only one of these should be present."
146  << endl;
147  }
148 
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>
154  (
155  IOobject::groupName(physicalProperties::typeName, phaseName_)
156  );
157  const scalarField hs
158  (
159  thermo.hs(scalarField(cells.size(), T), cells)
160  );
161 
162  forAll(cells, i)
163  {
164  eqn.source()[cells[i]] -=
165  mesh().V()[cells[i]]/set_.V()*massFlowRate*hs[i];
166  }
167  }
168  else
169  {
170  addGeneralSupType(eqn, fieldName);
171  }
172 }
173 
174 
175 template<class Type>
176 void Foam::fv::massSource::addSupType
177 (
178  const volScalarField& rho,
179  fvMatrix<Type>& eqn,
180  const word& fieldName
181 ) const
182 {
183  addSupType(eqn, fieldName);
184 }
185 
186 
187 template<class Type>
188 void Foam::fv::massSource::addSupType
189 (
190  const volScalarField& alpha,
191  const volScalarField& rho,
192  fvMatrix<Type>& eqn,
193  const word& fieldName
194 ) const
195 {
196  addSupType(eqn, fieldName);
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
201 
203 (
204  const word& name,
205  const word& modelType,
206  const dictionary& dict,
207  const fvMesh& mesh
208 )
209 :
210  fvModel(name, modelType, dict, mesh),
211  set_(coeffs(), mesh),
212  phaseName_(),
213  rhoName_(),
214  heName_(),
215  TName_(),
216  massFlowRate_(),
217  fieldValues_()
218 {
219  readCoeffs();
220 }
221 
222 
223 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 
225 bool Foam::fv::massSource::addsSupToField(const word& fieldName) const
226 {
227  const bool isThisPhase = IOobject::group(fieldName) == phaseName_;
228 
229  if
230  (
231  isThisPhase
232  && !(fieldName == rhoName_)
233  && !(fieldName == heName_ && fieldValues_.found(TName_))
234  && !fieldValues_.found(fieldName)
235  )
236  {
238  << "No value supplied for field " << fieldName << " in "
239  << type() << " fvModel " << name() << endl;
240 
241  return false;
242  }
243 
244  return isThisPhase;
245 }
246 
247 
249 {
250  wordList fieldNames = fieldValues_.toc();
251 
252  if (fieldValues_.found(TName_))
253  {
254  fieldNames[findIndex(fieldNames, TName_)] = heName_;
255  }
256 
257  return fieldNames;
258 }
259 
260 
262 
263 
265 
266 
268 
269 
271 {
272  set_.movePoints();
273  return true;
274 }
275 
276 
278 {
279  set_.topoChange(map);
280 }
281 
282 
284 {
285  set_.mapMesh(map);
286 }
287 
288 
290 {
291  set_.distribute(map);
292 }
293 
294 
296 {
297  if (fvModel::read(dict))
298  {
299  set_.read(coeffs());
300  readCoeffs();
301  return true;
302  }
303  else
304  {
305  return false;
306  }
307 }
308 
309 
310 // ************************************************************************* //
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: massSource.C:277
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: massSource.C:283
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
Definition: dictionary.H:156
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
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.
Definition: Ostream.H:251
Finite volume model abstract base class.
Definition: fvModel.H:57
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
static List< word > fieldNames
Definition: globalFoam.H:46
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
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.
Definition: massSource.C:203
static const word null
An empty word.
Definition: word.H:77
List< label > labelList
A List of labels.
Definition: labelList.H:56
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: massSource.C:248
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.
Definition: complex.C:47
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: massSource.C:225
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: massSource.C:289
#define WarningInFunction
Report a warning using Foam::Warning.
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Definition: fvModelM.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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.
Definition: massSource.C:270
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massSource.C:295
This fvModel applies a mass source to the continuity equation and to all field equations.
Definition: massSource.H:81
Namespace for OpenFOAM.
static autoPtr< Function1< scalar > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32