ReactingMultiphaseCloud.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-2018 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::ReactingMultiphaseCloud
26 
27 Description
28  Templated base class for multiphase reacting cloud
29 
30  - Adds to reacting cloud
31  - multiphase composition
32  - devolatilisatsion
33  - surface reactions
34 
35 SourceFiles
36  ReactingMultiphaseCloudI.H
37  ReactingMultiphaseCloud.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef ReactingMultiphaseCloud_H
42 #define ReactingMultiphaseCloud_H
43 
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of classes
52 
53 template<class CloudType>
55 
56 template<class CloudType>
58 
59 /*---------------------------------------------------------------------------*\
60  Class ReactingMultiphaseCloud Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 template<class CloudType>
65 :
66  public CloudType,
68 {
69 public:
70 
71  // Public typedefs
72 
73  //- Type of cloud this cloud was instantiated for
74  typedef CloudType cloudType;
75 
76  //- Type of parcel the cloud was instantiated for
77  typedef typename CloudType::particleType parcelType;
78 
79  //- Convenience typedef for this cloud type
81 
82 
83 private:
84 
85  // Private data
86 
87  //- Cloud copy pointer
89 
90 
91  // Private member functions
92 
93  //- Disallow default bitwise copy construct
95 
96  //- Disallow default bitwise assignment
97  void operator=(const ReactingMultiphaseCloud&);
98 
99 
100 protected:
101 
102  // Protected data
103 
104  //- Parcel constant properties
105  typename parcelType::constantProperties constProps_;
106 
107 
108  // References to the cloud sub-models
109 
110  //- Devolatilisation model
111  autoPtr
112  <
114  >
116 
117  //- Surface reaction model
118  autoPtr
119  <
121  >
123 
124 
125  // Check
126 
127  //- Total mass transferred to continuous phase via devolatilisation
128  scalar dMassDevolatilisation_;
129 
130  //- Total mass transferred to continuous phase via surface
131  // reactions
132  scalar dMassSurfaceReaction_;
133 
134 
135  // Protected Member Functions
136 
137  // Initialisation
138 
139  //- Set cloud sub-models
140  void setModels();
141 
142 
143  // Cloud evolution functions
144 
145  //- Reset state of cloud
147 
148 
149 public:
150 
151  // Constructors
152 
153  //- Construct given carrier gas fields
155  (
156  const word& cloudName,
157  const volScalarField& rho,
158  const volVectorField& U,
159  const dimensionedVector& g,
160  const SLGThermo& thermo,
161  bool readFields = true
162  );
163 
164 
165  //- Copy constructor with new name
167  (
169  const word& name
170  );
171 
172  //- Copy constructor with new name - creates bare cloud
174  (
175  const fvMesh& mesh,
176  const word& name,
178  );
179 
180  //- Construct and return clone based on (this) with new name
181  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
182  {
184  (
185  new ReactingMultiphaseCloud(*this, name)
186  );
187  }
188 
189  //- Construct and return bare clone based on (this) with new name
190  virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
191  {
193  (
194  new ReactingMultiphaseCloud(this->mesh(), name, *this)
195  );
196  }
197 
198 
199  //- Destructor
200  virtual ~ReactingMultiphaseCloud();
201 
202 
203  // Member Functions
204 
205  // Access
206 
207  //- Return a reference to the cloud copy
208  inline const ReactingMultiphaseCloud& cloudCopy() const;
209 
210  //- Return the constant properties
211  inline const typename parcelType::constantProperties&
212  constProps() const;
213 
214  //- Return access to the constant properties
215  inline typename parcelType::constantProperties& constProps();
216 
217 
218  // Sub-models
219 
220  //- Return const access to devolatilisation model
221  inline const DevolatilisationModel
222  <
224  >&
225  devolatilisation() const;
226 
227  //- Return reference to devolatilisation model
228  inline DevolatilisationModel
229  <
231  >&
233 
234  //- Return const access to reacting surface reaction model
235  inline const SurfaceReactionModel
236  <
238  >&
239  surfaceReaction() const;
240 
241  //- Return reference to reacting surface reaction model
242  inline SurfaceReactionModel
243  <
245  >&
246  surfaceReaction();
247 
248 
249  // Cloud evolution functions
250 
251  //- Set parcel thermo properties
253  (
254  parcelType& parcel,
255  const scalar lagrangianDt
256  );
257 
258  //- Check parcel properties
260  (
261  parcelType& parcel,
262  const scalar lagrangianDt,
263  const bool fullyDescribed
264  );
265 
266  //- Store the current cloud state
267  void storeState();
268 
269  //- Reset the current cloud to the previously stored state
270  void restoreState();
271 
272  //- Reset the cloud source terms
273  void resetSourceTerms();
274 
275  //- Evolve the cloud
276  void evolve();
277 
278 
279  // Mapping
280 
281  //- Remap the cells of particles corresponding to the
282  // mesh topology change with a default tracking data object
283  virtual void autoMap(const mapPolyMesh&);
284 
285 
286  // I-O
287 
288  //- Print cloud information
289  void info();
290 
291  //- Write the field data for the cloud
292  virtual void writeFields() const;
293 };
294 
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 } // End namespace Foam
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 #ifdef NoRepository
307  #include "ReactingMultiphaseCloud.C"
308 #endif
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 #endif
313 
314 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:105
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:270
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
const word & name() const
Return name.
Definition: IOobject.H:297
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
void setModels()
Set cloud sub-models.
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:34
parcelType::constantProperties constProps_
Parcel constant properties.
rhoReactionThermo & thermo
Definition: createFields.H:28
autoPtr< SurfaceReactionModel< ReactingMultiphaseCloud< CloudType > > > surfaceReactionModel_
Surface reaction model.
Templated base class for multiphase reacting cloud.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
void cloudReset(ReactingMultiphaseCloud< CloudType > &c)
Reset state of cloud.
void restoreState()
Reset the current cloud to the previously stored state.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
A class for handling words, derived from string.
Definition: word.H:59
CloudType cloudType
Type of cloud this cloud was instantiated for.
void resetSourceTerms()
Reset the cloud source terms.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
void storeState()
Store the current cloud state.
Virtual abstract base class for templated reactingMultiphaseCloud.
virtual ~ReactingMultiphaseCloud()
Destructor.
ReactingMultiphaseCloud< CloudType > reactingMultiphaseCloudType
Convenience typedef for this cloud type.
autoPtr< DevolatilisationModel< ReactingMultiphaseCloud< CloudType > > > devolatilisationModel_
Devolatilisation model.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
virtual void writeFields() const
Write the field data for the cloud.
U
Definition: pEqn.H:72
const ReactingMultiphaseCloud & cloudCopy() const
Return a reference to the cloud copy.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
scalar dMassDevolatilisation_
Total mass transferred to continuous phase via devolatilisation.
const SurfaceReactionModel< ReactingMultiphaseCloud< CloudType > > & surfaceReaction() const
Return const access to reacting surface reaction model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Templated devolatilisation model class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const dimensionedVector & g
void info()
Print cloud information.
scalar dMassSurfaceReaction_
Total mass transferred to continuous phase via surface.
Templated surface reaction model class.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
const DevolatilisationModel< ReactingMultiphaseCloud< CloudType > > & devolatilisation() const
Return const access to devolatilisation model.
Namespace for OpenFOAM.