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