ReactingHeterogeneousParcel.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) 2018-2024 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::ReactingHeterogeneousParcel
28 
29 Group
30  grpLagrangianIntermediateParcels
31 
32 Description
33  Reacting heterogeneous Parcel
34 
35 SourceFiles
36  ReactingHeterogeneousParcelI.H
37  ReactingHeterogeneousParcel.C
38  ReactingHeterogeneousParcelIO.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef ReactingHeterogeneousParcel_H
43 #define ReactingHeterogeneousParcel_H
44 
45 #include "demandDrivenEntry.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 template<class ParcelType>
54 
55 template<class ParcelType>
56 Ostream& operator<<
57 (
58  Ostream&,
60 );
61 
62 /*---------------------------------------------------------------------------*\
63  Class ReactingHeterogeneousParcel Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 template<class ParcelType>
68 :
69  public ParcelType
70 {
71 public:
72 
73  //- Size in bytes of the fields
74  static const std::size_t sizeofFields;
75 
76  //- Class to hold reacting particle constant properties
77  class constantProperties
78  :
79  public ParcelType::constantProperties
80  {
81  // Private data
82 
83  //- Fraction of enthalpy retained by parcel due to surface
84  // reactions
85  demandDrivenEntry<scalar> hRetentionCoeff_;
86 
87  public:
88 
89  // Constructors
90 
91  //- Null constructor
93 
94  //- Copy constructor
96 
97  //- Construct from dictionary
98  constantProperties(const dictionary& parentDict);
99 
100 
101  // Access
102 
103  //- Return const access to the fraction of enthalpy retained by
104  // parcel due to surface reactions
105  inline scalar hRetentionCoeff() const;
106  };
107 
108 
109  //- Use base tracking data
110  typedef typename ParcelType::trackingData trackingData;
111 
112 
113 private:
114 
115  // Private Member Functions
116 
117  //- Return the mixture effective specific heat capacity
118  template<class TrackCloudType>
119  scalar CpEff
120  (
121  TrackCloudType& cloud,
122  trackingData& td,
123  const scalar p,
124  const scalar T,
125  const label idS
126  ) const;
127 
128  //- Return the mixture effective sensible enthalpy
129  template<class TrackCloudType>
130  scalar HsEff
131  (
132  TrackCloudType& cloud,
133  trackingData& td,
134  const scalar p,
135  const scalar T,
136  const label idS
137  ) const;
138 
139  //- Return the mixture effective latent heat
140  template<class TrackCloudType>
141  scalar LEff
142  (
143  TrackCloudType& cloud,
144  trackingData& td,
145  const scalar p,
146  const scalar T,
147  const label idS
148  ) const;
149 
150 
151 protected:
152 
153  // Protected data
154 
155  // Parcel properties
156 
157  //- Progress variables for reactions
158  scalarField F_;
159 
160  //- Flag to identify if the particle can devolatilise and combust
161  // Combustion possible only after volatile content falls below
162  // threshold value. States include:
163  // 0 = can combust but can change
164  // 1 = can devolatilise, can combust
165  // -1 = cannot devolatilise or combust, and cannot change
166  label canCombust_;
167 
168 
169  // Protected Member Functions
170 
171  //- Return change of volume due to mass exchange
172  template<class TrackCloudType>
173  scalar updatedDeltaVolume
174  (
175  TrackCloudType& cloud,
176  const scalarField& dMass,
177  const scalar p,
178  const scalar T
179  );
180 
181 
182  //- Calculate surface reactions
183  template<class TrackCloudType>
185  (
186  TrackCloudType& cloud,
187  trackingData& td,
188  const scalar dt, // timestep
189  const scalar Res, // Re
190  const scalar nu, // nu
191  const scalar d, // diameter
192  const scalar T, // temperature
193  const scalar mass, // mass
194  const label canCombust, // 'can combust' flag
195  const scalar N, // flux of species emitted from particle
196  scalar& NCpW,
197  const scalarField& YSolid, // solid-phase mass fractions
198  scalarField& F, // progress of each reaction
199  scalarField& dMassSRSolid, // solid-phase mass transfer - local
200  scalarField& dMassSRCarrier, // carrier phase mass transfer
201  scalar& Sh, // explicit particle enthalpy source
202  scalar& dhsTrans // sensible enthalpy transfer to carrier
203  ) const;
204 
205 
206 public:
207 
208  // Static data members
209 
210  //- Runtime type information
211  TypeName("ReactingHeterogeneousParcel");
212 
213  //- String representation of properties
215  (
216  ParcelType,
217  + " nReactions(F1..FN)"
218  );
219 
220 
221  // Constructors
222 
223  //- Construct from mesh, position and topology
224  // Other properties initialised as null
226  (
227  const polyMesh& mesh,
228  const barycentric& coordinates,
229  const label celli,
230  const label tetFacei,
231  const label tetPti
232  );
233 
234  //- Construct from a position and a cell, searching for the rest of the
235  // required topology. Other properties are initialised as null.
237  (
238  const polyMesh& mesh,
239  const vector& position,
240  const label celli
241  );
242 
243  //- Construct from components
245  (
246  const polyMesh& mesh,
247  const barycentric& coordinates,
248  const label celli,
249  const label tetFacei,
250  const label tetPti,
251  const label typeId,
252  const scalar nParticle0,
253  const scalar d0,
254  const scalar dTarget0,
255  const vector& U0,
256  const vector& f0,
257  const vector& angularMomentum0,
258  const vector& torque0,
259  const scalarField& Y,
260  const scalarField& F,
261  const constantProperties& constProps
262  );
263 
264  //- Construct from Istream
266  (
267  const polyMesh& mesh,
268  Istream& is,
269  bool readFields = true,
270  bool newFormat = true
271  );
272 
273  //- Construct as a copy
275 
276  //- Construct as a copy
278  (
280  const polyMesh& mesh
281  );
282 
283  //- Return a (basic particle) clone
284  virtual autoPtr<particle> clone() const
285  {
286  return particle::Clone(*this);
287  }
288 
289  //- Return a (basic particle) clone
290  virtual autoPtr<particle> clone(const polyMesh& mesh) const
291  {
292  return particle::Clone(*this, mesh);
293  }
294 
295  //- Factory class to read-construct particles (for parallel transfer)
296  class iNew
297  {
298  const polyMesh& mesh_;
299 
300  public:
301 
302  iNew(const polyMesh& mesh)
303  :
304  mesh_(mesh)
305  {}
306 
307  autoPtr<ReactingHeterogeneousParcel<ParcelType>> operator()
308  (
309  Istream& is
310  ) const
311  {
312  return autoPtr<ReactingHeterogeneousParcel<ParcelType>>
313  (
314  new ReactingHeterogeneousParcel<ParcelType>
315  (mesh_, is, true)
316  );
317  }
318  };
319 
320 
321  // Member Functions
322 
323  // Access
324 
325  //- Return const access to F
326  inline const scalarField& F() const;
327 
328  //- Return const access to the canCombust flag
329  inline label canCombust() const;
330 
331 
332  // Edit
333 
334  //- Return access to F
335  inline scalarField& F();
336 
337  //- Return access to the canCombust flag
338  inline label& canCombust();
339 
341  // Main calculation loop
342 
343 
344  //- Update parcel properties over the time interval
345  template<class TrackCloudType>
346  void calc
347  (
348  TrackCloudType& cloud,
349  trackingData& td,
350  const scalar dt
351  );
352 
353 
354  // I-O
355 
356  //- Read - composition supplied
357  template<class CloudType, class CompositionType>
358  static void readFields
359  (
361  const CompositionType& compModel
362  );
363 
364  //- Read - no composition
365  template<class CloudType>
366  static void readFields(CloudType& c);
367 
368  //- Write - composition supplied
369  template<class CloudType, class CompositionType>
370  static void writeFields
371  (
372  const CloudType& c,
373  const CompositionType& compModel
374  );
375 
376  //- Read - no composition
377  template<class CloudType>
378  static void writeFields(const CloudType& c);
379 
380  //- Write individual parcel properties to stream
381  void writeProperties
382  (
383  Ostream& os,
384  const wordRes& filters,
385  const word& delim,
386  const bool namesOnly
387  ) const;
388 
389  //- Read particle fields as objects from the obr registry
390  // - no composition
391  template<class CloudType>
392  static void readObjects
393  (
394  CloudType& c,
395  const objectRegistry& obr
396  );
397 
398  //- Read particle fields as objects from the obr registry
399  template<class CloudType, class CompositionType>
400  static void readObjects
401  (
402  CloudType& c,
403  const CompositionType& compModel,
404  const objectRegistry& obr
405  );
406 
407  //- Write particle fields as objects into the obr registry
408  // - no composition
409  template<class CloudType>
410  static void writeObjects
411  (
412  const CloudType& c,
413  objectRegistry& obr
414  );
415 
416  //- Write particle fields as objects into the obr registry
417  template<class CloudType, class CompositionType>
418  static void writeObjects
419  (
420  const CloudType& c,
421  const CompositionType& compModel,
422  objectRegistry& obr
423  );
424 
425 
426  // Ostream Operator
427 
428  friend Ostream& operator<< <ParcelType>
429  (
430  Ostream&,
432  );
433 };
434 
435 
436 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
437 
438 } // End namespace Foam
439 
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441 
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 #ifdef NoRepository
448 #endif
449 
450 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
451 
452 #endif
453 
454 // ************************************************************************* //
TypeName("ReactingHeterogeneousParcel")
Runtime type information.
ParcelType::trackingData trackingData
Use base tracking data.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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
friend Ostream & operator(Ostream &, const ReactingHeterogeneousParcel< ParcelType > &)
void calcHeterogeneousReactions(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Res, const scalar nu, const scalar d, const scalar T, const scalar mass, const label canCombust, const scalar N, scalar &NCpW, const scalarField &YSolid, scalarField &F, scalarField &dMassSRSolid, scalarField &dMassSRCarrier, scalar &Sh, scalar &dhsTrans) const
Calculate surface reactions.
static const std::size_t sizeofFields
Size in bytes of the fields.
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
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
ReactingHeterogeneousParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
scalar hRetentionCoeff() const
Return const access to the fraction of enthalpy retained by.
scalar updatedDeltaVolume(TrackCloudType &cloud, const scalarField &dMass, const scalar p, const scalar T)
Return change of volume due to mass exchange.
const scalarField & F() const
Return const access to F.
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Class to hold reacting particle constant properties.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
label canCombust() const
Return const access to the canCombust flag.
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
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
label canCombust_
Flag to identify if the particle can devolatilise and combust.
const Vector< label > N(dict.get< Vector< label >>("N"))
AddToPropertyList(ParcelType,+" nReactions(F1..FN)")
String representation of properties.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
PtrList< coordinateSystem > coordinates(solidRegions.size())
PtrList< volScalarField > & Y
const dimensionedScalar c
Speed of light in a vacuum.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
scalarField F_
Progress variables for reactions.
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.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
volScalarField & nu
Namespace for OpenFOAM.