zeroDimensionalMassSource.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 
27 #include "basicThermo.H"
29 
30 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 Foam::fv::zeroDimensionalMassSourceBase::calcM0D() const
47 {
48  tmp<volScalarField> tm =
50  (
51  typedName("m0D"),
52  mesh(),
54  );
55 
56  HashTable<const basicThermo*> thermos(mesh().lookupClass<basicThermo>());
57 
58  forAllConstIter(HashTable<const basicThermo*>, thermos, thermoIter)
59  {
60  const basicThermo& thermo = *thermoIter();
61 
62  tmp<volScalarField> tRho = thermo.rho();
63  const volScalarField& rho = tRho();
64 
65  const word phaseName = thermo.phaseName();
66 
67  if (thermo.phaseName() != word::null)
68  {
69  const volScalarField& alpha =
71  (
72  IOobject::groupName("alpha", phaseName)
73  );
74 
75  tm.ref().ref() += alpha()*rho()*mesh().V();
76  }
77  else
78  {
79  tm.ref().ref() += rho()*mesh().V();
80  }
81  }
82 
83  return tm;
84 }
85 
86 
87 Foam::volScalarField& Foam::fv::zeroDimensionalMassSourceBase::initM0D() const
88 {
89  if (!mesh().foundObject<volScalarField>(typedName("m0D")))
90  {
91  volScalarField* mPtr =
92  new volScalarField
93  (
94  calcM0D()
95  );
96 
97  mPtr->store();
98  }
99 
100  return mesh().lookupObjectRef<volScalarField>(typedName("m0D"));
101 }
102 
103 
104 const Foam::volScalarField& Foam::fv::zeroDimensionalMassSourceBase::m() const
105 {
106  // If not registered, then read or create the mass field
107  if (!mesh().foundObject<volScalarField>(typedName("m")))
108  {
109  typeIOobject<volScalarField> mIo
110  (
111  typedName("m"),
112  mesh().time().timeName(),
113  mesh(),
116  );
117 
118  volScalarField* mPtr =
119  new volScalarField
120  (
121  mIo,
122  mesh(),
124  );
125 
126  mPtr->store();
127 
128  if (!mIo.headerOk())
129  {
130  *mPtr = m0D_;
131  }
132 
133  volScalarField* factorPtr =
134  new volScalarField
135  (
136  IOobject
137  (
138  typedName("factor"),
139  mesh().time().timeName(),
140  mesh(),
143  ),
144  *mPtr/m0D_
145  );
146 
147  factorPtr->store();
148  }
149 
150  volScalarField& m =
151  mesh().lookupObjectRef<volScalarField>(typedName("m"));
152 
153  volScalarField& factor =
154  mesh().lookupObjectRef<volScalarField>(typedName("factor"));
155 
156  // Update the mass if changes are available
157  if (mesh().foundObject<volScalarField>(typedName("deltaM")))
158  {
159  volScalarField& deltaM =
160  mesh().lookupObjectRef<volScalarField>(typedName("deltaM"));
161 
162  m = m.oldTime() + deltaM;
163 
164  factor = m/m0D_;
165 
166  deltaM.checkOut();
167  }
168 
169  return m;
170 }
171 
172 
173 Foam::scalar Foam::fv::zeroDimensionalMassSourceBase::massFlowRate() const
174 {
175  return zeroDimensionalMassFlowRate()*m0D_[0]/m()[0];
176 }
177 
178 
179 void Foam::fv::zeroDimensionalMassSource::readCoeffs()
180 {
181  readFieldValues();
182 
183  zeroDimensionalMassFlowRate_.reset
184  (
185  Function1<scalar>::New("massFlowRate", coeffs()).ptr()
186  );
187 }
188 
189 
190 Foam::scalar
191 Foam::fv::zeroDimensionalMassSource::zeroDimensionalMassFlowRate() const
192 {
193  return zeroDimensionalMassFlowRate_->value(mesh().time().userTimeValue());
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 
200 (
201  const word& name,
202  const word& modelType,
203  const fvMesh& mesh,
204  const dictionary& dict
205 )
206 :
207  massSourceBase(name, modelType, mesh, dict),
208  m0D_(initM0D())
209 {
210  if (mesh.nGeometricD() != 0)
211  {
213  << "Zero-dimensional fvModel applied to a "
214  << mesh.nGeometricD() << "-dimensional mesh"
215  << exit(FatalIOError);
216  }
217 }
218 
219 
221 (
222  const word& name,
223  const word& modelType,
224  const fvMesh& mesh,
225  const dictionary& dict
226 )
227 :
228  zeroDimensionalMassSourceBase(name, modelType, mesh, dict),
229  zeroDimensionalMassFlowRate_()
230 {
231  readCoeffs();
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
238 {
239  // Correct the zero-dimensional mass
240  m0D_ = calcM0D();
241 
242  // Create the mass change
243  if (!mesh().foundObject<volScalarField>(typedName("deltaM")))
244  {
245  volScalarField* dMPtr =
246  new volScalarField
247  (
248  IOobject
249  (
250  typedName("deltaM"),
251  mesh().time().timeName(),
252  mesh()
253  ),
254  mesh(),
256  );
257 
258  dMPtr->store();
259  }
260 
261  volScalarField& deltaM =
262  mesh().lookupObjectRef<volScalarField>(typedName("deltaM"));
263 
264  deltaM +=
265  mesh().time().deltaT()
266  *dimensionedScalar(dimMass/dimTime, zeroDimensionalMassFlowRate());
267 }
268 
269 
271 {
273  {
274  readCoeffs();
275  return true;
276  }
277  else
278  {
279  return false;
280  }
281 }
282 
283 
284 // ************************************************************************* //
#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
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:318
static word groupName(Name name, const word &group)
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massSource.C:372
zeroDimensionalMassSourceBase(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
This fvModel applies a mass source to the continuity equation and to all field equations,...
virtual bool read(const dictionary &dict)
Read source dictionary.
zeroDimensionalMassSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1055
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
word timeName
Definition: getTimeIndex.H:3
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
IOerror FatalIOError
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31