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-2017 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_(nullptr),
102  constProps_(this->particleProperties()),
103  compositionModel_(nullptr),
104  phaseChangeModel_(nullptr),
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  this->deleteLostParticles();
115  }
116  }
117 
118  // Set storage for mass source fields and initialise to zero
119  forAll(rhoTrans_, i)
120  {
121  const word& specieName = thermo.carrier().species()[i];
122  rhoTrans_.set
123  (
124  i,
126  (
127  IOobject
128  (
129  this->name() + ":rhoTrans_" + specieName,
130  this->db().time().timeName(),
131  this->db(),
132  IOobject::READ_IF_PRESENT,
133  IOobject::AUTO_WRITE
134  ),
135  this->mesh(),
136  dimensionedScalar("zero", dimMass, 0.0)
137  )
138  );
139  }
140 
141  if (this->solution().resetSourcesOnStartup())
142  {
143  resetSourceTerms();
144  }
145 }
146 
147 
148 template<class CloudType>
150 (
152  const word& name
153 )
154 :
155  CloudType(c, name),
156  reactingCloud(),
157  cloudCopyPtr_(nullptr),
158  constProps_(c.constProps_),
159  compositionModel_(c.compositionModel_->clone()),
160  phaseChangeModel_(c.phaseChangeModel_->clone()),
161  rhoTrans_(c.rhoTrans_.size())
162 {
163  forAll(c.rhoTrans_, i)
164  {
165  const word& specieName = this->thermo().carrier().species()[i];
166  rhoTrans_.set
167  (
168  i,
170  (
171  IOobject
172  (
173  this->name() + ":rhoTrans_" + specieName,
174  this->db().time().timeName(),
175  this->db(),
176  IOobject::NO_READ,
177  IOobject::NO_WRITE,
178  false
179  ),
180  c.rhoTrans_[i]
181  )
182  );
183  }
184 }
185 
186 
187 template<class CloudType>
189 (
190  const fvMesh& mesh,
191  const word& name,
192  const ReactingCloud<CloudType>& c
193 )
194 :
195  CloudType(mesh, name, c),
196  reactingCloud(),
197  cloudCopyPtr_(nullptr),
198  constProps_(),
199  compositionModel_(c.compositionModel_->clone()),
200 // compositionModel_(nullptr),
201  phaseChangeModel_(nullptr),
202  rhoTrans_(0)
203 {}
204 
205 
206 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
207 
208 template<class CloudType>
210 {}
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
215 template<class CloudType>
217 (
218  parcelType& parcel,
219  const scalar lagrangianDt
220 )
221 {
222  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
223 
224  parcel.pc() = this->thermo().thermo().p()[parcel.cell()];
225  parcel.Y() = composition().YMixture0();
226 }
227 
228 
229 template<class CloudType>
231 (
232  parcelType& parcel,
233  const scalar lagrangianDt,
234  const bool fullyDescribed
235 )
236 {
237  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
238 
239  if (fullyDescribed)
240  {
241  checkSuppliedComposition
242  (
243  parcel.Y(),
244  composition().YMixture0(),
245  "YMixture"
246  );
247  }
248 
249  // derived information - store initial mass
250  parcel.mass0() = parcel.mass();
251 }
252 
253 
254 template<class CloudType>
256 {
257  cloudCopyPtr_.reset
258  (
259  static_cast<ReactingCloud<CloudType>*>
260  (
261  clone(this->name() + "Copy").ptr()
262  )
263  );
264 }
265 
266 
267 template<class CloudType>
269 {
270  cloudReset(cloudCopyPtr_());
271  cloudCopyPtr_.clear();
272 }
273 
274 
275 template<class CloudType>
277 {
278  CloudType::resetSourceTerms();
279  forAll(rhoTrans_, i)
280  {
281  rhoTrans_[i].field() = 0.0;
282  }
283 }
284 
285 
286 template<class CloudType>
288 (
289  const ReactingCloud<CloudType>& cloudOldTime
290 )
291 {
292  CloudType::relaxSources(cloudOldTime);
293 
294  typedef volScalarField::Internal dsfType;
295 
296  forAll(rhoTrans_, fieldi)
297  {
298  dsfType& rhoT = rhoTrans_[fieldi];
299  const dsfType& rhoT0 = cloudOldTime.rhoTrans()[fieldi];
300  this->relax(rhoT, rhoT0, "rho");
301  }
302 }
303 
304 
305 template<class CloudType>
307 {
308  CloudType::scaleSources();
309 
310  typedef volScalarField::Internal dsfType;
311 
312  forAll(rhoTrans_, fieldi)
313  {
314  dsfType& rhoT = rhoTrans_[fieldi];
315  this->scale(rhoT, "rho");
316  }
317 }
318 
319 
320 template<class CloudType>
322 {
323  if (this->solution().canEvolve())
324  {
325  typename parcelType::template
326  TrackingData<ReactingCloud<CloudType>> td(*this);
327 
328  this->solve(td);
329  }
330 }
331 
332 
333 template<class CloudType>
335 {
337 
338  this->updateMesh();
339 }
340 
341 
342 template<class CloudType>
344 {
345  CloudType::info();
346 
347  this->phaseChange().info(Info);
348 }
349 
350 
351 template<class CloudType>
353 {
354  if (compositionModel_.valid())
355  {
356  CloudType::particleType::writeFields(*this, this->composition());
357  }
358 }
359 
360 
361 // ************************************************************************* //
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()
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
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:163
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void restoreState()
Reset the current cloud to the previously stored state.
const speciesTable & species() const
Return the table of species.
Virtual abstract base class for templated ReactingCloud.
Definition: reactingCloud.H:48
void setModels()
Set cloud sub-models.
Definition: ReactingCloud.C:34
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:158
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:218
volScalarField::Internal & 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.
const basicSpecieMixture & carrier() const
Return reference to the gaseous components.
Definition: SLGThermo.C:108
Templated base class for reacting cloud.
Definition: ReactingCloud.H:62
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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...
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
Definition: ReactingCloud.C:77
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
messageStream Info
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:92
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.