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 "basicSpecieMixture.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<basicSpecieMixture>(carrierThermo))
171  {
173  refCast<const basicSpecieMixture>(carrierThermo);
174 
176 
177  forAll(Y, i)
178  {
179  if (composition.solve(i))
180  {
181  fieldNames.append(Y[i].name());
182  }
183  }
184  }
185  }
186 
187  return fieldNames;
188 }
189 
190 
192 {
193  if (curTimeIndex_ == mesh().time().timeIndex())
194  {
195  return;
196  }
197 
198  if (!carrierHasThermo_)
199  {
200  tMu_.ref() = tRho_()*tCarrierViscosity_().nu();
201  }
202 
203  cloudsPtr_().evolve();
204 
205  curTimeIndex_ = mesh().time().timeIndex();
206 }
207 
208 
210 (
211  fvMatrix<scalar>& eqn,
212  const word& fieldName
213 ) const
214 {
215  if (debug)
216  {
217  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
218  }
219 
220  if (!carrierHasThermo_)
221  {
223  << "Applying source to compressible equation when carrier thermo "
224  << "is not available"
225  << exit(FatalError);
226  }
227 
228  if (fieldName == rhoName_)
229  {
230  eqn += cloudsPtr_().Srho(eqn.psi());
231  }
232 }
233 
234 
236 (
237  const volScalarField& rho,
238  fvMatrix<scalar>& eqn,
239  const word& fieldName
240 ) const
241 {
242  if (debug)
243  {
244  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
245  }
246 
247  if (!carrierHasThermo_)
248  {
250  << "Applying source to compressible equation when carrier thermo "
251  << "is not available"
252  << exit(FatalError);
253  }
254 
255  const fluidThermo& carrierThermo = tCarrierThermo_();
256 
257  if (fieldName == rhoName_)
258  {
259  eqn += cloudsPtr_().Srho(eqn.psi());
260  }
261  else if (fieldName == carrierThermo.he().name())
262  {
263  eqn += cloudsPtr_().Sh(eqn.psi());
264  }
265  else if
266  (
267  isA<basicSpecieMixture>(carrierThermo)
268  && refCast<const basicSpecieMixture>(carrierThermo).contains
269  (
270  eqn.psi().name()
271  )
272  )
273  {
274  eqn += cloudsPtr_().SYi
275  (
276  refCast<const basicSpecieMixture>(carrierThermo).index(eqn.psi()),
277  eqn.psi()
278  );
279  }
280 }
281 
282 
284 (
285  fvMatrix<vector>& eqn,
286  const word& fieldName
287 ) const
288 {
289  if (debug)
290  {
291  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
292  }
293 
294  if (carrierHasThermo_)
295  {
297  << "Applying source to incompressible equation when carrier thermo "
298  << "is available"
299  << exit(FatalError);
300  }
301 
302  if (fieldName == UName_)
303  {
304  eqn += cloudsPtr_().SU(eqn.psi())/tRho_();
305  }
306 }
307 
308 
310 (
311  const volScalarField& rho,
312  fvMatrix<vector>& eqn,
313  const word& fieldName
314 ) const
315 {
316  if (debug)
317  {
318  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
319  }
320 
321  if (!carrierHasThermo_)
322  {
324  << "Applying source to compressible equation when carrier thermo "
325  << "is not available"
326  << exit(FatalError);
327  }
328 
329  if (fieldName == UName_)
330  {
331  eqn += cloudsPtr_().SU(eqn.psi());
332  }
333 }
334 
335 
337 {
338  // Store the particle positions
339  cloudsPtr_().storeGlobalPositions();
340 }
341 
342 
344 {
345  cloudsPtr_().topoChange(map);
346 }
347 
348 
350 {
351  cloudsPtr_().mapMesh(map);
352 }
353 
354 
356 {
357  cloudsPtr_().distribute(map);
358 }
359 
360 
362 {
363  return true;
364 }
365 
366 
367 // ************************************************************************* //
#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
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
virtual volScalarField & he()=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:160
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:101
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:361
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:191
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add source to continuity equation.
Definition: clouds.C:210
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: clouds.C:343
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: clouds.C:355
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: clouds.C:336
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: clouds.C:349
clouds(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: clouds.C:52
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
An abstract 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:306
static List< word > fieldNames
Definition: globalFoam.H:46
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:251
messageStream Info
const dimensionSet dimAcceleration
const dimensionSet dimDensity
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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
basicSpecieMixture & composition
PtrList< volScalarField > & Y