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