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