ReactingHeterogeneousCloud.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) 2018-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30
31// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
32
33template<class CloudType>
35{
36 heterogeneousReactionModel_.reset
37 (
39 (
40 this->subModelProperties(),
41 *this
42 ).ptr()
43 );
44}
45
46
47template<class CloudType>
49(
51)
52{
53 CloudType::cloudReset(c);
54 heterogeneousReactionModel_.reset(c.heterogeneousReactionModel_.ptr());
55}
56
57
58// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59
60template<class CloudType>
62(
63 const word& cloudName,
64 const volScalarField& rho,
65 const volVectorField& U,
66 const dimensionedVector& g,
67 const SLGThermo& thermo,
68 bool readFields
69)
70:
71 CloudType(cloudName, rho, U, g, thermo, false),
73 cloudCopyPtr_(nullptr),
74 heterogeneousReactionModel_(nullptr)
75{
76 if (this->solution().active())
77 {
78 setModels();
79
80 if (readFields)
81 {
82 parcelType::readFields(*this, this->composition());
83 this->deleteLostParticles();
84 }
85 }
86}
87
88
89template<class CloudType>
91(
93 const word& name
94)
95:
96 CloudType(c, name),
98 cloudCopyPtr_(nullptr),
99 heterogeneousReactionModel_(c.heterogeneousReactionModel_->clone())
100{}
101
102
103template<class CloudType>
105(
106 const fvMesh& mesh,
107 const word& name,
109)
110:
111 CloudType(mesh, name, c),
113 cloudCopyPtr_(nullptr),
114 heterogeneousReactionModel_(c.heterogeneousReactionModel_->clone())
115{}
116
117
118// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119
120template<class CloudType>
122(
123 parcelType& parcel,
124 const scalar lagrangianDt
125)
126{
128 label idS = this->composition().idSolid();
129
130 // Set initial particle composition. Overwrite thermoProperties from
131 // ReactingCloud
132 parcel.Y() = this->composition().Y0(idS);
133
134 // Set initial progress to 0
135 parcel.F().setSize(heterogeneousReactionModel_->nF(), 0.0);
136
137 // Set the parcel to combust
138 parcel.canCombust() = 1;
139
140 // If rho0 was given in constProp use it. If not use the composition
141 // to set tho
142 if (this->constProps_.rho0() == -1)
143 {
144 const label idGas = this->composition().idGas();
145 const label idLiquid = this->composition().idLiquid();
146 const label idSolid = this->composition().idSolid();
147
148 const scalarField& Ygas = this->composition().Y0(idGas);
149 const scalarField& Yliq = this->composition().Y0(idLiquid);
150 const scalarField& Ysol = this->composition().Y0(idSolid);
151
152 const scalar p0 =
153 this->composition().thermo().thermo().p()[parcel.cell()];
154 const scalar T0 = this->constProps_.T0();
155
156 parcel.rho() = this->composition().rho(Ygas, Yliq, Ysol, T0, p0);
157 }
158}
159
160
161template<class CloudType>
163(
164 parcelType& parcel,
165 const scalar lagrangianDt,
166 const bool fullyDescribed
167)
168{
169 CloudType::checkParcelProperties(parcel, lagrangianDt, false);
170
171 const label solId = this->composition().idSolid();
172 const label liqId = this->composition().idLiquid();
173 const label gasId = this->composition().idGas();
174
175 // Check YMixture is pure solid
176 if
177 (
178 this->composition().YMixture0()[solId] != 1.0
179 || this->composition().YMixture0()[liqId] != 0.0
180 || this->composition().YMixture0()[gasId] != 0.0
181 )
182 {
184 << "The supplied composition must be : " << nl
185 << " YGasTot0 0 : " << nl
186 << " YLiquidTot0 0 : " << nl
187 << " YSolidTot0 1 : " << nl
188 << "This Cloud only works with pure solid particles."
189 << abort(FatalError);
190 }
191 if (this->composition().liquids().size() > 0)
192 {
194 << "The supplied composition has a liquid phase. " << nl
195 << "This Cloud only works with pure solid particles."
196 << abort(FatalError);
197 }
198}
199
200
201template<class CloudType>
203{
204 cloudCopyPtr_.reset
205 (
207 (
208 clone(this->name() + "Copy").ptr()
209 )
210 );
212
213
214template<class CloudType>
216{
217 cloudReset(cloudCopyPtr_());
218 cloudCopyPtr_.clear();
219}
220
221
222template<class CloudType>
224{
225 if (this->solution().canEvolve())
227 typename parcelType::trackingData td(*this);
228
229 this->solve(*this, td);
230 }
231}
233
234template<class CloudType>
236(
237 const mapPolyMesh& mapper
238)
241
242 this->updateMesh();
243}
244
246template<class CloudType>
250 heterogeneousReactionModel_->info(Info);
252
253
254template<class CloudType>
256{
258}
259
260
261template<class CloudType>
263readObjects(const objectRegistry& obr)
264{
266}
267
268
269template<class CloudType>
272{
274}
275
276
277// ************************************************************************* //
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
Base class for heterogeneous reacting models.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Templated base class for reacting heterogeneous cloud.
void storeState()
Store the current cloud state.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
virtual void readObjects(const objectRegistry &obr)
Read particle fields as objects from the obr registry.
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 cloudReset(ReactingHeterogeneousCloud< CloudType > &c)
Reset state of cloud.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
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
Allows specification of different writing frequency of objects registered to the database.
Definition: writeObjects.H:142
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
static void readObjects(Cloud< injectedParticle > &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
Virtual abstract base class for templated ReactingCloud.
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
CEqn solve()
scalar T0
Definition: createFields.H:22
const word cloudName(propsDict.get< word >("cloud"))