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-2024 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  //- Class to hold reacting multiphase particle constant properties
86  class constantProperties
87  :
88  public ParcelType::constantProperties
89  {
90  // Private Data
91 
92  //- Devolatilisation activation temperature [K]
94 
95  //- Latent heat of devolatilisation [J/kg]
97 
98  //- Fraction of enthalpy retained by parcel due to surface
99  // reactions
100  demandDrivenEntry<scalar> hRetentionCoeff_;
101 
102 
103  public:
104 
105  // Constructors
106 
107  //- Null constructor
109 
110  //- Copy constructor
112 
113  //- Construct from dictionary
114  constantProperties(const dictionary& parentDict);
115 
116 
117  // Access
118 
119  //- Return const access to the devolatilisation temperature
120  inline scalar TDevol() const;
121 
122  //- Return const access to the latent heat of devolatilisation
123  inline scalar LDevol() const;
124 
125  //- Return const access to the fraction of enthalpy retained by
126  // parcel due to surface reactions
127  inline scalar hRetentionCoeff() const;
128  };
129 
130 
131  //- Use base tracking data
132  typedef typename ParcelType::trackingData trackingData;
133 
134 
135 private:
136 
137  // Private Member Functions
138 
139  //- Return the mixture effective specific heat capacity
140  template<class TrackCloudType>
141  scalar CpEff
142  (
143  TrackCloudType& cloud,
144  trackingData& td,
145  const scalar p,
146  const scalar T,
147  const label idG,
148  const label idL,
149  const label idS
150  ) const;
151 
152  //- Return the mixture effective sensible enthalpy
153  template<class TrackCloudType>
154  scalar hsEff
155  (
156  TrackCloudType& cloud,
157  trackingData& td,
158  const scalar p,
159  const scalar T,
160  const label idG,
161  const label idL,
162  const label idS
163  ) const;
164 
165  //- Return the mixture effective latent heat
166  template<class TrackCloudType>
167  scalar LEff
168  (
169  TrackCloudType& cloud,
170  trackingData& td,
171  const scalar p,
172  const scalar T,
173  const label idG,
174  const label idL,
175  const label idS
176  ) const;
177 
178  //- Update the mass fractions (Y, YGas, YLiquid, YSolid)
179  scalar updateMassFractions
180  (
181  const scalar mass0,
182  const scalarField& dMassGas,
183  const scalarField& dMassLiquid,
184  const scalarField& dMassSolid,
185  const label idG,
186  const label idL,
187  const label idS
188  );
189 
190 
191 protected:
192 
193  // Protected data
194 
195  // Parcel properties
196 
197  //- Initial mass [kg]
198  scalar mass0_;
199 
200  //- Mass fractions of gases []
202 
203  //- Mass fractions of liquids []
205 
206  //- Mass fractions of solids []
208 
209  //- Flag to identify if the particle can devolatilise and combust
210  // Combustion possible only after volatile content falls below
211  // threshold value. States include:
212  // 0 = can devolatilise, cannot combust but can change
213  // 1 = can devolatilise, can combust
214  // -1 = cannot devolatilise or combust, and cannot change
216 
217 
218  // Protected Member Functions
219 
220  //- Calculate Devolatilisation
221  template<class TrackCloudType>
223  (
224  TrackCloudType& cloud,
225  trackingData& td,
226  const scalar dt, // timestep
227  const scalar Ts, // surface temperature
228  const scalar d, // diameter
229  const scalar T, // temperature
230  const scalar mass, // mass
231  const scalar mass0, // mass (initial on injection)
232  const scalarField& YGasEff,// gas component mass fractions
233  const scalarField& YLiquidEff,// liquid component mass fractions
234  const scalarField& YSolidEff,// solid component mass fractions
235  label& canCombust, // 'can combust' flag
236  scalarField& dMassDV, // mass transfer - local to particle
237  scalar& Sh, // explicit particle enthalpy source
238  scalar& N, // flux of species emitted from particle
239  scalar& NCpW, // sum of N*Cp*W of emission species
240  scalarField& Cs // carrier conc. of emission species
241  ) const;
242 
243  //- Calculate surface reactions
244  template<class TrackCloudType>
246  (
247  TrackCloudType& cloud,
248  trackingData& td,
249  const scalar dt, // timestep
250  const scalar d, // diameter
251  const scalar T, // temperature
252  const scalar mass, // mass
253  const label canCombust, // 'can combust' flag
254  const scalar N, // flux of species emitted from particle
255  const scalarField& YMix, // mixture mass fractions
256  const scalarField& YGas, // gas-phase mass fractions
257  const scalarField& YLiquid,// liquid-phase mass fractions
258  const scalarField& YSolid, // solid-phase mass fractions
259  scalarField& dMassSRGas, // gas-phase mass transfer - local
260  scalarField& dMassSRLiquid,// liquid-phase mass transfer - local
261  scalarField& dMassSRSolid, // solid-phase mass transfer - local
262  scalarField& dMassSRCarrier, // carrier phase mass transfer
263  scalar& Sh, // explicit particle enthalpy source
264  scalar& dhsTrans // sensible enthalpy transfer to carrier
265  ) const;
266 
267 
268 public:
269 
270  // Static Data Members
271 
272  //- String representation of properties
274  (
275  ParcelType,
276  " mass0"
277  + " nGas(Y1..YN)"
278  + " nLiquid(Y1..YN)"
279  + " nSolid(Y1..YN)"
280  );
281 
282 
283  // Constructors
284 
285  //- Construct from mesh, coordinates 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  const label facei
295  );
296 
297  //- Construct from a position and a cell, searching for the rest of the
298  // required topology. Other properties are initialised as null.
300  (
301  const polyMesh& mesh,
302  const vector& position,
303  const label celli,
304  label& nLocateBoundaryHits
305  );
306 
307  //- Construct from Istream
308  ReactingMultiphaseParcel(Istream& is, bool readFields = true);
309 
310  //- Construct as a copy
312 
313  //- Construct and return a clone
314  virtual autoPtr<particle> clone() const
315  {
316  return autoPtr<particle>(new ReactingMultiphaseParcel(*this));
317  }
318 
319  //- Construct from Istream and return
321  {
322  return
324  (
326  );
327  }
328 
329 
330  // Member Functions
331 
332  // Access
333 
334  //- Return const access to initial mass [kg]
335  inline scalar mass0() const;
336 
337  //- Return const access to mass fractions of gases
338  inline const scalarField& YGas() const;
339 
340  //- Return const access to mass fractions of liquids
341  inline const scalarField& YLiquid() const;
342 
343  //- Return const access to mass fractions of solids
344  inline const scalarField& YSolid() const;
345 
346  //- Return const access to the canCombust flag
347  inline label canCombust() const;
348 
349 
350  // Edit
351 
352  //- Return access to initial mass [kg]
353  inline scalar& mass0();
354 
355  //- Return access to mass fractions of gases
356  inline scalarField& YGas();
357 
358  //- Return access to mass fractions of liquids
359  inline scalarField& YLiquid();
360 
361  //- Return access to mass fractions of solids
362  inline scalarField& YSolid();
363 
364  //- Return access to the canCombust flag
365  inline label& canCombust();
366 
367 
368  // Main calculation loop
369 
370  //- Set cell values
371  template<class TrackCloudType>
372  void setCellValues(TrackCloudType& cloud, trackingData& td);
373 
374  //- Correct cell values using latest transfer information
375  template<class TrackCloudType>
377  (
378  TrackCloudType& cloud,
379  trackingData& td,
380  const scalar dt
381  );
382 
383  //- Update parcel properties over the time interval
384  template<class TrackCloudType>
385  void calc
386  (
387  TrackCloudType& cloud,
388  trackingData& td,
389  const scalar dt
390  );
391 
392 
393  // I-O
394 
395  //- Read
396  template<class CloudType, class CompositionType>
397  static void readFields
398  (
399  CloudType& c,
400  const CompositionType& compModel
401  );
402 
403  //- Read - no composition
404  template<class CloudType>
405  static void readFields(CloudType& c);
406 
407  //- Write
408  template<class CloudType, class CompositionType>
409  static void writeFields
410  (
411  const CloudType& c,
412  const CompositionType& compModel
413  );
414 
415  //- Read - composition supplied
416  template<class CloudType>
417  static void writeFields(const CloudType& c);
418 
419 
420  // Ostream Operator
421 
422  friend Ostream& operator<< <ParcelType>
423  (
424  Ostream&,
426  );
427 };
428 
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 } // End namespace Foam
433 
434 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
435 
437 
438 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439 
440 #ifdef NoRepository
441  #include "ReactingMultiphaseParcel.C"
442 #endif
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 #endif
447 
448 // ************************************************************************* //
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
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const dimensionedScalar c
Speed of light in a vacuum.
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
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)
volScalarField & p