ReactingCloud.H
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) 2011-2019 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 Class
25  Foam::ReactingCloud
26 
27 Description
28  Templated base class for reacting cloud
29 
30  - Adds to thermodynamic cloud
31  - Variable composition (single phase)
32  - Phase change
33 
34 SourceFiles
35  ReactingCloudI.H
36  ReactingCloud.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef ReactingCloud_H
41 #define ReactingCloud_H
42 
43 #include "reactingCloud.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // Forward declaration of classes
51 
52 template<class CloudType>
53 class CompositionModel;
54 
55 template<class CloudType>
56 class PhaseChangeModel;
57 
58 /*---------------------------------------------------------------------------*\
59  Class ReactingCloud Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 template<class CloudType>
63 class ReactingCloud
64 :
65  public CloudType,
66  public reactingCloud
67 {
68 public:
69 
70  // Public Typedefs
71 
72  //- Type of cloud this cloud was instantiated for
73  typedef CloudType cloudType;
74 
75  //- Type of parcel the cloud was instantiated for
76  typedef typename CloudType::particleType parcelType;
77 
78  //- Convenience typedef for this cloud type
80 
81 
82 private:
83 
84  // Private Data
85 
86  //- Cloud copy pointer
87  autoPtr<ReactingCloud<CloudType>> cloudCopyPtr_;
88 
89 
90 
91 protected:
92 
93  // Protected data
94 
95  //- Parcel constant properties
96  typename parcelType::constantProperties constProps_;
97 
98 
99  // References to the cloud sub-models
100 
101  //- Reacting composition model
104 
105  //- Reacting phase change model
108 
109 
110  // Sources
111 
112  //- Mass transfer fields - one per carrier phase specie
114 
115 
116  // Protected Member Functions
117 
118  // New parcel helper functions
119 
120  //- Check that size of a composition field is valid
122  (
123  const scalarField& YSupplied,
124  const scalarField& Y,
125  const word& YName
126  );
127 
128 
129  // Initialisation
130 
131  //- Set cloud sub-models
132  void setModels();
133 
134 
135  // Cloud evolution functions
136 
137  //- Reset state of cloud
139 
140 
141 public:
142 
143  // Constructors
144 
145  //- Construct given carrier gas fields
147  (
148  const word& cloudName,
149  const volScalarField& rho,
150  const volVectorField& U,
151  const dimensionedVector& g,
152  const SLGThermo& thermo,
153  bool readFields = true
154  );
155 
156  //- Copy constructor with new name
158 
159  //- Copy constructor with new name - creates bare cloud
161  (
162  const fvMesh& mesh,
163  const word& name,
165  );
166 
167  //- Disallow default bitwise copy construction
168  ReactingCloud(const ReactingCloud&) = delete;
169 
170  //- Construct and return clone based on (this) with new name
171  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
172  {
174  (
175  new ReactingCloud(*this, name)
176  );
177  }
178 
179  //- Construct and return bare clone based on (this) with new name
180  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
181  {
183  (
184  new ReactingCloud(this->mesh(), name, *this)
185  );
186  }
187 
188 
189  //- Destructor
190  virtual ~ReactingCloud();
191 
192 
193  // Member Functions
194 
195  // Access
196 
197  //- Return a reference to the cloud copy
198  inline const ReactingCloud& cloudCopy() const;
199 
200  //- Return the constant properties
201  inline const typename parcelType::constantProperties&
202  constProps() const;
203 
204  //- Return access to the constant properties
205  inline typename parcelType::constantProperties& constProps();
206 
207 
208  // Sub-models
209 
210  //- Return const access to reacting composition model
212  composition() const;
213 
214  //- Return const access to reacting phase change model
216  phaseChange() const;
217 
218  //- Return reference to reacting phase change model
220  phaseChange();
221 
222 
223  // Sources
224 
225  //- Mass
226 
227  //- Return reference to mass source for field i
229  rhoTrans(const label i);
230 
231  //- Return const access to mass source fields
233  rhoTrans() const;
234 
235  //- Return reference to mass source fields
237  rhoTrans();
238 
239  //- Return mass source term for specie i - specie eqn
240  inline tmp<fvScalarMatrix> SYi
241  (
242  const label i,
243  volScalarField& Yi
244  ) const;
245 
246  //- Return tmp mass source for field i - fully explicit
248  Srho(const label i) const;
249 
250  //- Return tmp total mass source for carrier phase
251  // - fully explicit
252  inline tmp<volScalarField::Internal> Srho() const;
253 
254  //- Return total mass source term [kg/m^3/s]
256 
257 
258  // Cloud evolution functions
259 
260  //- Set parcel thermo properties
262  (
263  parcelType& parcel,
264  const scalar lagrangianDt
265  );
266 
267  //- Check parcel properties
269  (
270  parcelType& parcel,
271  const scalar lagrangianDt,
272  const bool fullyDescribed
273  );
274 
275  //- Store the current cloud state
276  void storeState();
277 
278  //- Reset the current cloud to the previously stored state
279  void restoreState();
280 
281  //- Reset the cloud source terms
282  void resetSourceTerms();
283 
284  //- Apply relaxation to (steady state) cloud sources
285  void relaxSources(const ReactingCloud<CloudType>& cloudOldTime);
286 
287  //- Apply scaling to (transient) cloud sources
288  void scaleSources();
289 
290  //- Evolve the cloud
291  void evolve();
292 
293 
294  // Mapping
295 
296  //- Remap the cells of particles corresponding to the
297  // mesh topology change with a default tracking data object
298  virtual void autoMap(const mapPolyMesh&);
299 
300 
301  // I-O
302 
303  //- Print cloud information
304  void info();
305 
306  //- Write the field data for the cloud
307  virtual void writeFields() const;
308 
309 
310  // Member Operators
311 
312  //- Disallow default bitwise assignment
313  void operator=(const ReactingCloud&) = delete;
314 };
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace Foam
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 #include "ReactingCloudI.H"
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 #ifdef NoRepository
328  #include "ReactingCloud.C"
329 #endif
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 #endif
334 
335 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:105
virtual ~ReactingCloud()
Destructor.
void resetSourceTerms()
Reset the cloud source terms.
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:268
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
const word & name() const
Return name.
Definition: IOobject.H:295
Templated phase change model class.
Definition: ReactingCloud.H:55
void storeState()
Store the current cloud state.
parcelType::constantProperties constProps_
Parcel constant properties.
Definition: ReactingCloud.H:95
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void restoreState()
Reset the current cloud to the previously stored state.
Virtual abstract base class for templated ReactingCloud.
Definition: reactingCloud.H:48
void setModels()
Set cloud sub-models.
Definition: ReactingCloud.C:34
rhoReactionThermo & thermo
Definition: createFields.H:28
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
virtual void writeFields() const
Write the field data for the cloud.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
ReactingCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const SLGThermo &thermo, bool readFields=true)
Construct given carrier gas fields.
Definition: ReactingCloud.C:90
void evolve()
Evolve the cloud.
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
A class for handling words, derived from string.
Definition: word.H:59
autoPtr< PhaseChangeModel< ReactingCloud< CloudType > > > phaseChangeModel_
Reacting phase change model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
autoPtr< CompositionModel< ReactingCloud< CloudType > > > compositionModel_
Reacting composition model.
CloudType cloudType
Type of cloud this cloud was instantiated for.
Definition: ReactingCloud.H:72
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void scaleSources()
Apply scaling to (transient) cloud sources.
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
Templated base class for reacting cloud.
Definition: ReactingCloud.H:62
U
Definition: pEqn.H:72
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
const dimensionedScalar c
Speed of light in a vacuum.
void info()
Print cloud information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
Definition: ReactingCloud.C:77
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
Definition: ReactingCloud.C:58
void operator=(const ReactingCloud &)=delete
Disallow default bitwise assignment.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ReactingCloud.H:52
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition: ReactingCloud.H:75
ReactingCloud< CloudType > reactingCloudType
Convenience typedef for this cloud type.
Definition: ReactingCloud.H:78
Namespace for OpenFOAM.