ReactingMultiphaseParcel.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-2021 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 "fluidThermo.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 ReactingMultiphaseParcelName Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 
66 
67 /*---------------------------------------------------------------------------*\
68  Class ReactingMultiphaseParcel Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class ParcelType>
73 :
74  public ParcelType,
75  public ReactingMultiphaseParcelName
76 {
77  // Private Data
78 
79  //- Size in bytes of the fields
80  static const std::size_t sizeofFields_;
81 
82 
83 public:
84 
85  // IDs of phases in ReactingParcel phase list (Y)
86 
87  static const label GAS;
88  static const label LIQ;
89  static const label SLD;
90 
91 
92  //- Class to hold reacting multiphase particle constant properties
93  class constantProperties
94  :
95  public ParcelType::constantProperties
96  {
97  // Private Data
98 
99  //- Devolatilisation activation temperature [K]
101 
102  //- Latent heat of devolatilisation [J/kg]
104 
105  //- Fraction of enthalpy retained by parcel due to surface
106  // reactions
107  demandDrivenEntry<scalar> hRetentionCoeff_;
108 
109 
110  public:
111 
112  // Constructors
113 
114  //- Null constructor
116 
117  //- Copy constructor
119 
120  //- Construct from dictionary
121  constantProperties(const dictionary& parentDict);
122 
123 
124  // Access
125 
126  //- Return const access to the devolatilisation temperature
127  inline scalar TDevol() const;
128 
129  //- Return const access to the latent heat of devolatilisation
130  inline scalar LDevol() const;
131 
132  //- Return const access to the fraction of enthalpy retained by
133  // parcel due to surface reactions
134  inline scalar hRetentionCoeff() const;
135  };
136 
137 
138  //- Use base tracking data
139  typedef typename ParcelType::trackingData trackingData;
140 
141 
142 private:
143 
144  // Private Member Functions
145 
146  //- Return the mixture effective specific heat capacity
147  template<class TrackCloudType>
148  scalar CpEff
149  (
150  TrackCloudType& cloud,
151  trackingData& td,
152  const scalar p,
153  const scalar T,
154  const label idG,
155  const label idL,
156  const label idS
157  ) const;
158 
159  //- Return the mixture effective sensible enthalpy
160  template<class TrackCloudType>
161  scalar HsEff
162  (
163  TrackCloudType& cloud,
164  trackingData& td,
165  const scalar p,
166  const scalar T,
167  const label idG,
168  const label idL,
169  const label idS
170  ) const;
171 
172  //- Return the mixture effective latent heat
173  template<class TrackCloudType>
174  scalar LEff
175  (
176  TrackCloudType& cloud,
177  trackingData& td,
178  const scalar p,
179  const scalar T,
180  const label idG,
181  const label idL,
182  const label idS
183  ) const;
184 
185  //- Update the mass fractions (Y, YGas, YLiquid, YSolid)
186  scalar updateMassFractions
187  (
188  const scalar mass0,
189  const scalarField& dMassGas,
190  const scalarField& dMassLiquid,
191  const scalarField& dMassSolid
192  );
193 
194 
195 protected:
196 
197  // Protected data
198 
199  // Parcel properties
200 
201  //- Mass fractions of gases []
203 
204  //- Mass fractions of liquids []
206 
207  //- Mass fractions of solids []
209 
210  //- Flag to identify if the particle can devolatilise and combust
211  // Combustion possible only after volatile content falls below
212  // threshold value. States include:
213  // 0 = can devolatilise, cannot combust but can change
214  // 1 = can devolatilise, can combust
215  // -1 = cannot devolatilise or combust, and cannot change
217 
218 
219  // Protected Member Functions
220 
221  //- Calculate Devolatilisation
222  template<class TrackCloudType>
224  (
225  TrackCloudType& cloud,
226  trackingData& td,
227  const scalar dt, // timestep
228  const scalar Ts, // surface temperature
229  const scalar d, // diameter
230  const scalar T, // temperature
231  const scalar mass, // mass
232  const scalar mass0, // mass (initial on injection)
233  const scalarField& YGasEff,// gas component mass fractions
234  const scalarField& YLiquidEff,// liquid component mass fractions
235  const scalarField& YSolidEff,// solid component mass fractions
236  label& canCombust, // 'can combust' flag
237  scalarField& dMassDV, // mass transfer - local to particle
238  scalar& Sh, // explicit particle enthalpy source
239  scalar& N, // flux of species emitted from particle
240  scalar& NCpW, // sum of N*Cp*W of emission species
241  scalarField& Cs // carrier conc. of emission species
242  ) const;
243 
244  //- Calculate surface reactions
245  template<class TrackCloudType>
247  (
248  TrackCloudType& cloud,
249  trackingData& td,
250  const scalar dt, // timestep
251  const scalar d, // diameter
252  const scalar T, // temperature
253  const scalar mass, // mass
254  const label canCombust, // 'can combust' flag
255  const scalar N, // flux of species emitted from particle
256  const scalarField& YMix, // mixture mass fractions
257  const scalarField& YGas, // gas-phase mass fractions
258  const scalarField& YLiquid,// liquid-phase mass fractions
259  const scalarField& YSolid, // solid-phase mass fractions
260  scalarField& dMassSRGas, // gas-phase mass transfer - local
261  scalarField& dMassSRLiquid,// liquid-phase mass transfer - local
262  scalarField& dMassSRSolid, // solid-phase mass transfer - local
263  scalarField& dMassSRCarrier, // carrier phase mass transfer
264  scalar& Sh, // explicit particle enthalpy source
265  scalar& dhsTrans // sensible enthalpy transfer to carrier
266  ) const;
267 
268 
269 public:
270 
271  // Static Data Members
272 
273  //- String representation of properties
275  (
276  ParcelType,
277  " nGas(Y1..YN)"
278  + " nLiquid(Y1..YN)"
279  + " nSolid(Y1..YN)"
280  );
281 
282 
283  // Constructors
284 
285  //- Construct from mesh, position and topology
286  // Other properties initialised as null
288  (
289  const polyMesh& mesh,
290  const barycentric& coordinates,
291  const label celli,
292  const label tetFacei,
293  const label tetPti
294  );
295 
296  //- Construct from a position and a cell, searching for the rest of the
297  // required topology. Other properties are initialised as null.
299  (
300  const polyMesh& mesh,
301  const vector& position,
302  const label celli
303  );
304 
305  //- Construct from Istream
307  (
308  const polyMesh& mesh,
309  Istream& is,
310  bool readFields = true
311  );
312 
313  //- Construct as a copy
315 
316  //- Construct as a copy
318  (
320  const polyMesh& mesh
321  );
322 
323  //- Construct and return a (basic particle) clone
324  virtual autoPtr<particle> clone() const
325  {
326  return autoPtr<particle>(new ReactingMultiphaseParcel(*this));
327  }
328 
329  //- Construct and return a (basic particle) clone
330  virtual autoPtr<particle> clone(const polyMesh& mesh) const
331  {
332  return autoPtr<particle>(new ReactingMultiphaseParcel(*this, mesh));
333  }
334 
335  //- Factory class to read-construct particles used for
336  // parallel transfer
337  class iNew
338  {
339  const polyMesh& mesh_;
340 
341  public:
343  iNew(const polyMesh& mesh)
344  :
345  mesh_(mesh)
346  {}
347 
349  (
350  Istream& is
351  ) const
352  {
354  (
355  new ReactingMultiphaseParcel<ParcelType>(mesh_, is, true)
356  );
357  }
358  };
359 
360 
361  // Member Functions
362 
363  // Access
364 
365  //- Return const access to mass fractions of gases
366  inline const scalarField& YGas() const;
367 
368  //- Return const access to mass fractions of liquids
369  inline const scalarField& YLiquid() const;
370 
371  //- Return const access to mass fractions of solids
372  inline const scalarField& YSolid() const;
373 
374  //- Return const access to the canCombust flag
375  inline label canCombust() const;
376 
377 
378  // Edit
379 
380  //- Return access to mass fractions of gases
381  inline scalarField& YGas();
382 
383  //- Return access to mass fractions of liquids
384  inline scalarField& YLiquid();
385 
386  //- Return access to mass fractions of solids
387  inline scalarField& YSolid();
388 
389  //- Return access to the canCombust flag
390  inline label& canCombust();
391 
392 
393  // Main calculation loop
394 
395  //- Set cell values
396  template<class TrackCloudType>
397  void setCellValues(TrackCloudType& cloud, trackingData& td);
398 
399  //- Correct cell values using latest transfer information
400  template<class TrackCloudType>
402  (
403  TrackCloudType& cloud,
404  trackingData& td,
405  const scalar dt
406  );
407 
408  //- Update parcel properties over the time interval
409  template<class TrackCloudType>
410  void calc
411  (
412  TrackCloudType& cloud,
413  trackingData& td,
414  const scalar dt
415  );
416 
417 
418  // I-O
419 
420  //- Read
421  template<class CloudType, class CompositionType>
422  static void readFields
423  (
424  CloudType& c,
425  const CompositionType& compModel
426  );
427 
428  //- Read - no composition
429  template<class CloudType>
430  static void readFields(CloudType& c);
431 
432  //- Write
433  template<class CloudType, class CompositionType>
434  static void writeFields
435  (
436  const CloudType& c,
437  const CompositionType& compModel
438  );
439 
440  //- Read - composition supplied
441  template<class CloudType>
442  static void writeFields(const CloudType& c);
443 
444 
445  // Ostream Operator
446 
447  friend Ostream& operator<< <ParcelType>
448  (
449  Ostream&,
451  );
452 };
453 
454 
455 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
456 
457 } // End namespace Foam
458 
459 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
460 
462 
463 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
464 
465 #ifdef NoRepository
466  #include "ReactingMultiphaseParcel.C"
467 #endif
468 
469 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
470 
471 #endif
472 
473 // ************************************************************************* //
label canCombust() const
Return const access to the canCombust flag.
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:156
label canCombust_
Flag to identify if the particle can devolatilise and combust.
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:753
scalar LDevol() const
Return const access to the latent heat of devolatilisation.
fvMesh & mesh
const dimensionedScalar c
Speed of light in a vacuum.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
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.
ParcelType::trackingData trackingData
Use base tracking data.
scalarField YGas_
Mass fractions of gases [].
scalar hRetentionCoeff() const
Return const access to the fraction of enthalpy retained by.
const scalarField & YGas() const
Return const access to mass fractions of gases.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
scalar TDevol() const
Return const access to the devolatilisation temperature.
TemplateName(FvFaceCellWave)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
AddToPropertyList(ParcelType, " nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
String representation of properties.
void calcSurfaceReactions(TrackCloudType &cloud, trackingData &td, const scalar dt, 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.
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 [].
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:9
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:76
volScalarField & p
void calcDevolatilisation(TrackCloudType &cloud, trackingData &td, const scalar dt, 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.
Factory class to read-construct particles used for.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Class to hold reacting multiphase particle constant properties.
Namespace for OpenFOAM.
const scalarField & YSolid() const
Return const access to mass fractions of solids.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.