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-2023 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 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::fv::massSourceBase::readCoeffs()
47 {
48  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
49 
50  rhoName_ =
51  coeffs().lookupOrDefault<word>
52  (
53  "rho",
54  IOobject::groupName("rho", phaseName_)
55  );
56 
57  if
58  (
59  mesh().foundObject<basicThermo>
60  (
61  IOobject::groupName(physicalProperties::typeName, phaseName_)
62  )
63  )
64  {
65  const basicThermo& thermo =
66  mesh().lookupObject<basicThermo>
67  (
68  IOobject::groupName(physicalProperties::typeName, phaseName_)
69  );
70  heName_ = thermo.he().name();
71  TName_ = thermo.T().name();
72  }
73 }
74 
75 
76 template<class Type>
77 void Foam::fv::massSourceBase::addGeneralSupType
78 (
79  fvMatrix<Type>& eqn,
80  const word& fieldName
81 ) const
82 {
83  const labelUList cells = set_.cells();
84 
85  const scalar massFlowRate = this->massFlowRate();
86 
87  if (massFlowRate > 0)
88  {
89  const Type value =
90  fieldValues_[fieldName]->value<Type>(mesh().time().userTimeValue());
91 
92  forAll(cells, i)
93  {
94  eqn.source()[cells[i]] -=
95  mesh().V()[cells[i]]/set_.V()*massFlowRate*value;
96  }
97  }
98  else
99  {
100  forAll(cells, i)
101  {
102  eqn.diag()[cells[i]] +=
103  mesh().V()[cells[i]]/set_.V()*massFlowRate;
104  }
105  }
106 }
107 
108 
109 template<class Type>
110 void Foam::fv::massSourceBase::addSupType
111 (
112  fvMatrix<Type>& eqn,
113  const word& fieldName
114 ) const
115 {
116  addGeneralSupType(eqn, fieldName);
117 }
118 
119 
120 void Foam::fv::massSourceBase::addSupType
121 (
122  fvMatrix<scalar>& eqn,
123  const word& fieldName
124 ) const
125 {
126  const labelUList cells = set_.cells();
127 
128  if (fieldName == rhoName_)
129  {
130  const scalar massFlowRate = this->massFlowRate();
131 
132  forAll(cells, i)
133  {
134  eqn.source()[cells[i]] -=
135  mesh().V()[cells[i]]/set_.V()*massFlowRate;
136  }
137  }
138  else if (fieldName == heName_ && fieldValues_.found(TName_))
139  {
140  const scalar massFlowRate = this->massFlowRate();
141 
142  if (massFlowRate > 0)
143  {
144  if (fieldValues_.found(heName_))
145  {
147  << "Source " << name() << " defined for both field "
148  << heName_ << " and " << TName_
149  << ". Only one of these should be present." << endl;
150  }
151 
152  const basicThermo& thermo =
153  mesh().lookupObject<basicThermo>
154  (
156  (
157  physicalProperties::typeName,
158  phaseName_
159  )
160  );
161 
162  const scalar T =
163  fieldValues_[TName_]->value<scalar>
164  (
165  mesh().time().userTimeValue()
166  );
167 
168  const scalarField hs
169  (
170  thermo.hs(scalarField(cells.size(), T), cells)
171  );
172 
173  forAll(cells, i)
174  {
175  eqn.source()[cells[i]] -=
176  mesh().V()[cells[i]]/set_.V()*massFlowRate*hs[i];
177  }
178  }
179  else
180  {
181  forAll(cells, i)
182  {
183  eqn.diag()[cells[i]] +=
184  mesh().V()[cells[i]]/set_.V()*massFlowRate;
185  }
186  }
187  }
188  else
189  {
190  addGeneralSupType(eqn, fieldName);
191  }
192 }
193 
194 
195 template<class Type>
196 void Foam::fv::massSourceBase::addSupType
197 (
198  const volScalarField& rho,
199  fvMatrix<Type>& eqn,
200  const word& fieldName
201 ) const
202 {
203  addSupType(eqn, fieldName);
204 }
205 
206 
207 template<class Type>
208 void Foam::fv::massSourceBase::addSupType
209 (
210  const volScalarField& alpha,
211  const volScalarField& rho,
212  fvMatrix<Type>& eqn,
213  const word& fieldName
214 ) const
215 {
216  addSupType(eqn, fieldName);
217 }
218 
219 
220 void Foam::fv::massSource::readCoeffs()
221 {
222  readSet();
223 
224  readFieldValues();
225 
226  massFlowRate_.reset
227  (
228  Function1<scalar>::New("massFlowRate", coeffs()).ptr()
229  );
230 }
231 
232 
233 Foam::scalar Foam::fv::massSource::massFlowRate() const
234 {
235  return massFlowRate_->value(mesh().time().userTimeValue());
236 }
237 
238 
239 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
240 
242 {
243  set_.read(coeffs());
244 }
245 
246 
248 {
249  fieldValues_.clear();
250  const dictionary& fieldCoeffs = coeffs().subDict("fieldValues");
251  forAllConstIter(dictionary, fieldCoeffs, iter)
252  {
253  fieldValues_.set
254  (
255  iter().keyword(),
256  new unknownTypeFunction1(iter().keyword(), fieldCoeffs)
257  );
258  }
259 }
260 
261 
262 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
263 
265 (
266  const word& name,
267  const word& modelType,
268  const fvMesh& mesh,
269  const dictionary& dict
270 )
271 :
272  fvModel(name, modelType, mesh, dict),
273  phaseName_(),
274  rhoName_(),
275  heName_(),
276  TName_(),
277  set_(fvCellSet(mesh)),
278  fieldValues_()
279 {
280  readCoeffs();
281 }
282 
283 
285 (
286  const word& name,
287  const word& modelType,
288  const fvMesh& mesh,
289  const dictionary& dict
290 )
291 :
292  massSourceBase(name, modelType, mesh, dict),
293  massFlowRate_()
294 {
295  readCoeffs();
296 }
297 
298 
299 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
300 
301 bool Foam::fv::massSourceBase::addsSupToField(const word& fieldName) const
302 {
303  const bool isThisPhase = IOobject::group(fieldName) == phaseName_;
304 
305  if
306  (
307  isThisPhase
308  && massFlowRate() > 0
309  && !(fieldName == rhoName_)
310  && !(fieldName == heName_ && fieldValues_.found(TName_))
311  && !fieldValues_.found(fieldName)
312  )
313  {
315  << "No value supplied for field " << fieldName << " in "
316  << type() << " fvModel " << name() << endl;
317 
318  return false;
319  }
320 
321  return isThisPhase;
322 }
323 
324 
326 {
327  wordList fieldNames = fieldValues_.toc();
328 
329  if (fieldValues_.found(TName_))
330  {
331  fieldNames[findIndex(fieldNames, TName_)] = heName_;
332  }
333 
334  return fieldNames;
335 }
336 
337 
339 
340 
342 
343 
345 
346 
348 {
349  set_.movePoints();
350  return true;
351 }
352 
353 
355 {
356  set_.topoChange(map);
357 }
358 
359 
361 {
362  set_.mapMesh(map);
363 }
364 
365 
367 {
368  set_.distribute(map);
369 }
370 
371 
373 {
374  if (fvModel::read(dict))
375  {
376  readCoeffs();
377  return true;
378  }
379  else
380  {
381  return false;
382  }
383 }
384 
385 
387 {
389  {
390  readCoeffs();
391  return true;
392  }
393  else
394  {
395  return false;
396  }
397 }
398 
399 
400 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
static word groupName(Name name, const word &group)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
General run-time selected cell set selection class for fvMesh.
Definition: fvCellSet.H:88
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:187
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
virtual bool movePoints()
Update for mesh motion.
Definition: massSource.C:347
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: massSource.C:325
void readSet()
Read the set.
Definition: massSource.C:241
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: massSource.C:354
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: massSource.C:366
void readFieldValues()
Read the field values.
Definition: massSource.C:247
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massSource.C:372
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: massSource.C:360
massSourceBase(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: massSource.C:265
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: massSource.C:301
This fvModel applies a mass source to the continuity equation and to all field equations.
Definition: massSource.H:250
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massSource.C:386
massSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: massSource.C:285
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
bool foundObject(const word &name) const
Is the named Type in registry.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Wrapper around Function1 that constructs a function for an as yet unknown primitive type....
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations.
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Definition: fvModelM.H:51
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
static List< word > fieldNames
Definition: globalFoam.H:46
const cellShapeList & cells
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
UList< label > labelUList
Definition: UList.H:65
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31