ReactingMultiphaseParcel.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-2017 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::ReactingMultiphaseParcel
26 
27 Description
28  Multiphase variant of the reacting parcel class with one/two-way coupling
29  with the continuous phase.
30 
31 SourceFiles
32  ReactingMultiphaseParcelI.H
33  ReactingMultiphaseParcel.C
34  ReactingMultiphaseParcelIO.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef ReactingMultiphaseParcel_H
39 #define ReactingMultiphaseParcel_H
40 
41 #include "particle.H"
42 #include "SLGThermo.H"
43 #include "demandDrivenEntry.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 template<class ParcelType>
52 
53 template<class ParcelType>
54 Ostream& operator<<
55 (
56  Ostream&,
58 );
59 
60 /*---------------------------------------------------------------------------*\
61  Class ReactingMultiphaseParcel Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 template<class ParcelType>
66 :
67  public ParcelType
68 {
69  // Private data
70 
71  //- Size in bytes of the fields
72  static const std::size_t sizeofFields_;
73 
74 
75 public:
76 
77  // IDs of phases in ReacingParcel phase list (Y)
78 
79  static const label GAS;
80  static const label LIQ;
81  static const label SLD;
82 
83 
84  //- Class to hold reacting multiphase particle constant properties
85  class constantProperties
86  :
87  public ParcelType::constantProperties
88  {
89  // Private data
90 
91  //- Devolatilisation activation temperature [K]
93 
94  //- Latent heat of devolatilisation [J/kg]
96 
97  //- Fraction of enthalpy retained by parcel due to surface
98  // reactions
99  demandDrivenEntry<scalar> hRetentionCoeff_;
100 
101 
102  public:
103 
104  // Constructors
105 
106  //- Null constructor
108 
109  //- Copy constructor
111 
112  //- Construct from dictionary
113  constantProperties(const dictionary& parentDict);
114 
115 
116  // Access
117 
118  //- Return const access to the devolatilisation temperature
119  inline scalar TDevol() const;
120 
121  //- Return const access to the latent heat of devolatilisation
122  inline scalar LDevol() const;
123 
124  //- Return const access to the fraction of enthalpy retained by
125  // parcel due to surface reactions
126  inline scalar hRetentionCoeff() const;
127  };
128 
129 
130 private:
131 
132  // Private Member Functions
133 
134  //- Return the mixture effective specific heat capacity
135  template<class TrackData>
136  scalar CpEff
137  (
138  TrackData& td,
139  const scalar p,
140  const scalar T,
141  const label idG,
142  const label idL,
143  const label idS
144  ) const;
145 
146  //- Return the mixture effective sensible enthalpy
147  template<class TrackData>
148  scalar HsEff
149  (
150  TrackData& td,
151  const scalar p,
152  const scalar T,
153  const label idG,
154  const label idL,
155  const label idS
156  ) const;
157 
158  //- Return the mixture effective latent heat
159  template<class TrackData>
160  scalar LEff
161  (
162  TrackData& td,
163  const scalar p,
164  const scalar T,
165  const label idG,
166  const label idL,
167  const label idS
168  ) const;
169 
170  //- Update the mass fractions (Y, YGas, YLiquid, YSolid)
171  scalar updateMassFractions
172  (
173  const scalar mass0,
174  const scalarField& dMassGas,
175  const scalarField& dMassLiquid,
176  const scalarField& dMassSolid
177  );
178 
179 
180 protected:
181 
182  // Protected data
183 
184  // Parcel properties
185 
186  //- Mass fractions of gases []
188 
189  //- Mass fractions of liquids []
191 
192  //- Mass fractions of solids []
194 
195  //- Flag to identify if the particle can devolatilise and combust
196  // Combustion possible only after volatile content falls below
197  // threshold value. States include:
198  // 0 = can devolatilise, cannot combust but can change
199  // 1 = can devolatilise, can combust
200  // -1 = cannot devolatilise or combust, and cannot change
202 
203 
204  // Protected Member Functions
205 
206  //- Calculate Devolatilisation
207  template<class TrackData>
209  (
210  TrackData& td,
211  const scalar dt, // timestep
212  const scalar age, // age
213  const scalar Ts, // surface temperature
214  const scalar d, // diameter
215  const scalar T, // temperature
216  const scalar mass, // mass
217  const scalar mass0, // mass (initial on injection)
218  const scalarField& YGasEff,// gas component mass fractions
219  const scalarField& YLiquidEff,// liquid component mass fractions
220  const scalarField& YSolidEff,// solid component mass fractions
221  label& canCombust, // 'can combust' flag
222  scalarField& dMassDV, // mass transfer - local to particle
223  scalar& Sh, // explicit particle enthalpy source
224  scalar& N, // flux of species emitted from particle
225  scalar& NCpW, // sum of N*Cp*W of emission species
226  scalarField& Cs // carrier conc. of emission species
227  ) const;
228 
229  //- Calculate surface reactions
230  template<class TrackData>
232  (
233  TrackData& td,
234  const scalar dt, // timestep
235  const label celli, // owner cell
236  const scalar d, // diameter
237  const scalar T, // temperature
238  const scalar mass, // mass
239  const label canCombust, // 'can combust' flag
240  const scalar N, // flux of species emitted from particle
241  const scalarField& YMix, // mixture mass fractions
242  const scalarField& YGas, // gas-phase mass fractions
243  const scalarField& YLiquid,// liquid-phase mass fractions
244  const scalarField& YSolid, // solid-phase mass fractions
245  scalarField& dMassSRGas, // gas-phase mass transfer - local
246  scalarField& dMassSRLiquid,// liquid-phase mass transfer - local
247  scalarField& dMassSRSolid, // solid-phase mass transfer - local
248  scalarField& dMassSRCarrier, // carrier phase mass transfer
249  scalar& Sh, // explicit particle enthalpy source
250  scalar& dhsTrans // sensible enthalpy transfer to carrier
251  ) const;
252 
253 
254 public:
255 
256  // Static data members
257 
258  //- Runtime type information
259  TypeName("ReactingMultiphaseParcel");
260 
261  //- String representation of properties
263  (
264  ParcelType,
265  " nGas(Y1..YN)"
266  + " nLiquid(Y1..YN)"
267  + " nSolid(Y1..YN)"
268  );
269 
270 
271  // Constructors
272 
273  //- Construct from mesh, position and topology
274  // Other properties initialised as null
276  (
277  const polyMesh& mesh,
278  const barycentric& coordinates,
279  const label celli,
280  const label tetFacei,
281  const label tetPti
282  );
283 
284  //- Construct from a position and a cell, searching for the rest of the
285  // required topology. Other properties are initialised as null.
287  (
288  const polyMesh& mesh,
289  const vector& position,
290  const label celli
291  );
292 
293  //- Construct from components
295  (
296  const polyMesh& mesh,
297  const barycentric& coordinates,
298  const label celli,
299  const label tetFacei,
300  const label tetPti,
301  const label typeId,
302  const scalar nParticle0,
303  const scalar d0,
304  const scalar dTarget0,
305  const vector& U0,
306  const vector& f0,
307  const vector& angularMomentum0,
308  const vector& torque0,
309  const scalarField& Y0,
310  const scalarField& YGas0,
311  const scalarField& YLiquid0,
312  const scalarField& YSolid0,
313  const constantProperties& constProps
314  );
315 
316  //- Construct from Istream
318  (
319  const polyMesh& mesh,
320  Istream& is,
321  bool readFields = true
322  );
323 
324  //- Construct as a copy
326 
327  //- Construct as a copy
329  (
330  const ReactingMultiphaseParcel& p,
331  const polyMesh& mesh
332  );
333 
334  //- Construct and return a (basic particle) clone
335  virtual autoPtr<particle> clone() const
336  {
337  return autoPtr<particle>(new ReactingMultiphaseParcel(*this));
338  }
339 
340  //- Construct and return a (basic particle) clone
341  virtual autoPtr<particle> clone(const polyMesh& mesh) const
342  {
343  return autoPtr<particle>(new ReactingMultiphaseParcel(*this, mesh));
344  }
345 
346  //- Factory class to read-construct particles used for
347  // parallel transfer
348  class iNew
349  {
350  const polyMesh& mesh_;
351 
352  public:
354  iNew(const polyMesh& mesh)
355  :
356  mesh_(mesh)
357  {}
358 
360  (
361  Istream& is
362  ) const
363  {
365  (
366  new ReactingMultiphaseParcel<ParcelType>(mesh_, is, true)
367  );
368  }
369  };
370 
371 
372  // Member Functions
373 
374  // Access
375 
376  //- Return const access to mass fractions of gases
377  inline const scalarField& YGas() const;
378 
379  //- Return const access to mass fractions of liquids
380  inline const scalarField& YLiquid() const;
381 
382  //- Return const access to mass fractions of solids
383  inline const scalarField& YSolid() const;
384 
385  //- Return const access to the canCombust flag
386  inline label canCombust() const;
387 
388 
389  // Edit
390 
391  //- Return access to mass fractions of gases
392  inline scalarField& YGas();
393 
394  //- Return access to mass fractions of liquids
395  inline scalarField& YLiquid();
396 
397  //- Return access to mass fractions of solids
398  inline scalarField& YSolid();
399 
400  //- Return access to the canCombust flag
401  inline label& canCombust();
402 
403 
404  // Main calculation loop
405 
406  //- Set cell values
407  template<class TrackData>
408  void setCellValues
409  (
410  TrackData& td,
411  const scalar dt,
412  const label celli
413  );
414 
415  //- Correct cell values using latest transfer information
416  template<class TrackData>
418  (
419  TrackData& td,
420  const scalar dt,
421  const label celli
422  );
423 
424  //- Update parcel properties over the time interval
425  template<class TrackData>
426  void calc
427  (
428  TrackData& td,
429  const scalar dt,
430  const label celli
431  );
432 
433 
434  // I-O
435 
436  //- Read
437  template<class CloudType, class CompositionType>
438  static void readFields
439  (
440  CloudType& c,
441  const CompositionType& compModel
442  );
443 
444  //- Read - no composition
445  template<class CloudType>
446  static void readFields(CloudType& c);
447 
448  //- Write
449  template<class CloudType, class CompositionType>
450  static void writeFields
451  (
452  const CloudType& c,
453  const CompositionType& compModel
454  );
455 
456  //- Read - composition supplied
457  template<class CloudType>
458  static void writeFields(const CloudType& c);
459 
460 
461  // Ostream Operator
462 
463  friend Ostream& operator<< <ParcelType>
464  (
465  Ostream&,
467  );
468 };
469 
470 
471 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
472 
473 } // End namespace Foam
474 
475 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476 
478 
479 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
480 
481 #ifdef NoRepository
482  #include "ReactingMultiphaseParcel.C"
483 #endif
484 
485 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
486 
487 #endif
488 
489 // ************************************************************************* //
label canCombust() const
Return const access to the canCombust flag.
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
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
label canCombust_
Flag to identify if the particle can devolatilise and combust.
void setCellValues(TrackData &td, const scalar dt, const label celli)
Set cell values.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:739
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label celli)
Correct cell values using latest transfer information.
void calc(TrackData &td, const scalar dt, const label celli)
Update parcel properties over the time interval.
scalar LDevol() const
Return const access to the latent heat of devolatilisation.
scalarList Y0(nSpecie, 0.0)
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
scalarField YGas_
Mass fractions of gases [].
scalar hRetentionCoeff() const
Return const access to the fraction of enthalpy retained by.
dynamicFvMesh & mesh
const scalarField & YGas() const
Return const access to mass fractions of gases.
scalar TDevol() const
Return const access to the devolatilisation temperature.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
AddToPropertyList(ParcelType, " nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
String representation of properties.
scalarField YSolid_
Mass fractions of solids [].
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalarField YLiquid_
Mass fractions of liquids [].
PtrList< coordinateSystem > coordinates(solidRegions.size())
void calcDevolatilisation(TrackData &td, const scalar dt, const scalar age, const scalar Ts, const scalar d, const scalar T, const scalar mass, const scalar mass0, const scalarField &YGasEff, const scalarField &YLiquidEff, const scalarField &YSolidEff, label &canCombust, scalarField &dMassDV, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs) const
Calculate Devolatilisation.
const dimensionedScalar c
Speed of light in a vacuum.
friend Ostream & operator(Ostream &, const ReactingMultiphaseParcel< ParcelType > &)
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase...
label N
Definition: createFields.H:22
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
void calcSurfaceReactions(TrackData &td, const scalar dt, const label celli, const scalar d, const scalar T, const scalar mass, const label canCombust, const scalar N, const scalarField &YMix, const scalarField &YGas, const scalarField &YLiquid, const scalarField &YSolid, scalarField &dMassSRGas, scalarField &dMassSRLiquid, scalarField &dMassSRSolid, scalarField &dMassSRCarrier, scalar &Sh, scalar &dhsTrans) const
Calculate surface reactions.
TypeName("ReactingMultiphaseParcel")
Runtime type information.
Factory class to read-construct particles used for.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
Class to hold reacting multiphase particle constant properties.
Namespace for OpenFOAM.
const scalarField & YSolid() const
Return const access to mass fractions of solids.