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 | Copyright (C) 2011-2017 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 
115 
116  template<class CloudType>
117  class TrackingData
118  :
119  public ParcelType::template TrackingData<CloudType>
120  {
121  private:
122 
123  // Private data
124 
125  // Interpolators for continuous phase fields
126 
127  //- Interpolator for continuous phase pressure field
129 
130 
131  public:
132 
133  typedef typename ParcelType::template TrackingData<CloudType>::trackPart
134  trackPart;
135 
136  // Constructors
137 
138  //- Construct from components
139  inline TrackingData
140  (
141  CloudType& cloud,
142  trackPart part = ParcelType::template
144  );
145 
146 
147  // Member functions
148 
149  //- Return const access to the interpolator for continuous phase
150  // pressure field
151  inline const interpolation<scalar>& pInterp() const;
152  };
153 
154 
155 protected:
156 
157  // Protected data
158 
159  // Parcel properties
160 
161  //- Initial mass [kg]
162  scalar mass0_;
163 
164  //- Mass fractions of mixture []
165  scalarField Y_;
166 
167 
168  // Cell-based quantities
169 
170  //- Pressure [Pa]
171  scalar pc_;
172 
173 
174  // Protected Member Functions
175 
176  //- Calculate Phase change
177  template<class TrackData>
178  void calcPhaseChange
179  (
180  TrackData& td,
181  const scalar dt, // timestep
182  const label celli, // owner cell
183  const scalar Re, // Reynolds number
184  const scalar Pr, // Prandtl number
185  const scalar Ts, // Surface temperature
186  const scalar nus, // Surface kinematic viscosity
187  const scalar d, // diameter
188  const scalar T, // temperature
189  const scalar mass, // mass
190  const label idPhase, // id of phase involved in phase change
191  const scalar YPhase, // total mass fraction
192  const scalarField& YComponents, // component mass fractions
193  scalarField& dMassPC, // mass transfer - local to parcel
194  scalar& Sh, // explicit parcel enthalpy source
195  scalar& N, // flux of species emitted from parcel
196  scalar& NCpW, // sum of N*Cp*W of emission species
197  scalarField& Cs // carrier conc. of emission species
198  );
199 
200  //- Update mass fraction
201  scalar updateMassFraction
202  (
203  const scalar mass0,
204  const scalarField& dMass,
205  scalarField& Y
206  ) const;
207 
208 
209 public:
210 
211  // Static data members
212 
213  //- Runtime type information
214  TypeName("ReactingParcel");
215 
216  //- String representation of properties
218  (
219  ParcelType,
220  " mass0"
221  + " nPhases(Y1..YN)"
222  );
223 
224 
225  // Constructors
226 
227  //- Construct from mesh, coordinates and topology
228  // Other properties initialised as null
229  inline ReactingParcel
230  (
231  const polyMesh& mesh,
232  const barycentric& coordinates,
233  const label celli,
234  const label tetFacei,
235  const label tetPti
236  );
237 
238  //- Construct from a position and a cell, searching for the rest of the
239  // required topology. Other properties are initialised as null.
240  inline ReactingParcel
241  (
242  const polyMesh& mesh,
243  const vector& position,
244  const label celli
245  );
246 
247  //- Construct from components
248  inline ReactingParcel
249  (
250  const polyMesh& mesh,
251  const barycentric& coordinates,
252  const label celli,
253  const label tetFacei,
254  const label tetPti,
255  const label typeId,
256  const scalar nParticle0,
257  const scalar d0,
258  const scalar dTarget0,
259  const vector& U0,
260  const vector& f0,
261  const vector& angularMomentum0,
262  const vector& torque0,
263  const scalarField& Y0,
264  const constantProperties& constProps
265  );
266 
267  //- Construct from Istream
269  (
270  const polyMesh& mesh,
271  Istream& is,
272  bool readFields = true
273  );
274 
275  //- Construct as a copy
277  (
278  const ReactingParcel& p,
279  const polyMesh& mesh
280  );
281 
282  //- Construct as a copy
284 
285  //- Construct and return a (basic particle) clone
286  virtual autoPtr<particle> clone() const
287  {
289  }
290 
291  //- Construct and return a (basic particle) clone
292  virtual autoPtr<particle> clone(const polyMesh& mesh) const
293  {
294  return autoPtr<particle>
295  (
296  new ReactingParcel<ParcelType>(*this, mesh)
297  );
298  }
299 
300  //- Factory class to read-construct particles used for
301  // parallel transfer
302  class iNew
303  {
304  const polyMesh& mesh_;
305 
306  public:
308  iNew(const polyMesh& mesh)
309  :
310  mesh_(mesh)
311  {}
313  autoPtr<ReactingParcel<ParcelType>> operator()(Istream& is) const
314  {
316  (
317  new ReactingParcel<ParcelType>(mesh_, is, true)
318  );
319  }
320  };
321 
322 
323  // Member Functions
324 
325  // Access
326 
327  //- Return const access to initial mass [kg]
328  inline scalar mass0() const;
329 
330  //- Return const access to mass fractions of mixture []
331  inline const scalarField& Y() const;
332 
333  //- Return the owner cell pressure [Pa]
334  inline scalar pc() const;
335 
336  //- Return reference to the owner cell pressure [Pa]
337  inline scalar& pc();
338 
339 
340  // Edit
341 
342  //- Return access to initial mass [kg]
343  inline scalar& mass0();
344 
345  //- Return access to mass fractions of mixture []
346  inline scalarField& Y();
347 
348 
349  // Main calculation loop
350 
351  //- Set cell values
352  template<class TrackData>
353  void setCellValues
354  (
355  TrackData& td,
356  const scalar dt,
357  const label celli
358  );
359 
360  //- Correct cell values using latest transfer information
361  template<class TrackData>
363  (
364  TrackData& td,
365  const scalar dt,
366  const label celli
367  );
368 
369  //- Correct surface values due to emitted species
370  template<class TrackData>
372  (
373  TrackData& td,
374  const label celli,
375  const scalar T,
376  const scalarField& Cs,
377  scalar& rhos,
378  scalar& mus,
379  scalar& Prs,
380  scalar& kappas
381  );
382 
383  //- Update parcel properties over the time interval
384  template<class TrackData>
385  void calc
386  (
387  TrackData& td,
388  const scalar dt,
389  const label celli
390  );
391 
392 
393  // I-O
394 
395  //- Read
396  template<class CloudType, class CompositionType>
397  static void readFields
398  (
399  CloudType& c,
400  const CompositionType& compModel
401  );
402 
403  //- Read - no composition
404  template<class CloudType>
405  static void readFields(CloudType& c);
406 
407  //- Write
408  template<class CloudType, class CompositionType>
409  static void writeFields
410  (
411  const CloudType& c,
412  const CompositionType& compModel
413  );
414 
415  //- Write - composition supplied
416  template<class CloudType>
417  static void writeFields(const CloudType& c);
418 
419 
420  // Ostream Operator
421 
422  friend Ostream& operator<< <ParcelType>
423  (
424  Ostream&,
426  );
427 };
428 
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 } // End namespace Foam
433 
434 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
435 
436 #include "ReactingParcelI.H"
438 
439 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440 
441 #ifdef NoRepository
442  #include "ReactingParcel.C"
443 #endif
444 
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 
447 #endif
448 
449 // ************************************************************************* //
Class to hold reacting parcel constant properties.
scalar pc() const
Return the owner cell pressure [Pa].
ParcelType::template TrackingData< CloudType >::trackPart trackPart
dimensionedScalar Pr("Pr", dimless, laminarTransport)
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:137
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:739
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label celli)
Correct cell values using latest transfer information.
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.
scalarList Y0(nSpecie, 0.0)
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
void calc(TrackData &td, const scalar dt, const label celli)
Update parcel properties over the time interval.
dynamicFvMesh & mesh
scalar pMin() const
Return const access to the minimum pressure.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
scalar pc_
Pressure [Pa].
TypeName("ReactingParcel")
Runtime type information.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void setCellValues(TrackData &td, const scalar dt, const label celli)
Set cell values.
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 correctSurfaceValues(TrackData &td, const label celli, const scalar T, const scalarField &Cs, scalar &rhos, scalar &mus, scalar &Prs, scalar &kappas)
Correct surface values due to emitted species.
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 [].
void calcPhaseChange(TrackData &td, const scalar dt, const label celli, 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.
Namespace for OpenFOAM.