ReactingCloud.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "ReactingCloud.H"
27 
28 #include "CompositionModel.H"
29 #include "PhaseChangeModel.H"
30 
31 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
36  compositionModel_.reset
37  (
39  (
40  this->subModelProperties(),
41  *this
42  ).ptr()
43  );
44 
45  phaseChangeModel_.reset
46  (
48  (
49  this->subModelProperties(),
50  *this
51  ).ptr()
52  );
53 }
54 
55 
56 template<class CloudType>
58 (
59  const scalarField& YSupplied,
60  const scalarField& Y,
61  const word& YName
62 )
63 {
64  if (YSupplied.size() != Y.size())
65  {
67  << YName << " supplied, but size is not compatible with "
68  << "parcel composition: " << nl << " "
69  << YName << "(" << YSupplied.size() << ") vs required composition "
70  << YName << "(" << Y.size() << ")" << nl
71  << abort(FatalError);
72  }
73 }
74 
75 
76 template<class CloudType>
78 {
79  CloudType::cloudReset(c);
80 
81  compositionModel_.reset(c.compositionModel_.ptr());
82  phaseChangeModel_.reset(c.phaseChangeModel_.ptr());
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
88 template<class CloudType>
90 (
91  const word& cloudName,
92  const volScalarField& rho,
93  const volVectorField& U,
94  const dimensionedVector& g,
95  const SLGThermo& thermo,
96  bool readFields
97 )
98 :
99  CloudType(cloudName, rho, U, g, thermo, false),
100  reactingCloud(),
101  cloudCopyPtr_(NULL),
102  constProps_(this->particleProperties()),
103  compositionModel_(NULL),
104  phaseChangeModel_(NULL),
105  rhoTrans_(thermo.carrier().species().size())
106 {
107  if (this->solution().active())
108  {
109  setModels();
110 
111  if (readFields)
112  {
113  parcelType::readFields(*this, this->composition());
114  }
115  }
116 
117  // Set storage for mass source fields and initialise to zero
118  forAll(rhoTrans_, i)
119  {
120  const word& specieName = thermo.carrier().species()[i];
121  rhoTrans_.set
122  (
123  i,
125  (
126  IOobject
127  (
128  this->name() + ":rhoTrans_" + specieName,
129  this->db().time().timeName(),
130  this->db(),
131  IOobject::READ_IF_PRESENT,
132  IOobject::AUTO_WRITE
133  ),
134  this->mesh(),
135  dimensionedScalar("zero", dimMass, 0.0)
136  )
137  );
138  }
139 
140  if (this->solution().resetSourcesOnStartup())
141  {
142  resetSourceTerms();
143  }
144 }
145 
146 
147 template<class CloudType>
149 (
151  const word& name
152 )
153 :
154  CloudType(c, name),
155  reactingCloud(),
156  cloudCopyPtr_(NULL),
157  constProps_(c.constProps_),
158  compositionModel_(c.compositionModel_->clone()),
159  phaseChangeModel_(c.phaseChangeModel_->clone()),
160  rhoTrans_(c.rhoTrans_.size())
161 {
162  forAll(c.rhoTrans_, i)
163  {
164  const word& specieName = this->thermo().carrier().species()[i];
165  rhoTrans_.set
166  (
167  i,
169  (
170  IOobject
171  (
172  this->name() + ":rhoTrans_" + specieName,
173  this->db().time().timeName(),
174  this->db(),
175  IOobject::NO_READ,
176  IOobject::NO_WRITE,
177  false
178  ),
179  c.rhoTrans_[i]
180  )
181  );
182  }
183 }
184 
185 
186 template<class CloudType>
188 (
189  const fvMesh& mesh,
190  const word& name,
191  const ReactingCloud<CloudType>& c
192 )
193 :
194  CloudType(mesh, name, c),
195  reactingCloud(),
196  cloudCopyPtr_(NULL),
197  constProps_(),
198  compositionModel_(c.compositionModel_->clone()),
199 // compositionModel_(NULL),
200  phaseChangeModel_(NULL),
201  rhoTrans_(0)
202 {}
203 
204 
205 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
206 
207 template<class CloudType>
209 {}
210 
211 
212 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
213 
214 template<class CloudType>
216 (
217  parcelType& parcel,
218  const scalar lagrangianDt
219 )
220 {
221  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
222 
223  parcel.pc() = this->thermo().thermo().p()[parcel.cell()];
224  parcel.Y() = composition().YMixture0();
225 }
226 
227 
228 template<class CloudType>
230 (
231  parcelType& parcel,
232  const scalar lagrangianDt,
233  const bool fullyDescribed
234 )
235 {
236  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
237 
238  if (fullyDescribed)
239  {
240  checkSuppliedComposition
241  (
242  parcel.Y(),
243  composition().YMixture0(),
244  "YMixture"
245  );
246  }
247 
248  // derived information - store initial mass
249  parcel.mass0() = parcel.mass();
250 }
251 
252 
253 template<class CloudType>
255 {
256  cloudCopyPtr_.reset
257  (
258  static_cast<ReactingCloud<CloudType>*>
259  (
260  clone(this->name() + "Copy").ptr()
261  )
262  );
263 }
264 
265 
266 template<class CloudType>
268 {
269  cloudReset(cloudCopyPtr_());
270  cloudCopyPtr_.clear();
271 }
272 
273 
274 template<class CloudType>
276 {
277  CloudType::resetSourceTerms();
278  forAll(rhoTrans_, i)
279  {
280  rhoTrans_[i].field() = 0.0;
281  }
282 }
283 
284 
285 template<class CloudType>
287 (
288  const ReactingCloud<CloudType>& cloudOldTime
289 )
290 {
291  CloudType::relaxSources(cloudOldTime);
292 
293  typedef DimensionedField<scalar, volMesh> dsfType;
294 
295  forAll(rhoTrans_, fieldi)
296  {
297  dsfType& rhoT = rhoTrans_[fieldi];
298  const dsfType& rhoT0 = cloudOldTime.rhoTrans()[fieldi];
299  this->relax(rhoT, rhoT0, "rho");
300  }
301 }
302 
303 
304 template<class CloudType>
306 {
307  CloudType::scaleSources();
308 
309  typedef DimensionedField<scalar, volMesh> dsfType;
310 
311  forAll(rhoTrans_, fieldi)
312  {
313  dsfType& rhoT = rhoTrans_[fieldi];
314  this->scale(rhoT, "rho");
315  }
316 }
317 
318 
319 template<class CloudType>
321 {
322  if (this->solution().canEvolve())
323  {
324  typename parcelType::template
325  TrackingData<ReactingCloud<CloudType>> td(*this);
326 
327  this->solve(td);
328  }
329 }
330 
331 
332 template<class CloudType>
334 {
335  typedef typename particle::TrackingData<ReactingCloud<CloudType>> tdType;
336 
337  tdType td(*this);
338 
339  Cloud<parcelType>::template autoMap<tdType>(td, mapper);
340 
341  this->updateMesh();
342 }
343 
344 
345 template<class CloudType>
347 {
348  CloudType::info();
349 
350  this->phaseChange().info(Info);
351 }
352 
353 
354 template<class CloudType>
356 {
357  if (this->size())
358  {
359  CloudType::particleType::writeFields(*this, this->composition());
360  }
361 }
362 
363 
364 // ************************************************************************* //
virtual ~ReactingCloud()
Destructor.
void resetSourceTerms()
Reset the cloud source terms.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
DSMCCloud< dsmcParcel > CloudType
Templated phase change model class.
Definition: ReactingCloud.H:55
void storeState()
Store the current cloud state.
parcelType::constantProperties constProps_
Parcel constant properties.
UEqn relax()
virtual void writeFields() const
Write the field data for the cloud.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
basicMultiComponentMixture & composition
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
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:48
void setModels()
Set cloud sub-models.
Definition: ReactingCloud.C:34
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void evolve()
Evolve the cloud.
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
autoPtr< PhaseChangeModel< ReactingCloud< CloudType > > > phaseChangeModel_
Reacting phase change model.
rhoEqn solve()
word timeName
Definition: getTimeIndex.H:3
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Base cloud calls templated on particle type.
Definition: Cloud.H:52
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
DimensionedField< scalar, volMesh > & rhoTrans(const label i)
Mass.
static const char nl
Definition: Ostream.H:262
autoPtr< CompositionModel< ReactingCloud< CloudType > > > compositionModel_
Reacting composition model.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void scaleSources()
Apply scaling to (transient) cloud sources.
Templated base class for reacting cloud.
Definition: ReactingCloud.H:62
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< DimensionedField< scalar, volMesh > > rhoTrans_
Mass transfer fields - one per carrier phase specie.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void info()
Print cloud information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const basicSpecieMixture & carrier() const
Return reference to the gaseous components.
Definition: SLGThermo.C:108
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
Definition: ReactingCloud.C:77
messageStream Info
const speciesTable & species() const
Return the table of species.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
Definition: ReactingCloud.C:58
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ReactingCloud.H:52
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.