ReactingHeterogeneousParcel.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) 2018-2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::ReactingHeterogeneousParcel
28
29Group
30 grpLagrangianIntermediateParcels
31
32Description
33 Reacting heterogeneous Parcel
34
35SourceFiles
36 ReactingHeterogeneousParcelI.H
37 ReactingHeterogeneousParcel.C
38 ReactingHeterogeneousParcelIO.C
39
40\*---------------------------------------------------------------------------*/
41
42#ifndef ReactingHeterogeneousParcel_H
43#define ReactingHeterogeneousParcel_H
44
45#include "demandDrivenEntry.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52template<class ParcelType>
53class ReactingHeterogeneousParcel;
54
55template<class ParcelType>
56Ostream& operator<<
57(
58 Ostream&,
60);
61
62/*---------------------------------------------------------------------------*\
63 Class ReactingHeterogeneousParcel Declaration
64\*---------------------------------------------------------------------------*/
65
66template<class ParcelType>
68:
69 public ParcelType
70{
71public:
72
73 //- Size in bytes of the fields
74 static const std::size_t sizeofFields;
75
76 //- Class to hold reacting particle constant properties
78 :
79 public ParcelType::constantProperties
80 {
81 // Private data
82
83 //- Fraction of enthalpy retained by parcel due to surface
84 // reactions
85 demandDrivenEntry<scalar> hRetentionCoeff_;
86
87 public:
88
89 // Constructors
90
91 //- Null constructor
93
94 //- Copy constructor
96
97 //- Construct from dictionary
98 constantProperties(const dictionary& parentDict);
99
100
101 // Access
102
103 //- Return const access to the fraction of enthalpy retained by
104 // parcel due to surface reactions
105 inline scalar hRetentionCoeff() const;
106 };
107
108
109 //- Use base tracking data
110 typedef typename ParcelType::trackingData trackingData;
111
112
113private:
114
115 // Private Member Functions
116
117 //- Return the mixture effective specific heat capacity
118 template<class TrackCloudType>
119 scalar CpEff
120 (
121 TrackCloudType& cloud,
122 trackingData& td,
123 const scalar p,
124 const scalar T,
125 const label idS
126 ) const;
127
128 //- Return the mixture effective sensible enthalpy
129 template<class TrackCloudType>
130 scalar HsEff
131 (
132 TrackCloudType& cloud,
133 trackingData& td,
134 const scalar p,
135 const scalar T,
136 const label idS
137 ) const;
138
139 //- Return the mixture effective latent heat
140 template<class TrackCloudType>
141 scalar LEff
142 (
143 TrackCloudType& cloud,
144 trackingData& td,
145 const scalar p,
146 const scalar T,
147 const label idS
148 ) const;
149
150
151protected:
152
153 // Protected data
154
155 // Parcel properties
156
157 //- Progress variables for reactions
159
160 //- Flag to identify if the particle can devolatilise and combust
161 // Combustion possible only after volatile content falls below
162 // threshold value. States include:
163 // 0 = can combust but can change
164 // 1 = can devolatilise, can combust
165 // -1 = cannot devolatilise or combust, and cannot change
166 label canCombust_;
167
168
169 // Protected Member Functions
170
171 //- Return change of volume due to mass exchange
172 template<class TrackCloudType>
173 scalar updatedDeltaVolume
174 (
175 TrackCloudType& cloud,
176 const scalarField& dMass,
177 const scalar p,
178 const scalar T
179 );
180
181
182 //- Calculate surface reactions
183 template<class TrackCloudType>
185 (
186 TrackCloudType& cloud,
187 trackingData& td,
188 const scalar dt, // timestep
189 const scalar Res, // Re
190 const scalar nu, // nu
191 const scalar d, // diameter
192 const scalar T, // temperature
193 const scalar mass, // mass
194 const label canCombust, // 'can combust' flag
195 const scalar N, // flux of species emitted from particle
196 scalar& NCpW,
197 const scalarField& YSolid, // solid-phase mass fractions
198 scalarField& F, // progress of each reaction
199 scalarField& dMassSRSolid, // solid-phase mass transfer - local
200 scalarField& dMassSRCarrier, // carrier phase mass transfer
201 scalar& Sh, // explicit particle enthalpy source
202 scalar& dhsTrans // sensible enthalpy transfer to carrier
203 ) const;
204
205
206public:
207
208 // Static data members
209
210 //- Runtime type information
211 TypeName("ReactingHeterogeneousParcel");
212
213 //- String representation of properties
215 (
216 ParcelType,
217 + " nReactions(F1..FN)"
218 );
219
220
221 // Constructors
222
223 //- Construct from mesh, position and topology
224 // Other properties initialised as null
226 (
227 const polyMesh& mesh,
229 const label celli,
230 const label tetFacei,
231 const label tetPti
232 );
233
234 //- Construct from a position and a cell, searching for the rest of the
235 // required topology. Other properties are initialised as null.
237 (
238 const polyMesh& mesh,
239 const vector& position,
240 const label celli
241 );
242
243 //- Construct from components
245 (
246 const polyMesh& mesh,
248 const label celli,
249 const label tetFacei,
250 const label tetPti,
251 const label typeId,
252 const scalar nParticle0,
253 const scalar d0,
254 const scalar dTarget0,
255 const vector& U0,
256 const vector& f0,
257 const vector& angularMomentum0,
258 const vector& torque0,
259 const scalarField& Y,
260 const scalarField& F,
261 const constantProperties& constProps
262 );
263
264 //- Construct from Istream
266 (
267 const polyMesh& mesh,
268 Istream& is,
269 bool readFields = true,
270 bool newFormat = true
271 );
272
273 //- Construct as a copy
275
276 //- Construct as a copy
278 (
280 const polyMesh& mesh
281 );
282
283 //- Construct and return a (basic particle) clone
284 virtual autoPtr<particle> clone() const
285 {
287 }
288
289 //- Construct and return a (basic particle) clone
290 virtual autoPtr<particle> clone(const polyMesh& mesh) const
291 {
292 return autoPtr<particle>
293 (
295 );
296 }
297
298 //- Factory class to read-construct particles used for
299 // parallel transfer
300 class iNew
301 {
302 const polyMesh& mesh_;
303
304 public:
306 iNew(const polyMesh& mesh)
307 :
308 mesh_(mesh)
309 {}
312 (
313 Istream& is
314 ) const
315 {
317 (
319 (mesh_, is, true)
320 );
321 }
322 };
323
324
325 // Member Functions
326
327 // Access
328
329 //- Return const access to F
330 inline const scalarField& F() const;
331
332 //- Return const access to the canCombust flag
333 inline label canCombust() const;
334
335
336 // Edit
337
338 //- Return access to F
339 inline scalarField& F();
340
341 //- Return access to the canCombust flag
342 inline label& canCombust();
343
344
345 // Main calculation loop
346
347
348 //- Update parcel properties over the time interval
349 template<class TrackCloudType>
350 void calc
351 (
352 TrackCloudType& cloud,
353 trackingData& td,
354 const scalar dt
355 );
356
357
358 // I-O
359
360 //- Read - composition supplied
361 template<class CloudType, class CompositionType>
362 static void readFields
363 (
364 CloudType& c,
365 const CompositionType& compModel
366 );
367
368 //- Read - no composition
369 template<class CloudType>
370 static void readFields(CloudType& c);
371
372 //- Write - composition supplied
373 template<class CloudType, class CompositionType>
374 static void writeFields
375 (
376 const CloudType& c,
377 const CompositionType& compModel
378 );
379
380 //- Read - no composition
381 template<class CloudType>
382 static void writeFields(const CloudType& c);
383
384 //- Write individual parcel properties to stream
385 void writeProperties
386 (
387 Ostream& os,
388 const wordRes& filters,
389 const word& delim,
390 const bool namesOnly
391 ) const;
392
393 //- Read particle fields as objects from the obr registry
394 // - no composition
395 template<class CloudType>
396 static void readObjects
397 (
398 CloudType& c,
399 const objectRegistry& obr
400 );
401
402 //- Read particle fields as objects from the obr registry
403 template<class CloudType, class CompositionType>
404 static void readObjects
405 (
406 CloudType& c,
407 const CompositionType& compModel,
408 const objectRegistry& obr
409 );
410
411 //- Write particle fields as objects into the obr registry
412 // - no composition
413 template<class CloudType>
414 static void writeObjects
415 (
416 const CloudType& c,
417 objectRegistry& obr
418 );
419
420 //- Write particle fields as objects into the obr registry
421 template<class CloudType, class CompositionType>
422 static void writeObjects
423 (
424 const CloudType& c,
425 const CompositionType& compModel,
426 objectRegistry& obr
427 );
428
429
430 // Ostream Operator
432 friend Ostream& operator<< <ParcelType>
433 (
434 Ostream&,
436 );
437};
438
439
440// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441
442} // End namespace Foam
443
444// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445
447
448// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449
450#ifdef NoRepository
452#endif
453
454// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455
456#endif
457
458// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Class to hold reacting particle constant properties.
scalar hRetentionCoeff() const
Return const access to the fraction of enthalpy retained by.
Factory class to read-construct particles used for.
AddToPropertyList(ParcelType,+" nReactions(F1..FN)")
String representation of properties.
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
label canCombust_
Flag to identify if the particle can devolatilise and combust.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
friend Ostream & operator(Ostream &, const ReactingHeterogeneousParcel< ParcelType > &)
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
static const std::size_t sizeofFields
Size in bytes of the fields.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
TypeName("ReactingHeterogeneousParcel")
Runtime type information.
scalar updatedDeltaVolume(TrackCloudType &cloud, const scalarField &dMass, const scalar p, const scalar T)
Return change of volume due to mass exchange.
void calcHeterogeneousReactions(TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar Res, const scalar nu, const scalar d, const scalar T, const scalar mass, const label canCombust, const scalar N, scalar &NCpW, const scalarField &YSolid, scalarField &F, scalarField &dMassSRSolid, scalarField &dMassSRCarrier, scalar &Sh, scalar &dhsTrans) const
Calculate surface reactions.
label canCombust() const
Return const access to the canCombust flag.
scalarField F_
Progress variables for reactions.
ParcelType::trackingData trackingData
Use base tracking data.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
ReactingHeterogeneousParcel(const ReactingHeterogeneousParcel &p)
Construct as a copy.
const scalarField & F() const
Return const access to F.
ReactingHeterogeneousParcel(const ReactingHeterogeneousParcel &p, const polyMesh &mesh)
Construct as a copy.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
Class for demand-driven dictionary entries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
PtrList< volScalarField > & Y
const volScalarField & T
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
#define AddToPropertyList(ParcelType, str)
Add to existing static 'propertyList' for particle properties.
volScalarField & nu
const volScalarField & cp
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
const Vector< label > N(dict.get< Vector< label > >("N"))