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-2022 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  {
37  defineTypeNameAndDebug(clouds, 0);
38 
40  (
41  fvModel,
42  clouds,
43  dictionary
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& sourceName,
54  const word& modelType,
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
59  fvModel(sourceName, modelType, dict, mesh),
60  g_
61  (
62  IOobject
63  (
64  "g",
65  mesh.time().constant(),
66  mesh,
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().timeName(),
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  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
123  UName_(dict.lookupOrDefault<word>("U", "U")),
124  cloudsPtr_
125  (
126  carrierHasThermo_
127  ? new parcelCloudList
128  (
129  mesh.lookupObject<volScalarField>(rhoName_),
130  mesh.lookupObject<volVectorField>(UName_),
131  g_,
132  tCarrierThermo_()
133  )
134  : new parcelCloudList
135  (
136  tRho_(),
137  mesh.lookupObject<volVectorField>(UName_),
138  tMu_(),
139  g_
140  )
141  ),
142  curTimeIndex_(-1)
143 {}
144 
145 
146 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 
149 {
150  wordList fieldNames(1, UName_);
151 
152  if (carrierHasThermo_)
153  {
154  const fluidThermo& carrierThermo = tCarrierThermo_();
155 
156  fieldNames.append(rhoName_);
157 
158  fieldNames.append(carrierThermo.he().name());
159 
160  if (isA<basicSpecieMixture>(carrierThermo))
161  {
163  refCast<const basicSpecieMixture>(carrierThermo);
164 
165  const PtrList<volScalarField>& Y = composition.Y();
166 
167  forAll(Y, i)
168  {
169  if (composition.solve(i))
170  {
171  fieldNames.append(Y[i].name());
172  }
173  }
174  }
175  }
176 
177  return fieldNames;
178 }
179 
180 
182 {
183  if (curTimeIndex_ == mesh().time().timeIndex())
184  {
185  return;
186  }
187 
188  if (!carrierHasThermo_)
189  {
190  tMu_.ref() = tRho_()*tCarrierViscosity_().nu();
191  }
192 
193  cloudsPtr_().evolve();
194 
195  curTimeIndex_ = mesh().time().timeIndex();
196 }
197 
198 
200 (
201  fvMatrix<scalar>& eqn,
202  const word& fieldName
203 ) const
204 {
205  if (debug)
206  {
207  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
208  }
209 
210  if (!carrierHasThermo_)
211  {
213  << "Applying source to compressible equation when carrier thermo "
214  << "is not available"
215  << exit(FatalError);
216  }
217 
218  if (fieldName == rhoName_)
219  {
220  eqn += cloudsPtr_().Srho(eqn.psi());
221  }
222 }
223 
224 
226 (
227  const volScalarField& rho,
228  fvMatrix<scalar>& eqn,
229  const word& fieldName
230 ) const
231 {
232  if (debug)
233  {
234  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
235  }
236 
237  if (!carrierHasThermo_)
238  {
240  << "Applying source to compressible equation when carrier thermo "
241  << "is not available"
242  << exit(FatalError);
243  }
244 
245  const fluidThermo& carrierThermo = tCarrierThermo_();
246 
247  if (fieldName == rhoName_)
248  {
249  eqn += cloudsPtr_().Srho(eqn.psi());
250  }
251  else if (fieldName == carrierThermo.he().name())
252  {
253  eqn += cloudsPtr_().Sh(eqn.psi());
254  }
255  else if
256  (
257  isA<basicSpecieMixture>(carrierThermo)
258  && refCast<const basicSpecieMixture>(carrierThermo).contains
259  (
260  eqn.psi().name()
261  )
262  )
263  {
264  eqn += cloudsPtr_().SYi
265  (
266  refCast<const basicSpecieMixture>(carrierThermo).index(eqn.psi()),
267  eqn.psi()
268  );
269  }
270 }
271 
272 
274 (
275  fvMatrix<vector>& eqn,
276  const word& fieldName
277 ) const
278 {
279  if (debug)
280  {
281  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
282  }
283 
284  if (carrierHasThermo_)
285  {
287  << "Applying source to incompressible equation when carrier thermo "
288  << "is available"
289  << exit(FatalError);
290  }
291 
292  if (fieldName == UName_)
293  {
294  eqn += cloudsPtr_().SU(eqn.psi())/tRho_();
295  }
296 }
297 
298 
300 (
301  const volScalarField& rho,
302  fvMatrix<vector>& eqn,
303  const word& fieldName
304 ) const
305 {
306  if (debug)
307  {
308  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
309  }
310 
311  if (!carrierHasThermo_)
312  {
314  << "Applying source to compressible equation when carrier thermo "
315  << "is not available"
316  << exit(FatalError);
317  }
318 
319  if (fieldName == UName_)
320  {
321  eqn += cloudsPtr_().SU(eqn.psi());
322  }
323 }
324 
325 
327 {
328  // Store the particle positions
329  cloudsPtr_().storeGlobalPositions();
330 }
331 
332 
334 {}
335 
336 
338 {}
339 
340 
342 {
343  cloudsPtr_().distribute(map);
344 }
345 
346 
348 {
349  return true;
350 }
351 
352 
353 // ************************************************************************* //
clouds(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: clouds.C:52
bool solve(label speciei) const
Return true if the specie should be solved for.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: clouds.C:337
const word & name() const
Return name.
Definition: IOobject.H:315
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: clouds.C:333
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: clouds.C:326
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
basicSpecieMixture & composition
An abstract base class for physical properties.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool foundObject(const word &name) const
Is the named Type found?
Finite volume model abstract base class.
Definition: fvModel.H:57
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add source to continuity equation.
Definition: clouds.C:200
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
fvMesh & mesh
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
An abstract base class for Newtonian viscosity models.
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
static List< word > fieldNames
Definition: globalFoam.H:46
virtual void correct()
Solve the Lagrangian clouds and update the sources.
Definition: clouds.C:181
const dimensionSet dimAcceleration
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: clouds.C:148
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
const dimensionSet dimDensity
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
virtual bool movePoints()
Update for mesh motion.
Definition: clouds.C:347
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:52
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
label timeIndex
Definition: getTimeIndex.H:4
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: clouds.C:341
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Namespace for OpenFOAM.