ReactingMultiphaseCloud.C
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-2022 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 
27 
28 #include "DevolatilisationModel.H"
29 #include "SurfaceReactionModel.H"
30 
31 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
36  devolatilisationModel_.reset
37  (
39  (
40  this->subModelProperties(),
41  *this
42  ).ptr()
43  );
44 
45  surfaceReactionModel_.reset
46  (
48  (
49  this->subModelProperties(),
50  *this
51  ).ptr()
52  );
53 }
54 
55 
56 template<class CloudType>
58 (
60 )
61 {
62  CloudType::cloudReset(c);
63 
64  devolatilisationModel_.reset(c.devolatilisationModel_.ptr());
65  surfaceReactionModel_.reset(c.surfaceReactionModel_.ptr());
66 
67  dMassDevolatilisation_ = c.dMassDevolatilisation_;
68  dMassSurfaceReaction_ = c.dMassSurfaceReaction_;
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
74 template<class CloudType>
76 (
77  const word& cloudName,
78  const volScalarField& rho,
79  const volVectorField& U,
80  const dimensionedVector& g,
81  const fluidThermo& carrierThermo,
82  const bool readFields
83 )
84 :
85  CloudType(cloudName, rho, U, g, carrierThermo, false),
86  cloudCopyPtr_(nullptr),
87  constProps_(this->particleProperties()),
88  devolatilisationModel_(nullptr),
89  surfaceReactionModel_(nullptr),
90  dMassDevolatilisation_(0.0),
91  dMassSurfaceReaction_(0.0)
92 {
93  setModels();
94 
95  if (readFields)
96  {
97  parcelType::readFields(*this, this->composition());
98  this->deleteLostParticles();
99  }
100 
101  if (this->solution().resetSourcesOnStartup())
102  {
103  resetSourceTerms();
104  }
105 }
106 
107 
108 template<class CloudType>
110 (
112  const word& name
113 )
114 :
115  CloudType(c, name),
116  cloudCopyPtr_(nullptr),
117  constProps_(c.constProps_),
118  devolatilisationModel_(c.devolatilisationModel_->clone()),
119  surfaceReactionModel_(c.surfaceReactionModel_->clone()),
120  dMassDevolatilisation_(c.dMassDevolatilisation_),
121  dMassSurfaceReaction_(c.dMassSurfaceReaction_)
122 {}
123 
124 
125 template<class CloudType>
127 (
128  const fvMesh& mesh,
129  const word& name,
131 )
132 :
133  CloudType(mesh, name, c),
134  cloudCopyPtr_(nullptr),
135  constProps_(),
136  devolatilisationModel_(nullptr),
137  surfaceReactionModel_(nullptr),
138  dMassDevolatilisation_(0.0),
139  dMassSurfaceReaction_(0.0)
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
145 template<class CloudType>
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
152 template<class CloudType>
154 (
155  parcelType& parcel,
156  const scalar lagrangianDt
157 )
158 {
159  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
160 
161  label idGas = this->composition().idGas();
162  label idLiquid = this->composition().idLiquid();
163  label idSolid = this->composition().idSolid();
164 
165  parcel.YGas() = this->composition().Y0(idGas);
166  parcel.YLiquid() = this->composition().Y0(idLiquid);
167  parcel.YSolid() = this->composition().Y0(idSolid);
168 }
169 
170 
171 template<class CloudType>
173 (
174  parcelType& parcel,
175  const scalar lagrangianDt,
176  const bool fullyDescribed
177 )
178 {
179  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
180 
181  if (fullyDescribed)
182  {
183  label idGas = this->composition().idGas();
184  label idLiquid = this->composition().idLiquid();
185  label idSolid = this->composition().idSolid();
186 
187  this->checkSuppliedComposition
188  (
189  parcel.YGas(),
190  this->composition().Y0(idGas),
191  "YGas"
192  );
193  this->checkSuppliedComposition
194  (
195  parcel.YLiquid(),
196  this->composition().Y0(idLiquid),
197  "YLiquid"
198  );
199  this->checkSuppliedComposition
200  (
201  parcel.YSolid(),
202  this->composition().Y0(idSolid),
203  "YSolid"
204  );
205  }
206 }
207 
208 
209 template<class CloudType>
211 {
212  cloudCopyPtr_.reset
213  (
215  (
216  clone(this->name() + "Copy").ptr()
217  )
218  );
219 }
220 
221 
222 template<class CloudType>
224 {
225  cloudReset(cloudCopyPtr_());
226  cloudCopyPtr_.clear();
227 }
228 
229 
230 template<class CloudType>
232 {
233  CloudType::resetSourceTerms();
234 }
235 
236 
237 template<class CloudType>
239 {
240  if (this->solution().canEvolve())
241  {
242  typename parcelType::trackingData td(*this);
243 
244  this->solve(*this, td);
245  }
246 }
247 
248 
249 template<class CloudType>
251 (
252  const polyTopoChangeMap& mapper
253 )
254 {
256 
257  this->topoChange();
258 }
259 
260 
261 template<class CloudType>
263 {
264  CloudType::info();
265 
266  this->devolatilisation().info(Info);
267  this->surfaceReaction().info(Info);
268 }
269 
270 
271 template<class CloudType>
273 {
274  if (this->compositionModel_.valid())
275  {
276  CloudType::particleType::writeFields(*this, this->composition());
277  }
278  else
279  {
280  CloudType::particleType::writeFields(*this);
281  }
282 }
283 
284 
285 // ************************************************************************* //
virtual void autoMap(const polyTopoChangeMap &)
Remap the cells of particles corresponding to the.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
DSMCCloud< dsmcParcel > CloudType
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.
basicSpecieMixture & composition
parcelType::constantProperties constProps_
Parcel constant properties.
autoPtr< SurfaceReactionModel< ReactingMultiphaseCloud< CloudType > > > surfaceReactionModel_
Surface reaction model.
Templated base class for multiphase reacting cloud.
void cloudReset(ReactingMultiphaseCloud< CloudType > &c)
Reset state of cloud.
T clone(const T &t)
Definition: List.H:54
void restoreState()
Reset the current cloud to the previously stored state.
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
void resetSourceTerms()
Reset the cloud source terms.
Base cloud calls templated on particle type.
Definition: Cloud.H:52
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
void storeState()
Store the current cloud state.
virtual ~ReactingMultiphaseCloud()
Destructor.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
autoPtr< DevolatilisationModel< ReactingMultiphaseCloud< CloudType > > > devolatilisationModel_
Devolatilisation model.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
virtual void writeFields() const
Write the field data for the cloud.
rhoEqn solve()
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
scalar dMassDevolatilisation_
Total mass transferred to continuous phase via devolatilisation.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
messageStream Info
Templated devolatilisation model class.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
ReactingMultiphaseCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const fluidThermo &carrierThermo, const bool readFields=true)
Construct given carrier fields and thermo.
void info()
Print cloud information.
scalar dMassSurfaceReaction_
Total mass transferred to continuous phase via surface.
Templated surface reaction model class.