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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "ThermoCloud.H"
29 #include "ThermoParcel.H"
30 
31 #include "HeatTransferModel.H"
32 
33 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
34 
35 template<class CloudType>
37 {
38  heatTransferModel_.reset
39  (
41  (
42  this->subModelProperties(),
43  *this
44  ).ptr()
45  );
46 
47  TIntegrator_.reset
48  (
50  (
51  "T",
52  this->solution().integrationSchemes()
53  ).ptr()
54  );
55 
56  this->subModelProperties().readEntry("radiation", radiation_);
57 
58  if (radiation_)
59  {
60  radAreaP_.reset
61  (
63  (
64  IOobject
65  (
66  this->name() + ":radAreaP",
67  this->db().time().timeName(),
68  this->db(),
69  IOobject::READ_IF_PRESENT,
70  IOobject::AUTO_WRITE
71  ),
72  this->mesh(),
74  )
75  );
76 
77  radT4_.reset
78  (
80  (
81  IOobject
82  (
83  this->name() + ":radT4",
84  this->db().time().timeName(),
85  this->db(),
86  IOobject::READ_IF_PRESENT,
87  IOobject::AUTO_WRITE
88  ),
89  this->mesh(),
91  )
92  );
93 
94  radAreaPT4_.reset
95  (
97  (
98  IOobject
99  (
100  this->name() + ":radAreaPT4",
101  this->db().time().timeName(),
102  this->db(),
103  IOobject::READ_IF_PRESENT,
104  IOobject::AUTO_WRITE
105  ),
106  this->mesh(),
108  )
109  );
110  }
111 }
112 
113 
114 template<class CloudType>
116 {
117  CloudType::cloudReset(c);
118 
119  heatTransferModel_.reset(c.heatTransferModel_.ptr());
120  TIntegrator_.reset(c.TIntegrator_.ptr());
121 
122  radiation_ = c.radiation_;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
128 template<class CloudType>
130 (
131  const word& cloudName,
132  const volScalarField& rho,
133  const volVectorField& U,
134  const dimensionedVector& g,
135  const SLGThermo& thermo,
136  bool readFields
137 )
138 :
139  CloudType
140  (
141  cloudName,
142  rho,
143  U,
144  thermo.thermo().mu(),
145  g,
146  false
147  ),
148  thermoCloud(),
149  cloudCopyPtr_(nullptr),
150  constProps_(this->particleProperties()),
151  thermo_(thermo),
152  T_(thermo.thermo().T()),
153  p_(thermo.thermo().p()),
154  heatTransferModel_(nullptr),
155  TIntegrator_(nullptr),
156  radiation_(false),
157  radAreaP_(nullptr),
158  radT4_(nullptr),
159  radAreaPT4_(nullptr),
160  hsTrans_
161  (
163  (
164  IOobject
165  (
166  this->name() + ":hsTrans",
167  this->db().time().timeName(),
168  this->db(),
169  IOobject::READ_IF_PRESENT,
170  IOobject::AUTO_WRITE
171  ),
172  this->mesh(),
174  )
175  ),
176  hsCoeff_
177  (
179  (
180  IOobject
181  (
182  this->name() + ":hsCoeff",
183  this->db().time().timeName(),
184  this->db(),
185  IOobject::READ_IF_PRESENT,
186  IOobject::AUTO_WRITE
187  ),
188  this->mesh(),
190  )
191  )
192 {
193  if (this->solution().active())
194  {
195  setModels();
196 
197  if (readFields)
198  {
199  parcelType::readFields(*this);
200  this->deleteLostParticles();
201  }
202  }
203 
204  if (this->solution().resetSourcesOnStartup())
205  {
206  resetSourceTerms();
207  }
208 }
209 
210 
211 template<class CloudType>
213 (
215  const word& name
216 )
217 :
218  CloudType(c, name),
219  thermoCloud(),
220  cloudCopyPtr_(nullptr),
221  constProps_(c.constProps_),
222  thermo_(c.thermo_),
223  T_(c.T()),
224  p_(c.p()),
225  heatTransferModel_(c.heatTransferModel_->clone()),
226  TIntegrator_(c.TIntegrator_->clone()),
227  radiation_(c.radiation_),
228  radAreaP_(nullptr),
229  radT4_(nullptr),
230  radAreaPT4_(nullptr),
231  hsTrans_
232  (
234  (
235  IOobject
236  (
237  this->name() + ":hsTrans",
238  this->db().time().timeName(),
239  this->db(),
240  IOobject::NO_READ,
241  IOobject::NO_WRITE,
242  false
243  ),
244  c.hsTrans()
245  )
246  ),
247  hsCoeff_
248  (
250  (
251  IOobject
252  (
253  this->name() + ":hsCoeff",
254  this->db().time().timeName(),
255  this->db(),
256  IOobject::NO_READ,
257  IOobject::NO_WRITE,
258  false
259  ),
260  c.hsCoeff()
261  )
262  )
263 {
264  if (radiation_)
265  {
266  radAreaP_.reset
267  (
269  (
270  IOobject
271  (
272  this->name() + ":radAreaP",
273  this->db().time().timeName(),
274  this->db(),
275  IOobject::NO_READ,
276  IOobject::NO_WRITE,
277  false
278  ),
279  c.radAreaP()
280  )
281  );
282 
283  radT4_.reset
284  (
286  (
287  IOobject
288  (
289  this->name() + ":radT4",
290  this->db().time().timeName(),
291  this->db(),
292  IOobject::NO_READ,
293  IOobject::NO_WRITE,
294  false
295  ),
296  c.radT4()
297  )
298  );
299 
300  radAreaPT4_.reset
301  (
303  (
304  IOobject
305  (
306  this->name() + ":radAreaPT4",
307  this->db().time().timeName(),
308  this->db(),
309  IOobject::NO_READ,
310  IOobject::NO_WRITE,
311  false
312  ),
313  c.radAreaPT4()
314  )
315  );
316  }
317 }
318 
319 
320 template<class CloudType>
322 (
323  const fvMesh& mesh,
324  const word& name,
326 )
327 :
328  CloudType(mesh, name, c),
329  thermoCloud(),
330  cloudCopyPtr_(nullptr),
331  constProps_(),
332  thermo_(c.thermo()),
333  T_(c.T()),
334  p_(c.p()),
335  heatTransferModel_(nullptr),
336  TIntegrator_(nullptr),
337  radiation_(false),
338  radAreaP_(nullptr),
339  radT4_(nullptr),
340  radAreaPT4_(nullptr),
341  hsTrans_(nullptr),
342  hsCoeff_(nullptr)
343 {}
344 
345 
346 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
347 
348 template<class CloudType>
350 (
351  parcelType& parcel,
352  const scalar lagrangianDt
353 )
354 {
355  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
356 
357  parcel.T() = constProps_.T0();
358  parcel.Cp() = constProps_.Cp0();
359 }
360 
361 
362 template<class CloudType>
364 (
365  parcelType& parcel,
366  const scalar lagrangianDt,
367  const bool fullyDescribed
368 )
369 {
370  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
371 }
372 
373 
374 template<class CloudType>
376 {
377  cloudCopyPtr_.reset
378  (
379  static_cast<ThermoCloud<CloudType>*>
380  (
381  clone(this->name() + "Copy").ptr()
382  )
383  );
384 }
385 
386 
387 template<class CloudType>
389 {
390  cloudReset(cloudCopyPtr_());
391  cloudCopyPtr_.clear();
392 }
393 
394 
395 template<class CloudType>
397 {
398  CloudType::resetSourceTerms();
399  hsTrans_->field() = 0.0;
400  hsCoeff_->field() = 0.0;
401 
402  if (radiation_)
403  {
404  radAreaP_->field() = 0.0;
405  radT4_->field() = 0.0;
406  radAreaPT4_->field() = 0.0;
407  }
408 }
409 
410 
411 template<class CloudType>
413 (
414  const ThermoCloud<CloudType>& cloudOldTime
415 )
416 {
417  CloudType::relaxSources(cloudOldTime);
418 
419  this->relax(hsTrans_(), cloudOldTime.hsTrans(), "h");
420  this->relax(hsCoeff_(), cloudOldTime.hsCoeff(), "h");
421 
422  if (radiation_)
423  {
424  this->relax(radAreaP_(), cloudOldTime.radAreaP(), "radiation");
425  this->relax(radT4_(), cloudOldTime.radT4(), "radiation");
426  this->relax(radAreaPT4_(), cloudOldTime.radAreaPT4(), "radiation");
427  }
428 }
429 
430 
431 template<class CloudType>
433 {
434  CloudType::scaleSources();
435 
436  this->scale(hsTrans_(), "h");
437  this->scale(hsCoeff_(), "h");
438 
439  if (radiation_)
440  {
441  this->scale(radAreaP_(), "radiation");
442  this->scale(radT4_(), "radiation");
443  this->scale(radAreaPT4_(), "radiation");
444  }
445 }
446 
447 
448 template<class CloudType>
450 {
451  CloudType::preEvolve();
452 
453  this->pAmbient() = thermo_.thermo().p().average().value();
454 }
455 
456 
457 template<class CloudType>
459 {
460  if (this->solution().canEvolve())
461  {
462  typename parcelType::trackingData td(*this);
463 
464  this->solve(*this, td);
465  }
466 }
467 
468 
469 template<class CloudType>
471 {
473 
474  this->updateMesh();
475 }
476 
477 
478 template<class CloudType>
480 {
481  CloudType::info();
482 
483  Info<< " Temperature min/max = " << Tmin() << ", " << Tmax()
484  << endl;
485 }
486 
487 
488 // ************************************************************************* //
Foam::ThermoCloud::evolve
void evolve()
Evolve the cloud.
Definition: ThermoCloud.C:458
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::thermoCloud
Virtual abstract base class for templated ThermoCloud.
Definition: thermoCloud.H:50
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:50
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::ThermoCloud::scaleSources
void scaleSources()
Apply scaling to (transient) cloud sources.
Definition: ThermoCloud.C:432
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Cloud::parcelType
ParticleType parcelType
Parcels are just particles.
Definition: Cloud.H:116
rho
rho
Definition: readInitialConditions.H:96
Foam::dimensioned::T
dimensioned< Type > T() const
Return transpose.
Definition: dimensionedSphericalTensor.C:38
solve
CEqn solve()
Foam::HeatTransferModel
Templated heat transfer model class.
Definition: ThermoCloud.H:58
Foam::ThermoCloud::restoreState
void restoreState()
Reset the current cloud to the previously stored state.
Definition: ThermoCloud.C:388
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::ThermoCloud::cloudReset
void cloudReset(ThermoCloud< CloudType > &c)
Reset state of cloud.
Definition: ThermoCloud.C:115
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::ThermoCloud::radAreaPT4
volScalarField::Internal & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
Definition: ThermoCloudI.H:165
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::ThermoCloud::hsCoeff
volScalarField::Internal & hsCoeff()
Return coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:211
Foam::ThermoCloud::resetSourceTerms
void resetSourceTerms()
Reset the cloud source terms.
Definition: ThermoCloud.C:396
Foam::ThermoCloud::setParcelThermoProperties
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: ThermoCloud.C:350
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:38
timeName
word timeName
Definition: getTimeIndex.H:3
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< vector >
ThermoCloud.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam::ThermoCloud::autoMap
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: ThermoCloud.C:470
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::ThermoCloud::checkParcelProperties
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: ThermoCloud.C:364
Foam::ThermoCloud::setModels
void setModels()
Set cloud sub-models.
Definition: ThermoCloud.C:36
Foam::ThermoCloud::info
void info()
Print cloud information.
Definition: ThermoCloud.C:479
Foam::ThermoCloud
Templated base class for thermodynamic cloud.
Definition: ThermoCloud.H:65
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:54
HeatTransferModel.H
Foam::ThermoCloud::radT4
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
Definition: ThermoCloudI.H:135
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::ThermoCloud::relaxSources
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Definition: ThermoCloud.C:413
ThermoParcel.H
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::ThermoCloud::storeState
void storeState()
Store the current cloud state.
Definition: ThermoCloud.C:375
Foam::ThermoCloud::preEvolve
void preEvolve()
Pre-evolve.
Definition: ThermoCloud.C:449
Foam::ThermoCloud::radAreaP
volScalarField::Internal & radAreaP()
Radiation sum of parcel projected areas [m2].
Definition: ThermoCloudI.H:105
Foam::ThermoCloud::hsTrans
volScalarField::Internal & hsTrans()
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:195
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
relax
UEqn relax()