ThermoCloud.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-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 Class
25  Foam::ThermoCloud
26 
27 Description
28  Templated base class for thermodynamic cloud
29 
30  - Adds to momentum cloud
31  - Heat transfer
32 
33 SourceFiles
34  ThermoCloudI.H
35  ThermoCloud.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef ThermoCloud_H
40 #define ThermoCloud_H
41 
42 #include "volFieldsFwd.H"
43 #include "fvMatricesFwd.H"
44 #include "dimensionedTypes.H"
45 #include "fvMesh.H"
46 #include "parcelThermo.H"
47 #include "Cloud.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 // Forward declaration of classes
55 
56 class integrationScheme;
57 
58 template<class CloudType>
59 class HeatTransferModel;
60 
61 template<class CloudType>
62 class CompositionModel;
63 
64 
65 /*---------------------------------------------------------------------------*\
66  Class ThermoCloudName Declaration
67 \*---------------------------------------------------------------------------*/
68 
70 
71 
72 /*---------------------------------------------------------------------------*\
73  Class ThermoCloud Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 template<class CloudType>
77 class ThermoCloud
78 :
79  public CloudType,
80  public ThermoCloudName
81 {
82 public:
83 
84  // Public Typedefs
85 
86  //- Type of cloud this cloud was instantiated for
87  typedef CloudType cloudType;
88 
89  //- Type of parcel the cloud was instantiated for
90  typedef typename CloudType::particleType parcelType;
91 
92  //- Convenience typedef for this cloud type
94 
95 
96 private:
97 
98  // Private Data
99 
100  //- Cloud copy pointer
101  autoPtr<ThermoCloud<CloudType>> cloudCopyPtr_;
102 
103 
104 protected:
105 
106  // Protected data
107 
108  //- Thermo parcel constant properties
109  typename parcelType::constantProperties constProps_;
110 
111 
112  // References to the carrier gas fields
113 
114  //- Thermophysical properties of the carrier fluid
116 
117  //- SLG thermodynamics package
118  const parcelThermo thermo_;
119 
120  //- Temperature [K]
121  const volScalarField& T_;
122 
123  //- Pressure [Pa]
124  const volScalarField& p_;
125 
126 
127  // References to the cloud sub-models
128 
129  //- Heat transfer model
132 
133  //- Reacting composition model
136 
137 
138  // Reference to the particle integration schemes
139 
140  //- Temperature integration
142 
143 
144  // Modelling options
145 
146  //- Include radiation
148 
149  //- Radiation sum of parcel projected areas
151 
152  //- Radiation sum of parcel temperature^4
154 
155  //- Radiation sum of parcel projected areas * temperature^4
157 
158 
159  // Sources
160 
161  //- Sensible enthalpy transfer [J/kg]
163 
164  //- Coefficient for carrier phase hs equation [W/K]
166 
167 
168  // Protected Member Functions
169 
170  // Initialisation
171 
172  //- Set cloud sub-models
173  void setModels();
174 
175 
176  // Cloud evolution functions
177 
178  //- Reset state of cloud
180 
181 
182 public:
183 
184  // Constructors
185 
186  //- Construct given carrier fields and thermo
188  (
189  const word& cloudName,
190  const volScalarField& rho,
191  const volVectorField& U,
192  const dimensionedVector& g,
193  const fluidThermo& carrierThermo,
194  const bool readFields = true
195  );
196 
197  //- Copy constructor with new name
199 
200  //- Copy constructor with new name - creates bare cloud
202  (
203  const fvMesh& mesh,
204  const word& name,
206  );
207 
208  //- Disallow default bitwise copy construction
209  ThermoCloud(const ThermoCloud&) = delete;
210 
211  //- Construct and return clone based on (this) with new name
212  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
213  {
215  (
216  new ThermoCloud(*this, name)
217  );
218  }
219 
220  //- Construct and return bare clone based on (this) with new name
221  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
222  {
224  (
225  new ThermoCloud(this->mesh(), name, *this)
226  );
227  }
228 
229 
230  //- Destructor
231  virtual ~ThermoCloud();
232 
233 
234  // Member Functions
235 
236  // Access
237 
238  //- Return a reference to the cloud copy
239  inline const ThermoCloud& cloudCopy() const;
240 
241  //- Return the constant properties
242  inline const typename parcelType::constantProperties&
243  constProps() const;
244 
245  //- Return access to the constant properties
246  inline typename parcelType::constantProperties& constProps();
247 
248  //- Return const access to carrier thermo package
249  inline const fluidThermo& carrierThermo() const;
250 
251  //- Return const access to thermo package
252  inline const parcelThermo& thermo() const;
253 
254  //- Return const access to the carrier temperature field
255  inline const volScalarField& T() const;
256 
257  //- Return const access to the carrier pressure field
258  inline const volScalarField& p() const;
259 
260 
261  // Sub-models
262 
263  //- Return reference to heat transfer model
265  heatTransfer() const;
266 
267  //- Return const access to reacting composition model
269  composition() const;
270 
271 
272  // Integration schemes
273 
274  //-Return reference to velocity integration
275  inline const integrationScheme& TIntegrator() const;
276 
277 
278  // Modelling options
279 
280  //- Radiation flag
281  inline bool radiation() const;
282 
283  //- Radiation sum of parcel projected areas [m^2]
285 
286  //- Radiation sum of parcel projected areas [m^2]
287  inline const volScalarField::Internal&
288  radAreaP() const;
289 
290  //- Radiation sum of parcel temperature^4 [K4]
292 
293  //- Radiation sum of parcel temperature^4 [K4]
294  inline const volScalarField::Internal& radT4() const;
295 
296  //- Radiation sum of parcel projected area*temperature^4 [m2K4]
298 
299  //- Radiation sum of parcel temperature^4 [m2K4]
300  inline const volScalarField::Internal&
301  radAreaPT4() const;
302 
303 
304  // Sources
305 
306  // Enthalpy
307 
308  //- Return sensible enthalpy transfer [J/kg]
309  inline tmp<volScalarField::Internal> hsTrans() const;
310 
311  //- Access sensible enthalpy transfer [J/kg]
313 
314  //- Return coefficient for carrier phase hs equation
315  inline tmp<volScalarField::Internal> hsCoeff() const;
316 
317  //- Access coefficient for carrier phase hs equation
319 
320  //- Return sensible enthalpy source term [J/kg/m^3/s]
321  inline tmp<fvScalarMatrix> Sh
322  (
323  const volScalarField& hs
324  ) const;
325 
326 
327  // Radiation - overrides thermoCloud virtual abstract members
328 
329  //- Return tmp equivalent particulate emission
330  inline tmp<volScalarField> Ep() const;
331 
332  //- Return tmp equivalent particulate absorption
333  inline tmp<volScalarField> ap() const;
334 
335  //- Return tmp equivalent particulate scattering factor
336  inline tmp<volScalarField> sigmap() const;
337 
338 
339  // Check
340 
341  //- Maximum temperature
342  inline scalar Tmax() const;
343 
344  //- Minimum temperature
345  inline scalar Tmin() const;
346 
347 
348  // Cloud evolution functions
349 
350  //- Set parcel thermo properties
352  (
353  parcelType& parcel,
354  const scalar lagrangianDt
355  );
356 
357  //- Check parcel properties
359  (
360  parcelType& parcel,
361  const scalar lagrangianDt,
362  const bool fullyDescribed
363  );
364 
365  //- Store the current cloud state
366  void storeState();
367 
368  //- Reset the current cloud to the previously stored state
369  void restoreState();
370 
371  //- Reset the cloud source terms
372  void resetSourceTerms();
373 
374  //- Apply relaxation to (steady state) cloud sources
375  void relaxSources(const ThermoCloud<CloudType>& cloudOldTime);
376 
377  //- Apply scaling to (transient) cloud sources
378  void scaleSources();
379 
380  //- Pre-evolve
381  void preEvolve();
382 
383  //- Evolve the cloud
384  void evolve();
385 
386 
387  // Mapping
388 
389  //- Remap the cells of particles corresponding to the
390  // mesh topology change with a default tracking data object
391  virtual void autoMap(const polyTopoChangeMap&);
392 
393 
394  // I-O
395 
396  //- Print cloud information
397  void info();
398 
399  //- Write the field data for the cloud
400  virtual void writeFields() const;
401 
402 
403  // Member Operators
404 
405  //- Disallow default bitwise assignment
406  void operator=(const ThermoCloud&) = delete;
407 };
408 
409 
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
411 
412 } // End namespace Foam
413 
414 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
415 
416 #include "ThermoCloudI.H"
417 
418 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
419 
420 #ifdef NoRepository
421  #include "ThermoCloud.C"
422 #endif
423 
424 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
425 
426 #endif
427 
428 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:128
const volScalarField & T_
Temperature [K].
Definition: ThermoCloud.H:120
autoPtr< integrationScheme > TIntegrator_
Temperature integration.
Definition: ThermoCloud.H:140
const fluidThermo & carrierThermo_
Thermophysical properties of the carrier fluid.
Definition: ThermoCloud.H:114
CloudType cloudType
Type of cloud this cloud was instantiated for.
Definition: ThermoCloud.H:86
Switch radiation_
Include radiation.
Definition: ThermoCloud.H:146
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:288
const word & name() const
Return name.
Definition: IOobject.H:315
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: ThermoCloud.C:361
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
Definition: ThermoCloud.H:220
Templated heat transfer model class.
Definition: ThermoCloud.H:58
const parcelThermo thermo_
SLG thermodynamics package.
Definition: ThermoCloud.H:117
const ThermoCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: ThermoCloudI.H:35
U
Definition: pEqn.H:72
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
tmp< volScalarField > sigmap() const
Return tmp equivalent particulate scattering factor.
Definition: ThermoCloudI.H:347
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: parcelThermo.H:58
void restoreState()
Reset the current cloud to the previously stored state.
Definition: ThermoCloud.C:399
const volScalarField & p() const
Return const access to the carrier pressure field.
Definition: ThermoCloudI.H:80
autoPtr< volScalarField::Internal > hsTrans_
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloud.H:161
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
volScalarField::Internal & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
Definition: ThermoCloudI.H:179
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/kg/m^3/s].
Definition: ThermoCloudI.H:241
tmp< volScalarField::Internal > hsTrans() const
Return sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:209
void cloudReset(ThermoCloud< CloudType > &c)
Reset state of cloud.
Definition: ThermoCloud.C:125
void preEvolve()
Pre-evolve.
Definition: ThermoCloud.C:460
const dimensionedScalar c
Speed of light in a vacuum.
Templated base class for thermodynamic cloud.
Definition: ThermoCloud.H:76
ThermoCloud< CloudType > thermoCloudType
Convenience typedef for this cloud type.
Definition: ThermoCloud.H:92
parcelType::constantProperties constProps_
Thermo parcel constant properties.
Definition: ThermoCloud.H:108
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
Definition: ThermoCloudI.H:88
volScalarField::Internal & hsCoeffRef()
Access coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:233
volScalarField::Internal & radAreaP()
Radiation sum of parcel projected areas [m^2].
Definition: ThermoCloudI.H:119
autoPtr< volScalarField::Internal > radT4_
Radiation sum of parcel temperature^4.
Definition: ThermoCloud.H:152
void operator=(const ThermoCloud &)=delete
Disallow default bitwise assignment.
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
Definition: ThermoCloudI.H:149
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
Definition: ThermoCloudI.H:318
virtual void autoMap(const polyTopoChangeMap &)
Remap the cells of particles corresponding to the.
Definition: ThermoCloud.C:481
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const parcelType::constantProperties & constProps() const
Return the constant properties.
Definition: ThermoCloudI.H:43
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
ThermoCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const fluidThermo &carrierThermo, const bool readFields=true)
Construct given carrier fields and thermo.
Definition: ThermoCloud.C:142
TemplateName(FvFaceCellWave)
bool radiation() const
Radiation flag.
Definition: ThermoCloudI.H:111
const integrationScheme & TIntegrator() const
Return reference to velocity integration.
Definition: ThermoCloudI.H:104
const volScalarField & T() const
Return const access to the carrier temperature field.
Definition: ThermoCloudI.H:73
autoPtr< volScalarField::Internal > radAreaPT4_
Radiation sum of parcel projected areas * temperature^4.
Definition: ThermoCloud.H:155
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions...
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: ThermoCloud.C:375
void resetSourceTerms()
Reset the cloud source terms.
Definition: ThermoCloud.C:407
const volScalarField & p_
Pressure [Pa].
Definition: ThermoCloud.H:123
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
virtual void writeFields() const
Write the field data for the cloud.
Definition: ThermoCloud.C:500
Forward declarations of fvMatrix specialisations.
autoPtr< volScalarField::Internal > hsCoeff_
Coefficient for carrier phase hs equation [W/K].
Definition: ThermoCloud.H:164
scalar Tmax() const
Maximum temperature.
Definition: ThermoCloudI.H:376
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
Definition: ThermoCloudI.H:290
void evolve()
Evolve the cloud.
Definition: ThermoCloud.C:469
autoPtr< volScalarField::Internal > radAreaP_
Radiation sum of parcel projected areas.
Definition: ThermoCloud.H:149
volScalarField::Internal & hsTransRef()
Access sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:217
void info()
Print cloud information.
Definition: ThermoCloud.C:490
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void scaleSources()
Apply scaling to (transient) cloud sources.
Definition: ThermoCloud.C:443
const fluidThermo & carrierThermo() const
Return const access to carrier thermo package.
Definition: ThermoCloudI.H:59
tmp< volScalarField::Internal > hsCoeff() const
Return coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:225
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const parcelThermo & thermo() const
Return const access to thermo package.
Definition: ThermoCloudI.H:66
const CompositionModel< ThermoCloud< CloudType > > & composition() const
Return const access to reacting composition model.
Definition: ThermoCloudI.H:96
autoPtr< HeatTransferModel< ThermoCloud< CloudType > > > heatTransferModel_
Heat transfer model.
Definition: ThermoCloud.H:130
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: ThermoCloud.H:61
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
scalar Tmin() const
Minimum temperature.
Definition: ThermoCloudI.H:402
void storeState()
Store the current cloud state.
Definition: ThermoCloud.C:386
void setModels()
Set cloud sub-models.
Definition: ThermoCloud.C:34
Namespace for OpenFOAM.
autoPtr< CompositionModel< ThermoCloud< CloudType > > > compositionModel_
Reacting composition model.
Definition: ThermoCloud.H:134
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition: ThermoCloud.H:89
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: ThermoCloud.C:424
virtual ~ThermoCloud()
Destructor.
Definition: ThermoCloud.C:353