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-2021 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 "interpolation.H"
43 #include "fluidThermo.H"
44 #include "demandDrivenEntry.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 template<class ParcelType>
52 class ReactingParcel;
53 
54 template<class ParcelType>
55 Ostream& operator<<
56 (
57  Ostream&,
59 );
60 
61 
62 /*---------------------------------------------------------------------------*\
63  Class ReactingParcelName Declaration
64 \*---------------------------------------------------------------------------*/
65 
67 
68 
69 /*---------------------------------------------------------------------------*\
70  Class ReactingParcel Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 template<class ParcelType>
74 class ReactingParcel
75 :
76  public ParcelType,
77  public ReactingParcelName
78 {
79  // Private Data
80 
81  //- Size in bytes of the fields
82  static const std::size_t sizeofFields_;
83 
84 
85 public:
86 
87  //- Class to hold reacting parcel constant properties
88  class constantProperties
89  :
90  public ParcelType::constantProperties
91  {
92  // Private Data
93 
94  //- Minimum pressure [Pa]
96 
97  //- Constant volume flag - e.g. during mass transfer
98  demandDrivenEntry<bool> constantVolume_;
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 minimum pressure
118  inline scalar pMin() const;
119 
120  //- Return const access to the constant volume flag
121  inline bool constantVolume() const;
122  };
123 
124 
125  //- Use base tracking data
126  typedef typename ParcelType::trackingData trackingData;
127 
128 
129 protected:
130 
131  // Protected data
132 
133  // Parcel properties
134 
135  //- Initial mass [kg]
136  scalar mass0_;
137 
138  //- Mass fractions of mixture []
139  scalarField Y_;
140 
141 
142  // Protected Member Functions
143 
144  //- Calculate Phase change
145  template<class TrackCloudType>
146  void calcPhaseChange
147  (
148  TrackCloudType& cloud,
149  trackingData& td,
150  const scalar dt, // timestep
151  const scalar Re, // Reynolds number
152  const scalar Pr, // Prandtl number
153  const scalar Ts, // Surface temperature
154  const scalar nus, // Surface kinematic viscosity
155  const scalar d, // diameter
156  const scalar T, // temperature
157  const scalar mass, // mass
158  const label idPhase, // id of phase involved in phase change
159  const scalar YPhase, // total mass fraction
160  const scalarField& YComponents, // component mass fractions
161  scalarField& dMassPC, // mass transfer - local to parcel
162  scalar& Sh, // explicit parcel enthalpy source
163  scalar& N, // flux of species emitted from parcel
164  scalar& NCpW, // sum of N*Cp*W of emission species
165  scalarField& Cs // carrier conc. of emission species
166  );
167 
168  //- Update mass fraction
169  scalar updateMassFraction
170  (
171  const scalar mass0,
172  const scalarField& dMass,
173  scalarField& Y
174  ) const;
175 
176 
177 public:
178 
179  // Static Data Members
180 
181  //- String representation of properties
183  (
184  ParcelType,
185  " mass0"
186  + " nPhases(Y1..YN)"
187  );
188 
189 
190  // Constructors
191 
192  //- Construct from mesh, coordinates and topology
193  // Other properties initialised as null
194  inline ReactingParcel
195  (
196  const polyMesh& mesh,
197  const barycentric& coordinates,
198  const label celli,
199  const label tetFacei,
200  const label tetPti
201  );
202 
203  //- Construct from a position and a cell, searching for the rest of the
204  // required topology. Other properties are initialised as null.
205  inline ReactingParcel
206  (
207  const polyMesh& mesh,
208  const vector& position,
209  const label celli
210  );
211 
212  //- Construct from Istream
214  (
215  const polyMesh& mesh,
216  Istream& is,
217  bool readFields = true
218  );
219 
220  //- Construct as a copy
222  (
223  const ReactingParcel& p,
224  const polyMesh& mesh
225  );
226 
227  //- Construct as a copy
229 
230  //- Construct and return a (basic particle) clone
231  virtual autoPtr<particle> clone() const
232  {
234  }
235 
236  //- Construct and return a (basic particle) clone
237  virtual autoPtr<particle> clone(const polyMesh& mesh) const
238  {
239  return autoPtr<particle>
240  (
241  new ReactingParcel<ParcelType>(*this, mesh)
242  );
243  }
244 
245  //- Factory class to read-construct particles used for
246  // parallel transfer
247  class iNew
248  {
249  const polyMesh& mesh_;
250 
251  public:
253  iNew(const polyMesh& mesh)
254  :
255  mesh_(mesh)
256  {}
258  autoPtr<ReactingParcel<ParcelType>> operator()(Istream& is) const
259  {
261  (
262  new ReactingParcel<ParcelType>(mesh_, is, true)
263  );
264  }
265  };
266 
267 
268  // Member Functions
269 
270  // Access
271 
272  //- Return const access to initial mass [kg]
273  inline scalar mass0() const;
274 
275  //- Return const access to mass fractions of mixture []
276  inline const scalarField& Y() const;
277 
278 
279  // Edit
280 
281  //- Return access to initial mass [kg]
282  inline scalar& mass0();
283 
284  //- Return access to mass fractions of mixture []
285  inline scalarField& Y();
286 
287 
288  // Main calculation loop
289 
290  //- Set cell values
291  template<class TrackCloudType>
292  void setCellValues(TrackCloudType& cloud, trackingData& td);
293 
294  //- Correct cell values using latest transfer information
295  template<class TrackCloudType>
297  (
298  TrackCloudType& cloud,
299  trackingData& td,
300  const scalar dt
301  );
302 
303  //- Correct surface values due to emitted species
304  template<class TrackCloudType>
306  (
307  TrackCloudType& cloud,
308  trackingData& td,
309  const scalar T,
310  const scalarField& Cs,
311  scalar& rhos,
312  scalar& mus,
313  scalar& Prs,
314  scalar& kappas
315  );
316 
317  //- Update parcel properties over the time interval
318  template<class TrackCloudType>
319  void calc
320  (
321  TrackCloudType& cloud,
322  trackingData& td,
323  const scalar dt
324  );
325 
326 
327  // I-O
328 
329  //- Read
330  template<class CloudType, class CompositionType>
331  static void readFields
332  (
333  CloudType& c,
334  const CompositionType& compModel
335  );
336 
337  //- Read - no composition
338  template<class CloudType>
339  static void readFields(CloudType& c);
340 
341  //- Write
342  template<class CloudType, class CompositionType>
343  static void writeFields
344  (
345  const CloudType& c,
346  const CompositionType& compModel
347  );
348 
349  //- Write - composition supplied
350  template<class CloudType>
351  static void writeFields(const CloudType& c);
352 
353 
354  // Ostream Operator
355 
356  friend Ostream& operator<< <ParcelType>
357  (
358  Ostream&,
360  );
361 };
362 
363 
364 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 
366 } // End namespace Foam
367 
368 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 
370 #include "ReactingParcelI.H"
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 
374 #ifdef NoRepository
375  #include "ReactingParcel.C"
376 #endif
377 
378 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
379 
380 #endif
381 
382 // ************************************************************************* //
Class to hold reacting parcel constant properties.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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.
fvMesh & mesh
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
const dimensionedScalar c
Speed of light in a vacuum.
ParcelType::trackingData trackingData
Use base tracking data.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
scalar pMin() const
Return const access to the minimum pressure.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
TemplateName(FvFaceCellWave)
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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.
scalar mass0() const
Return const access to initial mass [kg].
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
label N
Definition: createFields.H:9
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:76
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:75
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.