clouds.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 "clouds.H"
27 #include "multicomponentThermo.H"
28 #include "fvMatrix.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace fv
36  {
38 
40  (
41  fvModel,
42  clouds,
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& sourceName,
54  const word& modelType,
55  const fvMesh& mesh,
56  const dictionary& dict
57 )
58 :
59  fvModel(sourceName, modelType, mesh, dict),
60  g_
61  (
62  IOobject
63  (
64  "g",
65  mesh.time().constant(),
66  mesh,
67  IOobject::READ_IF_PRESENT,
68  IOobject::NO_WRITE
69  ),
71  ),
72  carrierHasThermo_
73  (
74  mesh.foundObject<fluidThermo>(physicalProperties::typeName)
75  ),
76  tCarrierThermo_
77  (
78  carrierHasThermo_
80  (
81  mesh.lookupObject<fluidThermo>(physicalProperties::typeName)
82  )
83  : tmpNrc<fluidThermo>(nullptr)
84  ),
85  tCarrierViscosity_
86  (
87  carrierHasThermo_
88  ? tmpNrc<viscosityModel>(nullptr)
90  (
91  mesh.lookupObject<viscosityModel>(physicalProperties::typeName)
92  )
93  ),
94  tRho_
95  (
96  carrierHasThermo_
97  ? tmp<volScalarField>(nullptr)
99  (
100  new volScalarField
101  (
102  IOobject
103  (
104  "rho",
105  mesh.time().name(),
106  mesh
107  ),
108  mesh,
109  dimensionedScalar("rho", dimDensity, tCarrierViscosity_())
110  )
111  )
112  ),
113  tMu_
114  (
115  carrierHasThermo_
116  ? tmp<volScalarField>(nullptr)
118  (
119  new volScalarField("mu", tRho_()*tCarrierViscosity_().nu())
120  )
121  ),
122  cloudNames_
123  (
124  dict.lookupOrDefault<wordList>
125  (
126  "clouds",
127  parcelCloudList::defaultCloudNames
128  )
129  ),
130  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
131  UName_(dict.lookupOrDefault<word>("U", "U")),
132  cloudsPtr_
133  (
134  carrierHasThermo_
135  ? new parcelCloudList
136  (
137  cloudNames_,
138  mesh.lookupObject<volScalarField>(rhoName_),
139  mesh.lookupObject<volVectorField>(UName_),
140  g_,
141  tCarrierThermo_()
142  )
143  : new parcelCloudList
144  (
145  cloudNames_,
146  tRho_(),
147  mesh.lookupObject<volVectorField>(UName_),
148  tMu_(),
149  g_
150  )
151  ),
152  curTimeIndex_(-1)
153 {}
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
159 {
160  wordList fieldNames(1, UName_);
161 
162  if (carrierHasThermo_)
163  {
164  const fluidThermo& carrierThermo = tCarrierThermo_();
165 
166  fieldNames.append(rhoName_);
167 
168  fieldNames.append(carrierThermo.he().name());
169 
170  if (isA<multicomponentThermo>(carrierThermo))
171  {
172  const multicomponentThermo& carrierMcThermo =
173  refCast<const multicomponentThermo>(carrierThermo);
174 
175  forAll(carrierMcThermo.Y(), i)
176  {
177  if (carrierMcThermo.solveSpecie(i))
178  {
179  fieldNames.append(carrierMcThermo.Y()[i].name());
180  }
181  }
182  }
183  }
184 
185  return fieldNames;
186 }
187 
188 
190 {
191  if (curTimeIndex_ == mesh().time().timeIndex())
192  {
193  return;
194  }
195 
196  if (!carrierHasThermo_)
197  {
198  tMu_.ref() = tRho_()*tCarrierViscosity_().nu();
199  }
200 
201  cloudsPtr_().evolve();
202 
203  curTimeIndex_ = mesh().time().timeIndex();
204 }
205 
206 
208 (
209  const volScalarField& rho,
210  fvMatrix<scalar>& eqn
211 ) const
212 {
213  if (debug)
214  {
215  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
216  }
217 
218  if (!carrierHasThermo_)
219  {
221  << "Applying source to compressible equation when carrier thermo "
222  << "is not available"
223  << exit(FatalError);
224  }
225 
226  if (rho.name() == rhoName_)
227  {
228  eqn += cloudsPtr_().Srho(eqn.psi());
229  }
230  else
231  {
233  << "Support for field " << rho.name() << " is not implemented"
234  << exit(FatalError);
235  }
236 }
237 
238 
240 (
241  const volScalarField& rho,
242  const volScalarField& field,
243  fvMatrix<scalar>& eqn
244 ) const
245 {
246  if (debug)
247  {
248  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
249  }
250 
251  if (!carrierHasThermo_)
252  {
254  << "Applying source to compressible equation when carrier thermo "
255  << "is not available"
256  << exit(FatalError);
257  }
258 
259  const fluidThermo& carrierThermo = tCarrierThermo_();
260 
261  if (&field == &carrierThermo.he())
262  {
263  eqn += cloudsPtr_().Sh(eqn.psi());
264  return;
265  }
266 
267  if (isA<multicomponentThermo>(carrierThermo))
268  {
269  const multicomponentThermo& carrierMcThermo =
270  refCast<const multicomponentThermo>(carrierThermo);
271 
272  if (carrierMcThermo.containsSpecie(field.name()))
273  {
274  eqn +=
275  cloudsPtr_().SYi
276  (
277  carrierMcThermo.specieIndex(field),
278  field
279  );
280  return;
281  }
282  }
283 
284  {
286  << "Support for field " << field.name() << " is not implemented"
287  << exit(FatalError);
288  }
289 }
290 
291 
293 (
294  const volVectorField& U,
295  fvMatrix<vector>& eqn
296 ) const
297 {
298  if (debug)
299  {
300  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
301  }
302 
303  if (carrierHasThermo_)
304  {
306  << "Applying source to incompressible equation when carrier thermo "
307  << "is available"
308  << exit(FatalError);
309  }
310 
311  if (U.name() == UName_)
312  {
313  eqn += cloudsPtr_().SU(eqn.psi())/tRho_();
314  }
315  else
316  {
318  << "Support for field " << U.name() << " is not implemented"
319  << exit(FatalError);
320  }
321 }
322 
323 
325 (
326  const volScalarField& rho,
327  const volVectorField& U,
328  fvMatrix<vector>& eqn
329 ) const
330 {
331  if (debug)
332  {
333  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
334  }
335 
336  if (!carrierHasThermo_)
337  {
339  << "Applying source to compressible equation when carrier thermo "
340  << "is not available"
341  << exit(FatalError);
342  }
343 
344  if (U.name() == UName_)
345  {
346  eqn += cloudsPtr_().SU(eqn.psi());
347  }
348  else
349  {
351  << "Support for field " << U.name() << " is not implemented"
352  << exit(FatalError);
353  }
354 }
355 
356 
358 {
359  // Store the particle positions
360  cloudsPtr_().storeGlobalPositions();
361 }
362 
363 
365 {
366  cloudsPtr_().topoChange(map);
367 }
368 
369 
371 {
372  cloudsPtr_().mapMesh(map);
373 }
374 
375 
377 {
378  cloudsPtr_().distribute(map);
379 }
380 
381 
383 {
384  return true;
385 }
386 
387 
388 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
This fvModel adds any number of Lagrangian clouds to any single-phase solver. The particles are track...
Definition: clouds.H:115
virtual bool movePoints()
Update for mesh motion.
Definition: clouds.C:382
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: clouds.C:158
virtual void correct()
Solve the Lagrangian clouds and update the sources.
Definition: clouds.C:189
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn) const
Add source to continuity equation.
Definition: clouds.C:208
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: clouds.C:364
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: clouds.C:376
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: clouds.C:357
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: clouds.C:370
clouds(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: clouds.C:52
Base-class for multi-component thermodynamic properties.
bool solveSpecie(const label speciei) const
Should the given specie be solved for? I.e., is it active and.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
bool containsSpecie(const word &specieName) const
Does the mixture include this specie?
label specieIndex(const volScalarField &Yi) const
Access the specie index of the given mass-fraction field.
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
A base class for physical properties.
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.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:53
A class for managing temporary objects.
Definition: tmp.H:55
An abstract base class for Newtonian viscosity models.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
static List< word > fieldNames
Definition: globalFoam.H:46
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
const dimensionSet dimAcceleration
const dimensionSet dimDensity
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label timeIndex
Definition: getTimeIndex.H:4
labelList fv(nPoints)
dictionary dict