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 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  carrierThermo_
61  (
63  ),
64  clouds_
65  (
66  mesh.lookupObject<volScalarField>("rho"),
67  mesh.lookupObject<volVectorField>("U"),
69  carrierThermo_
70  ),
71  curTimeIndex_(-1)
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 {
79  wordList fieldNames({"rho", "U", carrierThermo_.he().name()});
80 
81  if (isA<basicSpecieMixture>(carrierThermo_))
82  {
84  refCast<const basicSpecieMixture>(carrierThermo_);
85 
86  const PtrList<volScalarField>& Y = composition.Y();
87 
88  forAll(Y, i)
89  {
90  if (composition.solve(i))
91  {
92  fieldNames.append(Y[i].name());
93  }
94  }
95  }
96 
97  return fieldNames;
98 }
99 
100 
102 {
103  if (curTimeIndex_ == mesh().time().timeIndex())
104  {
105  return;
106  }
107 
108  clouds_.evolve();
109 
110  curTimeIndex_ = mesh().time().timeIndex();
111 }
112 
113 
115 (
116  fvMatrix<scalar>& eqn,
117  const word& fieldName
118 ) const
119 {
120  if (debug)
121  {
122  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
123  }
124 
125  if (fieldName == "rho")
126  {
127  eqn += clouds_.Srho(eqn.psi());
128  }
129  else
130  {
132  << "Support for field " << fieldName << " is not implemented"
133  << exit(FatalError);
134  }
135 }
136 
137 
139 (
140  const volScalarField& rho,
141  fvMatrix<scalar>& eqn,
142  const word& fieldName
143 ) const
144 {
145  if (debug)
146  {
147  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
148  }
149 
150  if (fieldName == "rho")
151  {
152  eqn += clouds_.Srho(eqn.psi());
153  }
154  else if (fieldName == carrierThermo_.he().name())
155  {
156  eqn += clouds_.Sh(eqn.psi());
157  }
158  else if
159  (
160  isA<basicSpecieMixture>(carrierThermo_)
161  && refCast<const basicSpecieMixture>(carrierThermo_).contains
162  (
163  eqn.psi().name()
164  )
165  )
166  {
167  eqn += clouds_.SYi
168  (
169  refCast<const basicSpecieMixture>(carrierThermo_).index(eqn.psi()),
170  eqn.psi()
171  );
172  }
173  else
174  {
176  << "Support for field " << fieldName << " is not implemented"
177  << exit(FatalError);
178  }
179 }
180 
181 
183 (
184  const volScalarField& rho,
185  fvMatrix<vector>& eqn,
186  const word& fieldName
187 ) const
188 {
189  if (debug)
190  {
191  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
192  }
193 
194  if (fieldName == "U")
195  {
196  eqn += clouds_.SU(eqn.psi());
197  }
198  else
199  {
201  << "Support for field " << fieldName << " is not implemented"
202  << exit(FatalError);
203  }
204 }
205 
206 
208 {
209  // Store the particle positions
210  clouds_.storeGlobalPositions();
211 }
212 
213 
214 // ************************************************************************* //
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
const word & name() const
Return name.
Definition: IOobject.H:303
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: clouds.C:207
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:323
basicSpecieMixture & composition
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Finite volume model abstract base class.
Definition: fvModel.H:55
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add source to continuity equation.
Definition: clouds.C:115
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
static List< word > fieldNames
Definition: globalFoam.H:46
virtual void correct()
Solve the Lagrangian clouds and update the sources.
Definition: clouds.C:101
dynamicFvMesh & mesh
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: clouds.C:77
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
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
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
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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:78
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
messageStream Info
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Namespace for OpenFOAM.