ReactingParcel.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::ReactingParcel
29 
30 Group
31  grpLagrangianIntermediateParcels
32 
33 Description
34  Reacting parcel class with one/two-way coupling with the continuous
35  phase.
36 
37 SourceFiles
38  ReactingParcelI.H
39  ReactingParcel.C
40  ReactingParcelIO.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef ReactingParcel_H
45 #define ReactingParcel_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>
57 class ReactingParcel;
58 
59 template<class ParcelType>
60 Ostream& operator<<
61 (
62  Ostream&,
64 );
65 
66 
67 /*---------------------------------------------------------------------------*\
68  Class ReactingParcel Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class ParcelType>
72 class ReactingParcel
73 :
74  public ParcelType
75 {
76 public:
77 
78  //- Size in bytes of the fields
79  static const std::size_t sizeofFields;
80 
81 
82 
83  //- Class to hold reacting parcel constant properties
84  class constantProperties
85  :
86  public ParcelType::constantProperties
87  {
88  public:
89 
90  //- Type of activity coefficient models
91  enum volumeUpdateType
92  {
93  mConstRho,
97  };
98 
99  private:
100 
101  // Private Data
102 
103  //- Minimum pressure [Pa]
105 
106  //- Method to update particle rho and diameter
107  //- 0: constant rho
108  //- 1: constant volume
109  //- 2: recalculation of rho and volume based on the lost mass
110  demandDrivenEntry<bool> constantVolume_;
111 
112  //- Method to update vol and rho
113  demandDrivenEntry<label> volUpdateType_;
114 
115 
116  public:
117 
118  // Constructors
119 
120  //- Null constructor
122 
123  //- Copy constructor
125 
126  //- Construct from dictionary
127  constantProperties(const dictionary& parentDict);
128 
129 
130  // Access
131 
132  //- Return const access to the minimum pressure
133  inline scalar pMin() const;
134 
135  //- Return const access to the constant volume flag
136  inline bool constantVolume() const;
137 
138  //- Return const access to the constant volume flag
139  inline label volUpdateType() const;
140  };
141 
142 
143  class trackingData
144  :
145  public ParcelType::trackingData
146  {
147  private:
148 
149  // Private data
150 
151  // Interpolators for continuous phase fields
152 
153  //- Interpolator for continuous phase pressure field
155 
156 
157  // Cached continuous phase properties
158 
159  //- Pressure [Pa]
160  scalar pc_;
161 
163  public:
164 
165  typedef typename ParcelType::trackingData::trackPart trackPart;
166 
167  // Constructors
168 
169  //- Construct from components
170  template<class TrackCloudType>
171  inline trackingData
172  (
173  const TrackCloudType& cloud,
174  trackPart part = ParcelType::trackingData::tpLinearTrack
175  );
176 
177 
178  // Member functions
179 
180  //- Return const access to the interpolator for continuous phase
181  // pressure field
182  inline const interpolation<scalar>& pInterp() const;
183 
184  //- Return the continuous phase pressure
185  inline scalar pc() const;
186 
187  //- Access the continuous phase pressure
188  inline scalar& pc();
189  };
190 
191 
192 protected:
193 
194  // Protected data
195 
196  // Parcel properties
197 
198  //- Initial mass [kg]
199  scalar mass0_;
200 
201  //- Mass fractions of mixture []
202  scalarField Y_;
203 
204 
205  // Protected Member Functions
206 
207  //- Return change of volume due to mass exchange
208  template<class TrackCloudType>
209  scalar updatedDeltaVolume
210  (
211  TrackCloudType& cloud,
212  const scalarField& dMass,
213  const scalar p,
214  const scalar T
215  );
216 
217 
218  //- Calculate Phase change
219  template<class TrackCloudType>
220  void calcPhaseChange
221  (
222  TrackCloudType& cloud,
223  trackingData& td,
224  const scalar dt, // timestep
225  const scalar Re, // Reynolds number
226  const scalar Pr, // Prandtl number
227  const scalar Ts, // Surface temperature
228  const scalar nus, // Surface kinematic viscosity
229  const scalar d, // diameter
230  const scalar T, // temperature
231  const scalar mass, // mass
232  const scalar rho, // density
233  const label idPhase, // id of phase involved in phase change
234  const scalar YPhase, // total mass fraction
235  const scalarField& YLiq, // liquid component mass fractions
236  const scalarField& YSol, // solid component mass fractions
237  scalarField& dMassPC, // mass transfer - local to parcel
238  scalar& Sh, // explicit parcel enthalpy source
239  scalar& N, // flux of species emitted from parcel
240  scalar& NCpW, // sum of N*Cp*W of emission species
241  scalarField& Cs // carrier conc. of emission species
242  );
243 
244  //- Update mass fraction
245  scalar updateMassFraction
246  (
247  const scalar mass0,
248  const scalarField& dMass,
249  scalarField& Y
250  ) const;
251 
252 
253 public:
254 
255  // Static data members
256 
257  //- Runtime type information
258  TypeName("ReactingParcel");
259 
260  //- String representation of properties
262  (
263  ParcelType,
264  " mass0"
265  + " nPhases(Y1..YN)"
266  );
267 
268 
269  // Constructors
270 
271  //- Construct from mesh, coordinates and topology
272  // Other properties initialised as null
273  inline ReactingParcel
274  (
275  const polyMesh& mesh,
276  const barycentric& coordinates,
277  const label celli,
278  const label tetFacei,
279  const label tetPti
280  );
281 
282  //- Construct from a position and a cell, searching for the rest of the
283  // required topology. Other properties are initialised as null.
284  inline ReactingParcel
285  (
286  const polyMesh& mesh,
287  const vector& position,
288  const label celli
289  );
290 
291  //- Construct from components
292  inline ReactingParcel
293  (
294  const polyMesh& mesh,
295  const barycentric& coordinates,
296  const label celli,
297  const label tetFacei,
298  const label tetPti,
299  const label typeId,
300  const scalar nParticle0,
301  const scalar d0,
302  const scalar dTarget0,
303  const vector& U0,
304  const vector& f0,
305  const vector& angularMomentum0,
306  const vector& torque0,
307  const scalarField& Y0,
308  const constantProperties& constProps
309  );
310 
311  //- Construct from Istream
313  (
314  const polyMesh& mesh,
315  Istream& is,
316  bool readFields = true,
317  bool newFormat = true
318  );
319 
320  //- Construct as a copy
322  (
323  const ReactingParcel& p,
324  const polyMesh& mesh
325  );
326 
327  //- Construct as a copy
329 
330  //- Return a (basic particle) clone
331  virtual autoPtr<particle> clone() const
332  {
333  return particle::Clone(*this);
334  }
335 
336  //- Return a (basic particle) clone
337  virtual autoPtr<particle> clone(const polyMesh& mesh) const
338  {
339  return particle::Clone(*this, mesh);
340  }
341 
342  //- Factory class to read-construct particles (for parallel transfer)
343  class iNew
344  {
345  const polyMesh& mesh_;
346 
347  public:
348 
349  iNew(const polyMesh& mesh)
350  :
351  mesh_(mesh)
352  {}
353 
354  autoPtr<ReactingParcel<ParcelType>> operator()(Istream& is) const
355  {
356  return autoPtr<ReactingParcel<ParcelType>>
357  (
358  new ReactingParcel<ParcelType>(mesh_, is, true)
359  );
360  }
361  };
362 
363 
364  // Member Functions
365 
366  // Access
367 
368  //- Return const access to initial mass [kg]
369  inline scalar mass0() const;
370 
371  //- Return const access to mass fractions of mixture []
372  inline const scalarField& Y() const;
373 
374  //- Return const access to mass fractions of gases
375  // Note: for compatibilty only - returns Y()
376  inline const scalarField& YGas() const;
377 
378  //- Return const access to mass fractions of liquids
379  // Note: for compatibilty only - returns Y()
380  inline const scalarField& YLiquid() const;
381 
382  //- Return const access to mass fractions of solids
383  // Note: for compatibilty only - returns Y()
384  inline const scalarField& YSolid() const;
385 
386 
387  // Edit
388 
389  //- Return access to initial mass [kg]
390  inline scalar& mass0();
391 
392  //- Return access to mass fractions of mixture []
393  inline scalarField& Y();
394 
395 
396  // Main calculation loop
397 
398  //- Set cell values
399  template<class TrackCloudType>
400  void setCellValues(TrackCloudType& cloud, trackingData& td);
402  //- Correct cell values using latest transfer information
403  template<class TrackCloudType>
405  (
406  TrackCloudType& cloud,
407  trackingData& td,
408  const scalar dt
409  );
410 
411  //- Correct surface values due to emitted species
412  template<class TrackCloudType>
414  (
415  TrackCloudType& cloud,
416  trackingData& td,
417  const scalar T,
418  const scalarField& Cs,
419  scalar& rhos,
420  scalar& mus,
421  scalar& Prs,
422  scalar& kappas
423  );
424 
425  //- Update parcel properties over the time interval
426  template<class TrackCloudType>
427  void calc
428  (
429  TrackCloudType& cloud,
430  trackingData& td,
431  const scalar dt
432  );
433 
434 
435  // I-O
436 
437  //- Read - composition supplied
438  template<class CloudType, class CompositionType>
439  static void readFields
440  (
441  CloudType& c,
442  const CompositionType& compModel
443  );
444 
445  //- Read - no composition
446  template<class CloudType>
447  static void readFields(CloudType& c);
448 
449  //- Write - composition supplied
450  template<class CloudType, class CompositionType>
451  static void writeFields
452  (
453  const CloudType& c,
454  const CompositionType& compModel
455  );
456 
457  //- Write - no composition
458  template<class CloudType>
459  static void writeFields(const CloudType& c);
460 
461  //- Write individual parcel properties to stream
462  void writeProperties
463  (
464  Ostream& os,
465  const wordRes& filters,
466  const word& delim,
467  const bool namesOnly = false
468  ) const;
469 
470  //- Read particle fields as objects from the obr registry
471  // - no composition
472  template<class CloudType>
473  static void readObjects
474  (
475  CloudType& c,
476  const objectRegistry& obr
477  );
478 
479  //- Read particle fields as objects from the obr registry
480  template<class CloudType, class CompositionType>
481  static void readObjects
482  (
483  CloudType& c,
484  const CompositionType& compModel,
485  const objectRegistry& obr
486  );
487 
488  //- Write particle fields as objects into the obr registry
489  // - no composition
490  template<class CloudType>
491  static void writeObjects
492  (
493  const CloudType& c,
494  objectRegistry& obr
495  );
496 
497  //- Write particle fields as objects into the obr registry
498  template<class CloudType, class CompositionType>
499  static void writeObjects
500  (
501  const CloudType& c,
502  const CompositionType& compModel,
503  objectRegistry& obr
504  );
505 
506 
507  // Ostream Operator
508 
509  friend Ostream& operator<< <ParcelType>
510  (
511  Ostream&,
513  );
514 };
515 
516 
517 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
518 
519 } // End namespace Foam
520 
521 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
522 
523 #include "ReactingParcelI.H"
525 
526 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
527 
528 #ifdef NoRepository
529  #include "ReactingParcel.C"
530 #endif
531 
532 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
533 
534 #endif
535 
536 // ************************************************************************* //
Class to hold reacting parcel constant properties.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
volumeUpdateType
Type of activity coefficient models.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
scalar updatedDeltaVolume(TrackCloudType &cloud, const scalarField &dMass, const scalar p, const scalar T)
Return change of volume due to mass exchange.
scalarList Y0(nSpecie, Zero)
autoPtr< ReactingParcel< ParcelType > > operator()(Istream &is) const
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)
scalar pc() const
Return the continuous phase pressure.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar updateMassFraction(const scalar mass0, const scalarField &dMass, scalarField &Y) const
Update mass fraction.
scalarField Y_
Mass fractions of mixture [].
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
const interpolation< scalar > & pInterp() const
Return const access to the interpolator for continuous phase.
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
bool constantVolume() const
Return const access to the constant volume flag.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
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
scalar pMin() const
Return const access to the minimum pressure.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
label volUpdateType() const
Return const access to the constant volume flag.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
ParcelType::trackingData::trackPart trackPart
void correctSurfaceValues(TrackCloudType &cloud, trackingData &td, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas)
Correct surface values due to emitted species.
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
TypeName("ReactingParcel")
Runtime type information.
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)
AddToPropertyList(ParcelType, " mass0"+" nPhases(Y1..YN)")
String representation of properties.
iNew(const polyMesh &mesh)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static const std::size_t sizeofFields
Size in bytes of the fields.
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
const Vector< label > N(dict.get< Vector< label >>("N"))
const scalarField & YSolid() const
Return const access to mass fractions of solids.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
PtrList< coordinateSystem > coordinates(solidRegions.size())
scalar mass0() const
Return const access to initial mass [kg].
const dimensionedScalar c
Speed of light in a vacuum.
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
const scalarField & YGas() const
Return const access to mass fractions of gases.
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
Reacting parcel class with one/two-way coupling with the continuous phase.
volScalarField & p
Registry of regIOobjects.
scalar mass0_
Initial mass [kg].
trackingData(const TrackCloudType &cloud, trackPart part=ParcelType::trackingData::tpLinearTrack)
Construct from components.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
const scalarField & Y() const
Return const access to mass fractions of mixture [].
const scalarField & Cs
Namespace for OpenFOAM.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
void calcPhaseChange(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Re, const scalar Pr, const scalar Ts, const scalar nus, const scalar d, const scalar T, const scalar mass, const scalar rho, const label idPhase, const scalar YPhase, const scalarField &YLiq, const scalarField &YSol, scalarField &dMassPC, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs)
Calculate Phase change.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.