Loading...
Searching...
No Matches
ThermoCloud.C
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-2023 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
27\*---------------------------------------------------------------------------*/
28
29#include "ThermoCloud.H"
31
32#include "HeatTransferModel.H"
33
34// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
35
36template<class CloudType>
38{
40 (
41 HeatTransferModel<ThermoCloud<CloudType>>::New
42 (
43 this->subModelProperties(),
44 *this
45 ).ptr()
46 );
47
48 TIntegrator_.reset
49 (
51 (
52 "T",
54 ).ptr()
55 );
56
57 this->subModelProperties().readEntry("radiation", radiation_);
58
59 if (radiation_)
60 {
61 radAreaP_.reset
62 (
64 (
66 (
67 IOobject::scopedName(this->name(), "radAreaP"),
68 this->db().time().timeName(),
69 this->db(),
73 ),
74 this->mesh(),
76 )
77 );
78
79 radT4_.reset
80 (
82 (
84 (
85 IOobject::scopedName(this->name(), "radT4"),
86 this->db().time().timeName(),
87 this->db(),
91 ),
92 this->mesh(),
94 )
95 );
96
97 radAreaPT4_.reset
98 (
100 (
102 (
103 IOobject::scopedName(this->name(), "radAreaPT4"),
104 this->db().time().timeName(),
105 this->db(),
109 ),
110 this->mesh(),
113 );
114 }
115}
116
117
118template<class CloudType>
119void Foam::ThermoCloud<CloudType>::cloudReset(ThermoCloud<CloudType>& c)
120{
121 CloudType::cloudReset(c);
122
123 heatTransferModel_.reset(c.heatTransferModel_.ptr());
124 TIntegrator_.reset(c.TIntegrator_.ptr());
125
126 radiation_ = c.radiation_;
127}
128
129
130// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
131
132template<class CloudType>
133Foam::ThermoCloud<CloudType>::ThermoCloud
134(
135 const word& cloudName,
136 const volScalarField& rho,
137 const volVectorField& U,
138 const dimensionedVector& g,
139 const SLGThermo& thermo,
140 bool readFields
141)
142:
144 (
145 cloudName,
146 rho,
147 U,
148 thermo.thermo().mu(),
149 g,
150 false
151 ),
152 thermoCloud(),
153 cloudCopyPtr_(nullptr),
154 constProps_(this->particleProperties()),
155 thermo_(thermo),
156 T_(thermo.thermo().T()),
157 p_(thermo.thermo().p()),
158 heatTransferModel_(nullptr),
159 TIntegrator_(nullptr),
160 radiation_(false),
161 radAreaP_(nullptr),
162 radT4_(nullptr),
163 radAreaPT4_(nullptr),
164 hsTrans_
165 (
166 new volScalarField::Internal
167 (
169 (
170 IOobject::scopedName(this->name(), "hsTrans"),
171 this->db().time().timeName(),
172 this->db(),
173 IOobject::READ_IF_PRESENT,
174 IOobject::AUTO_WRITE,
175 IOobject::REGISTER
176 ),
177 this->mesh(),
179 )
180 ),
181 hsCoeff_
182 (
183 new volScalarField::Internal
184 (
186 (
187 IOobject::scopedName(this->name(), "hsCoeff"),
188 this->db().time().timeName(),
189 this->db(),
190 IOobject::READ_IF_PRESENT,
191 IOobject::AUTO_WRITE,
192 IOobject::REGISTER
193 ),
194 this->mesh(),
196 )
197 )
199 if (this->solution().active())
200 {
201 setModels();
202
203 if (readFields)
204 {
205 parcelType::readFields(*this);
207 }
208 }
209
210 if (this->solution().resetSourcesOnStartup())
213 }
214}
215
216
217template<class CloudType>
218Foam::ThermoCloud<CloudType>::ThermoCloud
219(
220 ThermoCloud<CloudType>& c,
221 const word& name
222)
223:
224 CloudType(c, name),
225 thermoCloud(),
226 cloudCopyPtr_(nullptr),
227 constProps_(c.constProps_),
228 thermo_(c.thermo_),
229 T_(c.T()),
230 p_(c.p()),
231 heatTransferModel_(c.heatTransferModel_->clone()),
232 TIntegrator_(c.TIntegrator_->clone()),
233 radiation_(c.radiation_),
234 radAreaP_(nullptr),
235 radT4_(nullptr),
236 radAreaPT4_(nullptr),
237 hsTrans_
238 (
239 new volScalarField::Internal
240 (
241 IOobject
242 (
243 IOobject::scopedName(this->name(), "hsTrans"),
244 this->db().time().timeName(),
245 this->db(),
246 IOobject::NO_READ,
247 IOobject::NO_WRITE,
248 IOobject::NO_REGISTER
249 ),
250 c.hsTrans()
251 )
252 ),
253 hsCoeff_
254 (
255 new volScalarField::Internal
256 (
257 IOobject
258 (
259 IOobject::scopedName(this->name(), "hsCoeff"),
260 this->db().time().timeName(),
261 this->db(),
262 IOobject::NO_READ,
263 IOobject::NO_WRITE,
264 IOobject::NO_REGISTER
265 ),
266 c.hsCoeff()
267 )
268 )
269{
270 if (radiation_)
271 {
272 radAreaP_.reset
273 (
275 (
276 IOobject
277 (
278 IOobject::scopedName(this->name(), "radAreaP"),
279 this->db().time().timeName(),
280 this->db(),
284 ),
285 c.radAreaP()
286 )
287 );
288
289 radT4_.reset
290 (
292 (
293 IOobject
294 (
295 IOobject::scopedName(this->name(), "radT4"),
296 this->db().time().timeName(),
297 this->db(),
301 ),
302 c.radT4()
303 )
304 );
305
306 radAreaPT4_.reset
307 (
309 (
310 IOobject
311 (
312 IOobject::scopedName(this->name(), "radAreaPT4"),
313 this->db().time().timeName(),
314 this->db(),
318 ),
319 c.radAreaPT4()
321 );
322 }
323}
324
325
326template<class CloudType>
327Foam::ThermoCloud<CloudType>::ThermoCloud
328(
329 const fvMesh& mesh,
330 const word& name,
331 const ThermoCloud<CloudType>& c
332)
333:
334 CloudType(mesh, name, c),
335 thermoCloud(),
336 cloudCopyPtr_(nullptr),
337 constProps_(),
338 thermo_(c.thermo()),
339 T_(c.T()),
340 p_(c.p()),
341 heatTransferModel_(nullptr),
342 TIntegrator_(nullptr),
343 radiation_(false),
344 radAreaP_(nullptr),
345 radT4_(nullptr),
346 radAreaPT4_(nullptr),
347 hsTrans_(nullptr),
348 hsCoeff_(nullptr)
349{}
350
351
352// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
353
354template<class CloudType>
356(
357 parcelType& parcel,
358 const scalar lagrangianDt
359)
360{
361 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
363 parcel.T() = constProps_.T0();
364 parcel.Cp() = constProps_.Cp0();
365}
366
367
368template<class CloudType>
370(
371 parcelType& parcel,
372 const scalar lagrangianDt,
373 const bool fullyDescribed
375{
376 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
377}
378
379
380template<class CloudType>
382{
383 cloudCopyPtr_.reset
384 (
385 static_cast<ThermoCloud<CloudType>*>
386 (
387 clone(this->name() + "Copy").ptr()
388 )
389 );
390}
391
392
393template<class CloudType>
396 cloudReset(cloudCopyPtr_());
397 cloudCopyPtr_.clear();
398}
399
400
401template<class CloudType>
403{
404 CloudType::resetSourceTerms();
405 hsTrans_->field() = 0.0;
406 hsCoeff_->field() = 0.0;
407
408 if (radiation_)
409 {
410 radAreaP_->field() = 0.0;
411 radT4_->field() = 0.0;
412 radAreaPT4_->field() = 0.0;
413 }
414}
415
416
417template<class CloudType>
419(
420 const ThermoCloud<CloudType>& cloudOldTime
421)
422{
423 CloudType::relaxSources(cloudOldTime);
424
425 this->relax(hsTrans_(), cloudOldTime.hsTrans(), "h");
426 this->relax(hsCoeff_(), cloudOldTime.hsCoeff(), "h");
427
428 if (radiation_)
429 {
430 this->relax(radAreaP_(), cloudOldTime.radAreaP(), "radiation");
431 this->relax(radT4_(), cloudOldTime.radT4(), "radiation");
432 this->relax(radAreaPT4_(), cloudOldTime.radAreaPT4(), "radiation");
433 }
434}
435
436
437template<class CloudType>
439{
440 CloudType::scaleSources();
442 this->scale(hsTrans_(), "h");
443 this->scale(hsCoeff_(), "h");
444
445 if (radiation_)
446 {
447 this->scale(radAreaP_(), "radiation");
448 this->scale(radT4_(), "radiation");
449 this->scale(radAreaPT4_(), "radiation");
451}
452
453
454template<class CloudType>
456(
457 const typename parcelType::trackingData& td
458)
459{
460 CloudType::preEvolve(td);
461
462 this->pAmbient() = thermo_.thermo().p().average().value();
463}
464
466template<class CloudType>
468{
469 if (this->solution().canEvolve())
471 typename parcelType::trackingData td(*this);
473 this->solve(*this, td);
474 }
476
477
478template<class CloudType>
482
483 this->updateMesh();
484}
486
487template<class CloudType>
489{
491
492 Log_<< " Temperature min/max = " << Tmin() << ", " << Tmax()
493 << endl;
494}
495
496
497// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition Cloud.C:108
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition Cloud.C:297
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition Cloud.C:125
const word & cloudName() const
const IOdictionary & particleProperties() const
const fvMesh & mesh() const
DimensionedField< scalar, volMesh > Internal
Templated class to calculate the fluid-particle heat transfer coefficients based on a specified Nusse...
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition SLGThermo.H:63
Templated base class for thermodynamic cloud.
Definition ThermoCloud.H:66
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].
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.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
const volScalarField & p_
Pressure [Pa].
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
autoPtr< volScalarField::Internal > hsCoeff_
Coefficient for carrier phase hs equation [W/K].
autoPtr< volScalarField::Internal > radT4_
Radiation sum of parcel temperature^4.
void scaleSources()
Apply scaling to (transient) cloud sources.
autoPtr< integrationScheme > TIntegrator_
Temperature integration.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition ThermoCloud.H:79
Switch radiation_
Include radiation.
const SLGThermo & thermo() const
Return const access to thermo package.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
void cloudReset(ThermoCloud< CloudType > &c)
Reset state of cloud.
void evolve()
Evolve the cloud.
autoPtr< HeatTransferModel< ThermoCloud< CloudType > > > heatTransferModel_
Heat transfer model.
parcelType::constantProperties constProps_
Thermo parcel constant properties.
autoPtr< volScalarField::Internal > radAreaPT4_
Radiation sum of parcel projected areas * temperature^4.
volScalarField::Internal & hsTrans()
Sensible enthalpy transfer [J/kg].
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.
void resetSourceTerms()
Reset the cloud source terms.
const SLGThermo & thermo_
SLG thermodynamics package.
const volScalarField & T_
Temperature [K].
const volScalarField & p() const
Return const access to the carrier pressure field.
autoPtr< volScalarField::Internal > hsTrans_
Sensible enthalpy transfer [J/kg].
autoPtr< volScalarField::Internal > radAreaP_
Radiation sum of parcel projected areas.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static autoPtr< integrationScheme > New(const word &phiName, const dictionary &dict)
Select an integration scheme.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const Time & time() const noexcept
Return time registry.
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
Virtual abstract base class for templated ThermoCloud.
Definition thermoCloud.H:47
thermoCloud()=default
Null constructor.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
const volScalarField & T
UEqn relax()
dynamicFvMesh & mesh
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
word timeName
Definition getTimeIndex.H:3
#define Log_
Report write to Foam::Info if the class log switch is true.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
CEqn solve()