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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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
30
31using namespace Foam::constant;
32
33// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34
35template<class CloudType>
38{
39 return *cloudCopyPtr_;
40}
41
42
43template<class CloudType>
46{
47 return constProps_;
48}
49
50
51template<class CloudType>
54{
55 return constProps_;
56}
57
58
59template<class CloudType>
61{
62 return thermo_;
63}
64
65
66template<class CloudType>
68{
69 return T_;
70}
71
72
73template<class CloudType>
75{
76 return p_;
77}
78
79
80template<class CloudType>
83{
84 return *heatTransferModel_;
85}
86
87
88template<class CloudType>
89inline const Foam::integrationScheme&
91{
92 return *TIntegrator_;
93}
94
95
96template<class CloudType>
98{
99 return radiation_;
100}
101
102
103template<class CloudType>
106{
107 if (!radiation_)
108 {
110 << "Radiation field requested, but radiation model not active"
111 << abort(FatalError);
112 }
113
114 return *radAreaP_;
115}
116
117
118template<class CloudType>
121{
122 if (!radiation_)
123 {
125 << "Radiation field requested, but radiation model not active"
126 << abort(FatalError);
127 }
128
129 return *radAreaP_;
130}
131
132
133template<class CloudType>
136{
137 if (!radiation_)
138 {
140 << "Radiation field requested, but radiation model not active"
141 << abort(FatalError);
142 }
143
144 return *radT4_;
145}
146
147
148template<class CloudType>
151{
152 if (!radiation_)
153 {
155 << "Radiation field requested, but radiation model not active"
156 << abort(FatalError);
157 }
158
159 return *radT4_;
160}
161
162
163template<class CloudType>
166{
167 if (!radiation_)
168 {
170 << "Radiation field requested, but radiation model not active"
171 << abort(FatalError);
172 }
173
174 return *radAreaPT4_;
175}
176
177
178template<class CloudType>
181{
182 if (!radiation_)
183 {
185 << "Radiation field requested, but radiation model not active"
186 << abort(FatalError);
187 }
188
189 return *radAreaPT4_;
190}
191
192
193template<class CloudType>
196{
197 return *hsTrans_;
198}
199
200
201template<class CloudType>
204{
205 return *hsTrans_;
206}
207
208
209template<class CloudType>
212{
213 return *hsCoeff_;
214}
215
216
217template<class CloudType>
220{
221 return *hsCoeff_;
222}
223
224
225template<class CloudType>
228{
229 if (debug)
230 {
231 Info<< "hsTrans min/max = " << min(hsTrans()).value() << ", "
232 << max(hsTrans()).value() << nl
233 << "hsCoeff min/max = " << min(hsCoeff()).value() << ", "
234 << max(hsCoeff()).value() << endl;
235 }
236
237 if (this->solution().coupled())
238 {
239 if (this->solution().semiImplicit("h"))
241 const volScalarField Cp(thermo_.thermo().Cp());
243 Vdt(this->mesh().V()*this->db().time().deltaT());
244
245 return
246 hsTrans()/Vdt
247 - fvm::SuSp(hsCoeff()/(Cp*Vdt), hs)
248 + hsCoeff()/(Cp*Vdt)*hs;
249 }
250 else
251 {
253 fvScalarMatrix& fvm = tfvm.ref();
254
255 fvm.source() = -hsTrans()/(this->db().time().deltaT());
257 return tfvm;
258 }
259 }
260
261 return tmp<fvScalarMatrix>::New(hs, dimEnergy/dimTime);
263
264
265template<class CloudType>
267{
271 (
273 (
274 this->name() + ":radiation:Ep",
275 this->db().time().timeName(),
276 this->db(),
277 IOobject::NO_READ,
278 IOobject::NO_WRITE,
279 false
280 ),
281 this->mesh(),
282 dimensionedScalar(dimMass/dimLength/pow3(dimTime), Zero)
283 )
284 );
285
286 if (radiation_)
287 {
288 scalarField& Ep = tEp.ref().primitiveFieldRef();
289 const scalar dt = this->db().time().deltaTValue();
290 const scalarField& V = this->mesh().V();
291 const scalar epsilon = constProps_.epsilon0();
292 const scalarField& sumAreaPT4 = radAreaPT4_->field();
293
294 Ep = sumAreaPT4*epsilon*physicoChemical::sigma.value()/V/dt;
295 }
296
297 return tEp;
298}
299
300
301template<class CloudType>
303{
305 (
307 (
309 (
310 this->name() + ":radiation:ap",
311 this->db().time().timeName(),
312 this->db(),
313 IOobject::NO_READ,
314 IOobject::NO_WRITE,
315 false
316 ),
317 this->mesh(),
318 dimensionedScalar(dimless/dimLength, Zero)
319 )
320 );
321
322 if (radiation_)
323 {
324 scalarField& ap = tap.ref().primitiveFieldRef();
325 const scalar dt = this->db().time().deltaTValue();
326 const scalarField& V = this->mesh().V();
327 const scalar epsilon = constProps_.epsilon0();
328 const scalarField& sumAreaP = radAreaP_->field();
329
330 ap = sumAreaP*epsilon/V/dt;
331 }
332
333 return tap;
334}
335
336
337template<class CloudType>
340{
341 tmp<volScalarField> tsigmap
342 (
344 (
346 (
347 this->name() + ":radiation:sigmap",
348 this->db().time().timeName(),
349 this->db(),
350 IOobject::NO_READ,
351 IOobject::NO_WRITE,
352 false
353 ),
354 this->mesh(),
356 )
357 );
358
359 if (radiation_)
360 {
361 scalarField& sigmap = tsigmap.ref().primitiveFieldRef();
362 const scalar dt = this->db().time().deltaTValue();
363 const scalarField& V = this->mesh().V();
364 const scalar epsilon = constProps_.epsilon0();
365 const scalar f = constProps_.f0();
366 const scalarField& sumAreaP = radAreaP_->field();
367
368 sigmap = sumAreaP*(1.0 - f)*(1.0 - epsilon)/V/dt;
369 }
370
371 return tsigmap;
372}
373
374
375template<class CloudType>
376inline Foam::scalar Foam::ThermoCloud<CloudType>::Tmax() const
377{
378 scalar val = -GREAT;
379 bool nonEmpty = false;
380
381 for (const parcelType& p : *this)
382 {
383 val = max(val, p.T());
384 nonEmpty = true;
385 }
386
387 if (returnReduce(nonEmpty, orOp<bool>()))
388 {
389 return returnReduce(val, maxOp<scalar>());
390 }
391
392 return 0;
393}
394
395
396template<class CloudType>
397inline Foam::scalar Foam::ThermoCloud<CloudType>::Tmin() const
398{
399 scalar val = GREAT;
400 bool nonEmpty = false;
401
402 for (const parcelType& p : *this)
403 {
404 val = min(val, p.T());
405 nonEmpty = true;
406 }
407
408 if (returnReduce(nonEmpty, orOp<bool>()))
409 {
410 return returnReduce(val, minOp<scalar>());
411 }
412
413 return 0;
414}
415
416
417// ************************************************************************* //
Y[inertIndex] max(0.0)
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:82
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Templated class to calculate the fluid-particle heat transfer coefficients based on a specified Nusse...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
scalar Sh() const
Sherwood number.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
Templated base class for thermodynamic cloud.
Definition: ThermoCloud.H:70
const parcelType::constantProperties & constProps() const
Return the constant properties.
Definition: ThermoCloudI.H:45
volScalarField::Internal & radAreaP()
Radiation sum of parcel projected areas [m2].
Definition: ThermoCloudI.H:105
volScalarField::Internal & hsCoeff()
Return coefficient for carrier phase hs equation.
Definition: ThermoCloudI.H:211
volScalarField::Internal & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
Definition: ThermoCloudI.H:165
const volScalarField & T() const
Return const access to the carrier temperature field.
Definition: ThermoCloudI.H:67
bool radiation() const
Radiation flag.
Definition: ThermoCloudI.H:97
scalar Tmin() const
Minimum temperature.
Definition: ThermoCloudI.H:397
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
Definition: ThermoCloudI.H:135
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition: ThermoCloud.H:79
const SLGThermo & thermo() const
Return const access to thermo package.
Definition: ThermoCloudI.H:60
const integrationScheme & TIntegrator() const
Return reference to velocity integration.
Definition: ThermoCloudI.H:90
const ThermoCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: ThermoCloudI.H:37
scalar Tmax() const
Maximum temperature.
Definition: ThermoCloudI.H:376
tmp< volScalarField > sigmap() const
Return tmp equivalent particulate scattering factor.
Definition: ThermoCloudI.H:339
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
Definition: ThermoCloudI.H:82
volScalarField::Internal & hsTrans()
Sensible enthalpy transfer [J/kg].
Definition: ThermoCloudI.H:195
const volScalarField & p() const
Return const access to the carrier pressure field.
Definition: ThermoCloudI.H:74
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
Definition: ThermoCloudI.H:266
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
Definition: ThermoCloudI.H:302
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Field< Type > & source() noexcept
Definition: fvMatrix.H:458
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
volScalarField & p
const volScalarField & Cp
Definition: EEqn.H:7
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
scalar epsilon
word timeName
Definition: getTimeIndex.H:3
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)