Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::ThermoCloud
29
30Group
31 grpLagrangianIntermediateClouds
32
33Description
34 Templated base class for thermodynamic cloud
35
36 - Adds to kinematic cloud
37 - Heat transfer
38
39SourceFiles
40 ThermoCloudI.H
41 ThermoCloud.C
42
43\*---------------------------------------------------------------------------*/
44
45#ifndef ThermoCloud_H
46#define ThermoCloud_H
47
48#include "KinematicCloud.H"
49#include "thermoCloud.H"
50#include "SLGThermo.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57// Forward declaration of classes
58
59template<class CloudType>
61
62/*---------------------------------------------------------------------------*\
63 Class ThermoCloud Declaration
64\*---------------------------------------------------------------------------*/
65
66template<class CloudType>
67class ThermoCloud
68:
69 public CloudType,
70 public thermoCloud
71{
72public:
73
74 // Public typedefs
75
76 //- Type of cloud this cloud was instantiated for
77 typedef CloudType cloudType;
78
79 //- Type of parcel the cloud was instantiated for
80 typedef typename CloudType::particleType parcelType;
81
82 //- Convenience typedef for this cloud type
83 typedef ThermoCloud<CloudType> thermoCloudType;
85
86private:
87
88 // Private data
89
90 //- Cloud copy pointer
91 autoPtr<ThermoCloud<CloudType>> cloudCopyPtr_;
92
93
94 // Private member functions
95
96 //- No copy construct
97 ThermoCloud(const ThermoCloud&) = delete;
98
99 //- No copy assignment
100 void operator=(const ThermoCloud&) = delete;
101
102
103protected:
104
105 // Protected data
106
107 //- Thermo parcel constant properties
108 typename parcelType::constantProperties constProps_;
109
110
111 // References to the carrier gas fields
112
113 //- SLG thermodynamics package
114 const SLGThermo& thermo_;
115
116 //- Temperature [K]
118
119 //- Pressure [Pa]
120 const volScalarField& p_;
121
122
123 // References to the cloud sub-models
124
125 //- Heat transfer model
128
129
130 // Reference to the particle integration schemes
131
132 //- Temperature integration
134
136 // Modelling options
137
138 //- Include radiation
140
141 //- Radiation sum of parcel projected areas
143
144 //- Radiation sum of parcel temperature^4
146
147 //- Radiation sum of parcel projected areas * temperature^4
149
150
151 // Sources
153 //- Sensible enthalpy transfer [J/kg]
155
156 //- Coefficient for carrier phase hs equation [W/K]
158
159
160 // Protected Member Functions
161
162 // Initialisation
163
164 //- Set cloud sub-models
165 void setModels();
166
167
168 // Cloud evolution functions
169
170 //- Reset state of cloud
171 void cloudReset(ThermoCloud<CloudType>& c);
172
173
174public:
176 // Constructors
177
178 //- Construct given carrier gas fields
179 ThermoCloud
180 (
181 const word& cloudName,
182 const volScalarField& rho,
184 const dimensionedVector& g,
185 const SLGThermo& thermo,
186 bool readFields = true
187 );
189 //- Copy constructor with new name
190 ThermoCloud(ThermoCloud<CloudType>& c, const word& name);
191
192 //- Copy constructor with new name - creates bare cloud
193 ThermoCloud
194 (
195 const fvMesh& mesh,
196 const word& name,
197 const ThermoCloud<CloudType>& c
198 );
199
200 //- Construct and return clone based on (this) with new name
202 {
204 (
205 new ThermoCloud(*this, name)
206 );
207 }
208
209 //- Construct and return bare clone based on (this) with new name
210 virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
211 {
213 (
214 new ThermoCloud(this->mesh(), name, *this)
215 );
217
218
219 //- Destructor
220 virtual ~ThermoCloud() = default;
221
222
223 // Member Functions
224
225 // Access
226
227 //- Return a reference to the cloud copy
228 inline const ThermoCloud& cloudCopy() const;
230 //- Return the constant properties
231 inline const typename parcelType::constantProperties&
232 constProps() const;
233
234 //- Return access to the constant properties
235 inline typename parcelType::constantProperties& constProps();
236
237 //- Return const access to thermo package
238 inline const SLGThermo& thermo() const;
239
240 //- Return const access to the carrier temperature field
241 inline const volScalarField& T() const;
242
243 //- Return const access to the carrier pressure field
244 inline const volScalarField& p() const;
245
246
247 // Sub-models
248
249 //- Return reference to heat transfer model
251 heatTransfer() const;
252
253
254 // Integration schemes
256 //-Return reference to velocity integration
257 inline const integrationScheme& TIntegrator() const;
258
259
260 // Modelling options
261
262 //- Radiation flag
263 inline bool radiation() const;
264
265 //- Radiation sum of parcel projected areas [m2]
268 //- Radiation sum of parcel projected areas [m2]
269 inline const volScalarField::Internal&
270 radAreaP() const;
271
272 //- Radiation sum of parcel temperature^4 [K4]
274
275 //- Radiation sum of parcel temperature^4 [K4]
276 inline const volScalarField::Internal& radT4() const;
278 //- Radiation sum of parcel projected area*temperature^4 [m2K4]
280
281 //- Radiation sum of parcel temperature^4 [m2K4]
282 inline const volScalarField::Internal&
283 radAreaPT4() const;
284
285
286 // Sources
287
288 //- Transfer the effect of parcel to the carrier phase
289 inline void transferToCarrier
290 (
291 const parcelType& p,
292 const typename parcelType::trackingData& td
293 );
294
295
296 // Enthalpy
297
298 //- Sensible enthalpy transfer [J/kg]
300
301 //- Sensible enthalpy transfer [J/kg]
302 inline const volScalarField::Internal&
303 hsTrans() const;
304
305 //- Return coefficient for carrier phase hs equation
307
308 //- Return const coefficient for carrier phase hs equation
309 inline const volScalarField::Internal&
310 hsCoeff() const;
311
312 //- Return sensible enthalpy source term [J/kg/m3/s]
313 inline tmp<fvScalarMatrix> Sh(volScalarField& hs) const;
314
315
316 // Radiation - overrides thermoCloud virtual abstract members
317
318 //- Return tmp equivalent particulate emission
319 inline tmp<volScalarField> Ep() const;
321 //- Return tmp equivalent particulate absorption
322 inline tmp<volScalarField> ap() const;
323
324 //- Return tmp equivalent particulate scattering factor
325 inline tmp<volScalarField> sigmap() const;
326
327
328 // Check
329
330 //- Maximum temperature
331 inline scalar Tmax() const;
332
333 //- Minimum temperature
334 inline scalar Tmin() const;
335
336
337 // Cloud evolution functions
338
339 //- Set parcel thermo properties
341 (
342 parcelType& parcel,
343 const scalar lagrangianDt
344 );
345
346 //- Check parcel properties
348 (
349 parcelType& parcel,
350 const scalar lagrangianDt,
351 const bool fullyDescribed
352 );
353
354 //- Store the current cloud state
355 void storeState();
356
357 //- Reset the current cloud to the previously stored state
358 void restoreState();
359
360 //- Reset the cloud source terms
361 void resetSourceTerms();
362
363 //- Apply relaxation to (steady state) cloud sources
364 void relaxSources(const ThermoCloud<CloudType>& cloudOldTime);
365
366 //- Apply scaling to (transient) cloud sources
367 void scaleSources();
369 //- Pre-evolve
370 void preEvolve(const typename parcelType::trackingData& td);
371
372 //- Evolve the cloud
373 void evolve();
374
375
376 // Mapping
377
378 //- Remap the cells of particles corresponding to the
379 // mesh topology change with a default tracking data object
380 virtual void autoMap(const mapPolyMesh&);
381
382
383 // I-O
384
385 //- Print cloud information
386 void info();
387};
388
389
390// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
392} // End namespace Foam
393
394// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395
396#include "ThermoCloudI.H"
398// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399
400#ifdef NoRepository
401 #include "ThermoCloud.C"
402#endif
403
404// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
405
406#endif
407
408// ************************************************************************* //
ParticleType particleType
Definition Cloud.H:130
const word & cloudName() const
const fvMesh & mesh() const
DimensionedField< scalar, volMesh > Internal
Templated class to calculate the fluid-particle heat transfer coefficients based on a specified Nusse...
autoPtr< IOobject > clone() const
Clone.
Definition IOobject.H:641
const volVectorField & U() const
Return carrier gas velocity.
const volScalarField & rho() const
Return carrier gas density.
const dimensionedVector & g() const
Gravity.
const fvMesh & mesh() const
Return reference to the mesh.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition SLGThermo.H:63
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Templated base class for thermodynamic cloud.
Definition ThermoCloud.H:66
const parcelType::constantProperties & constProps() const
Return the constant properties.
volScalarField::Internal & radAreaP()
Radiation sum of parcel projected areas [m2].
void setModels()
Set cloud sub-models.
Definition ThermoCloud.C:30
volScalarField::Internal & hsCoeff()
Return coefficient for carrier phase hs equation.
volScalarField::Internal & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
virtual ~ThermoCloud()=default
Destructor.
void storeState()
Store the current cloud state.
const volScalarField & T() const
Return const access to the carrier temperature field.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
bool radiation() const
Radiation flag.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
ThermoCloud(ThermoCloud< CloudType > &c, const word &name)
Copy constructor with new name.
scalar Tmin() const
Minimum temperature.
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
const volScalarField::Internal & hsTrans() const
Sensible enthalpy transfer [J/kg].
void scaleSources()
Apply scaling to (transient) cloud sources.
KinematicCloud< Cloud< basicThermoParcel > >::particleType parcelType
Definition ThermoCloud.H:79
const integrationScheme & TIntegrator() const
Return reference to velocity integration.
KinematicCloud< Cloud< basicThermoParcel > > cloudType
Definition ThermoCloud.H:74
const ThermoCloud & cloudCopy() const
Return a reference to the cloud copy.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
scalar Tmax() const
Maximum temperature.
tmp< volScalarField > sigmap() const
Return tmp equivalent particulate scattering factor.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
void cloudReset(ThermoCloud< CloudType > &c)
Reset state of cloud.
const volScalarField::Internal & radAreaP() const
Radiation sum of parcel projected areas [m2].
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
tmp< fvScalarMatrix > Sh(volScalarField &hs) const
Return sensible enthalpy source term [J/kg/m3/s].
void evolve()
Evolve the cloud.
autoPtr< HeatTransferModel< ThermoCloud< KinematicCloud< Cloud< basicThermoParcel > > > > > heatTransferModel_
ThermoCloud< KinematicCloud< Cloud< basicThermoParcel > > > thermoCloudType
Definition ThermoCloud.H:84
const volScalarField::Internal & hsCoeff() const
Return const coefficient for carrier phase hs equation.
volScalarField::Internal & hsTrans()
Sensible enthalpy transfer [J/kg].
parcelType::constantProperties & constProps()
Return access to the constant properties.
void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
ThermoCloud(const fvMesh &mesh, const word &name, const ThermoCloud< CloudType > &c)
Copy constructor with new name - creates bare cloud.
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
ThermoCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const SLGThermo &thermo, bool readFields=true)
Construct given carrier gas fields.
const volScalarField::Internal & radT4() const
Radiation sum of parcel temperature^4 [K4].
void resetSourceTerms()
Reset the cloud source terms.
const volScalarField & p() const
Return const access to the carrier pressure field.
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
const volScalarField::Internal & radAreaPT4() const
Radiation sum of parcel temperature^4 [m2K4].
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
thermoCloud()=default
Null constructor.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.