ThermoCloud.C
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-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "ThermoCloud.H"
27 #include "ThermoParcel.H"
28 
29 #include "HeatTransferModel.H"
30 
31 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
36  heatTransferModel_.reset
37  (
39  (
40  this->subModelProperties(),
41  *this
42  ).ptr()
43  );
44 
45  TIntegrator_.reset
46  (
48  (
49  "T",
50  this->solution().integrationSchemes()
51  ).ptr()
52  );
53 
54  if (this->solution().coupled())
55  {
56  this->subModelProperties().lookup("radiation") >> radiation_;
57  }
58 
59  if (radiation_)
60  {
61  radAreaP_.reset
62  (
64  (
65  IOobject
66  (
67  this->name() + ":radAreaP",
68  this->db().time().timeName(),
69  this->db(),
70  IOobject::READ_IF_PRESENT,
71  IOobject::AUTO_WRITE
72  ),
73  this->mesh(),
74  dimensionedScalar("zero", dimArea, 0.0)
75  )
76  );
77 
78  radT4_.reset
79  (
81  (
82  IOobject
83  (
84  this->name() + ":radT4",
85  this->db().time().timeName(),
86  this->db(),
87  IOobject::READ_IF_PRESENT,
88  IOobject::AUTO_WRITE
89  ),
90  this->mesh(),
92  )
93  );
94 
95  radAreaPT4_.reset
96  (
98  (
99  IOobject
100  (
101  this->name() + ":radAreaPT4",
102  this->db().time().timeName(),
103  this->db(),
104  IOobject::READ_IF_PRESENT,
105  IOobject::AUTO_WRITE
106  ),
107  this->mesh(),
109  (
110  "zero",
112  0.0
113  )
114  )
115  );
116  }
117 }
118 
119 
120 template<class CloudType>
122 {
123  CloudType::cloudReset(c);
124 
125  heatTransferModel_.reset(c.heatTransferModel_.ptr());
126  TIntegrator_.reset(c.TIntegrator_.ptr());
127 
128  radiation_ = c.radiation_;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
134 template<class CloudType>
136 (
137  const word& cloudName,
138  const volScalarField& rho,
139  const volVectorField& U,
140  const dimensionedVector& g,
141  const SLGThermo& thermo,
142  bool readFields
143 )
144 :
145  CloudType
146  (
147  cloudName,
148  rho,
149  U,
150  thermo.thermo().mu(),
151  g,
152  false
153  ),
154  thermoCloud(),
155  cloudCopyPtr_(nullptr),
156  constProps_(this->particleProperties()),
157  thermo_(thermo),
158  T_(thermo.thermo().T()),
159  p_(thermo.thermo().p()),
160  heatTransferModel_(nullptr),
161  TIntegrator_(nullptr),
162  radiation_(false),
163  radAreaP_(nullptr),
164  radT4_(nullptr),
165  radAreaPT4_(nullptr),
166  hsTrans_
167  (
169  (
170  IOobject
171  (
172  this->name() + ":hsTrans",
173  this->db().time().timeName(),
174  this->db(),
175  IOobject::READ_IF_PRESENT,
176  IOobject::AUTO_WRITE
177  ),
178  this->mesh(),
179  dimensionedScalar("zero", dimEnergy, 0.0)
180  )
181  ),
182  hsCoeff_
183  (
185  (
186  IOobject
187  (
188  this->name() + ":hsCoeff",
189  this->db().time().timeName(),
190  this->db(),
191  IOobject::READ_IF_PRESENT,
192  IOobject::AUTO_WRITE
193  ),
194  this->mesh(),
196  )
197  )
198 {
199  if (this->solution().active())
200  {
201  setModels();
202 
203  if (readFields)
204  {
205  parcelType::readFields(*this);
206  this->deleteLostParticles();
207  }
208  }
209 
210  if (this->solution().resetSourcesOnStartup())
211  {
212  resetSourceTerms();
213  }
214 }
215 
216 
217 template<class CloudType>
219 (
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  (
240  (
241  IOobject
242  (
243  this->name() + ":hsTrans",
244  this->db().time().timeName(),
245  this->db(),
246  IOobject::NO_READ,
247  IOobject::NO_WRITE,
248  false
249  ),
250  c.hsTrans()
251  )
252  ),
253  hsCoeff_
254  (
256  (
257  IOobject
258  (
259  this->name() + ":hsCoeff",
260  this->db().time().timeName(),
261  this->db(),
262  IOobject::NO_READ,
263  IOobject::NO_WRITE,
264  false
265  ),
266  c.hsCoeff()
267  )
268  )
269 {
270  if (radiation_)
271  {
272  radAreaP_.reset
273  (
275  (
276  IOobject
277  (
278  this->name() + ":radAreaP",
279  this->db().time().timeName(),
280  this->db(),
281  IOobject::NO_READ,
282  IOobject::NO_WRITE,
283  false
284  ),
285  c.radAreaP()
286  )
287  );
288 
289  radT4_.reset
290  (
292  (
293  IOobject
294  (
295  this->name() + ":radT4",
296  this->db().time().timeName(),
297  this->db(),
298  IOobject::NO_READ,
299  IOobject::NO_WRITE,
300  false
301  ),
302  c.radT4()
303  )
304  );
305 
306  radAreaPT4_.reset
307  (
309  (
310  IOobject
311  (
312  this->name() + ":radAreaPT4",
313  this->db().time().timeName(),
314  this->db(),
315  IOobject::NO_READ,
316  IOobject::NO_WRITE,
317  false
318  ),
319  c.radAreaPT4()
320  )
321  );
322  }
323 }
324 
325 
326 template<class CloudType>
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 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
353 
354 template<class CloudType>
356 {}
357 
358 
359 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
360 
361 template<class CloudType>
363 (
364  parcelType& parcel,
365  const scalar lagrangianDt
366 )
367 {
368  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
369 
370  parcel.T() = constProps_.T0();
371  parcel.Cp() = constProps_.Cp0();
372 }
373 
374 
375 template<class CloudType>
377 (
378  parcelType& parcel,
379  const scalar lagrangianDt,
380  const bool fullyDescribed
381 )
382 {
383  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
384 }
385 
386 
387 template<class CloudType>
389 {
390  cloudCopyPtr_.reset
391  (
392  static_cast<ThermoCloud<CloudType>*>
393  (
394  clone(this->name() + "Copy").ptr()
395  )
396  );
397 }
398 
399 
400 template<class CloudType>
402 {
403  cloudReset(cloudCopyPtr_());
404  cloudCopyPtr_.clear();
405 }
406 
407 
408 template<class CloudType>
410 {
411  CloudType::resetSourceTerms();
412  hsTrans_->field() = 0.0;
413  hsCoeff_->field() = 0.0;
414 
415  if (radiation_)
416  {
417  radAreaP_->field() = 0.0;
418  radT4_->field() = 0.0;
419  radAreaPT4_->field() = 0.0;
420  }
421 }
422 
423 
424 template<class CloudType>
426 (
427  const ThermoCloud<CloudType>& cloudOldTime
428 )
429 {
430  CloudType::relaxSources(cloudOldTime);
431 
432  this->relax(hsTrans_(), cloudOldTime.hsTrans(), "h");
433  this->relax(hsCoeff_(), cloudOldTime.hsCoeff(), "h");
434 
435  if (radiation_)
436  {
437  this->relax(radAreaP_(), cloudOldTime.radAreaP(), "radiation");
438  this->relax(radT4_(), cloudOldTime.radT4(), "radiation");
439  this->relax(radAreaPT4_(), cloudOldTime.radAreaPT4(), "radiation");
440  }
441 }
442 
443 
444 template<class CloudType>
446 {
447  CloudType::scaleSources();
448 
449  this->scale(hsTrans_(), "h");
450  this->scale(hsCoeff_(), "h");
451 
452  if (radiation_)
453  {
454  this->scale(radAreaP_(), "radiation");
455  this->scale(radT4_(), "radiation");
456  this->scale(radAreaPT4_(), "radiation");
457  }
458 }
459 
460 
461 template<class CloudType>
463 {
464  CloudType::preEvolve();
465 
466  this->pAmbient() = thermo_.thermo().p().average().value();
467 }
468 
469 
470 template<class CloudType>
472 {
473  if (this->solution().canEvolve())
474  {
475  typename parcelType::trackingData td(*this);
476 
477  this->solve(*this, td);
478  }
479 }
480 
481 
482 template<class CloudType>
484 {
486 
487  this->updateMesh();
488 }
489 
490 
491 template<class CloudType>
493 {
494  CloudType::info();
495 
496  Info<< " Temperature min/max = " << Tmin() << ", " << Tmax()
497  << endl;
498 }
499 
500 
501 // ************************************************************************* //
autoPtr< integrationScheme > TIntegrator_
Temperature integration.
Definition: ThermoCloud.H:126
Switch radiation_
Include radiation.
Definition: ThermoCloud.H:132
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
DSMCCloud< dsmcParcel > CloudType
const SLGThermo & thermo_
SLG thermodynamics package.
Definition: ThermoCloud.H:107
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: ThermoCloud.C:363
Templated heat transfer model class.
Definition: ThermoCloud.H:53
volScalarField::Internal & hsTrans()
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:192
UEqn relax()
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void restoreState()
Reset the current cloud to the previously stored state.
Definition: ThermoCloud.C:401
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const volScalarField & p() const
Return const access to the carrier prressure field.
Definition: ThermoCloudI.H:71
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:499
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:121
void preEvolve()
Pre-evolve.
Definition: ThermoCloud.C:462
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
parcelType::constantProperties constProps_
Thermo parcel constant properties.
Definition: ThermoCloud.H:101
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:487
volScalarField::Internal & radAreaP()
Radiation sum of parcel projected areas [m2].
Definition: ThermoCloudI.H:102
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
Definition: ThermoCloudI.H:132
dynamicFvMesh & mesh
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: ThermoCloud.C:483
A class for handling words, derived from string.
Definition: word.H:59
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
rhoEqn solve()
word timeName
Definition: getTimeIndex.H:3
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
Base cloud calls templated on particle type.
Definition: Cloud.H:52
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:218
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: ThermoCloud.C:377
void resetSourceTerms()
Reset the cloud source terms.
Definition: ThermoCloud.C:409
void evolve()
Evolve the cloud.
Definition: ThermoCloud.C:471
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void info()
Print cloud information.
Definition: ThermoCloud.C:492
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...
dimensionedScalar pow4(const dimensionedScalar &ds)
void scaleSources()
Apply scaling to (transient) cloud sources.
Definition: ThermoCloud.C:445
messageStream Info
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
autoPtr< HeatTransferModel< ThermoCloud< CloudType > > > heatTransferModel_
Heat transfer model.
Definition: ThermoCloud.H:120
const dimensionedVector & g
volScalarField::Internal & hsCoeff()
Return coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:208
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void storeState()
Store the current cloud state.
Definition: ThermoCloud.C:388
void setModels()
Set cloud sub-models.
Definition: ThermoCloud.C:34
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
const fluidThermo & thermo() const
Return reference to the thermo database.
Definition: SLGThermo.C:102
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: ThermoCloud.C:426
virtual ~ThermoCloud()
Destructor.
Definition: ThermoCloud.C:355