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-2016 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 owner, position, and cloud owner
228  // Other properties initialised as null
229  inline ReactingParcel
230  (
231  const polyMesh& mesh,
232  const vector& position,
233  const label celli,
234  const label tetFacei,
235  const label tetPtI
236  );
237 
238  //- Construct from components
239  inline ReactingParcel
240  (
241  const polyMesh& mesh,
242  const vector& position,
243  const label celli,
244  const label tetFacei,
245  const label tetPtI,
246  const label typeId,
247  const scalar nParticle0,
248  const scalar d0,
249  const scalar dTarget0,
250  const vector& U0,
251  const vector& f0,
252  const vector& angularMomentum0,
253  const vector& torque0,
254  const scalarField& Y0,
255  const constantProperties& constProps
256  );
257 
258  //- Construct from Istream
260  (
261  const polyMesh& mesh,
262  Istream& is,
263  bool readFields = true
264  );
265 
266  //- Construct as a copy
268  (
269  const ReactingParcel& p,
270  const polyMesh& mesh
271  );
272 
273  //- Construct as a copy
275 
276  //- Construct and return a (basic particle) clone
277  virtual autoPtr<particle> clone() const
278  {
280  }
281 
282  //- Construct and return a (basic particle) clone
283  virtual autoPtr<particle> clone(const polyMesh& mesh) const
284  {
285  return autoPtr<particle>
286  (
287  new ReactingParcel<ParcelType>(*this, mesh)
288  );
289  }
290 
291  //- Factory class to read-construct particles used for
292  // parallel transfer
293  class iNew
294  {
295  const polyMesh& mesh_;
296 
297  public:
299  iNew(const polyMesh& mesh)
300  :
301  mesh_(mesh)
302  {}
304  autoPtr<ReactingParcel<ParcelType>> operator()(Istream& is) const
305  {
307  (
308  new ReactingParcel<ParcelType>(mesh_, is, true)
309  );
310  }
311  };
312 
313 
314  // Member Functions
315 
316  // Access
317 
318  //- Return const access to initial mass [kg]
319  inline scalar mass0() const;
320 
321  //- Return const access to mass fractions of mixture []
322  inline const scalarField& Y() const;
323 
324  //- Return the owner cell pressure [Pa]
325  inline scalar pc() const;
326 
327  //- Return reference to the owner cell pressure [Pa]
328  inline scalar& pc();
329 
330 
331  // Edit
332 
333  //- Return access to initial mass [kg]
334  inline scalar& mass0();
335 
336  //- Return access to mass fractions of mixture []
337  inline scalarField& Y();
338 
339 
340  // Main calculation loop
341 
342  //- Set cell values
343  template<class TrackData>
344  void setCellValues
345  (
346  TrackData& td,
347  const scalar dt,
348  const label celli
349  );
350 
351  //- Correct cell values using latest transfer information
352  template<class TrackData>
354  (
355  TrackData& td,
356  const scalar dt,
357  const label celli
358  );
359 
360  //- Correct surface values due to emitted species
361  template<class TrackData>
363  (
364  TrackData& td,
365  const label celli,
366  const scalar T,
367  const scalarField& Cs,
368  scalar& rhos,
369  scalar& mus,
370  scalar& Prs,
371  scalar& kappas
372  );
373 
374  //- Update parcel properties over the time interval
375  template<class TrackData>
376  void calc
377  (
378  TrackData& td,
379  const scalar dt,
380  const label celli
381  );
382 
383 
384  // I-O
385 
386  //- Read
387  template<class CloudType, class CompositionType>
388  static void readFields
389  (
390  CloudType& c,
391  const CompositionType& compModel
392  );
393 
394  //- Read - no composition
395  template<class CloudType>
396  static void readFields(CloudType& c);
397 
398  //- Write
399  template<class CloudType, class CompositionType>
400  static void writeFields
401  (
402  const CloudType& c,
403  const CompositionType& compModel
404  );
405 
406  //- Write - composition supplied
407  template<class CloudType>
408  static void writeFields(const CloudType& c);
409 
410 
411  // Ostream Operator
412 
413  friend Ostream& operator<< <ParcelType>
414  (
415  Ostream&,
417  );
418 };
419 
420 
421 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
422 
423 } // End namespace Foam
424 
425 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
426 
427 #include "ReactingParcelI.H"
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 #ifdef NoRepository
433  #include "ReactingParcel.C"
434 #endif
435 
436 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
437 
438 #endif
439 
440 // ************************************************************************* //
Class to hold reacting parcel constant properties.
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
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
scalar Sh
Definition: solveChemistry.H:2
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
scalarField Y_
Mass fractions of mixture [].
void cellValueSourceCorrection(TrackData &td, const scalar dt, const label celli)
Correct cell values using latest transfer information.
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:620
const scalarField & Y() const
Return const access to mass fractions of mixture [].
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
bool constantVolume() const
Return const access to the constant volume flag.
scalar mass0() const
Return const access to initial mass [kg].
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
AddToPropertyList(ParcelType," mass0"+" nPhases(Y1..YN)")
String representation of properties.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
scalar pMin() const
Return const access to the minimum pressure.
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.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
ReactingParcel(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from owner, position, and cloud owner.
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.
scalar pc() const
Return the owner cell pressure [Pa].
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:53
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:68
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
scalar updateMassFraction(const scalar mass0, const scalarField &dMass, scalarField &Y) const
Update mass fraction.
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.