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-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::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 owner, position, and cloud owner
274  // Other properties initialised as null
276  (
277  const polyMesh& mesh,
278  const vector& position,
279  const label celli,
280  const label tetFacei,
281  const label tetPtI
282  );
283 
284 
285  //- Construct from components
287  (
288  const polyMesh& mesh,
289  const vector& position,
290  const label celli,
291  const label tetFacei,
292  const label tetPtI,
293  const label typeId,
294  const scalar nParticle0,
295  const scalar d0,
296  const scalar dTarget0,
297  const vector& U0,
298  const vector& f0,
299  const vector& angularMomentum0,
300  const vector& torque0,
301  const scalarField& Y0,
302  const scalarField& YGas0,
303  const scalarField& YLiquid0,
304  const scalarField& YSolid0,
305  const constantProperties& constProps
306  );
307 
308  //- Construct from Istream
310  (
311  const polyMesh& mesh,
312  Istream& is,
313  bool readFields = true
314  );
315 
316  //- Construct as a copy
318 
319  //- Construct as a copy
321  (
322  const ReactingMultiphaseParcel& p,
323  const polyMesh& mesh
324  );
325 
326  //- Construct and return a (basic particle) clone
327  virtual autoPtr<particle> clone() const
328  {
329  return autoPtr<particle>(new ReactingMultiphaseParcel(*this));
330  }
331 
332  //- Construct and return a (basic particle) clone
333  virtual autoPtr<particle> clone(const polyMesh& mesh) const
334  {
335  return autoPtr<particle>(new ReactingMultiphaseParcel(*this, mesh));
336  }
337 
338  //- Factory class to read-construct particles used for
339  // parallel transfer
340  class iNew
341  {
342  const polyMesh& mesh_;
343 
344  public:
346  iNew(const polyMesh& mesh)
347  :
348  mesh_(mesh)
349  {}
350 
352  (
353  Istream& is
354  ) const
355  {
357  (
358  new ReactingMultiphaseParcel<ParcelType>(mesh_, is, true)
359  );
360  }
361  };
362 
363 
364  // Member Functions
365 
366  // Access
367 
368  //- Return const access to mass fractions of gases
369  inline const scalarField& YGas() const;
370 
371  //- Return const access to mass fractions of liquids
372  inline const scalarField& YLiquid() const;
373 
374  //- Return const access to mass fractions of solids
375  inline const scalarField& YSolid() const;
376 
377  //- Return const access to the canCombust flag
378  inline label canCombust() const;
379 
380 
381  // Edit
382 
383  //- Return access to mass fractions of gases
384  inline scalarField& YGas();
385 
386  //- Return access to mass fractions of liquids
387  inline scalarField& YLiquid();
388 
389  //- Return access to mass fractions of solids
390  inline scalarField& YSolid();
391 
392  //- Return access to the canCombust flag
393  inline label& canCombust();
394 
395 
396  // Main calculation loop
397 
398  //- Set cell values
399  template<class TrackData>
400  void setCellValues
401  (
402  TrackData& td,
403  const scalar dt,
404  const label celli
405  );
406 
407  //- Correct cell values using latest transfer information
408  template<class TrackData>
410  (
411  TrackData& td,
412  const scalar dt,
413  const label celli
414  );
415 
416  //- Update parcel properties over the time interval
417  template<class TrackData>
418  void calc
419  (
420  TrackData& td,
421  const scalar dt,
422  const label celli
423  );
424 
425 
426  // I-O
427 
428  //- Read
429  template<class CloudType, class CompositionType>
430  static void readFields
431  (
432  CloudType& c,
433  const CompositionType& compModel
434  );
435 
436  //- Read - no composition
437  template<class CloudType>
438  static void readFields(CloudType& c);
439 
440  //- Write
441  template<class CloudType, class CompositionType>
442  static void writeFields
443  (
444  const CloudType& c,
445  const CompositionType& compModel
446  );
447 
448  //- Read - composition supplied
449  template<class CloudType>
450  static void writeFields(const CloudType& c);
451 
452 
453  // Ostream Operator
454 
455  friend Ostream& operator<< <ParcelType>
456  (
457  Ostream&,
459  );
460 };
461 
462 
463 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
464 
465 } // End namespace Foam
466 
467 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
468 
470 
471 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
472 
473 #ifdef NoRepository
474  #include "ReactingMultiphaseParcel.C"
475 #endif
476 
477 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
478 
479 #endif
480 
481 // ************************************************************************* //
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
AddToPropertyList(ParcelType," nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
String representation of properties.
label canCombust_
Flag to identify if the particle can devolatilise and combust.
scalar hRetentionCoeff() const
Return const access to the fraction of enthalpy retained by.
scalar Sh
Definition: solveChemistry.H:2
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
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label celli)
Correct cell values using latest transfer information.
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:620
void calc(TrackData &td, const scalar dt, const label celli)
Update parcel properties over the time interval.
scalar TDevol() const
Return const access to the devolatilisation temperature.
scalarList Y0(nSpecie, 0.0)
label canCombust() const
Return const access to the canCombust flag.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
scalarField YGas_
Mass fractions of gases [].
dynamicFvMesh & mesh
ReactingMultiphaseParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
const scalarField & YGas() const
Return const access to mass fractions of gases.
scalarField YSolid_
Mass fractions of solids [].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalar LDevol() const
Return const access to the latent heat of devolatilisation.
scalarField YLiquid_
Mass fractions of liquids [].
const scalarField & YSolid() const
Return const access to mass fractions of solids.
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
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.
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.
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
TypeName("ReactingMultiphaseParcel")
Runtime type information.
Factory class to read-construct particles used for.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Class to hold reacting multiphase particle constant properties.
Namespace for OpenFOAM.