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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::ReactingMultiphaseParcel
29 
30 Group
31  grpLagrangianIntermediateParcels
32 
33 Description
34  Multiphase variant of the reacting parcel class with one/two-way coupling
35  with the continuous phase.
36 
37 SourceFiles
38  ReactingMultiphaseParcelI.H
39  ReactingMultiphaseParcel.C
40  ReactingMultiphaseParcelIO.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef ReactingMultiphaseParcel_H
45 #define ReactingMultiphaseParcel_H
46 
47 #include "particle.H"
48 #include "SLGThermo.H"
49 #include "demandDrivenEntry.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 template<class ParcelType> class ReactingMultiphaseParcel;
57 
58 template<class ParcelType>
59 Ostream& operator<<(Ostream&, const ReactingMultiphaseParcel<ParcelType>&);
60 
61 /*---------------------------------------------------------------------------*\
62  Class ReactingMultiphaseParcel Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class ParcelType>
67 :
68  public ParcelType
69 {
70 public:
71 
72  //- Size in bytes of the fields
73  static const std::size_t sizeofFields;
74 
75 
76  // IDs of phases in ReactingParcel phase list (Y)
77 
78  static const label GAS;
79  static const label LIQ;
80  static const label SLD;
81 
82 
83  //- Class to hold reacting multiphase particle constant properties
84  class constantProperties
85  :
86  public ParcelType::constantProperties
87  {
88  // Private data
89 
90  //- Devolatilisation activation temperature [K]
92 
93  //- Latent heat of devolatilisation [J/kg]
95 
96  //- Fraction of enthalpy retained by parcel due to surface
97  // reactions
98  demandDrivenEntry<scalar> hRetentionCoeff_;
99 
100 
101  public:
102 
103  // Constructors
104 
105  //- Null constructor
107 
108  //- Copy constructor
110 
111  //- Construct from dictionary
112  constantProperties(const dictionary& parentDict);
113 
114 
115  // Access
116 
117  //- Return const access to the devolatilisation temperature
118  inline scalar TDevol() const;
119 
120  //- Return const access to the latent heat of devolatilisation
121  inline scalar LDevol() const;
122 
123  //- Return const access to the fraction of enthalpy retained by
124  // parcel due to surface reactions
125  inline scalar hRetentionCoeff() const;
126  };
127 
128 
129  //- Use base tracking data
130  typedef typename ParcelType::trackingData trackingData;
131 
132 
133 private:
134 
135  // Private Member Functions
136 
137  //- Return the mixture effective specific heat capacity
138  template<class TrackCloudType>
139  scalar CpEff
140  (
141  TrackCloudType& cloud,
142  trackingData& td,
143  const scalar p,
144  const scalar T,
145  const label idG,
146  const label idL,
147  const label idS
148  ) const;
149 
150  //- Return the mixture effective sensible enthalpy
151  template<class TrackCloudType>
152  scalar HsEff
153  (
154  TrackCloudType& cloud,
155  trackingData& td,
156  const scalar p,
157  const scalar T,
158  const label idG,
159  const label idL,
160  const label idS
161  ) const;
162 
163  //- Return the mixture effective latent heat
164  template<class TrackCloudType>
165  scalar LEff
166  (
167  TrackCloudType& cloud,
168  trackingData& td,
169  const scalar p,
170  const scalar T,
171  const label idG,
172  const label idL,
173  const label idS
174  ) const;
175 
176  //- Update the mass fractions (Y, YGas, YLiquid, YSolid)
177  scalar updateMassFractions
178  (
179  const scalar mass0,
180  const scalarField& dMassGas,
181  const scalarField& dMassLiquid,
182  const scalarField& dMassSolid
183  );
184 
185 
186 protected:
187 
188  // Protected data
189 
190  // Parcel properties
191 
192  //- Mass fractions of gases []
194 
195  //- Mass fractions of liquids []
197 
198  //- Mass fractions of solids []
200 
201  //- Flag to identify if the particle can devolatilise and combust
202  // Combustion possible only after volatile content falls below
203  // threshold value. States include:
204  // 0 = can devolatilise, cannot combust but can change
205  // 1 = can devolatilise, can combust
206  // -1 = cannot devolatilise or combust, and cannot change
207  label canCombust_;
208 
209 
210  // Protected Member Functions
211 
212 
213  //- Return change of volume due to mass exchange
214  template<class TrackCloudType>
215  scalar updatedDeltaVolume
216  (
217  TrackCloudType& cloud,
218  const scalarField& dMassGas,
219  const scalarField& dMassLiquid,
220  const scalarField& dMassSolid,
221  const label idG,
222  const label idL,
223  const label idS,
224  const scalar p,
225  const scalar T
226  );
227 
228  //- Calculate Devolatilisation
229  template<class TrackCloudType>
231  (
232  TrackCloudType& cloud,
233  trackingData& td,
234  const scalar dt, // timestep
235  const scalar age, // age
236  const scalar Ts, // surface temperature
237  const scalar d, // diameter
238  const scalar T, // temperature
239  const scalar mass, // mass
240  const scalar mass0, // mass (initial on injection)
241  const scalarField& YGasEff,// gas component mass fractions
242  const scalarField& YLiquidEff,// liquid component mass fractions
243  const scalarField& YSolidEff,// solid component mass fractions
244  label& canCombust, // 'can combust' flag
245  scalarField& dMassDV, // mass transfer - local to particle
246  scalar& Sh, // explicit particle enthalpy source
247  scalar& N, // flux of species emitted from particle
248  scalar& NCpW, // sum of N*Cp*W of emission species
249  scalarField& Cs // carrier conc. of emission species
250  ) const;
251 
252  //- Calculate surface reactions
253  template<class TrackCloudType>
255  (
256  TrackCloudType& cloud,
257  trackingData& td,
258  const scalar dt, // timestep
259  const scalar d, // diameter
260  const scalar Re, // Re
261  const scalar nu, // nu
262  const scalar T, // temperature
263  const scalar mass, // mass
264  const label canCombust, // 'can combust' flag
265  const scalar N, // flux of species emitted from particle
266  const scalarField& YMix, // mixture mass fractions
267  const scalarField& YGas, // gas-phase mass fractions
268  const scalarField& YLiquid,// liquid-phase mass fractions
269  const scalarField& YSolid, // solid-phase mass fractions
270  scalarField& dMassSRGas, // gas-phase mass transfer - local
271  scalarField& dMassSRLiquid,// liquid-phase mass transfer - local
272  scalarField& dMassSRSolid, // solid-phase mass transfer - local
273  scalarField& dMassSRCarrier, // carrier phase mass transfer
274  scalar& Sh, // explicit particle enthalpy source
275  scalar& dhsTrans // sensible enthalpy transfer to carrier
276  ) const;
277 
278 
279 public:
280 
281  // Static data members
282 
283  //- Runtime type information
284  TypeName("ReactingMultiphaseParcel");
285 
286  //- String representation of properties
288  (
289  ParcelType,
290  " nGas(Y1..YN)"
291  + " nLiquid(Y1..YN)"
292  + " nSolid(Y1..YN)"
293  );
294 
295 
296  // Constructors
297 
298  //- Construct from mesh, position and topology
299  // Other properties initialised as null
301  (
302  const polyMesh& mesh,
303  const barycentric& coordinates,
304  const label celli,
305  const label tetFacei,
306  const label tetPti
307  );
308 
309  //- Construct from a position and a cell, searching for the rest of the
310  // required topology. Other properties are initialised as null.
312  (
313  const polyMesh& mesh,
314  const vector& position,
315  const label celli
316  );
317 
318  //- Construct from components
320  (
321  const polyMesh& mesh,
322  const barycentric& coordinates,
323  const label celli,
324  const label tetFacei,
325  const label tetPti,
326  const label typeId,
327  const scalar nParticle0,
328  const scalar d0,
329  const scalar dTarget0,
330  const vector& U0,
331  const vector& f0,
332  const vector& angularMomentum0,
333  const vector& torque0,
334  const scalarField& Y0,
335  const scalarField& YGas0,
336  const scalarField& YLiquid0,
337  const scalarField& YSolid0,
338  const constantProperties& constProps
339  );
340 
341  //- Construct from Istream
343  (
344  const polyMesh& mesh,
345  Istream& is,
346  bool readFields = true,
347  bool newFormat = true
348  );
349 
350  //- Construct as a copy
352 
353  //- Construct as a copy
355  (
357  const polyMesh& mesh
358  );
359 
360  //- Return a (basic particle) clone
361  virtual autoPtr<particle> clone() const
362  {
363  return particle::Clone(*this);
364  }
365 
366  //- Return a (basic particle) clone
367  virtual autoPtr<particle> clone(const polyMesh& mesh) const
368  {
369  return particle::Clone(*this, mesh);
370  }
371 
372  //- Factory class to read-construct particles (for parallel transfer)
373  class iNew
374  {
375  const polyMesh& mesh_;
376 
377  public:
378 
379  iNew(const polyMesh& mesh)
380  :
381  mesh_(mesh)
382  {}
383 
384  autoPtr<ReactingMultiphaseParcel<ParcelType>> operator()
385  (
386  Istream& is
387  ) const
388  {
389  return autoPtr<ReactingMultiphaseParcel<ParcelType>>
390  (
391  new ReactingMultiphaseParcel<ParcelType>(mesh_, is, true)
392  );
393  }
394  };
395 
396 
397  // Member Functions
398 
399  // Access
400 
401  //- Return const access to mass fractions of gases
402  inline const scalarField& YGas() const;
403 
404  //- Return const access to mass fractions of liquids
405  inline const scalarField& YLiquid() const;
406 
407  //- Return const access to mass fractions of solids
408  inline const scalarField& YSolid() const;
409 
410  //- Return const access to the canCombust flag
411  inline label canCombust() const;
412 
413 
414  // Edit
415 
416  //- Return access to mass fractions of gases
417  inline scalarField& YGas();
418 
419  //- Return access to mass fractions of liquids
420  inline scalarField& YLiquid();
421 
422  //- Return access to mass fractions of solids
423  inline scalarField& YSolid();
424 
425  //- Return access to the canCombust flag
426  inline label& canCombust();
427 
428 
429  // Main calculation loop
430 
431  //- Set cell values
432  template<class TrackCloudType>
433  void setCellValues(TrackCloudType& cloud, trackingData& td);
434 
435  //- Correct cell values using latest transfer information
436  template<class TrackCloudType>
438  (
439  TrackCloudType& cloud,
440  trackingData& td,
441  const scalar dt
442  );
443 
444  //- Update parcel properties over the time interval
445  template<class TrackCloudType>
446  void calc
447  (
448  TrackCloudType& cloud,
449  trackingData& td,
450  const scalar dt
451  );
452 
454  // I-O
455 
456  //- Read - composition supplied
457  template<class CloudType, class CompositionType>
458  static void readFields
459  (
460  CloudType& c,
461  const CompositionType& compModel
462  );
463 
464  //- Read - no composition
465  template<class CloudType>
466  static void readFields(CloudType& c);
467 
468  //- Write - composition supplied
469  template<class CloudType, class CompositionType>
470  static void writeFields
471  (
472  const CloudType& c,
473  const CompositionType& compModel
474  );
475 
476  //- Read - no composition
477  template<class CloudType>
478  static void writeFields(const CloudType& c);
479 
480  //- Write individual parcel properties to stream
481  void writeProperties
482  (
483  Ostream& os,
484  const wordRes& filters,
485  const word& delim,
486  const bool namesOnly = false
487  ) const;
488 
489  //- Read particle fields as objects from the obr registry
490  // - no composition
491  template<class CloudType>
492  static void readObjects
493  (
494  CloudType& c,
495  const objectRegistry& obr
496  );
497 
498  //- Read particle fields as objects from the obr registry
499  template<class CloudType, class CompositionType>
500  static void readObjects
501  (
502  CloudType& c,
503  const CompositionType& compModel,
504  const objectRegistry& obr
505  );
506 
507  //- Write particle fields as objects into the obr registry
508  // - no composition
509  template<class CloudType>
510  static void writeObjects
511  (
512  const CloudType& c,
513  objectRegistry& obr
514  );
515 
516  //- Write particle fields as objects into the obr registry
517  template<class CloudType, class CompositionType>
518  static void writeObjects
519  (
520  const CloudType& c,
521  const CompositionType& compModel,
522  objectRegistry& obr
523  );
524 
525 
526  // Ostream Operator
527 
528  friend Ostream& operator<< <ParcelType>
529  (
530  Ostream&,
532  );
533 };
534 
535 
536 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
537 
538 } // End namespace Foam
539 
540 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
541 
543 
544 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
545 
546 #ifdef NoRepository
547  #include "ReactingMultiphaseParcel.C"
548 #endif
549 
550 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
551 
552 #endif
553 
554 // ************************************************************************* //
label canCombust() const
Return const access to the canCombust flag.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
scalarList Y0(nSpecie, Zero)
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label canCombust_
Flag to identify if the particle can devolatilise and combust.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void calcDevolatilisation(TrackCloudType &cloud, trackingData &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.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: POSIX.C:1063
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
scalar LDevol() const
Return const access to the latent heat of devolatilisation.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
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.
dynamicFvMesh & mesh
const scalarField & YGas() const
Return const access to mass fractions of gases.
void calcSurfaceReactions(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar d, const scalar Re, const scalar nu, 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.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
scalar TDevol() const
Return const access to the devolatilisation temperature.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
AddToPropertyList(ParcelType, " nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
String representation of properties.
scalarField YSolid_
Mass fractions of solids [].
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
OBJstream os(runTime.globalPath()/outputName)
static const std::size_t sizeofFields
Size in bytes of the fields.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
scalarField YLiquid_
Mass fractions of liquids [].
const Vector< label > N(dict.get< Vector< label >>("N"))
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
PtrList< coordinateSystem > coordinates(solidRegions.size())
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...
scalar updatedDeltaVolume(TrackCloudType &cloud, const scalarField &dMassGas, const scalarField &dMassLiquid, const scalarField &dMassSolid, const label idG, const label idL, const label idS, const scalar p, const scalar T)
Return change of volume due to mass exchange.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
volScalarField & p
Registry of regIOobjects.
TypeName("ReactingMultiphaseParcel")
Runtime type information.
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:67
volScalarField & nu
const scalarField & Cs
Class to hold reacting multiphase particle constant properties.
Namespace for OpenFOAM.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
const scalarField & YSolid() const
Return const access to mass fractions of solids.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.