ReactingParcel.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-2019 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::ReactingParcel
26 
27 Description
28  Reacting parcel class with one/two-way coupling with the continuous
29  phase.
30 
31 SourceFiles
32  ReactingParcelI.H
33  ReactingParcel.C
34  ReactingParcelIO.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef ReactingParcel_H
39 #define ReactingParcel_H
40 
41 #include "particle.H"
42 #include "SLGThermo.H"
43 #include "demandDrivenEntry.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 template<class ParcelType>
51 class ReactingParcel;
52 
53 template<class ParcelType>
54 Ostream& operator<<
55 (
56  Ostream&,
58 );
59 
60 
61 /*---------------------------------------------------------------------------*\
62  Class ReactingParcel Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class ParcelType>
66 class ReactingParcel
67 :
68  public ParcelType
69 {
70  // Private Data
71 
72  //- Size in bytes of the fields
73  static const std::size_t sizeofFields_;
74 
75 
76 public:
77 
78  //- Class to hold reacting parcel constant properties
79  class constantProperties
80  :
81  public ParcelType::constantProperties
82  {
83  // Private Data
84 
85  //- Minimum pressure [Pa]
87 
88  //- Constant volume flag - e.g. during mass transfer
89  demandDrivenEntry<bool> constantVolume_;
90 
91 
92  public:
93 
94  // Constructors
95 
96  //- Null constructor
98 
99  //- Copy constructor
101 
102  //- Construct from dictionary
103  constantProperties(const dictionary& parentDict);
104 
105 
106  // Access
107 
108  //- Return const access to the minimum pressure
109  inline scalar pMin() const;
110 
111  //- Return const access to the constant volume flag
112  inline bool constantVolume() const;
113  };
114 
116  class trackingData
117  :
118  public ParcelType::trackingData
119  {
120  private:
121 
122  // Private Data
123 
124  // Interpolators for continuous phase fields
125 
126  //- Interpolator for continuous phase pressure field
128 
129 
130  // Cached continuous phase properties
131 
132  //- Pressure [Pa]
133  scalar pc_;
134 
135 
136  public:
138  typedef typename ParcelType::trackingData::trackPart trackPart;
139 
140  // Constructors
141 
142  //- Construct from components
143  template<class TrackCloudType>
144  inline trackingData
145  (
146  const TrackCloudType& cloud,
147  trackPart part = ParcelType::trackingData::tpLinearTrack
148  );
149 
150 
151  // Member Functions
152 
153  //- Return const access to the interpolator for continuous phase
154  // pressure field
155  inline const interpolation<scalar>& pInterp() const;
156 
157  //- Return the continuous phase pressure
158  inline scalar pc() const;
159 
160  //- Access the continuous phase pressure
161  inline scalar& pc();
162  };
163 
164 
165 protected:
166 
167  // Protected data
168 
169  // Parcel properties
170 
171  //- Initial mass [kg]
172  scalar mass0_;
173 
174  //- Mass fractions of mixture []
175  scalarField Y_;
176 
177 
178  // Protected Member Functions
179 
180  //- Calculate Phase change
181  template<class TrackCloudType>
182  void calcPhaseChange
183  (
184  TrackCloudType& cloud,
185  trackingData& td,
186  const scalar dt, // timestep
187  const scalar Re, // Reynolds number
188  const scalar Pr, // Prandtl number
189  const scalar Ts, // Surface temperature
190  const scalar nus, // Surface kinematic viscosity
191  const scalar d, // diameter
192  const scalar T, // temperature
193  const scalar mass, // mass
194  const label idPhase, // id of phase involved in phase change
195  const scalar YPhase, // total mass fraction
196  const scalarField& YComponents, // component mass fractions
197  scalarField& dMassPC, // mass transfer - local to parcel
198  scalar& Sh, // explicit parcel enthalpy source
199  scalar& N, // flux of species emitted from parcel
200  scalar& NCpW, // sum of N*Cp*W of emission species
201  scalarField& Cs // carrier conc. of emission species
202  );
203 
204  //- Update mass fraction
205  scalar updateMassFraction
206  (
207  const scalar mass0,
208  const scalarField& dMass,
209  scalarField& Y
210  ) const;
211 
212 
213 public:
214 
215  // Static Data Members
216 
217  //- Runtime type information
218  TypeName("ReactingParcel");
219 
220  //- String representation of properties
222  (
223  ParcelType,
224  " mass0"
225  + " nPhases(Y1..YN)"
226  );
227 
228 
229  // Constructors
230 
231  //- Construct from mesh, coordinates and topology
232  // Other properties initialised as null
233  inline ReactingParcel
234  (
235  const polyMesh& mesh,
236  const barycentric& coordinates,
237  const label celli,
238  const label tetFacei,
239  const label tetPti
240  );
241 
242  //- Construct from a position and a cell, searching for the rest of the
243  // required topology. Other properties are initialised as null.
244  inline ReactingParcel
245  (
246  const polyMesh& mesh,
247  const vector& position,
248  const label celli
249  );
250 
251  //- Construct from components
252  inline ReactingParcel
253  (
254  const polyMesh& mesh,
255  const barycentric& coordinates,
256  const label celli,
257  const label tetFacei,
258  const label tetPti,
259  const label typeId,
260  const scalar nParticle0,
261  const scalar d0,
262  const scalar dTarget0,
263  const vector& U0,
264  const vector& f0,
265  const vector& angularMomentum0,
266  const vector& torque0,
267  const scalarField& Y0,
268  const constantProperties& constProps
269  );
270 
271  //- Construct from Istream
273  (
274  const polyMesh& mesh,
275  Istream& is,
276  bool readFields = true
277  );
278 
279  //- Construct as a copy
281  (
282  const ReactingParcel& p,
283  const polyMesh& mesh
284  );
285 
286  //- Construct as a copy
288 
289  //- Construct and return a (basic particle) clone
290  virtual autoPtr<particle> clone() const
291  {
293  }
294 
295  //- Construct and return a (basic particle) clone
296  virtual autoPtr<particle> clone(const polyMesh& mesh) const
297  {
298  return autoPtr<particle>
299  (
300  new ReactingParcel<ParcelType>(*this, mesh)
301  );
302  }
303 
304  //- Factory class to read-construct particles used for
305  // parallel transfer
306  class iNew
307  {
308  const polyMesh& mesh_;
309 
310  public:
312  iNew(const polyMesh& mesh)
313  :
314  mesh_(mesh)
315  {}
317  autoPtr<ReactingParcel<ParcelType>> operator()(Istream& is) const
318  {
320  (
321  new ReactingParcel<ParcelType>(mesh_, is, true)
322  );
323  }
324  };
325 
326 
327  // Member Functions
328 
329  // Access
330 
331  //- Return const access to initial mass [kg]
332  inline scalar mass0() const;
333 
334  //- Return const access to mass fractions of mixture []
335  inline const scalarField& Y() const;
336 
337 
338  // Edit
339 
340  //- Return access to initial mass [kg]
341  inline scalar& mass0();
342 
343  //- Return access to mass fractions of mixture []
344  inline scalarField& Y();
345 
346 
347  // Main calculation loop
348 
349  //- Set cell values
350  template<class TrackCloudType>
351  void setCellValues(TrackCloudType& cloud, trackingData& td);
352 
353  //- Correct cell values using latest transfer information
354  template<class TrackCloudType>
356  (
357  TrackCloudType& cloud,
358  trackingData& td,
359  const scalar dt
360  );
361 
362  //- Correct surface values due to emitted species
363  template<class TrackCloudType>
365  (
366  TrackCloudType& cloud,
367  trackingData& td,
368  const scalar T,
369  const scalarField& Cs,
370  scalar& rhos,
371  scalar& mus,
372  scalar& Prs,
373  scalar& kappas
374  );
375 
376  //- Update parcel properties over the time interval
377  template<class TrackCloudType>
378  void calc
379  (
380  TrackCloudType& cloud,
381  trackingData& td,
382  const scalar dt
383  );
384 
385 
386  // I-O
387 
388  //- Read
389  template<class CloudType, class CompositionType>
390  static void readFields
391  (
392  CloudType& c,
393  const CompositionType& compModel
394  );
395 
396  //- Read - no composition
397  template<class CloudType>
398  static void readFields(CloudType& c);
399 
400  //- Write
401  template<class CloudType, class CompositionType>
402  static void writeFields
403  (
404  const CloudType& c,
405  const CompositionType& compModel
406  );
407 
408  //- Write - composition supplied
409  template<class CloudType>
410  static void writeFields(const CloudType& c);
411 
412 
413  // Ostream Operator
414 
415  friend Ostream& operator<< <ParcelType>
416  (
417  Ostream&,
419  );
420 };
421 
422 
423 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
424 
425 } // End namespace Foam
426 
427 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
428 
429 #include "ReactingParcelI.H"
431 
432 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
433 
434 #ifdef NoRepository
435  #include "ReactingParcel.C"
436 #endif
437 
438 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439 
440 #endif
441 
442 // ************************************************************************* //
Class to hold reacting parcel constant properties.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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 label idPhase, const scalar YPhase, const scalarField &YComponents, scalarField &dMassPC, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs)
Calculate Phase change.
Factory class to read-construct particles used for.
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, recursively if necessary, the source to the destination.
Definition: POSIX.C:753
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.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
scalarList Y0(nSpecie, 0.0)
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
dynamicFvMesh & mesh
scalar pMin() const
Return const access to the minimum pressure.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
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.
TypeName("ReactingParcel")
Runtime type information.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
AddToPropertyList(ParcelType, " mass0"+" nPhases(Y1..YN)")
String representation of properties.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
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.
label N
Definition: createFields.H:22
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Reacting parcel class with one/two-way coupling with the continuous phase.
volScalarField & p
scalar mass0_
Initial mass [kg].
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
const scalarField & Y() const
Return const access to mass fractions of mixture [].
Namespace for OpenFOAM.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.