Loading...
Searching...
No Matches
ReactingCloud.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::ReactingCloud
29
30Group
31 grpLagrangianIntermediateClouds
32
33Description
34 Templated base class for reacting cloud
35
36 - Adds to thermodynamic cloud
37 - Variable composition (single phase)
38 - Phase change
39
40SourceFiles
41 ReactingCloudI.H
42 ReactingCloud.C
43
44\*---------------------------------------------------------------------------*/
45
46#ifndef ReactingCloud_H
47#define ReactingCloud_H
48
49#include "reactingCloud.H"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
55
56// Forward declaration of classes
57
58template<class CloudType>
60
61template<class CloudType>
63
64/*---------------------------------------------------------------------------*\
65 Class ReactingCloud Declaration
66\*---------------------------------------------------------------------------*/
67
68template<class CloudType>
69class ReactingCloud
70:
71 public CloudType,
72 public reactingCloud
73{
74public:
75
76 // Public typedefs
77
78 //- Type of cloud this cloud was instantiated for
79 typedef CloudType cloudType;
80
81 //- Type of parcel the cloud was instantiated for
82 typedef typename CloudType::particleType parcelType;
83
84 //- Convenience typedef for this cloud type
85 typedef ReactingCloud<CloudType> reactingCloudType;
87
88private:
89
90 // Private data
91
92 //- Cloud copy pointer
94
95
96 // Private member functions
97
98 //- No copy construct
99 ReactingCloud(const ReactingCloud&) = delete;
100
101 //- No copy assignment
102 void operator=(const ReactingCloud&) = delete;
103
104
105protected:
106
107 // Protected data
108
109 //- Parcel constant properties
110 typename parcelType::constantProperties constProps_;
111
112
113 // References to the cloud sub-models
114
115 //- Reacting composition model
118
119 //- Reacting phase change model
122
123
124 // Sources
125
126 //- Mass transfer fields - one per carrier phase specie
129
130 // Protected Member Functions
131
132 // New parcel helper functions
133
134 //- Check that size of a composition field is valid
136 (
137 const scalarField& YSupplied,
138 const scalarField& Y,
139 const word& YName
140 );
141
143 // Initialisation
144
145 //- Set cloud sub-models
146 void setModels();
147
148
149 // Cloud evolution functions
150
151 //- Reset state of cloud
152 void cloudReset(ReactingCloud<CloudType>& c);
153
154
155public:
156
157 // Constructors
158
159 //- Construct given carrier gas fields
160 ReactingCloud
161 (
162 const word& cloudName,
163 const volScalarField& rho,
164 const volVectorField& U,
166 const SLGThermo& thermo,
167 bool readFields = true
168 );
169
170 //- Copy constructor with new name
171 ReactingCloud(ReactingCloud<CloudType>& c, const word& name);
172
173 //- Copy constructor with new name - creates bare cloud
174 ReactingCloud
175 (
176 const fvMesh& mesh,
177 const word& name,
178 const ReactingCloud<CloudType>& c
179 );
180
181 //- Construct and return clone based on (this) with new name
185 (
186 new ReactingCloud(*this, name)
187 );
188 }
189
190 //- Construct and return bare clone based on (this) with new name
191 virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
192 {
194 (
195 new ReactingCloud(this->mesh(), name, *this)
196 );
197 }
198
199
200 //- Destructor
201 virtual ~ReactingCloud() = default;
202
203
204 // Member Functions
205
206 // Access
207
208 //- Return a reference to the cloud copy
209 inline const ReactingCloud& cloudCopy() const;
210
211 //- Return the constant properties
212 inline const typename parcelType::constantProperties&
213 constProps() const;
214
215 //- Return access to the constant properties
216 inline typename parcelType::constantProperties& constProps();
217
218
219 // Sub-models
220
221 //- Return const access to reacting composition model
223 composition() const;
224
225 //- Return const access to reacting phase change model
227 phaseChange() const;
228
229 //- Return reference to reacting phase change model
231 phaseChange();
232
233
234 // Sources
235
236 //- Transfer the effect of parcel to the carrier phase
237 inline void transferToCarrier
238 (
239 const parcelType& p,
240 const typename parcelType::trackingData& td
241 );
242
243
244 // Mass
245
246 //- Return reference to mass source for field i
248 rhoTrans(const label i);
249
250 //- Return const access to mass source fields
252 rhoTrans() const;
253
254 //- Return reference to mass source fields
256 rhoTrans();
257
258 //- Return mass source term for specie i - specie eqn
260 (
261 const label i,
263 ) const;
265 //- Return tmp mass source for field i - fully explicit
267 Srho(const label i) const;
268
269 //- Return tmp total mass source for carrier phase
270 // - fully explicit
271 inline tmp<volScalarField::Internal> Srho() const;
272
273 //- Return total mass source term [kg/m3/s]
275
277 // Cloud evolution functions
278
279 //- Set parcel thermo properties
281 (
282 parcelType& parcel,
283 const scalar lagrangianDt
284 );
285
286 //- Check parcel properties
288 (
289 parcelType& parcel,
290 const scalar lagrangianDt,
291 const bool fullyDescribed
292 );
293
294 //- Store the current cloud state
295 void storeState();
296
297 //- Reset the current cloud to the previously stored state
298 void restoreState();
299
300 //- Reset the cloud source terms
301 void resetSourceTerms();
302
303 //- Apply relaxation to (steady state) cloud sources
304 void relaxSources(const ReactingCloud<CloudType>& cloudOldTime);
305
306 //- Apply scaling to (transient) cloud sources
307 void scaleSources();
308
309 //- Evolve the cloud
310 void evolve();
311
312
313 // Mapping
315 //- Remap the cells of particles corresponding to the
316 // mesh topology change with a default tracking data object
317 virtual void autoMap(const mapPolyMesh&);
318
319
320 // I-O
321
322 //- Print cloud information
323 void info();
325 //- Write the field data for the cloud
326 virtual void writeFields() const;
327
328 //- Write particle fields as objects into the obr registry
329 virtual void writeObjects(objectRegistry& obr) const;
330};
332
333// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334
335} // End namespace Foam
337// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338
339#include "ReactingCloudI.H"
340
341// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
342
343#ifdef NoRepository
345#endif
346
347// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348
349#endif
350
351// ************************************************************************* //
const uniformDimensionedVectorField & g
ParticleType particleType
Definition Cloud.H:130
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
const word & cloudName() const
const fvMesh & mesh() const
DimensionedField< scalar, volMesh > Internal
autoPtr< IOobject > clone() const
Clone.
Definition IOobject.H:641
Templated phase change model class.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Templated base class for reacting cloud.
virtual ~ReactingCloud()=default
Destructor.
autoPtr< PhaseChangeModel< ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicReactingParcel > > > > > > phaseChangeModel_
const parcelType::constantProperties & constProps() const
Return the constant properties.
ReactingCloud(ReactingCloud< CloudType > &c, const word &name)
Copy constructor with new name.
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
void setModels()
Set cloud sub-models.
ReactingCloud(const fvMesh &mesh, const word &name, const ReactingCloud< CloudType > &c)
Copy constructor with new name - creates bare cloud.
void storeState()
Store the current cloud state.
ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicReactingParcel > > > > reactingCloudType
ReactingCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const SLGThermo &thermo, bool readFields=true)
Construct given carrier gas fields.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
tmp< volScalarField::Internal > Srho(const label i) const
Return tmp mass source for field i - fully explicit.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
void scaleSources()
Apply scaling to (transient) cloud sources.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
tmp< fvScalarMatrix > Srho(volScalarField &rho) const
Return total mass source term [kg/m3/s].
ThermoCloud< KinematicCloud< Cloud< basicReactingParcel > > >::particleType parcelType
PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange()
Return reference to reacting phase change model.
volScalarField::Internal & rhoTrans(const label i)
Return reference to mass source for field i.
virtual void writeFields() const
Write the field data for the cloud.
ThermoCloud< KinematicCloud< Cloud< basicReactingParcel > > > cloudType
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
void evolve()
Evolve the cloud.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
autoPtr< CompositionModel< ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicReactingParcel > > > > > > compositionModel_
parcelType::constantProperties & constProps()
Return access to the constant properties.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
PtrList< volScalarField::Internal > & rhoTrans()
Return reference to mass source fields.
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
void resetSourceTerms()
Reset the cloud source terms.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition SLGThermo.H:63
const SLGThermo & thermo() const
Return const access to thermo package.
const volScalarField & p() const
Return const access to the carrier pressure field.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
reactingCloud()=default
Null constructor.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
PtrList< volScalarField > & Y
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.