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-2023 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>
51 class ReactingMultiphaseParcel;
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  //- Initial mass [kg]
202  scalar mass0_;
203 
204  //- Mass fractions of gases []
206 
207  //- Mass fractions of liquids []
209 
210  //- Mass fractions of solids []
212 
213  //- Flag to identify if the particle can devolatilise and combust
214  // Combustion possible only after volatile content falls below
215  // threshold value. States include:
216  // 0 = can devolatilise, cannot combust but can change
217  // 1 = can devolatilise, can combust
218  // -1 = cannot devolatilise or combust, and cannot change
220 
221 
222  // Protected Member Functions
223 
224  //- Calculate Devolatilisation
225  template<class TrackCloudType>
227  (
228  TrackCloudType& cloud,
229  trackingData& td,
230  const scalar dt, // timestep
231  const scalar Ts, // surface temperature
232  const scalar d, // diameter
233  const scalar T, // temperature
234  const scalar mass, // mass
235  const scalar mass0, // mass (initial on injection)
236  const scalarField& YGasEff,// gas component mass fractions
237  const scalarField& YLiquidEff,// liquid component mass fractions
238  const scalarField& YSolidEff,// solid component mass fractions
239  label& canCombust, // 'can combust' flag
240  scalarField& dMassDV, // mass transfer - local to particle
241  scalar& Sh, // explicit particle enthalpy source
242  scalar& N, // flux of species emitted from particle
243  scalar& NCpW, // sum of N*Cp*W of emission species
244  scalarField& Cs // carrier conc. of emission species
245  ) const;
246 
247  //- Calculate surface reactions
248  template<class TrackCloudType>
250  (
251  TrackCloudType& cloud,
252  trackingData& td,
253  const scalar dt, // timestep
254  const scalar d, // diameter
255  const scalar T, // temperature
256  const scalar mass, // mass
257  const label canCombust, // 'can combust' flag
258  const scalar N, // flux of species emitted from particle
259  const scalarField& YMix, // mixture mass fractions
260  const scalarField& YGas, // gas-phase mass fractions
261  const scalarField& YLiquid,// liquid-phase mass fractions
262  const scalarField& YSolid, // solid-phase mass fractions
263  scalarField& dMassSRGas, // gas-phase mass transfer - local
264  scalarField& dMassSRLiquid,// liquid-phase mass transfer - local
265  scalarField& dMassSRSolid, // solid-phase mass transfer - local
266  scalarField& dMassSRCarrier, // carrier phase mass transfer
267  scalar& Sh, // explicit particle enthalpy source
268  scalar& dhsTrans // sensible enthalpy transfer to carrier
269  ) const;
270 
271 
272 public:
273 
274  // Static Data Members
275 
276  //- String representation of properties
278  (
279  ParcelType,
280  " mass0"
281  + " nGas(Y1..YN)"
282  + " nLiquid(Y1..YN)"
283  + " nSolid(Y1..YN)"
284  );
285 
286 
287  // Constructors
288 
289  //- Construct from mesh, coordinates and topology
290  // Other properties initialised as null
292  (
293  const polyMesh& mesh,
294  const barycentric& coordinates,
295  const label celli,
296  const label tetFacei,
297  const label tetPti,
298  const label facei
299  );
300 
301  //- Construct from a position and a cell, searching for the rest of the
302  // required topology. Other properties are initialised as null.
304  (
305  const polyMesh& mesh,
306  const vector& position,
307  const label celli,
308  label& nLocateBoundaryHits
309  );
310 
311  //- Construct from Istream
312  ReactingMultiphaseParcel(Istream& is, bool readFields = true);
313 
314  //- Construct as a copy
316 
317  //- Construct and return a clone
318  virtual autoPtr<particle> clone() const
319  {
320  return autoPtr<particle>(new ReactingMultiphaseParcel(*this));
321  }
322 
323  //- Construct from Istream and return
325  {
326  return
328  (
330  );
331  }
332 
333 
334  // Member Functions
335 
336  // Access
337 
338  //- Return const access to initial mass [kg]
339  inline scalar mass0() const;
340 
341  //- Return const access to mass fractions of gases
342  inline const scalarField& YGas() const;
343 
344  //- Return const access to mass fractions of liquids
345  inline const scalarField& YLiquid() const;
346 
347  //- Return const access to mass fractions of solids
348  inline const scalarField& YSolid() const;
349 
350  //- Return const access to the canCombust flag
351  inline label canCombust() const;
352 
353 
354  // Edit
355 
356  //- Return access to initial mass [kg]
357  inline scalar& mass0();
358 
359  //- Return access to mass fractions of gases
360  inline scalarField& YGas();
361 
362  //- Return access to mass fractions of liquids
363  inline scalarField& YLiquid();
364 
365  //- Return access to mass fractions of solids
366  inline scalarField& YSolid();
367 
368  //- Return access to the canCombust flag
369  inline label& canCombust();
370 
371 
372  // Main calculation loop
373 
374  //- Set cell values
375  template<class TrackCloudType>
376  void setCellValues(TrackCloudType& cloud, trackingData& td);
377 
378  //- Correct cell values using latest transfer information
379  template<class TrackCloudType>
381  (
382  TrackCloudType& cloud,
383  trackingData& td,
384  const scalar dt
385  );
386 
387  //- Update parcel properties over the time interval
388  template<class TrackCloudType>
389  void calc
390  (
391  TrackCloudType& cloud,
392  trackingData& td,
393  const scalar dt
394  );
395 
396 
397  // I-O
398 
399  //- Read
400  template<class CloudType, class CompositionType>
401  static void readFields
402  (
403  CloudType& c,
404  const CompositionType& compModel
405  );
406 
407  //- Read - no composition
408  template<class CloudType>
409  static void readFields(CloudType& c);
410 
411  //- Write
412  template<class CloudType, class CompositionType>
413  static void writeFields
414  (
415  const CloudType& c,
416  const CompositionType& compModel
417  );
418 
419  //- Read - composition supplied
420  template<class CloudType>
421  static void writeFields(const CloudType& c);
422 
423 
424  // Ostream Operator
425 
426  friend Ostream& operator<< <ParcelType>
427  (
428  Ostream&,
430  );
431 };
432 
433 
434 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
435 
436 } // End namespace Foam
437 
438 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439 
441 
442 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
443 
444 #ifdef NoRepository
445  #include "ReactingMultiphaseParcel.C"
446 #endif
447 
448 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449 
450 #endif
451 
452 // ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to hold reacting multiphase particle constant properties.
scalar hRetentionCoeff() const
Return const access to the fraction of enthalpy retained by.
scalar TDevol() const
Return const access to the devolatilisation temperature.
scalar LDevol() const
Return const access to the latent heat of devolatilisation.
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.
scalarField YLiquid_
Mass fractions of liquids [].
scalarField YSolid_
Mass fractions of solids [].
label canCombust_
Flag to identify if the particle can devolatilise and combust.
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.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
AddToPropertyList(ParcelType, " mass0"+" nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
String representation of properties.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
scalarField YGas_
Mass fractions of gases [].
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.
label canCombust() const
Return const access to the canCombust flag.
ParcelType::trackingData trackingData
Use base tracking data.
const scalarField & YGas() const
Return const access to mass fractions of gases.
virtual autoPtr< particle > clone() const
Construct and return a clone.
scalar mass0() const
Return const access to initial mass [kg].
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
static autoPtr< ReactingMultiphaseParcel > New(Istream &is)
Construct from Istream and return.
const scalarField & YSolid() const
Return const access to mass fractions of solids.
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
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
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
TemplateName(FvFaceCellWave)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p