Loading...
Searching...
No Matches
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
30#include "CompositionModel.H"
31#include "PhaseChangeModel.H"
32
33// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
34
35template<class CloudType>
37{
39 (
40 CompositionModel<ReactingCloud<CloudType>>::New
41 (
42 this->subModelProperties(),
43 *this
44 ).ptr()
45 );
46
48 (
49 PhaseChangeModel<ReactingCloud<CloudType>>::New
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>
79void Foam::ReactingCloud<CloudType>::cloudReset(ReactingCloud<CloudType>& c)
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>
91Foam::ReactingCloud<CloudType>::ReactingCloud
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 fieldName
124 (
126 (
127 this->name(),
128 "rhoTrans_" + thermo.carrier().species()[i]
129 )
130 );
131
132 rhoTrans_.set
133 (
134 i,
136 (
137 IOobject
138 (
139 fieldName,
140 this->db().time().timeName(),
141 this->db(),
145 ),
146 this->mesh(),
148 )
149 );
150 }
151
152 if (this->solution().resetSourcesOnStartup())
155 }
156}
157
158
159template<class CloudType>
160Foam::ReactingCloud<CloudType>::ReactingCloud
161(
162 ReactingCloud<CloudType>& c,
163 const word& name
164)
166 CloudType(c, name),
168 cloudCopyPtr_(nullptr),
174 forAll(c.rhoTrans_, i)
175 {
176 const word fieldName
177 (
179 (
180 this->name(),
181 "rhoTrans_" + this->thermo().carrier().species()[i]
182 )
183 );
184
185 rhoTrans_.set
186 (
187 i,
189 (
191 (
192 fieldName,
193 this->db().time().timeName(),
194 this->db(),
198 ),
199 c.rhoTrans_[i]
201 );
202 }
203}
204
205
206template<class CloudType>
207Foam::ReactingCloud<CloudType>::ReactingCloud
208(
209 const fvMesh& mesh,
210 const word& name,
211 const ReactingCloud<CloudType>& c
212)
213:
214 CloudType(mesh, name, c),
215 reactingCloud(),
216 cloudCopyPtr_(nullptr),
217 constProps_(),
218 compositionModel_(c.compositionModel_->clone()),
219 phaseChangeModel_(nullptr),
221{}
222
223
224// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225
226template<class CloudType>
228(
229 parcelType& parcel,
230 const scalar lagrangianDt
231)
232{
233 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
234
235 parcel.Y() = composition().YMixture0();
236}
237
238
239template<class CloudType>
241(
242 parcelType& parcel,
243 const scalar lagrangianDt,
244 const bool fullyDescribed
245)
246{
247 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
248
249 if (fullyDescribed)
250 {
251 checkSuppliedComposition
252 (
253 parcel.Y(),
254 composition().YMixture0(),
255 "YMixture"
256 );
257 }
259 // derived information - store initial mass
260 parcel.mass0() = parcel.mass();
261}
262
263
264template<class CloudType>
266{
267 cloudCopyPtr_.reset
268 (
269 static_cast<ReactingCloud<CloudType>*>
270 (
271 clone(this->name() + "Copy").ptr()
272 )
273 );
274}
275
276
277template<class CloudType>
280 cloudReset(cloudCopyPtr_());
281 cloudCopyPtr_.clear();
282}
283
284
285template<class CloudType>
287{
288 CloudType::resetSourceTerms();
289 forAll(rhoTrans_, i)
291 rhoTrans_[i].field() = 0.0;
292 }
293}
294
295
296template<class CloudType>
298(
299 const ReactingCloud<CloudType>& cloudOldTime
300)
301{
302 CloudType::relaxSources(cloudOldTime);
303
304 typedef volScalarField::Internal dsfType;
305
306 forAll(rhoTrans_, fieldi)
307 {
308 dsfType& rhoT = rhoTrans_[fieldi];
309 const dsfType& rhoT0 = cloudOldTime.rhoTrans()[fieldi];
310 this->relax(rhoT, rhoT0, "rho");
311 }
312}
313
314
315template<class CloudType>
317{
318 CloudType::scaleSources();
319
320 typedef volScalarField::Internal dsfType;
321
322 forAll(rhoTrans_, fieldi)
323 {
324 dsfType& rhoT = rhoTrans_[fieldi];
325 this->scale(rhoT, "rho");
326 }
327}
328
329
330template<class CloudType>
332{
333 if (this->solution().canEvolve())
334 {
335 typename parcelType::trackingData td(*this);
337 this->solve(*this, td);
338 }
339}
340
341
342template<class CloudType>
343void Foam::ReactingCloud<CloudType>::autoMap(const mapPolyMesh& mapper)
346
347 this->updateMesh();
348}
349
350
351template<class CloudType>
355
356 this->phaseChange().info();
357}
358
359
360template<class CloudType>
362{
366 }
367}
369
370template<class CloudType>
372{
374}
375
376
377// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition Cloud.C:108
virtual void writeFields() const
Write the field data for the cloud of particles Dummy at.
Definition CloudIO.C:285
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition Cloud.C:297
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition Cloud.C:125
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
const word & cloudName() const
const IOdictionary & particleProperties() const
const fvMesh & mesh() const
DimensionedField< scalar, volMesh > Internal
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Templated phase change model class.
Templated base class for reacting cloud.
autoPtr< PhaseChangeModel< ReactingCloud< CloudType > > > phaseChangeModel_
Reacting phase change model.
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
void setModels()
Set cloud sub-models.
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.
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.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
volScalarField::Internal & rhoTrans(const label i)
Return reference to mass source for field i.
virtual void writeFields() const
Write the field data for the cloud.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
void evolve()
Evolve the cloud.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
parcelType::constantProperties constProps_
Parcel constant properties.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
autoPtr< CompositionModel< ReactingCloud< CloudType > > > compositionModel_
Reacting composition model.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void resetSourceTerms()
Reset the cloud source terms.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition SLGThermo.H:63
const basicSpecieMixture & carrier() const
Return reference to the gaseous components.
Definition SLGThermo.C:104
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const speciesTable & species() const
Return the table of species.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Definition cloud.C:91
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
Virtual abstract base class for templated ReactingCloud.
reactingCloud()=default
Null constructor.
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
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:66
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:600
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
word timeName
Definition getTimeIndex.H:3
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
CEqn solve()
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299