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 | Copyright (C) 2011-2016 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 
148 
149  template<class CloudType>
150  class TrackingData
151  :
152  public ParcelType::template TrackingData<CloudType>
153  {
154  private:
155 
156  // Private data
157 
158  //- Local copy of carrier specific heat field
159  // Cp not stored on carrier thermo, but returned as tmp<...>
160  const volScalarField Cp_;
161 
162  //- Local copy of carrier thermal conductivity field
163  // kappa not stored on carrier thermo, but returned as tmp<...>
164  const volScalarField kappa_;
165 
166 
167  // Interpolators for continuous phase fields
168 
169  //- Temperature field interpolator
171 
172  //- Specific heat capacity field interpolator
174 
175  //- Thermal conductivity field interpolator
176  autoPtr<interpolation<scalar>> kappaInterp_;
177 
178  //- Radiation field interpolator
180 
181 
182  public:
183 
184  typedef typename ParcelType::template TrackingData<CloudType>::trackPart
185  trackPart;
186 
187  // Constructors
188 
189  //- Construct from components
190  inline TrackingData
191  (
192  CloudType& cloud,
193  trackPart part = ParcelType::template
195  );
196 
197 
198  // Member functions
199 
200  //- Return access to the locally stored carrier Cp field
201  inline const volScalarField& Cp() const;
202 
203  //- Return access to the locally stored carrier kappa field
204  inline const volScalarField& kappa() const;
205 
206  //- Return const access to the interpolator for continuous
207  // phase temperature field
208  inline const interpolation<scalar>& TInterp() const;
209 
210  //- Return const access to the interpolator for continuous
211  // phase specific heat capacity field
212  inline const interpolation<scalar>& CpInterp() const;
213 
214  //- Return const access to the interpolator for continuous
215  // phase thermal conductivity field
216  inline const interpolation<scalar>& kappaInterp() const;
217 
218  //- Return const access to the interpolator for continuous
219  // radiation field
220  inline const interpolation<scalar>& GInterp() const;
221  };
222 
223 
224 protected:
225 
226  // Protected data
227 
228  // Parcel properties
229 
230  //- Temperature [K]
231  scalar T_;
232 
233  //- Specific heat capacity [J/(kg.K)]
234  scalar Cp_;
235 
236 
237  // Cell-based quantities
238 
239  //- Temperature [K]
240  scalar Tc_;
241 
242  //- Specific heat capacity [J/(kg.K)]
243  scalar Cpc_;
244 
245 
246  // Protected Member Functions
247 
248  //- Calculate new particle temperature
249  template<class TrackData>
250  scalar calcHeatTransfer
251  (
252  TrackData& td,
253  const scalar dt, // timestep
254  const label celli, // owner cell
255  const scalar Re, // Reynolds number
256  const scalar Pr, // Prandtl number - surface
257  const scalar kappa, // Thermal conductivity - surface
258  const scalar NCpW, // Sum of N*Cp*W of emission species
259  const scalar Sh, // explicit particle enthalpy source
260  scalar& dhsTrans, // sensible enthalpy transfer to carrier
261  scalar& Sph // linearised heat transfer coefficient
262  );
263 
264 
265 public:
266 
267  // Static data members
268 
269  //- Runtime type information
270  TypeName("ThermoParcel");
271 
272  //- String representation of properties
274  (
275  ParcelType,
276  " T"
277  + " Cp"
278  );
279 
280 
281  // Constructors
282 
283  //- Construct from owner, position, and cloud owner
284  // Other properties initialised as null
285  inline ThermoParcel
286  (
287  const polyMesh& mesh,
288  const vector& position,
289  const label celli,
290  const label tetFacei,
291  const label tetPtI
292  );
293 
294  //- Construct from components
295  inline ThermoParcel
296  (
297  const polyMesh& mesh,
298  const vector& position,
299  const label celli,
300  const label tetFacei,
301  const label tetPtI,
302  const label typeId,
303  const scalar nParticle0,
304  const scalar d0,
305  const scalar dTarget0,
306  const vector& U0,
307  const vector& f0,
308  const vector& angularMomentum0,
309  const vector& torque0,
310  const constantProperties& constProps
311  );
312 
313  //- Construct from Istream
315  (
316  const polyMesh& mesh,
317  Istream& is,
318  bool readFields = true
319  );
320 
321  //- Construct as a copy
322  ThermoParcel(const ThermoParcel& p);
323 
324  //- Construct as a copy
325  ThermoParcel(const ThermoParcel& p, const polyMesh& mesh);
326 
327  //- Construct and return a (basic particle) clone
328  virtual autoPtr<particle> clone() const
329  {
330  return autoPtr<particle>(new ThermoParcel(*this));
331  }
332 
333  //- Construct and return a (basic particle) clone
334  virtual autoPtr<particle> clone(const polyMesh& mesh) const
335  {
336  return autoPtr<particle>(new ThermoParcel(*this, mesh));
337  }
338 
339  //- Factory class to read-construct particles used for
340  // parallel transfer
341  class iNew
342  {
343  const polyMesh& mesh_;
344 
345  public:
347  iNew(const polyMesh& mesh)
348  :
349  mesh_(mesh)
350  {}
352  autoPtr<ThermoParcel<ParcelType>> operator()(Istream& is) const
353  {
355  (
356  new ThermoParcel<ParcelType>(mesh_, is, true)
357  );
358  }
359  };
360 
361 
362  // Member Functions
363 
364  // Access
365 
366  //- Return const access to temperature
367  inline scalar T() const;
368 
369  //- Return const access to specific heat capacity
370  inline scalar Cp() const;
371 
372  //- Return the parcel sensible enthalpy
373  inline scalar hs() const;
374 
375  //- Return const access to carrier temperature
376  inline scalar Tc() const;
377 
378  //- Return const access to carrier specific heat capacity
379  inline scalar Cpc() const;
380 
381 
382  // Edit
383 
384  //- Return access to temperature
385  inline scalar& T();
386 
387  //- Return access to specific heat capacity
388  inline scalar& Cp();
389 
390 
391  // Main calculation loop
392 
393  //- Set cell values
394  template<class TrackData>
395  void setCellValues
396  (
397  TrackData& td,
398  const scalar dt,
399  const label celli
400  );
401 
402  //- Correct cell values using latest transfer information
403  template<class TrackData>
405  (
406  TrackData& td,
407  const scalar dt,
408  const label celli
409  );
410 
411  //- Calculate surface thermo properties
412  template<class TrackData>
413  void calcSurfaceValues
414  (
415  TrackData& td,
416  const label celli,
417  const scalar T,
418  scalar& Ts,
419  scalar& rhos,
420  scalar& mus,
421  scalar& Pr,
422  scalar& kappas
423  ) const;
424 
425  //- Update parcel properties over the time interval
426  template<class TrackData>
427  void calc
428  (
429  TrackData& td,
430  const scalar dt,
431  const label celli
432  );
433 
434 
435  // I-O
436 
437  //- Read
438  template<class CloudType>
439  static void readFields(CloudType& c);
440 
441  //- Write
442  template<class CloudType>
443  static void writeFields(const CloudType& c);
444 
445 
446  // Ostream Operator
447 
448  friend Ostream& operator<< <ParcelType>
449  (
450  Ostream&,
452  );
453 };
454 
455 
456 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
457 
458 } // End namespace Foam
459 
460 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
461 
462 #include "ThermoParcelI.H"
464 
465 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
466 
467 #ifdef NoRepository
468  #include "ThermoParcel.C"
469 #endif
470 
471 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
472 
473 #endif
474 
475 // ************************************************************************* //
scalar Cpc() const
Return const access to carrier specific heat capacity.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
void calcSurfaceValues(TrackData &td, const label celli, const scalar T, scalar &Ts, scalar &rhos, scalar &mus, scalar &Pr, scalar &kappas) const
Calculate surface thermo properties.
Definition: ThermoParcel.C:95
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
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label celli)
Correct cell values using latest transfer information.
Definition: ThermoParcel.C:67
scalar f0() const
Return const access to the particle scattering factor [].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
scalar Cp_
Specific heat capacity [J/(kg.K)].
Definition: ThermoParcel.H:233
Class to hold thermo particle constant properties.
Definition: ThermoParcel.H:78
scalar calcHeatTransfer(TrackData &td, const scalar dt, const label celli, const scalar Re, const scalar Pr, const scalar kappa, const scalar NCpW, const scalar Sh, scalar &dhsTrans, scalar &Sph)
Calculate new particle temperature.
void calc(TrackData &td, const scalar dt, const label celli)
Update parcel properties over the time interval.
Definition: ThermoParcel.C:138
scalar Sh
Definition: solveChemistry.H:2
scalar TMax() const
Return const access to maximum temperature [K].
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
Definition: ThermoParcel.H:327
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].
scalar T() const
Return const access to temperature.
scalar TMin() const
Return const access to minimum temperature [K].
scalar T0() const
Return const access to the particle initial temperature [K].
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:620
void setCellValues(TrackData &td, const scalar dt, const label celli)
Set cell values.
Definition: ThermoParcel.C:36
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
ParcelType::template TrackingData< CloudType >::trackPart trackPart
Definition: ThermoParcel.H:184
Factory class to read-construct particles used for.
Definition: ThermoParcel.H:340
TypeName("ThermoParcel")
Runtime type information.
scalar T_
Temperature [K].
Definition: ThermoParcel.H:230
dynamicFvMesh & mesh
static void readFields(CloudType &c)
Read.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
ThermoParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: ThermoParcelI.H:75
scalar Cpc_
Specific heat capacity [J/(kg.K)].
Definition: ThermoParcel.H:242
static void writeFields(const CloudType &c)
Write.
scalar Tc() const
Return const access to carrier temperature.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const dimensionedScalar c
Speed of light in a vacuum.
scalar Cp0() const
Return const access to the particle specific heat capacity.
scalar hs() const
Return the parcel sensible enthalpy.
scalar Cp() const
Return const access to specific heat capacity.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
scalar epsilon0() const
Return const access to the particle emissivity [].
Namespace for OpenFOAM.
Thermodynamic parcel class with one/two-way coupling with the continuous phase. Includes Kinematic pa...
Definition: ThermoParcel.H:51
AddToPropertyList(ParcelType," T"+" Cp")
String representation of properties.
scalar Tc_
Temperature [K].
Definition: ThermoParcel.H:239