Loading...
Searching...
No Matches
ThermoParcel.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) 2016-2024 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
27Class
28 Foam::ThermoParcel
29
30Group
31 grpLagrangianIntermediateParcels
32
33Description
34 Thermodynamic parcel class with one/two-way coupling with the continuous
35 phase. Includes Kinematic parcel sub-models, plus:
36 - heat transfer
37
38SourceFiles
39 ThermoParcelI.H
40 ThermoParcel.C
41 ThermoParcelIO.C
42
43\*---------------------------------------------------------------------------*/
44
45#ifndef ThermoParcel_H
46#define ThermoParcel_H
47
48#include "particle.H"
49#include "SLGThermo.H"
50#include "demandDrivenEntry.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
57template<class ParcelType>
58class ThermoParcel;
59
60template<class ParcelType>
61Ostream& operator<<
62(
63 Ostream&,
65);
66
67/*---------------------------------------------------------------------------*\
68 Class ThermoParcel Declaration
69\*---------------------------------------------------------------------------*/
70
71template<class ParcelType>
72class ThermoParcel
73:
74 public ParcelType
75{
76public:
77
78 //- Size in bytes of the fields
79 static const std::size_t sizeofFields;
80
81
82 //- Class to hold thermo particle constant properties
84 :
85 public ParcelType::constantProperties
86 {
87
88 // Private data
89
90 //- Particle initial temperature [K]
92
93 //- Minimum temperature [K]
95
96 //- Maximum temperature [K]
98
99 //- Particle specific heat capacity [J/(kg.K)]
101
102 //- Particle emissivity [] (radiation)
104
105 //- Particle scattering factor [] (radiation)
107
108
109 public:
110
111 // Constructors
112
113 //- Null constructor
115
116 //- Copy constructor
118
119 //- Construct from dictionary
120 constantProperties(const dictionary& parentDict);
121
122
123 // Member functions
124
125 // Access
126
127 //- Return const access to the particle initial temperature [K]
128 inline scalar T0() const;
129
130 //- Return const access to minimum temperature [K]
131 inline scalar TMin() const;
132
133 //- Return const access to maximum temperature [K]
134 inline scalar TMax() const;
135
136 //- Set the maximum temperature [K]
137 inline void setTMax(const scalar TMax);
138
139 //- Return const access to the particle specific heat capacity
140 // [J/(kg.K)]
141 inline scalar Cp0() const;
142
143 //- Return const access to the particle emissivity []
144 // Active for radiation only
145 inline scalar epsilon0() const;
146
147 //- Return const access to the particle scattering factor []
148 // Active for radiation only
149 inline scalar f0() const;
150 };
151
152
153 class trackingData
154 :
155 public ParcelType::trackingData
156 {
157 private:
158
159 // Private data
160
161 //- Local copy of carrier specific heat field
162 // Cp not stored on carrier thermo, but returned as tmp<...>
163 const volScalarField Cp_;
164
165 //- Local copy of carrier thermal conductivity field
166 // kappa not stored on carrier thermo, but returned as tmp<...>
167 const volScalarField kappa_;
168
169
170 // Interpolators for continuous phase fields
171
172 //- Temperature field interpolator
174
175 //- Specific heat capacity field interpolator
177
178 //- Thermal conductivity field interpolator
179 autoPtr<interpolation<scalar>> kappaInterp_;
180
181 //- Radiation field interpolator
183
184
185 // Cached continuous phase properties
186
187 //- Temperature [K]
188 scalar Tc_;
189
190 //- Specific heat capacity [J/(kg.K)]
191 scalar Cpc_;
192
193
194 public:
195
196 typedef typename ParcelType::trackingData::trackPart trackPart;
197
198 // Constructors
199
200 //- Construct from components
201 template <class TrackCloudType>
202 inline trackingData
203 (
204 const TrackCloudType& cloud,
205 trackPart part = ParcelType::trackingData::tpLinearTrack
206 );
207
208
209 // Member functions
210
211 //- Return access to the locally stored carrier Cp field
212 inline const volScalarField& Cp() const;
213
214 //- Return access to the locally stored carrier kappa field
215 inline const volScalarField& kappa() const;
216
217 //- Return const access to the interpolator for continuous
218 // phase temperature field
219 inline const interpolation<scalar>& TInterp() const;
220
221 //- Return const access to the interpolator for continuous
222 // phase specific heat capacity field
223 inline const interpolation<scalar>& CpInterp() const;
224
225 //- Return const access to the interpolator for continuous
226 // phase thermal conductivity field
227 inline const interpolation<scalar>& kappaInterp() const;
228
229 //- Return const access to the interpolator for continuous
230 // radiation field
231 inline const interpolation<scalar>& GInterp() const;
232
233 //- Return the continuous phase temperature
234 inline scalar Tc() const;
235
236 //- Access the continuous phase temperature
237 inline scalar& Tc();
238
239 //- Return the continuous phase specific heat capacity
240 inline scalar Cpc() const;
241
242 //- Access the continuous phase specific heat capacity
243 inline scalar& Cpc();
244 };
245
246
247protected:
249 // Protected data
250
251 // Parcel properties
252
253 //- Temperature [K]
254 scalar T_;
255
256 //- Specific heat capacity [J/(kg.K)]
257 scalar Cp_;
258
259
260 // Protected Member Functions
261
262 //- Calculate new particle temperature
263 template<class TrackCloudType>
264 scalar calcHeatTransfer
265 (
266 TrackCloudType& cloud,
268 const scalar dt, // timestep
269 const scalar Re, // Reynolds number
270 const scalar Pr, // Prandtl number - surface
271 const scalar kappa, // Thermal conductivity - surface
272 const scalar NCpW, // Sum of N*Cp*W of emission species
273 const scalar Sh, // explicit particle enthalpy source
274 scalar& dhsTrans, // sensible enthalpy transfer to carrier
275 scalar& Sph // linearised heat transfer coefficient
276 );
277
278
279public:
280
281 // Static data members
282
283 //- Runtime type information
284 TypeName("ThermoParcel");
285
286 //- String representation of properties
288 (
289 ParcelType,
290 " T"
291 + " Cp"
292 );
293
294
295 // Constructors
296
297 //- Construct from mesh, coordinates and topology
298 // Other properties initialised as null
299 inline ThermoParcel
300 (
301 const polyMesh& mesh,
303 const label celli,
304 const label tetFacei,
305 const label tetPti
306 );
307
308 //- Construct from a position and a cell, searching for the rest of the
309 // required topology. Other properties are initialised as null.
310 inline ThermoParcel
311 (
312 const polyMesh& mesh,
313 const vector& position,
314 const label celli
315 );
316
317 //- Construct from components
318 inline ThermoParcel
319 (
320 const polyMesh& mesh,
322 const label celli,
323 const label tetFacei,
324 const label tetPti,
325 const label typeId,
326 const scalar nParticle0,
327 const scalar d0,
328 const scalar dTarget0,
329 const vector& U0,
330 const vector& f0,
331 const vector& angularMomentum0,
332 const vector& torque0,
333 const constantProperties& constProps
334 );
335
336 //- Construct from Istream
338 (
340 Istream& is,
341 bool readFields = true,
342 bool newFormat = true
343 );
344
345 //- Construct as a copy
347
348 //- Construct as a copy
349 ThermoParcel(const ThermoParcel& p, const polyMesh& mesh);
350
351 //- Return a (basic particle) clone
352 virtual autoPtr<particle> clone() const
353 {
354 return particle::Clone(*this);
355 }
356
357 //- Return a (basic particle) clone
358 virtual autoPtr<particle> clone(const polyMesh& mesh) const
359 {
360 return particle::Clone(*this, mesh);
361 }
362
363 //- Factory class to read-construct particles (for parallel transfer)
364 class iNew
365 {
366 const polyMesh& mesh_;
367
368 public:
369
371 :
372 mesh_(mesh)
373 {}
374
376 {
378 (
379 new ThermoParcel<ParcelType>(mesh_, is, true)
380 );
381 }
382 };
383
384
385 // Member Functions
386
387 // Access
388
389 //- Return const access to temperature
390 inline scalar T() const;
391
392 //- Return const access to specific heat capacity
393 inline scalar Cp() const;
394
395 //- Return the parcel sensible enthalpy
396 inline scalar hs() const;
397
398
399 // Edit
400
401 //- Return access to temperature
402 inline scalar& T();
403
404 //- Return access to specific heat capacity
405 inline scalar& Cp();
406
407
408 // Main calculation loop
409
410 //- Set cell values
411 template<class TrackCloudType>
412 void setCellValues(TrackCloudType& cloud, trackingData& td);
413
414 //- Correct cell values using latest transfer information
415 template<class TrackCloudType>
417 (
418 TrackCloudType& cloud,
419 trackingData& td,
420 const scalar dt
421 );
422
423 //- Calculate surface thermo properties
424 template<class TrackCloudType>
426 (
427 TrackCloudType& cloud,
428 trackingData& td,
429 const scalar T,
430 scalar& Ts,
431 scalar& rhos,
432 scalar& mus,
433 scalar& Pr,
434 scalar& kappas
435 ) const;
436
437 //- Update parcel properties over the time interval
438 template<class TrackCloudType>
439 void calc
440 (
441 TrackCloudType& cloud,
442 trackingData& td,
443 const scalar dt
444 );
445
447 // I-O
448
449 //- Read
450 template<class CloudType>
451 static void readFields(CloudType& c);
452
453 //- Write
454 template<class CloudType>
455 static void writeFields(const CloudType& c);
457 //- Write individual parcel properties to stream
458 void writeProperties
459 (
460 Ostream& os,
461 const wordRes& filters,
462 const word& delim,
463 const bool namesOnly = false
464 ) const;
465
466 //- Read particle fields as objects from the obr registry
467 template<class CloudType>
468 static void readObjects(CloudType& c, const objectRegistry& obr);
469
470 //- Write particle fields as objects into the obr registry
471 template<class CloudType>
472 static void writeObjects(const CloudType& c, objectRegistry& obr);
473
474
475 // Ostream Operator
476
477 friend Ostream& operator<< <ParcelType>
479 Ostream&,
481 );
482};
484
485// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
486
487} // End namespace Foam
488
489// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
490
491#include "ThermoParcelI.H"
493
494// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
495
496#ifdef NoRepository
497 #include "ThermoParcel.C"
498#endif
499
500// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
501
502#endif
503
504// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
scalar Re(const trackingData &td) const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to hold thermo particle constant properties.
void setTMax(const scalar TMax)
Set the maximum temperature [K].
scalar f0() const
Return const access to the particle scattering factor [].
scalar TMin() const
Return const access to minimum temperature [K].
scalar epsilon0() const
Return const access to the particle emissivity [].
scalar Cp0() const
Return const access to the particle specific heat capacity.
scalar TMax() const
Return const access to maximum temperature [K].
scalar T0() const
Return const access to the particle initial temperature [K].
iNew(const polyMesh &mesh)
autoPtr< ThermoParcel< ParcelType > > operator()(Istream &is) const
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
const volScalarField & kappa() const
Return access to the locally stored carrier kappa field.
const interpolation< scalar > & CpInterp() const
Return const access to the interpolator for continuous.
scalar Cpc() const
Return the continuous phase specific heat capacity.
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
ParcelType::trackingData::trackPart trackPart
const volScalarField & Cp() const
Return access to the locally stored carrier Cp field.
trackingData(const TrackCloudType &cloud, trackPart part=ParcelType::trackingData::tpLinearTrack)
Construct from components.
scalar Tc() const
Return the continuous phase temperature.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
scalar Cp() const
Return const access to specific heat capacity.
scalar calcHeatTransfer(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
scalar T() const
Return const access to temperature.
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Return a (basic particle) clone.
AddToPropertyList(ParcelType, " T"+" Cp")
String representation of properties.
scalar hs() const
Return the parcel sensible enthalpy.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
scalar & T()
Return access to temperature.
ThermoParcel(const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
Construct from Istream.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label typeId, const scalar nParticle0, const scalar d0, const scalar dTarget0, const vector &U0, const vector &f0, const vector &angularMomentum0, const vector &torque0, const constantProperties &constProps)
Construct from components.
ThermoParcel(const ThermoParcel &p, const polyMesh &mesh)
Construct as a copy.
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
ThermoParcel(const ThermoParcel &p)
Construct as a copy.
ThermoParcel(const polyMesh &mesh, const vector &position, const label celli)
Construct from a position and a cell, searching for the rest of the.
static void writeFields(const CloudType &c)
Write.
void calcSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
scalar & Cp()
Return access to specific heat capacity.
TypeName("ThermoParcel")
Runtime type information.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
Class for demand-driven dictionary entries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for volume field interpolation.
Registry of regIOobjects.
vector position() const
Return current particle position.
Definition particleI.H:283
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition particle.H:552
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
OBJstream os(runTime.globalPath()/outputName)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
DSMCCloud< dsmcParcel > CloudType
GeometricField< scalar, fvPatchField, volMesh > volScalarField
scalarField Re(const UList< complex > &cmplx)
Extract real component.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68