ThermoParcel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 Class
25  Foam::ThermoParcel
26 
27 Description
28  Thermodynamic parcel class with one/two-way coupling with the continuous
29  phase
30 
31 SourceFiles
32  ThermoParcelI.H
33  ThermoParcel.C
34  ThermoParcelIO.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef ThermoParcel_H
39 #define ThermoParcel_H
40 
41 #include "particle.H"
42 #include "interpolation.H"
43 #include "fluidThermo.H"
44 #include "demandDrivenEntry.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 template<class ParcelType>
52 class ThermoParcel;
53 
54 template<class ParcelType>
55 Ostream& operator<<
56 (
57  Ostream&,
59 );
60 
61 
62 /*---------------------------------------------------------------------------*\
63  Class ThermoParcelName Declaration
64 \*---------------------------------------------------------------------------*/
65 
67 
68 
69 /*---------------------------------------------------------------------------*\
70  Class ThermoParcel Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 template<class ParcelType>
74 class ThermoParcel
75 :
76  public ParcelType,
77  public ThermoParcelName
78 {
79  // Private Data
80 
81  //- Size in bytes of the fields
82  static const std::size_t sizeofFields_;
83 
84 
85 public:
86 
87  //- Class to hold thermo particle constant properties
88  class constantProperties
89  :
90  public ParcelType::constantProperties
91  {
92 
93  // Private Data
94 
95  //- Particle initial temperature [K]
97 
98  //- Minimum temperature [K]
100 
101  //- Maximum temperature [K]
103 
104  //- Particle specific heat capacity [J/kg/K]
106 
107  //- Particle emissivity [] (radiation)
108  demandDrivenEntry<scalar> epsilon0_;
109 
110  //- Particle scattering factor [] (radiation)
112 
113 
114  public:
115 
116  // Constructors
117 
118  //- Null constructor
120 
121  //- Copy constructor
123 
124  //- Construct from dictionary
125  constantProperties(const dictionary& parentDict);
126 
127 
128  // Member Functions
129 
130  // Access
131 
132  //- Return const access to the particle initial temperature [K]
133  inline scalar T0() const;
134 
135  //- Return const access to minimum temperature [K]
136  inline scalar TMin() const;
137 
138  //- Return const access to maximum temperature [K]
139  inline scalar TMax() const;
140 
141  //- Set the maximum temperature [K]
142  inline void setTMax(const scalar TMax);
143 
144  //- Return const access to the particle specific heat capacity
145  // [J/kg/K]
146  inline scalar Cp0() const;
147 
148  //- Return const access to the particle emissivity []
149  // Active for radiation only
150  inline scalar epsilon0() const;
151 
152  //- Return const access to the particle scattering factor []
153  // Active for radiation only
154  inline scalar f0() const;
155  };
156 
157 
158  class trackingData
159  :
160  public ParcelType::trackingData
161  {
162  private:
163 
164  // Private Data
165 
166  //- Local copy of carrier specific heat field
167  // Cp not stored on carrier thermo, but returned as tmp<...>
168  const volScalarField Cp_;
169 
170  //- Local copy of carrier thermal conductivity field
171  // kappa not stored on carrier thermo, but returned as tmp<...>
172  const volScalarField kappa_;
173 
174 
175  // Interpolators for continuous phase fields
176 
177  //- Interpolator for continuous phase pressure field
179 
180  //- Temperature field interpolator
182 
183  //- Specific heat capacity field interpolator
185 
186  //- Thermal conductivity field interpolator
187  autoPtr<interpolation<scalar>> kappaInterp_;
188 
189  //- Radiation field interpolator
191 
192 
193  // Cached continuous phase properties
194 
195  //- Pressure [Pa]
196  scalar pc_;
197 
198  //- Temperature [K]
199  scalar Tc_;
200 
201  //- Specific heat capacity [J/kg/K]
202  scalar Cpc_;
203 
204 
205  public:
206 
207  // Constructors
208 
209  //- Construct from components
210  template <class TrackCloudType>
211  inline trackingData(const TrackCloudType& cloud);
212 
213 
214  // Member Functions
215 
216  //- Return access to the locally stored carrier Cp field
217  inline const volScalarField& Cp() const;
218 
219  //- Return access to the locally stored carrier kappa field
220  inline const volScalarField& kappa() const;
221 
222  //- Return const access to the interpolator for continuous phase
223  // pressure field
224  inline const interpolation<scalar>& pInterp() const;
225 
226  //- Return const access to the interpolator for continuous
227  // phase temperature field
228  inline const interpolation<scalar>& TInterp() const;
229 
230  //- Return const access to the interpolator for continuous
231  // phase specific heat capacity field
232  inline const interpolation<scalar>& CpInterp() const;
233 
234  //- Return const access to the interpolator for continuous
235  // phase thermal conductivity field
236  inline const interpolation<scalar>& kappaInterp() const;
237 
238  //- Return const access to the interpolator for continuous
239  // radiation field
240  inline const interpolation<scalar>& GInterp() const;
241 
242  //- Return the continuous phase pressure
243  inline scalar pc() const;
244 
245  //- Access the continuous phase pressure
246  inline scalar& pc();
247 
248  //- Return the continuous phase temperature
249  inline scalar Tc() const;
250 
251  //- Access the continuous phase temperature
252  inline scalar& Tc();
253 
254  //- Return the continuous phase specific heat capacity
255  inline scalar Cpc() const;
256 
257  //- Access the continuous phase specific heat capacity
258  inline scalar& Cpc();
259  };
260 
261 
262 protected:
263 
264  // Protected data
265 
266  // Parcel properties
267 
268  //- Temperature [K]
269  scalar T_;
270 
271  //- Specific heat capacity [J/kg/K]
272  scalar Cp_;
273 
274 
275  // Protected Member Functions
276 
277  //- Calculate new particle temperature
278  template<class TrackCloudType>
279  scalar calcHeatTransfer
280  (
281  TrackCloudType& cloud,
282  trackingData& td,
283  const scalar dt, // timestep
284  const scalar Re, // Reynolds number
285  const scalar Pr, // Prandtl number - surface
286  const scalar kappa, // Thermal conductivity - surface
287  const scalar NCpW, // Sum of N*Cp*W of emission species
288  const scalar Sh, // explicit particle enthalpy source
289  scalar& dhsTrans, // sensible enthalpy transfer to carrier
290  scalar& Sph // linearised heat transfer coefficient
291  );
292 
293 
294 public:
295 
296  // Static Data Members
297 
298  //- String representation of properties
300  (
301  ParcelType,
302  " T"
303  + " Cp"
304  );
305 
306 
307  // Constructors
308 
309  //- Construct from mesh, coordinates and topology
310  // Other properties initialised as null
311  inline ThermoParcel
312  (
313  const polyMesh& mesh,
314  const barycentric& coordinates,
315  const label celli,
316  const label tetFacei,
317  const label tetPti,
318  const label facei
319  );
320 
321  //- Construct from a position and a cell, searching for the rest of the
322  // required topology. Other properties are initialised as null.
323  inline ThermoParcel
324  (
325  const polyMesh& mesh,
326  const vector& position,
327  const label celli
328  );
329 
330  //- Construct from Istream
331  ThermoParcel(Istream& is, bool readFields = true);
332 
333  //- Construct as a copy
334  ThermoParcel(const ThermoParcel& p);
335 
336  //- Construct and return a clone
337  virtual autoPtr<particle> clone() const
338  {
339  return autoPtr<particle>(new ThermoParcel(*this));
340  }
341 
342  //- Construct from Istream and return
343  static autoPtr<ThermoParcel> New(Istream& is)
344  {
345  return autoPtr<ThermoParcel>(new ThermoParcel(is));
346  }
347 
348 
349  // Member Functions
350 
351  // Access
352 
353  //- Return const access to temperature
354  inline scalar T() const;
355 
356  //- Return const access to specific heat capacity
357  inline scalar Cp() const;
358 
359 
360  // Edit
361 
362  //- Return access to temperature
363  inline scalar& T();
364 
365  //- Return access to specific heat capacity
366  inline scalar& Cp();
367 
368 
369  // Main calculation loop
370 
371  //- Set cell values
372  template<class TrackCloudType>
373  void setCellValues(TrackCloudType& cloud, trackingData& td);
374 
375  //- Correct cell values using latest transfer information
376  template<class TrackCloudType>
378  (
379  TrackCloudType& cloud,
380  trackingData& td,
381  const scalar dt
382  );
383 
384  //- Calculate surface thermo properties
385  template<class TrackCloudType>
386  void calcSurfaceValues
387  (
388  const TrackCloudType& cloud,
389  const trackingData& td,
390  const scalar T,
391  scalar& Ts,
392  scalar& rhos,
393  scalar& mus,
394  scalar& Pr,
395  scalar& kappas
396  ) const;
397 
398  //- Update parcel properties over the time interval
399  template<class TrackCloudType>
400  void calc
401  (
402  TrackCloudType& cloud,
403  trackingData& td,
404  const scalar dt
405  );
406 
407 
408  // I-O
409 
410  //- Read
411  template<class CloudType>
412  static void readFields(CloudType& c);
413 
414  //- Write
415  template<class CloudType>
416  static void writeFields(const CloudType& c);
417 
418  //- Write
419  template<class CloudType, class CompositionType>
420  static void writeFields
421  (
422  const CloudType& c,
423  const CompositionType& compModel
424  );
425 
426 
427  // Ostream Operator
428 
429  friend Ostream& operator<< <ParcelType>
430  (
431  Ostream&,
433  );
434 };
435 
436 
437 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
438 
439 } // End namespace Foam
440 
441 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
442 
443 #include "ThermoParcelI.H"
445 
446 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
447 
448 #ifdef NoRepository
449  #include "ThermoParcel.C"
450 #endif
451 
452 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
453 
454 #endif
455 
456 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:79
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to hold thermo particle constant properties.
Definition: ThermoParcel.H:90
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].
const interpolation< scalar > & GInterp() const
Return const access to the interpolator for continuous.
trackingData(const TrackCloudType &cloud)
Construct from components.
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.
scalar pc() const
Return the continuous phase pressure.
const interpolation< scalar > & TInterp() const
Return const access to the interpolator for continuous.
const interpolation< scalar > & pInterp() const
Return const access to the interpolator for continuous phase.
const interpolation< scalar > & kappaInterp() const
Return const access to the interpolator for continuous.
const volScalarField & Cp() const
Return access to the locally stored carrier Cp field.
scalar Tc() const
Return the continuous phase temperature.
Thermodynamic parcel class with one/two-way coupling with the continuous phase.
Definition: ThermoParcel.H:77
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.
AddToPropertyList(ParcelType, " T"+" Cp")
String representation of properties.
scalar Cp_
Specific heat capacity [J/kg/K].
Definition: ThermoParcel.H:271
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Definition: ThermoParcel.C:67
static autoPtr< ThermoParcel > New(Istream &is)
Construct from Istream and return.
Definition: ThermoParcel.H:342
virtual autoPtr< particle > clone() const
Construct and return a clone.
Definition: ThermoParcel.H:336
static void writeFields(const CloudType &c)
Write.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
Definition: ThermoParcel.C:37
scalar T_
Temperature [K].
Definition: ThermoParcel.H:268
static void readFields(CloudType &c)
Read.
void calcSurfaceValues(const TrackCloudType &cloud, const trackingData &td, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
Definition: ThermoParcel.C:95
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Definition: ThermoParcel.C:138
ThermoParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
Definition: ThermoParcelI.H:77
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base class for interpolation.
Definition: interpolation.H:55
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
TemplateName(FvFaceCellWave)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
volScalarField & p