ReactingMultiphaseCloud.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
30
33
34// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
35
36template<class CloudType>
38{
39 devolatilisationModel_.reset
40 (
42 (
43 this->subModelProperties(),
44 *this
45 ).ptr()
46 );
47
48 surfaceReactionModel_.reset
49 (
51 (
52 this->subModelProperties(),
53 *this
54 ).ptr()
55 );
56}
57
58
59template<class CloudType>
61(
63)
64{
65 CloudType::cloudReset(c);
66
67 devolatilisationModel_.reset(c.devolatilisationModel_.ptr());
68 surfaceReactionModel_.reset(c.surfaceReactionModel_.ptr());
69
70 dMassDevolatilisation_ = c.dMassDevolatilisation_;
71 dMassSurfaceReaction_ = c.dMassSurfaceReaction_;
72}
73
74
75// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76
77template<class CloudType>
79(
80 const word& cloudName,
81 const volScalarField& rho,
82 const volVectorField& U,
83 const dimensionedVector& g,
84 const SLGThermo& thermo,
85 bool readFields
86)
87:
88 CloudType(cloudName, rho, U, g, thermo, false),
90 cloudCopyPtr_(nullptr),
91 constProps_(this->particleProperties()),
92 devolatilisationModel_(nullptr),
93 surfaceReactionModel_(nullptr),
94 dMassDevolatilisation_(0.0),
95 dMassSurfaceReaction_(0.0)
96{
97 if (this->solution().active())
98 {
99 setModels();
100
101 if (readFields)
102 {
103 parcelType::readFields(*this, this->composition());
104 this->deleteLostParticles();
105 }
106 }
107
108 if (this->solution().resetSourcesOnStartup())
109 {
111 }
112}
113
114
115template<class CloudType>
117(
119 const word& name
120)
121:
122 CloudType(c, name),
124 cloudCopyPtr_(nullptr),
125 constProps_(c.constProps_),
126 devolatilisationModel_(c.devolatilisationModel_->clone()),
127 surfaceReactionModel_(c.surfaceReactionModel_->clone()),
128 dMassDevolatilisation_(c.dMassDevolatilisation_),
129 dMassSurfaceReaction_(c.dMassSurfaceReaction_)
130{}
131
132
133template<class CloudType>
135(
136 const fvMesh& mesh,
137 const word& name,
139)
140:
141 CloudType(mesh, name, c),
143 cloudCopyPtr_(nullptr),
144 constProps_(),
145 devolatilisationModel_(nullptr),
146 surfaceReactionModel_(nullptr),
147 dMassDevolatilisation_(0.0),
148 dMassSurfaceReaction_(0.0)
149{}
151
152// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153
154template<class CloudType>
156{}
157
159// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160
161template<class CloudType>
163(
164 parcelType& parcel,
165 const scalar lagrangianDt
166)
167{
168 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
169
170 label idGas = this->composition().idGas();
171 label idLiquid = this->composition().idLiquid();
172 label idSolid = this->composition().idSolid();
173
174 parcel.YGas() = this->composition().Y0(idGas);
175 parcel.YLiquid() = this->composition().Y0(idLiquid);
176 parcel.YSolid() = this->composition().Y0(idSolid);
178
179 // If rho0 was given in constProp use it. If not use the composition
180 // to set tho
181 if (constProps_.rho0() == -1)
182 {
183 const scalarField& Ygas = this->composition().Y0(idGas);
184 const scalarField& Yliq = this->composition().Y0(idLiquid);
185 const scalarField& Ysol = this->composition().Y0(idSolid);
186 const scalar p0 =
187 this->composition().thermo().thermo().p()[parcel.cell()];
188 const scalar T0 = constProps_.T0();
189
190 parcel.rho() = this->composition().rho(Ygas, Yliq, Ysol, T0, p0);
191 }
192}
193
194
195template<class CloudType>
197(
198 parcelType& parcel,
199 const scalar lagrangianDt,
200 const bool fullyDescribed
201)
202{
203 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
205 if (fullyDescribed)
206 {
207 label idGas = this->composition().idGas();
208 label idLiquid = this->composition().idLiquid();
209 label idSolid = this->composition().idSolid();
210
211 this->checkSuppliedComposition
212 (
213 parcel.YGas(),
214 this->composition().Y0(idGas),
215 "YGas"
216 );
217 this->checkSuppliedComposition
218 (
219 parcel.YLiquid(),
220 this->composition().Y0(idLiquid),
221 "YLiquid"
222 );
223 this->checkSuppliedComposition
224 (
225 parcel.YSolid(),
226 this->composition().Y0(idSolid),
227 "YSolid"
228 );
229 }
230}
231
232
233template<class CloudType>
235{
236 cloudCopyPtr_.reset
237 (
239 (
240 clone(this->name() + "Copy").ptr()
241 )
242 );
243}
244
245
246template<class CloudType>
248{
249 cloudReset(cloudCopyPtr_());
250 cloudCopyPtr_.clear();
251}
252
253
254template<class CloudType>
258}
259
260
261template<class CloudType>
264 if (this->solution().canEvolve())
265 {
266 typename parcelType::trackingData td(*this);
267
268 this->solve(*this, td);
269 }
270}
272
273template<class CloudType>
275(
276 const mapPolyMesh& mapper
278{
281 this->updateMesh();
282}
283
284
285template<class CloudType>
289
290 this->devolatilisation().info(Info);
291 this->surfaceReaction().info(Info);
292}
294
295template<class CloudType>
297{
298 if (this->compositionModel_)
299 {
301 }
302}
303
304
305// ************************************************************************* //
const uniformDimensionedVectorField & g
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:118
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:308
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
void info() const
Print cloud information.
Definition: DSMCCloud.C:940
Templated devolatilisation model class.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
void resetSourceTerms()
Reset the cloud source terms.
Templated base class for multiphase reacting cloud.
void setModels()
Set cloud sub-models.
void storeState()
Store the current cloud state.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
virtual void writeFields() const
Write the field data for the cloud.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void cloudReset(ReactingMultiphaseCloud< CloudType > &c)
Reset state of cloud.
void resetSourceTerms()
Reset the cloud source terms.
virtual ~ReactingMultiphaseCloud()
Destructor.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
Templated surface reaction model class.
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m3].
Class used to pass data into container.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Virtual abstract base class for templated reactingMultiphaseCloud.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
basicSpecieMixture & composition
const volScalarField & p0
Definition: EEqn.H:36
dynamicFvMesh & mesh
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
CEqn solve()
scalar T0
Definition: createFields.H:22
const word cloudName(propsDict.get< word >("cloud"))