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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-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
29#include "ReactingCloud.H"
30#include "CompositionModel.H"
31#include "PhaseChangeModel.H"
32
33// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
34
35template<class CloudType>
37{
38 compositionModel_.reset
39 (
41 (
42 this->subModelProperties(),
43 *this
44 ).ptr()
45 );
46
47 phaseChangeModel_.reset
48 (
50 (
51 this->subModelProperties(),
52 *this
53 ).ptr()
54 );
55}
56
57
58template<class CloudType>
60(
61 const scalarField& YSupplied,
62 const scalarField& Y,
63 const word& YName
64)
65{
66 if (YSupplied.size() != Y.size())
67 {
69 << YName << " supplied, but size is not compatible with "
70 << "parcel composition: " << nl << " "
71 << YName << "(" << YSupplied.size() << ") vs required composition "
72 << YName << "(" << Y.size() << ")" << nl
73 << abort(FatalError);
74 }
75}
76
77
78template<class CloudType>
80{
81 CloudType::cloudReset(c);
82
83 compositionModel_.reset(c.compositionModel_.ptr());
84 phaseChangeModel_.reset(c.phaseChangeModel_.ptr());
85}
86
87
88// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89
90template<class CloudType>
92(
93 const word& cloudName,
94 const volScalarField& rho,
95 const volVectorField& U,
96 const dimensionedVector& g,
97 const SLGThermo& thermo,
98 bool readFields
99)
100:
101 CloudType(cloudName, rho, U, g, thermo, false),
103 cloudCopyPtr_(nullptr),
104 constProps_(this->particleProperties()),
105 compositionModel_(nullptr),
106 phaseChangeModel_(nullptr),
107 rhoTrans_(thermo.carrier().species().size())
108{
109 if (this->solution().active())
110 {
111 setModels();
112
113 if (readFields)
114 {
115 parcelType::readFields(*this, this->composition());
116 this->deleteLostParticles();
117 }
118 }
119
120 // Set storage for mass source fields and initialise to zero
121 forAll(rhoTrans_, i)
122 {
123 const word& specieName = thermo.carrier().species()[i];
124 rhoTrans_.set
125 (
126 i,
128 (
130 (
131 this->name() + ":rhoTrans_" + specieName,
132 this->db().time().timeName(),
133 this->db(),
136 ),
137 this->mesh(),
139 )
140 );
141 }
142
143 if (this->solution().resetSourcesOnStartup())
144 {
146 }
147}
148
149
150template<class CloudType>
152(
154 const word& name
155)
156:
157 CloudType(c, name),
159 cloudCopyPtr_(nullptr),
160 constProps_(c.constProps_),
161 compositionModel_(c.compositionModel_->clone()),
162 phaseChangeModel_(c.phaseChangeModel_->clone()),
163 rhoTrans_(c.rhoTrans_.size())
164{
165 forAll(c.rhoTrans_, i)
166 {
167 const word& specieName = this->thermo().carrier().species()[i];
168 rhoTrans_.set
169 (
170 i,
172 (
174 (
175 this->name() + ":rhoTrans_" + specieName,
176 this->db().time().timeName(),
177 this->db(),
180 false
181 ),
182 c.rhoTrans_[i]
183 )
184 );
185 }
186}
187
188
189template<class CloudType>
191(
192 const fvMesh& mesh,
193 const word& name,
195)
196:
197 CloudType(mesh, name, c),
199 cloudCopyPtr_(nullptr),
200 constProps_(),
201 compositionModel_(c.compositionModel_->clone()),
202 phaseChangeModel_(nullptr),
203 rhoTrans_(0)
204{}
205
206
207// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
208
209template<class CloudType>
211(
212 parcelType& parcel,
213 const scalar lagrangianDt
214)
215{
216 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
217
218 parcel.Y() = composition().YMixture0();
219}
220
221
222template<class CloudType>
224(
225 parcelType& parcel,
226 const scalar lagrangianDt,
227 const bool fullyDescribed
228)
229{
230 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
231
232 if (fullyDescribed)
233 {
234 checkSuppliedComposition
235 (
236 parcel.Y(),
237 composition().YMixture0(),
238 "YMixture"
239 );
240 }
241
242 // derived information - store initial mass
243 parcel.mass0() = parcel.mass();
244}
245
246
247template<class CloudType>
249{
250 cloudCopyPtr_.reset
251 (
252 static_cast<ReactingCloud<CloudType>*>
253 (
254 clone(this->name() + "Copy").ptr()
255 )
256 );
257}
258
259
260template<class CloudType>
262{
263 cloudReset(cloudCopyPtr_());
264 cloudCopyPtr_.clear();
265}
266
267
268template<class CloudType>
270{
272 forAll(rhoTrans_, i)
273 {
274 rhoTrans_[i].field() = 0.0;
275 }
276}
277
279template<class CloudType>
281(
282 const ReactingCloud<CloudType>& cloudOldTime
283)
284{
285 CloudType::relaxSources(cloudOldTime);
287 typedef volScalarField::Internal dsfType;
288
289 forAll(rhoTrans_, fieldi)
290 {
291 dsfType& rhoT = rhoTrans_[fieldi];
292 const dsfType& rhoT0 = cloudOldTime.rhoTrans()[fieldi];
293 this->relax(rhoT, rhoT0, "rho");
294 }
296
297
298template<class CloudType>
300{
302
303 typedef volScalarField::Internal dsfType;
304
305 forAll(rhoTrans_, fieldi)
306 {
307 dsfType& rhoT = rhoTrans_[fieldi];
308 this->scale(rhoT, "rho");
309 }
310}
311
312
313template<class CloudType>
315{
316 if (this->solution().canEvolve())
318 typename parcelType::trackingData td(*this);
319
320 this->solve(*this, td);
321 }
322}
323
324
325template<class CloudType>
327{
329
330 this->updateMesh();
331}
332
333
334template<class CloudType>
336{
338
339 this->phaseChange().info(Info);
340}
341
342
343template<class CloudType>
345{
346 if (compositionModel_)
347 {
349 }
350}
351
352
353template<class CloudType>
355{
357}
358
359
360// ************************************************************************* //
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 reacting parcel composition model class Consists of carrier species (via thermo package),...
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
void info() const
Print cloud information.
Definition: DSMCCloud.C:940
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
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 scaleSources()
Apply scaling to (transient) cloud sources.
void resetSourceTerms()
Reset the cloud source terms.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
Templated phase change model class.
Templated base class for reacting cloud.
Definition: ReactingCloud.H:72
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
void setModels()
Set cloud sub-models.
Definition: ReactingCloud.C:36
void storeState()
Store the current cloud state.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
Definition: ReactingCloud.C:79
void scaleSources()
Apply scaling to (transient) cloud sources.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
Definition: ReactingCloud.C:60
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition: ReactingCloud.H:81
volScalarField::Internal & rhoTrans(const label i)
Mass.
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 evolve()
Evolve the cloud.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void resetSourceTerms()
Reset the cloud source terms.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
const basicSpecieMixture & carrier() const
Return reference to the gaseous components.
Definition: SLGThermo.C:111
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const speciesTable & species() const
Return the table of species.
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
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
Virtual abstract base class for templated ReactingCloud.
Definition: reactingCloud.H:51
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
PtrList< volScalarField > & Y
UEqn relax()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
word timeName
Definition: getTimeIndex.H:3
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const word cloudName(propsDict.get< word >("cloud"))