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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::ReactingCloud
29 
30 Group
31  grpLagrangianIntermediateClouds
32 
33 Description
34  Templated base class for reacting cloud
35 
36  - Adds to thermodynamic cloud
37  - Variable composition (single phase)
38  - Phase change
39 
40 SourceFiles
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 
53 namespace Foam
54 {
55 
56 // Forward declaration of classes
57 
58 template<class CloudType>
59 class CompositionModel;
60 
61 template<class CloudType>
62 class PhaseChangeModel;
63 
64 /*---------------------------------------------------------------------------*\
65  Class ReactingCloud Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 template<class CloudType>
69 class ReactingCloud
70 :
71  public CloudType,
72  public reactingCloud
73 {
74 public:
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
86 
87 
88 private:
89 
90  // Private data
91 
92  //- Cloud copy pointer
93  autoPtr<ReactingCloud<CloudType>> cloudCopyPtr_;
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 
105 protected:
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
153 
154 
155 public:
156 
157  // Constructors
158 
159  //- Construct given carrier gas fields
161  (
162  const word& cloudName,
163  const volScalarField& rho,
164  const volVectorField& U,
165  const dimensionedVector& g,
166  const SLGThermo& thermo,
167  bool readFields = true
168  );
169 
170  //- Copy constructor with new name
172 
173  //- Copy constructor with new name - creates bare cloud
175  (
176  const fvMesh& mesh,
177  const word& name,
179  );
180 
181  //- Construct and return clone based on (this) with new name
182  virtual autoPtr<Cloud<parcelType>> clone(const word& name)
183  {
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  //- Mass
237 
238  //- Return reference to mass source for field i
240  rhoTrans(const label i);
241 
242  //- Return const access to mass source fields
244  rhoTrans() const;
245 
246  //- Return reference to mass source fields
248  rhoTrans();
249 
250  //- Return mass source term for specie i - specie eqn
251  inline tmp<fvScalarMatrix> SYi
252  (
253  const label i,
254  volScalarField& Yi
255  ) const;
256 
257  //- Return tmp mass source for field i - fully explicit
259  Srho(const label i) const;
260 
261  //- Return tmp total mass source for carrier phase
262  // - fully explicit
263  inline tmp<volScalarField::Internal> Srho() const;
264 
265  //- Return total mass source term [kg/m3/s]
267 
268 
269  // Cloud evolution functions
270 
271  //- Set parcel thermo properties
273  (
274  parcelType& parcel,
275  const scalar lagrangianDt
276  );
277 
278  //- Check parcel properties
280  (
281  parcelType& parcel,
282  const scalar lagrangianDt,
283  const bool fullyDescribed
284  );
285 
286  //- Store the current cloud state
287  void storeState();
288 
289  //- Reset the current cloud to the previously stored state
290  void restoreState();
291 
292  //- Reset the cloud source terms
293  void resetSourceTerms();
294 
295  //- Apply relaxation to (steady state) cloud sources
296  void relaxSources(const ReactingCloud<CloudType>& cloudOldTime);
297 
298  //- Apply scaling to (transient) cloud sources
299  void scaleSources();
300 
301  //- Evolve the cloud
302  void evolve();
303 
304 
305  // Mapping
306 
307  //- Remap the cells of particles corresponding to the
308  // mesh topology change with a default tracking data object
309  virtual void autoMap(const mapPolyMesh&);
310 
311 
312  // I-O
313 
314  //- Print cloud information
315  void info();
316 
317  //- Write the field data for the cloud
318  virtual void writeFields() const;
319 
320  //- Write particle fields as objects into the obr registry
321  virtual void writeObjects(objectRegistry& obr) const;
322 };
323 
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 } // End namespace Foam
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
331 #include "ReactingCloudI.H"
332 
333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334 
335 #ifdef NoRepository
336  #include "ReactingCloud.C"
337 #endif
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 #endif
342 
343 // ************************************************************************* //
ParcelType particleType
Definition: Cloud.H:130
void resetSourceTerms()
Reset the cloud source terms.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:570
Templated phase change model class.
Definition: ReactingCloud.H:57
void storeState()
Store the current cloud state.
parcelType::constantProperties constProps_
Parcel constant properties.
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:30
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void restoreState()
Reset the current cloud to the previously stored state.
Virtual abstract base class for templated ReactingCloud.
Definition: reactingCloud.H:46
void setModels()
Set cloud sub-models.
Definition: ReactingCloud.C:29
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
virtual void writeFields() const
Write the field data for the cloud.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
void evolve()
Evolve the cloud.
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
A class for handling words, derived from Foam::string.
Definition: word.H:63
autoPtr< PhaseChangeModel< ReactingCloud< CloudType > > > phaseChangeModel_
Reacting phase change model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:60
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
const uniformDimensionedVectorField & g
autoPtr< CompositionModel< ReactingCloud< CloudType > > > compositionModel_
Reacting composition model.
CloudType cloudType
Type of cloud this cloud was instantiated for.
Definition: ReactingCloud.H:76
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void scaleSources()
Apply scaling to (transient) cloud sources.
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:37
Templated base class for reacting cloud.
Definition: ReactingCloud.H:64
U
Definition: pEqn.H:72
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const dimensionedScalar c
Speed of light in a vacuum.
virtual ~ReactingCloud()=default
Destructor.
void info()
Print cloud information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
Definition: ReactingCloud.C:72
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
Definition: ReactingCloud.C:53
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ReactingCloud.H:54
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
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.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition: ReactingCloud.H:81
ReactingCloud< CloudType > reactingCloudType
Convenience typedef for this cloud type.
Definition: ReactingCloud.H:86
Namespace for OpenFOAM.