All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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:
191 
192  // Constructors
193 
194  //- Construct from components
195  template <class TrackCloudType>
196  inline trackingData(const TrackCloudType& cloud);
197 
198 
199  // Member Functions
200 
201  //- Return access to the locally stored carrier Cp field
202  inline const volScalarField& Cp() const;
203 
204  //- Return access to the locally stored carrier kappa field
205  inline const volScalarField& kappa() const;
206 
207  //- Return const access to the interpolator for continuous
208  // phase temperature field
209  inline const interpolation<scalar>& TInterp() const;
210 
211  //- Return const access to the interpolator for continuous
212  // phase specific heat capacity field
213  inline const interpolation<scalar>& CpInterp() const;
214 
215  //- Return const access to the interpolator for continuous
216  // phase thermal conductivity field
217  inline const interpolation<scalar>& kappaInterp() const;
218 
219  //- Return const access to the interpolator for continuous
220  // radiation field
221  inline const interpolation<scalar>& GInterp() const;
222 
223  //- Return the continuous phase temperature
224  inline scalar Tc() const;
225 
226  //- Access the continuous phase temperature
227  inline scalar& Tc();
228 
229  //- Return the continuous phase specific heat capacity
230  inline scalar Cpc() const;
231 
232  //- Access the continuous phase specific heat capacity
233  inline scalar& Cpc();
234  };
235 
236 
237 protected:
238 
239  // Protected data
240 
241  // Parcel properties
242 
243  //- Temperature [K]
244  scalar T_;
245 
246  //- Specific heat capacity [J/kg/K]
247  scalar Cp_;
248 
249 
250  // Protected Member Functions
251 
252  //- Calculate new particle temperature
253  template<class TrackCloudType>
254  scalar calcHeatTransfer
255  (
256  TrackCloudType& cloud,
257  trackingData& td,
258  const scalar dt, // timestep
259  const scalar Re, // Reynolds number
260  const scalar Pr, // Prandtl number - surface
261  const scalar kappa, // Thermal conductivity - surface
262  const scalar NCpW, // Sum of N*Cp*W of emission species
263  const scalar Sh, // explicit particle enthalpy source
264  scalar& dhsTrans, // sensible enthalpy transfer to carrier
265  scalar& Sph // linearised heat transfer coefficient
266  );
267 
268 
269 public:
270 
271  // Static Data Members
272 
273  //- Runtime type information
274  TypeName("ThermoParcel");
275 
276  //- String representation of properties
278  (
279  ParcelType,
280  " T"
281  + " Cp"
282  );
283 
284 
285  // Constructors
286 
287  //- Construct from mesh, coordinates and topology
288  // Other properties initialised as null
289  inline ThermoParcel
290  (
291  const polyMesh& mesh,
292  const barycentric& coordinates,
293  const label celli,
294  const label tetFacei,
295  const label tetPti
296  );
297 
298  //- Construct from a position and a cell, searching for the rest of the
299  // required topology. Other properties are initialised as null.
300  inline ThermoParcel
301  (
302  const polyMesh& mesh,
303  const vector& position,
304  const label celli
305  );
306 
307  //- Construct from components
308  inline ThermoParcel
309  (
310  const polyMesh& mesh,
311  const barycentric& coordinates,
312  const label celli,
313  const label tetFacei,
314  const label tetPti,
315  const label typeId,
316  const scalar nParticle0,
317  const scalar d0,
318  const scalar dTarget0,
319  const vector& U0,
320  const vector& f0,
321  const vector& angularMomentum0,
322  const vector& torque0,
323  const constantProperties& constProps
324  );
325 
326  //- Construct from Istream
328  (
329  const polyMesh& mesh,
330  Istream& is,
331  bool readFields = true
332  );
333 
334  //- Construct as a copy
335  ThermoParcel(const ThermoParcel& p);
336 
337  //- Construct as a copy
338  ThermoParcel(const ThermoParcel& p, const polyMesh& mesh);
339 
340  //- Construct and return a (basic particle) clone
341  virtual autoPtr<particle> clone() const
342  {
343  return autoPtr<particle>(new ThermoParcel(*this));
344  }
345 
346  //- Construct and return a (basic particle) clone
347  virtual autoPtr<particle> clone(const polyMesh& mesh) const
348  {
349  return autoPtr<particle>(new ThermoParcel(*this, mesh));
350  }
351 
352  //- Factory class to read-construct particles used for
353  // parallel transfer
354  class iNew
355  {
356  const polyMesh& mesh_;
357 
358  public:
360  iNew(const polyMesh& mesh)
361  :
362  mesh_(mesh)
363  {}
365  autoPtr<ThermoParcel<ParcelType>> operator()(Istream& is) const
366  {
368  (
369  new ThermoParcel<ParcelType>(mesh_, is, true)
370  );
371  }
372  };
373 
374 
375  // Member Functions
376 
377  // Access
378 
379  //- Return const access to temperature
380  inline scalar T() const;
381 
382  //- Return const access to specific heat capacity
383  inline scalar Cp() const;
384 
385  //- Return the parcel sensible enthalpy
386  inline scalar hs() const;
387 
388 
389  // Edit
390 
391  //- Return access to temperature
392  inline scalar& T();
393 
394  //- Return access to specific heat capacity
395  inline scalar& Cp();
396 
397 
398  // Main calculation loop
399 
400  //- Set cell values
401  template<class TrackCloudType>
402  void setCellValues(TrackCloudType& cloud, trackingData& td);
403 
404  //- Correct cell values using latest transfer information
405  template<class TrackCloudType>
407  (
408  TrackCloudType& cloud,
409  trackingData& td,
410  const scalar dt
411  );
412 
413  //- Calculate surface thermo properties
414  template<class TrackCloudType>
415  void calcSurfaceValues
416  (
417  TrackCloudType& cloud,
418  trackingData& td,
419  const scalar T,
420  scalar& Ts,
421  scalar& rhos,
422  scalar& mus,
423  scalar& Pr,
424  scalar& kappas
425  ) const;
426 
427  //- Update parcel properties over the time interval
428  template<class TrackCloudType>
429  void calc
430  (
431  TrackCloudType& cloud,
432  trackingData& td,
433  const scalar dt
434  );
435 
436 
437  // I-O
438 
439  //- Read
440  template<class CloudType>
441  static void readFields(CloudType& c);
442 
443  //- Write
444  template<class CloudType>
445  static void writeFields(const CloudType& c);
446 
447 
448  // Ostream Operator
449 
450  friend Ostream& operator<< <ParcelType>
451  (
452  Ostream&,
454  );
455 };
456 
457 
458 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
459 
460 } // End namespace Foam
461 
462 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
463 
464 #include "ThermoParcelI.H"
466 
467 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
468 
469 #ifdef NoRepository
470  #include "ThermoParcel.C"
471 #endif
472 
473 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
474 
475 #endif
476 
477 // ************************************************************************* //
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
const dimensionedScalar & kappa
Coulomb constant: default SI units: [N.m2/C2].
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:246
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 & c
Speed of light in a vacuum.
Factory class to read-construct particles used for.
Definition: ThermoParcel.H:353
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: ThermoParcel.H:340
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:243
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:54
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.
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 [].