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 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  (
61  )
62  )
63  {
64  const basicThermo& thermo =
65  mesh().lookupObject<basicThermo>
66  (
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(iter().keyword(), nullptr);
78  }
79 
80  #define callReadFieldValues(Type, nullArg) readFieldValues<Type>();
82  #undef callReadFieldValues
83 
84  massFlowRate_.reset(Function1<scalar>::New("massFlowRate", coeffs()).ptr());
85 }
86 
87 
88 template<class Type>
89 void Foam::fv::massSource::readFieldValues()
90 {
91  const dictionary& fieldCoeffs = coeffs().subDict("fieldValues");
92 
93  forAllConstIter(dictionary, fieldCoeffs, iter)
94  {
95  const word& fieldName = iter().keyword();
96 
97  if (mesh().foundObject<VolField<Type>>(fieldName))
98  {
99  fieldValues_.set
100  (
101  fieldName,
102  objectFunction1::New<VolField>
103  (
104  fieldName,
105  fieldCoeffs,
106  fieldName,
107  mesh()
108  ).ptr()
109  );
110  }
111  }
112 }
113 
114 
115 template<class Type>
116 void Foam::fv::massSource::addGeneralSupType
117 (
118  fvMatrix<Type>& eqn,
119  const word& fieldName
120 ) const
121 {
122  const scalar t = mesh().time().value();
123  const scalar massFlowRate = massFlowRate_->value(t);
124  const Type value = fieldValues_[fieldName]->value<Type>(t);
125 
126  const labelList& cells = set_.cells();
127 
128  forAll(cells, i)
129  {
130  eqn.source()[cells[i]] -=
131  mesh().V()[cells[i]]/set_.V()*massFlowRate*value;
132  }
133 }
134 
135 
136 template<class Type>
137 void Foam::fv::massSource::addSupType
138 (
139  fvMatrix<Type>& eqn,
140  const word& fieldName
141 ) const
142 {
143  addGeneralSupType(eqn, fieldName);
144 }
145 
146 
147 void Foam::fv::massSource::addSupType
148 (
149  fvMatrix<scalar>& eqn,
150  const word& fieldName
151 ) const
152 {
153  const labelList& cells = set_.cells();
154 
155  if (fieldName == rhoName_)
156  {
157  const scalar t = mesh().time().value();
158  const scalar massFlowRate = massFlowRate_->value(t);
159 
160  forAll(cells, i)
161  {
162  eqn.source()[cells[i]] -=
163  mesh().V()[cells[i]]/set_.V()*massFlowRate;
164  }
165  }
166  else if (fieldName == heName_ && fieldValues_.found(TName_))
167  {
168  if (fieldValues_.found(heName_))
169  {
171  << "Source " << name() << " defined for both field " << heName_
172  << " and " << TName_ << ". Only one of these should be present."
173  << endl;
174  }
175 
176  const scalar t = mesh().time().value();
177  const scalar massFlowRate = massFlowRate_->value(t);
178  const scalar T = fieldValues_[TName_]->value<scalar>(t);
179  const basicThermo& thermo =
180  mesh().lookupObject<basicThermo>
181  (
183  );
184  const scalarField hs
185  (
186  thermo.hs(scalarField(cells.size(), T), cells)
187  );
188 
189  forAll(cells, i)
190  {
191  eqn.source()[cells[i]] -=
192  mesh().V()[cells[i]]/set_.V()*massFlowRate*hs[i];
193  }
194  }
195  else
196  {
197  addGeneralSupType(eqn, fieldName);
198  }
199 }
200 
201 
202 template<class Type>
203 void Foam::fv::massSource::addSupType
204 (
205  const volScalarField& rho,
206  fvMatrix<Type>& eqn,
207  const word& fieldName
208 ) const
209 {
210  addSupType(eqn, fieldName);
211 }
212 
213 
214 template<class Type>
215 void Foam::fv::massSource::addSupType
216 (
217  const volScalarField& alpha,
218  const volScalarField& rho,
219  fvMatrix<Type>& eqn,
220  const word& fieldName
221 ) const
222 {
223  addSupType(eqn, fieldName);
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
228 
230 (
231  const word& name,
232  const word& modelType,
233  const dictionary& dict,
234  const fvMesh& mesh
235 )
236 :
237  fvModel(name, modelType, dict, mesh),
238  set_(coeffs(), mesh),
239  phaseName_(),
240  rhoName_(),
241  heName_(),
242  TName_(),
243  massFlowRate_(),
244  fieldValues_()
245 {
246  readCoeffs();
247 }
248 
249 
250 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
251 
252 bool Foam::fv::massSource::addsSupToField(const word& fieldName) const
253 {
254  const bool isThisPhase = IOobject::group(fieldName) == phaseName_;
255 
256  if
257  (
258  isThisPhase
259  && !(fieldName == rhoName_)
260  && !(fieldName == heName_ && fieldValues_.found(TName_))
261  && !fieldValues_.found(fieldName)
262  )
263  {
265  << "No value supplied for field " << fieldName << " in "
266  << type() << " fvModel " << name() << endl;
267 
268  return false;
269  }
270 
271  return isThisPhase;
272 }
273 
274 
276 {
277  wordList fieldNames = fieldValues_.toc();
278 
279  if (fieldValues_.found(TName_))
280  {
281  fieldNames[findIndex(fieldNames, TName_)] = heName_;
282  }
283 
284  return fieldNames;
285 }
286 
287 
289 
290 
292 
293 
295 
296 
298 {
299  set_.updateMesh(mpm);
300 }
301 
302 
304 {
305  if (fvModel::read(dict))
306  {
307  set_.read(coeffs());
308  readCoeffs();
309  return true;
310  }
311  else
312  {
313  return false;
314  }
315 }
316 
317 
318 // ************************************************************************* //
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:158
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:55
word group() const
Return group (extension part of name)
Definition: IOobject.C:330
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
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:57
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
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:230
static const word null
An empty word.
Definition: word.H:77
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: massSource.C:275
#define callReadFieldValues(Type, nullArg)
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)
FOR_ALL_FIELD_TYPES(IMPLEMENT_FV_MODEL_ADD_SUP, fv::massSource)
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: massSource.C:252
#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:78
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.
virtual void updateMesh(const mapPolyMesh &)
Update for mesh changes.
Definition: massSource.C:297
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massSource.C:303
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