ThermoCloud.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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  // Private member functions
89 
90  //- Disallow default bitwise copy construct
91  ThermoCloud(const ThermoCloud&);
92 
93  //- Disallow default bitwise assignment
94  void operator=(const ThermoCloud&);
95 
96 
97 protected:
98 
99  // Protected data
100 
101  //- Thermo parcel constant properties
102  typename parcelType::constantProperties constProps_;
103 
104 
105  // References to the carrier gas fields
106 
107  //- SLG thermodynamics package
108  const SLGThermo& thermo_;
109 
110  //- Temperature [K]
111  const volScalarField& T_;
112 
113  //- Pressure [Pa]
114  const volScalarField& p_;
115 
116 
117  // References to the cloud sub-models
118 
119  //- Heat transfer model
122 
123 
124  // Reference to the particle integration schemes
125 
126  //- Temperature integration
128 
129 
130  // Modelling options
131 
132  //- Include radiation
134 
135  //- Radiation sum of parcel projected areas
137 
138  //- Radiation sum of parcel temperature^4
140 
141  //- Radiation sum of parcel projected areas * temperature^4
143 
144 
145  // Sources
146 
147  //- Sensible enthalpy transfer [J/kg]
149 
150  //- Coefficient for carrier phase hs equation [W/K]
152 
153 
154  // Protected Member Functions
155 
156  // Initialisation
157 
158  //- Set cloud sub-models
159  void setModels();
160 
161 
162  // Cloud evolution functions
163 
164  //- Reset state of cloud
166 
167 
168 public:
169 
170  // Constructors
171 
172  //- Construct given carrier gas fields
174  (
175  const word& cloudName,
176  const volScalarField& rho,
177  const volVectorField& U,
178  const dimensionedVector& g,
179  const SLGThermo& thermo,
180  bool readFields = true
181  );
182 
183  //- Copy constructor with new name
185 
186  //- Copy constructor with new name - creates bare cloud
188  (
189  const fvMesh& mesh,
190  const word& name,
192  );
193 
194  //- Construct and return clone based on (this) with new name
195  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
196  {
198  (
199  new ThermoCloud(*this, name)
200  );
201  }
202 
203  //- Construct and return bare clone based on (this) with new name
204  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
205  {
207  (
208  new ThermoCloud(this->mesh(), name, *this)
209  );
210  }
211 
212 
213  //- Destructor
214  virtual ~ThermoCloud();
215 
216 
217  // Member Functions
218 
219  // Access
220 
221  //- Return a reference to the cloud copy
222  inline const ThermoCloud& cloudCopy() const;
223 
224  //- Return the constant properties
225  inline const typename parcelType::constantProperties&
226  constProps() const;
227 
228  //- Return access to the constant properties
229  inline typename parcelType::constantProperties& constProps();
230 
231  //- Return const access to thermo package
232  inline const SLGThermo& thermo() const;
233 
234  //- Return const access to the carrier temperature field
235  inline const volScalarField& T() const;
236 
237  //- Return const access to the carrier prressure field
238  inline const volScalarField& p() const;
239 
240 
241  // Sub-models
242 
243  //- Return reference to heat transfer model
245  heatTransfer() const;
246 
247 
248  // Integration schemes
249 
250  //-Return reference to velocity integration
251  inline const scalarIntegrationScheme& TIntegrator() const;
252 
253 
254  // Modelling options
255 
256  //- Radiation flag
257  inline bool radiation() const;
258 
259  //- Radiation sum of parcel projected areas [m2]
261 
262  //- Radiation sum of parcel projected areas [m2]
264  radAreaP() const;
265 
266  //- Radiation sum of parcel temperature^4 [K4]
268 
269  //- Radiation sum of parcel temperature^4 [K4]
270  inline const DimensionedField<scalar, volMesh>& radT4() const;
271 
272  //- Radiation sum of parcel projected area*temperature^4 [m2K4]
274 
275  //- Radiation sum of parcel temperature^4 [m2K4]
277  radAreaPT4() const;
278 
279 
280  // Sources
281 
282  // Enthalpy
283 
284  //- Sensible enthalpy transfer [J/kg]
286 
287  //- Sensible enthalpy transfer [J/kg]
289  hsTrans() const;
290 
291  //- Return coefficient for carrier phase hs equation
293 
294  //- Return const coefficient for carrier phase hs equation
296  hsCoeff() const;
297 
298  //- Return sensible enthalpy source term [J/kg/m3/s]
299  inline tmp<fvScalarMatrix> Sh(volScalarField& hs) const;
300 
301 
302  // Radiation - overrides thermoCloud virtual abstract members
303 
304  //- Return tmp equivalent particulate emission
305  inline tmp<volScalarField> Ep() const;
306 
307  //- Return tmp equivalent particulate absorption
308  inline tmp<volScalarField> ap() const;
309 
310  //- Return tmp equivalent particulate scattering factor
311  inline tmp<volScalarField> sigmap() const;
312 
313 
314  // Check
315 
316  //- Maximum temperature
317  inline scalar Tmax() const;
318 
319  //- Minimum temperature
320  inline scalar Tmin() const;
321 
322 
323  // Cloud evolution functions
324 
325  //- Set parcel thermo properties
327  (
328  parcelType& parcel,
329  const scalar lagrangianDt
330  );
331 
332  //- Check parcel properties
334  (
335  parcelType& parcel,
336  const scalar lagrangianDt,
337  const bool fullyDescribed
338  );
339 
340  //- Store the current cloud state
341  void storeState();
342 
343  //- Reset the current cloud to the previously stored state
344  void restoreState();
345 
346  //- Reset the cloud source terms
347  void resetSourceTerms();
348 
349  //- Apply relaxation to (steady state) cloud sources
350  void relaxSources(const ThermoCloud<CloudType>& cloudOldTime);
351 
352  //- Apply scaling to (transient) cloud sources
353  void scaleSources();
354 
355  //- Pre-evolve
356  void preEvolve();
357 
358  //- Evolve the cloud
359  void evolve();
360 
361 
362  // Mapping
363 
364  //- Remap the cells of particles corresponding to the
365  // mesh topology change with a default tracking data object
366  virtual void autoMap(const mapPolyMesh&);
367 
368 
369  // I-O
370 
371  //- Print cloud information
372  void info();
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:114
const volScalarField & T_
Temperature [K].
Definition: ThermoCloud.H:110
CloudType cloudType
Type of cloud this cloud was instantiated for.
Definition: ThermoCloud.H:70
Switch radiation_
Include radiation.
Definition: ThermoCloud.H:132
autoPtr< DimensionedField< scalar, volMesh > > radT4_
Radiation sum of parcel temperature^4.
Definition: ThermoCloud.H:138
U
Definition: pEqn.H:83
const SLGThermo & thermo_
SLG thermodynamics package.
Definition: ThermoCloud.H:107
tmp< volScalarField > sigmap() const
Return tmp equivalent particulate scattering factor.
Definition: ThermoCloudI.H:336
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: ThermoCloud.C:362
Templated heat transfer model class.
Definition: ThermoCloud.H:53
DimensionedField< scalar, volMesh > & hsTrans()
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:192
scalar Tmin() const
Minimum temperature.
Definition: ThermoCloudI.H:399
const parcelType::constantProperties & constProps() const
Return the constant properties.
Definition: ThermoCloudI.H:42
autoPtr< DimensionedField< scalar, volMesh > > hsTrans_
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloud.H:147
void restoreState()
Reset the current cloud to the previously stored state.
Definition: ThermoCloud.C:400
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
Definition: ThermoCloud.H:203
const SLGThermo & thermo() const
Return const access to thermo package.
Definition: ThermoCloudI.H:57
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.
Definition: Switch.H:60
DimensionedField< scalar, volMesh > & 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:121
void preEvolve()
Pre-evolve.
Definition: ThermoCloud.C:461
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:101
DimensionedField< scalar, volMesh > & radAreaP()
Radiation sum of parcel projected areas [m2].
Definition: ThermoCloudI.H:102
DimensionedField< scalar, volMesh > & radT4()
Radiation sum of parcel temperature^4 [K4].
Definition: ThermoCloudI.H:132
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: ThermoCloud.C:483
autoPtr< scalarIntegrationScheme > TIntegrator_
Temperature integration.
Definition: ThermoCloud.H:126
A class for handling words, derived from string.
Definition: word.H:59
scalar Tmax() const
Maximum temperature.
Definition: ThermoCloudI.H:373
Foam::autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:239
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
Definition: ThermoCloudI.H:79
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const dimensionedVector & g
bool radiation() const
Radiation flag.
Definition: ThermoCloudI.H:94
Virtual abstract base class for templated ThermoCloud.
Definition: thermoCloud.H:48
virtual void readFields()
Read the field data for the cloud of particles. Dummy at.
const volScalarField & T() const
Return const access to the carrier temperature field.
Definition: ThermoCloudI.H:64
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: ThermoCloud.C:376
void resetSourceTerms()
Reset the cloud source terms.
Definition: ThermoCloud.C:408
const volScalarField & p_
Pressure [Pa].
Definition: ThermoCloud.H:113
const ThermoCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: ThermoCloudI.H:34
void evolve()
Evolve the cloud.
Definition: ThermoCloud.C:470
const fvMesh & mesh() const
Return refernce to the mesh.
Definition: DSMCCloudI.H:41
autoPtr< DimensionedField< scalar, volMesh > > radAreaPT4_
Radiation sum of parcel projected areas * temperature^4.
Definition: ThermoCloud.H:141
void info()
Print cloud information.
Definition: ThermoCloud.C:496
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
Top level model for Integration schemes.
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:444
autoPtr< DimensionedField< scalar, volMesh > > hsCoeff_
Coefficient for carrier phase hs equation [W/K].
Definition: ThermoCloud.H:150
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
const volScalarField & p() const
Return const access to the carrier prressure field.
Definition: ThermoCloudI.H:71
tmp< fvScalarMatrix > Sh(volScalarField &hs) const
Return sensible enthalpy source term [J/kg/m3/s].
Definition: ThermoCloudI.H:224
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
autoPtr< HeatTransferModel< ThermoCloud< CloudType > > > heatTransferModel_
Heat transfer model.
Definition: ThermoCloud.H:120
A class for managing temporary objects.
Definition: PtrList.H:54
autoPtr< DimensionedField< scalar, volMesh > > radAreaP_
Radiation sum of parcel projected areas.
Definition: ThermoCloud.H:135
DimensionedField< scalar, volMesh > & hsCoeff()
Return coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:208
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
Definition: ThermoCloudI.H:299
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
const scalarIntegrationScheme & TIntegrator() const
Return reference to velocity integration.
Definition: ThermoCloudI.H:87
const word & name() const
Return name.
Definition: IOobject.H:260
void storeState()
Store the current cloud state.
Definition: ThermoCloud.C:387
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:425
virtual ~ThermoCloud()
Destructor.
Definition: ThermoCloud.C:354
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
Definition: ThermoCloudI.H:263