ThermoCloudI.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 | Copyright (C) 2011-2015 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 
27 
28 using namespace Foam::constant;
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class CloudType>
33 inline const Foam::ThermoCloud<CloudType>&
35 {
36  return cloudCopyPtr_();
37 }
38 
39 
40 template<class CloudType>
41 inline const typename CloudType::particleType::constantProperties&
43 {
44  return constProps_;
45 }
46 
47 
48 template<class CloudType>
49 inline typename CloudType::particleType::constantProperties&
51 {
52  return constProps_;
53 }
54 
55 
56 template<class CloudType>
58 {
59  return thermo_;
60 }
61 
62 
63 template<class CloudType>
65 {
66  return T_;
67 }
68 
69 
70 template<class CloudType>
72 {
73  return p_;
74 }
75 
76 
77 template<class CloudType>
80 {
81  return heatTransferModel_;
82 }
83 
84 
85 template<class CloudType>
88 {
89  return TIntegrator_;
90 }
91 
92 
93 template<class CloudType>
95 {
96  return radiation_;
97 }
98 
99 
100 template<class CloudType>
103 {
104  if (!radiation_)
105  {
107  << "Radiation field requested, but radiation model not active"
108  << abort(FatalError);
109  }
110 
111  return radAreaP_();
112 }
113 
114 
115 template<class CloudType>
118 {
119  if (!radiation_)
120  {
122  << "Radiation field requested, but radiation model not active"
123  << abort(FatalError);
124  }
125 
126  return radAreaP_();
127 }
128 
129 
130 template<class CloudType>
133 {
134  if (!radiation_)
135  {
137  << "Radiation field requested, but radiation model not active"
138  << abort(FatalError);
139  }
140 
141  return radT4_();
142 }
143 
144 
145 template<class CloudType>
148 {
149  if (!radiation_)
150  {
152  << "Radiation field requested, but radiation model not active"
153  << abort(FatalError);
154  }
155 
156  return radT4_();
157 }
158 
159 
160 template<class CloudType>
163 {
164  if (!radiation_)
165  {
167  << "Radiation field requested, but radiation model not active"
168  << abort(FatalError);
169  }
170 
171  return radAreaPT4_();
172 }
173 
174 
175 template<class CloudType>
178 {
179  if (!radiation_)
180  {
182  << "Radiation field requested, but radiation model not active"
183  << abort(FatalError);
184  }
185 
186  return radAreaPT4_();
187 }
188 
189 
190 template<class CloudType>
193 {
194  return hsTrans_();
195 }
196 
197 
198 template<class CloudType>
201 {
202  return hsTrans_();
203 }
204 
205 
206 template<class CloudType>
209 {
210  return hsCoeff_();
211 }
212 
213 
214 template<class CloudType>
217 {
218  return hsCoeff_();
219 }
220 
221 
222 template<class CloudType>
225 {
226  if (debug)
227  {
228  Info<< "hsTrans min/max = " << min(hsTrans()).value() << ", "
229  << max(hsTrans()).value() << nl
230  << "hsCoeff min/max = " << min(hsCoeff()).value() << ", "
231  << max(hsCoeff()).value() << endl;
232  }
233 
234  if (this->solution().coupled())
235  {
236  if (this->solution().semiImplicit("h"))
237  {
238  const volScalarField Cp(thermo_.thermo().Cp());
240  Vdt(this->mesh().V()*this->db().time().deltaT());
241 
242  return
243  hsTrans()/Vdt
244  - fvm::SuSp(hsCoeff()/(Cp*Vdt), hs)
245  + hsCoeff()/(Cp*Vdt)*hs;
246  }
247  else
248  {
250  fvScalarMatrix& fvm = tfvm();
251 
252  fvm.source() = -hsTrans()/(this->db().time().deltaT());
253 
254  return tfvm;
255  }
256  }
257 
259 }
260 
261 
262 template<class CloudType>
264 {
266  (
267  new volScalarField
268  (
269  IOobject
270  (
271  this->name() + ":radiation:Ep",
272  this->db().time().timeName(),
273  this->db(),
276  false
277  ),
278  this->mesh(),
280  )
281  );
282 
283  if (radiation_)
284  {
285  scalarField& Ep = tEp().internalField();
286  const scalar dt = this->db().time().deltaTValue();
287  const scalarField& V = this->mesh().V();
288  const scalar epsilon = constProps_.epsilon0();
289  const scalarField& sumAreaPT4 = radAreaPT4_->field();
290 
291  Ep = sumAreaPT4*epsilon*physicoChemical::sigma.value()/V/dt;
292  }
293 
294  return tEp;
295 }
296 
297 
298 template<class CloudType>
300 {
302  (
303  new volScalarField
304  (
305  IOobject
306  (
307  this->name() + ":radiation:ap",
308  this->db().time().timeName(),
309  this->db(),
312  false
313  ),
314  this->mesh(),
315  dimensionedScalar("zero", dimless/dimLength, 0.0)
316  )
317  );
318 
319  if (radiation_)
320  {
321  scalarField& ap = tap().internalField();
322  const scalar dt = this->db().time().deltaTValue();
323  const scalarField& V = this->mesh().V();
324  const scalar epsilon = constProps_.epsilon0();
325  const scalarField& sumAreaP = radAreaP_->field();
326 
327  ap = sumAreaP*epsilon/V/dt;
328  }
329 
330  return tap;
331 }
332 
333 
334 template<class CloudType>
337 {
338  tmp<volScalarField> tsigmap
339  (
340  new volScalarField
341  (
342  IOobject
343  (
344  this->name() + ":radiation:sigmap",
345  this->db().time().timeName(),
346  this->db(),
349  false
350  ),
351  this->mesh(),
352  dimensionedScalar("zero", dimless/dimLength, 0.0)
353  )
354  );
355 
356  if (radiation_)
357  {
358  scalarField& sigmap = tsigmap().internalField();
359  const scalar dt = this->db().time().deltaTValue();
360  const scalarField& V = this->mesh().V();
361  const scalar epsilon = constProps_.epsilon0();
362  const scalar f = constProps_.f0();
363  const scalarField& sumAreaP = radAreaP_->field();
364 
365  sigmap *= sumAreaP*(1.0 - f)*(1.0 - epsilon)/V/dt;
366  }
367 
368  return tsigmap;
369 }
370 
371 
372 template<class CloudType>
373 inline Foam::scalar Foam::ThermoCloud<CloudType>::Tmax() const
374 {
375  scalar T = -GREAT;
376  scalar n = 0;
377  forAllConstIter(typename ThermoCloud<CloudType>, *this, iter)
378  {
379  const parcelType& p = iter();
380  T = max(T, p.T());
381  n++;
382  }
383 
384  reduce(T, maxOp<scalar>());
385  reduce(n, sumOp<label>());
386 
387  if (n > 0)
388  {
389  return T;
390  }
391  else
392  {
393  return 0.0;
394  }
395 }
396 
397 
398 template<class CloudType>
399 inline Foam::scalar Foam::ThermoCloud<CloudType>::Tmin() const
400 {
401  scalar T = GREAT;
402  scalar n = 0;
403  forAllConstIter(typename ThermoCloud<CloudType>, *this, iter)
404  {
405  const parcelType& p = iter();
406  T = min(T, p.T());
407  n++;
408  }
409 
410  reduce(T, minOp<scalar>());
411  reduce(n, sumOp<label>());
412 
413  if (n > 0)
414  {
415  return T;
416  }
417  else
418  {
419  return 0.0;
420  }
421 }
422 
423 
424 // ************************************************************************* //
Foam::ThermoCloud::p
const volScalarField & p() const
Return const access to the carrier prressure field.
Definition: ThermoCloudI.H:71
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
p
p
Definition: pEqn.H:62
Foam::ThermoCloud::Tmin
scalar Tmin() const
Minimum temperature.
Definition: ThermoCloudI.H:399
Foam::ThermoCloud::Sh
tmp< fvScalarMatrix > Sh(volScalarField &hs) const
Return sensible enthalpy source term [J/kg/m3/s].
Definition: ThermoCloudI.H:224
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
Foam::ThermoCloud::TIntegrator
const scalarIntegrationScheme & TIntegrator() const
Return reference to velocity integration.
Definition: ThermoCloudI.H:87
Foam::minOp
Definition: ops.H:173
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::ThermoCloud::constProps
const parcelType::constantProperties & constProps() const
Return the constant properties.
Definition: ThermoCloudI.H:42
Foam::ThermoCloud::sigmap
tmp< volScalarField > sigmap() const
Return tmp equivalent particulate scattering factor.
Definition: ThermoCloudI.H:336
Foam::HeatTransferModel
Templated heat transfer model class.
Definition: ThermoCloud.H:53
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::IntegrationScheme
Top level model for Integration schemes.
Definition: IntegrationScheme.H:51
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::ThermoCloud::ap
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
Definition: ThermoCloudI.H:299
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::ThermoCloud::radAreaP
DimensionedField< scalar, volMesh > & radAreaP()
Radiation sum of parcel projected areas [m2].
Definition: ThermoCloudI.H:102
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::ThermoCloud::hsTrans
DimensionedField< scalar, volMesh > & hsTrans()
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:192
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::ThermoCloud::thermo
const SLGThermo & thermo() const
Return const access to thermo package.
Definition: ThermoCloudI.H:57
Foam::ThermoCloud::radT4
DimensionedField< scalar, volMesh > & radT4()
Radiation sum of parcel temperature^4 [K4].
Definition: ThermoCloudI.H:132
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::ThermoCloud::Ep
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
Definition: ThermoCloudI.H:263
Foam::FatalError
error FatalError
Foam::ThermoCloud::radAreaPT4
DimensionedField< scalar, volMesh > & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
Definition: ThermoCloudI.H:162
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::ThermoCloud::cloudCopy
const ThermoCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: ThermoCloudI.H:34
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
epsilon
epsilon
Definition: createTDFields.H:56
Foam::sumOp
Definition: ops.H:162
f
labelList f(nPoints)
physicoChemicalConstants.H
Foam::ThermoCloud
Templated base class for thermodynamic cloud.
Definition: ThermoCloud.H:60
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
Foam::ThermoCloud::heatTransfer
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
Definition: ThermoCloudI.H:79
Foam::ThermoCloud::Tmax
scalar Tmax() const
Maximum temperature.
Definition: ThermoCloudI.H:373
Foam::ThermoCloud::radiation
bool radiation() const
Radiation flag.
Definition: ThermoCloudI.H:94
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
Foam::ThermoCloud::T
const volScalarField & T() const
Return const access to the carrier temperature field.
Definition: ThermoCloudI.H:64
Foam::ThermoCloud::hsCoeff
DimensionedField< scalar, volMesh > & hsCoeff()
Return coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:208
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51