zeroDimensionalMassSourceBase.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-2024 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 "fvCellSet.H"
28 #include "basicThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 Foam::fv::zeroDimensionalMassSourceBase::calcM0D() const
46 {
47  tmp<volScalarField> tm =
49  (
50  typedName("m0D"),
51  mesh(),
53  );
54 
55  HashTable<const basicThermo*> thermos(mesh().lookupClass<basicThermo>());
56 
57  forAllConstIter(HashTable<const basicThermo*>, thermos, thermoIter)
58  {
59  const basicThermo& thermo = *thermoIter();
60 
61  tmp<volScalarField> tRho = thermo.rho();
62  const volScalarField& rho = tRho();
63 
64  const word phaseName = thermo.phaseName();
65 
66  if (thermo.phaseName() != word::null)
67  {
68  const volScalarField& alpha =
70  (
72  );
73 
74  tm.ref().internalFieldRef() += alpha()*rho()*mesh().V();
75  }
76  else
77  {
78  tm.ref().internalFieldRef() += rho()*mesh().V();
79  }
80  }
81 
82  return tm;
83 }
84 
85 
86 Foam::volScalarField& Foam::fv::zeroDimensionalMassSourceBase::initM0D() const
87 {
88  if (!mesh().foundObject<volScalarField>(typedName("m0D")))
89  {
90  volScalarField* mPtr =
91  new volScalarField
92  (
93  calcM0D()
94  );
95 
96  mPtr->store();
97  }
98 
99  return mesh().lookupObjectRef<volScalarField>(typedName("m0D"));
100 }
101 
102 
103 const Foam::volScalarField& Foam::fv::zeroDimensionalMassSourceBase::m() const
104 {
105  // If not registered, then read or create the mass field
106  if (!mesh().foundObject<volScalarField>(typedName("m")))
107  {
108  typeIOobject<volScalarField> mIo
109  (
110  typedName("m"),
111  mesh().time().name(),
112  mesh(),
115  );
116 
117  volScalarField* mPtr =
118  new volScalarField
119  (
120  mIo,
121  mesh(),
123  );
124 
125  mPtr->store();
126 
127  if (!mIo.headerOk())
128  {
129  *mPtr = m0D_;
130  }
131 
132  volScalarField* factorPtr =
133  new volScalarField
134  (
135  IOobject
136  (
137  typedName("factor"),
138  mesh().time().name(),
139  mesh(),
142  ),
143  *mPtr/m0D_
144  );
145 
146  factorPtr->store();
147  }
148 
149  volScalarField& m =
150  mesh().lookupObjectRef<volScalarField>(typedName("m"));
151 
152  volScalarField& factor =
153  mesh().lookupObjectRef<volScalarField>(typedName("factor"));
154 
155  // Update the mass if changes are available
156  if (mesh().foundObject<volScalarField>(typedName("deltaM")))
157  {
158  volScalarField& deltaM =
159  mesh().lookupObjectRef<volScalarField>(typedName("deltaM"));
160 
161  m = m.oldTime() + deltaM;
162 
163  factor = m/m0D_;
164 
165  deltaM.checkOut();
166  }
167 
168  return m;
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173 
175 (
176  const word& name,
177  const word& modelType,
178  const fvMesh& mesh,
179  const dictionary& dict
180 )
181 :
182  massSourceBase(name, modelType, mesh, dict),
183  m0D_(initM0D())
184 {
185  if (mesh.nGeometricD() != 0)
186  {
188  << "Zero-dimensional fvModel applied to a "
189  << mesh.nGeometricD() << "-dimensional mesh"
190  << exit(FatalIOError);
191  }
192 }
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
198 {
199  static labelList zero(1, Zero);
200  return labelUList(zero);
201 }
202 
203 
205 {
206  return 1;
207 }
208 
209 
211 {
212  return mesh().V()[0];
213 }
214 
215 
217 {
218  return
220  (
222  massFlowRate()*m0D_[0]/m()[0]
223  );
224 }
225 
226 
228 {
229  return true;
230 }
231 
232 
234 (
235  const polyTopoChangeMap& map
236 )
237 {}
238 
239 
241 (
242  const polyMeshMap& map
243 )
244 {}
245 
246 
248 (
249  const polyDistributionMap& map
250 )
251 {}
252 
253 
255 {
256  // Correct the zero-dimensional mass
257  m0D_ = calcM0D();
258 
259  // Create the mass change
260  if (!mesh().foundObject<volScalarField>(typedName("deltaM")))
261  {
262  volScalarField* dMPtr =
263  new volScalarField
264  (
265  IOobject
266  (
267  typedName("deltaM"),
268  mesh().time().name(),
269  mesh()
270  ),
271  mesh(),
273  );
274 
275  dMPtr->store();
276  }
277 
278  volScalarField& deltaM =
279  mesh().lookupObjectRef<volScalarField>(typedName("deltaM"));
280 
281  deltaM +=
282  mesh().time().deltaT()
283  *dimensionedScalar(dimMass/dimTime, massFlowRate());
284 }
285 
286 
287 // ************************************************************************* //
#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.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
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:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
const word & phaseName() const
Return the phase name.
Base class for mass source models.
Base class for zero-dimensional mass source models.
virtual bool movePoints()
Update for mesh motion.
virtual scalar V() const
Return the volume of cells that the source applies to.
zeroDimensionalMassSourceBase(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual labelUList cells() const
Return the cells that the source applies to.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual label nCells() const
Return the number of cells that the source applies to.
virtual dimensionedScalar S() const
Return the source value.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1000
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
IOerror FatalIOError
const dimensionSet dimMass
UList< label > labelUList
Definition: UList.H:65
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31