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 | 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 \*---------------------------------------------------------------------------*/
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_(NULL),
156  constProps_(this->particleProperties()),
157  thermo_(thermo),
158  T_(thermo.thermo().T()),
159  p_(thermo.thermo().p()),
160  heatTransferModel_(NULL),
161  TIntegrator_(NULL),
162  radiation_(false),
163  radAreaP_(NULL),
164  radT4_(NULL),
165  radAreaPT4_(NULL),
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  }
207  }
208 
209  if (this->solution().resetSourcesOnStartup())
210  {
211  resetSourceTerms();
212  }
213 }
214 
215 
216 template<class CloudType>
218 (
220  const word& name
221 )
222 :
223  CloudType(c, name),
224  thermoCloud(),
225  cloudCopyPtr_(NULL),
226  constProps_(c.constProps_),
227  thermo_(c.thermo_),
228  T_(c.T()),
229  p_(c.p()),
230  heatTransferModel_(c.heatTransferModel_->clone()),
231  TIntegrator_(c.TIntegrator_->clone()),
232  radiation_(c.radiation_),
233  radAreaP_(NULL),
234  radT4_(NULL),
235  radAreaPT4_(NULL),
236  hsTrans_
237  (
239  (
240  IOobject
241  (
242  this->name() + ":hsTrans",
243  this->db().time().timeName(),
244  this->db(),
245  IOobject::NO_READ,
246  IOobject::NO_WRITE,
247  false
248  ),
249  c.hsTrans()
250  )
251  ),
252  hsCoeff_
253  (
255  (
256  IOobject
257  (
258  this->name() + ":hsCoeff",
259  this->db().time().timeName(),
260  this->db(),
261  IOobject::NO_READ,
262  IOobject::NO_WRITE,
263  false
264  ),
265  c.hsCoeff()
266  )
267  )
268 {
269  if (radiation_)
270  {
271  radAreaP_.reset
272  (
274  (
275  IOobject
276  (
277  this->name() + ":radAreaP",
278  this->db().time().timeName(),
279  this->db(),
280  IOobject::NO_READ,
281  IOobject::NO_WRITE,
282  false
283  ),
284  c.radAreaP()
285  )
286  );
287 
288  radT4_.reset
289  (
291  (
292  IOobject
293  (
294  this->name() + ":radT4",
295  this->db().time().timeName(),
296  this->db(),
297  IOobject::NO_READ,
298  IOobject::NO_WRITE,
299  false
300  ),
301  c.radT4()
302  )
303  );
304 
305  radAreaPT4_.reset
306  (
308  (
309  IOobject
310  (
311  this->name() + ":radAreaPT4",
312  this->db().time().timeName(),
313  this->db(),
314  IOobject::NO_READ,
315  IOobject::NO_WRITE,
316  false
317  ),
318  c.radAreaPT4()
319  )
320  );
321  }
322 }
323 
324 
325 template<class CloudType>
327 (
328  const fvMesh& mesh,
329  const word& name,
330  const ThermoCloud<CloudType>& c
331 )
332 :
333  CloudType(mesh, name, c),
334  thermoCloud(),
335  cloudCopyPtr_(NULL),
336  constProps_(),
337  thermo_(c.thermo()),
338  T_(c.T()),
339  p_(c.p()),
340  heatTransferModel_(NULL),
341  TIntegrator_(NULL),
342  radiation_(false),
343  radAreaP_(NULL),
344  radT4_(NULL),
345  radAreaPT4_(NULL),
346  hsTrans_(NULL),
347  hsCoeff_(NULL)
348 {}
349 
350 
351 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
352 
353 template<class CloudType>
355 {}
356 
357 
358 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
359 
360 template<class CloudType>
362 (
363  parcelType& parcel,
364  const scalar lagrangianDt
365 )
366 {
367  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
368 
369  parcel.T() = constProps_.T0();
370  parcel.Cp() = constProps_.Cp0();
371 }
372 
373 
374 template<class CloudType>
376 (
377  parcelType& parcel,
378  const scalar lagrangianDt,
379  const bool fullyDescribed
380 )
381 {
382  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
383 }
384 
385 
386 template<class CloudType>
388 {
389  cloudCopyPtr_.reset
390  (
391  static_cast<ThermoCloud<CloudType>*>
392  (
393  clone(this->name() + "Copy").ptr()
394  )
395  );
396 }
397 
398 
399 template<class CloudType>
401 {
402  cloudReset(cloudCopyPtr_());
403  cloudCopyPtr_.clear();
404 }
405 
406 
407 template<class CloudType>
409 {
410  CloudType::resetSourceTerms();
411  hsTrans_->field() = 0.0;
412  hsCoeff_->field() = 0.0;
413 
414  if (radiation_)
415  {
416  radAreaP_->field() = 0.0;
417  radT4_->field() = 0.0;
418  radAreaPT4_->field() = 0.0;
419  }
420 }
421 
422 
423 template<class CloudType>
425 (
426  const ThermoCloud<CloudType>& cloudOldTime
427 )
428 {
429  CloudType::relaxSources(cloudOldTime);
430 
431  this->relax(hsTrans_(), cloudOldTime.hsTrans(), "h");
432  this->relax(hsCoeff_(), cloudOldTime.hsCoeff(), "h");
433 
434  if (radiation_)
435  {
436  this->relax(radAreaP_(), cloudOldTime.radAreaP(), "radiation");
437  this->relax(radT4_(), cloudOldTime.radT4(), "radiation");
438  this->relax(radAreaPT4_(), cloudOldTime.radAreaPT4(), "radiation");
439  }
440 }
441 
442 
443 template<class CloudType>
445 {
446  CloudType::scaleSources();
447 
448  this->scale(hsTrans_(), "h");
449  this->scale(hsCoeff_(), "h");
450 
451  if (radiation_)
452  {
453  this->scale(radAreaP_(), "radiation");
454  this->scale(radT4_(), "radiation");
455  this->scale(radAreaPT4_(), "radiation");
456  }
457 }
458 
459 
460 template<class CloudType>
462 {
463  CloudType::preEvolve();
464 
465  this->pAmbient() = thermo_.thermo().p().average().value();
466 }
467 
468 
469 template<class CloudType>
471 {
472  if (this->solution().canEvolve())
473  {
474  typename parcelType::template
475  TrackingData<ThermoCloud<CloudType>> td(*this);
476 
477  this->solve(td);
478  }
479 }
480 
481 
482 template<class CloudType>
484 {
485  typedef typename particle::TrackingData<ThermoCloud<CloudType>> tdType;
486 
487  tdType td(*this);
488 
489  Cloud<parcelType>::template autoMap<tdType>(td, mapper);
490 
491  this->updateMesh();
492 }
493 
494 
495 template<class CloudType>
497 {
498  CloudType::info();
499 
500  Info<< " Temperature min/max = " << Tmin() << ", " << Tmax()
501  << endl;
502 }
503 
504 
505 // ************************************************************************* //
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:362
Templated heat transfer model class.
Definition: ThermoCloud.H:53
DimensionedField< scalar, volMesh > & hsTrans()
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:192
UEqn relax()
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void restoreState()
Reset the current cloud to the previously stored state.
Definition: ThermoCloud.C:400
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const SLGThermo & thermo() const
Return const access to thermo package.
Definition: ThermoCloudI.H:57
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
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< 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:466
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
dynamicFvMesh & mesh
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
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:217
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:478
const dimensionedVector & g
Virtual abstract base class for templated ThermoCloud.
Definition: thermoCloud.H:48
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
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
void evolve()
Evolve the cloud.
Definition: ThermoCloud.C:470
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:496
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:444
messageStream Info
const volScalarField & p() const
Return const access to the carrier prressure field.
Definition: ThermoCloudI.H:71
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
DimensionedField< scalar, volMesh > & 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:91
const fluidThermo & thermo() const
Return reference to the thermo database.
Definition: SLGThermo.C:102
void storeState()
Store the current cloud state.
Definition: ThermoCloud.C:387
void setModels()
Set cloud sub-models.
Definition: ThermoCloud.C:34
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: ThermoCloud.C:425
virtual ~ThermoCloud()
Destructor.
Definition: ThermoCloud.C:354