reactingOneDim.H
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) 2016 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
27Class
28 Foam::regionModels::pyrolysisModels::reactingOneDim
29
30Description
31 Reacting, 1-D pyrolysis model
32 NOTE: The moveMesh option can only be applied to solid reaction such as:
33 PMMA -> gas at the moment.
34
35SourceFiles
36 reactingOneDim.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef reactingOneDim_H
41#define reactingOneDim_H
42
43#include "pyrolysisModel.H"
45#include "radiationModel.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51namespace regionModels
52{
53namespace pyrolysisModels
54{
55
56
57/*---------------------------------------------------------------------------*\
58 Class reactingOneDim Declaration
59\*---------------------------------------------------------------------------*/
62:
63 public pyrolysisModel
64{
65private:
66
67 // Private member functions
68
69 //- No copy construct
70 reactingOneDim(const reactingOneDim&) = delete;
71
72 //- No copy assignment
73 void operator=(const reactingOneDim&) = delete;
74
75 //- Read model controls
76 void readReactingOneDimControls();
77
78
79protected:
80
81 // Protected data
82
83 //- Reference to solid thermo
85
86 //- Reference to the solid chemistry model
88
89 //- Pointer to radiation model
91
92
93 // Reference to solid thermo properties
94
95 //- Density [kg/m3]
97
98 //- List of solid components
100
101 // Non-const access to enthalpy
103
104
105 // Solution parameters
106
107 //- Number of non-orthogonal correctors
108 label nNonOrthCorr_;
109
110 //- Maximum diffusivity
111 scalar maxDiff_;
112
113 //- Minimum delta for combustion
114 scalar minimumDelta_;
115
116
117 // Fields
118
119 //- Total gas mass flux to the primary region [kg/m2/s]
121
122 //- Sensible enthalpy gas flux [J/m2/s]
124
125 //- Heat release rate [J/s/m3]
127
128
129 // Source term fields
130
131 //- Coupled region radiative heat flux [W/m2]
132 // Requires user to input mapping info for coupled patches
133 //volScalarField qrCoupled_;
134
135 //- In depth radiative heat flux [W/m2]
137
138
139 // Checks
140
141 //- Cumulative lost mass of the condensed phase [kg]
143
144 //- Cumulative mass generation of the gas phase [kg]
146
147 //- Total mass gas flux at the pyrolysing walls [kg/s]
148 scalar totalGasMassFlux_;
149
150 //- Total heat release rate [J/s]
152
153
154 // Options
155
156 //- Add gas enthalpy source term
157 bool gasHSource_;
158
159 //- Add in depth radiation source term
160 bool qrHSource_;
161
162 //- Use chemistry solvers (ode or sequential)
164
165
166 // Protected member functions
167
168 //- Read control parameters from dictionary
169 bool read();
170
171 //- Read control parameters from dict
172 bool read(const dictionary& dict);
173
174 //- Update submodels
175 void updateFields();
176
177 //- Update/move mesh based on change in mass
178 void updateMesh(const scalarField& mass0);
179
180 //- Update radiative flux in pyrolysis region
181 void updateqr();
182
183 //- Update enthalpy flux for pyrolysis gases
184 void updatePhiGas();
185
186 //- Mass check
188
189
190 // Equations
191
192 //- Solve continuity equation
193 void solveContinuity();
194
195 //- Solve energy
196 void solveEnergy();
197
198 //- Solve solid species mass conservation
199 void solveSpeciesMass();
200
201
202public:
203
204 //- Runtime type information
205 TypeName("reactingOneDim");
206
207
208 // Constructors
209
210 //- Construct from type name and mesh
212 (
213 const word& modelType,
214 const fvMesh& mesh,
215 const word& regionType
216 );
217
218 //- Construct from type name, mesh and dictionary
220 (
221 const word& modelType,
222 const fvMesh& mesh,
223 const dictionary& dict,
224 const word& regionType
225 );
226
227
228
229
230 //- Destructor
231 virtual ~reactingOneDim();
232
233
234 // Member Functions
235
236 // Access
237
238 //- Fields
239
240 //- Return const density [Kg/m3]
241 const volScalarField& rho() const;
242
243 //- Return const temperature [K]
244 virtual const volScalarField& T() const;
245
246 //- Return specific heat capacity [J/kg/K]
247 virtual const tmp<volScalarField> Cp() const;
248
249 //- Return the region absorptivity [1/m]
250 virtual tmp<volScalarField> kappaRad() const;
251
252 //- Return the region thermal conductivity [W/m/k]
253 virtual tmp<volScalarField> kappa() const;
254
255 //- Return the total gas mass flux to primary region [kg/m2/s]
256 virtual const surfaceScalarField& phiGas() const;
257
258
259 // Solution parameters
260
261 //- Return the number of non-orthogonal correctors
262 inline label nNonOrthCorr() const;
263
264 //- Return max diffusivity allowed in the solid
265 virtual scalar maxDiff() const;
266
267
268 // Helper functions
269
270 //- External hook to add mass to the primary region
271 virtual scalar addMassSources
272 (
273 const label patchi, // patchi on primary region
274 const label facei // facei of patchi
275 );
276
277 //- Mean diffusion number of the solid region
278 virtual scalar solidRegionDiffNo() const;
279
280
281 // Evolution
282
283 //- Pre-evolve region
284 virtual void preEvolveRegion();
285
286 //- Evolve the pyrolysis equations
287 virtual void evolveRegion();
288
289
290 // I-O
291
292 //- Provide some feedback
293 virtual void info();
294};
295
296
297// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298
299} // End namespace pyrolysisModels
300} // End namespace regionModels
301} // End namespace Foam
302
303// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304
305#include "reactingOneDimI.H"
306
307// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308
309#endif
310
311// ************************************************************************* //
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Reacting, 1-D pyrolysis model NOTE: The moveMesh option can only be applied to solid reaction such as...
void solveSpeciesMass()
Solve solid species mass conservation.
bool qrHSource_
Add in depth radiation source term.
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
void solveContinuity()
Solve continuity equation.
const volScalarField & rho() const
Fields.
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
dimensionedScalar lostSolidMass_
Cumulative lost mass of the condensed phase [kg].
void updateqr()
Update radiative flux in pyrolysis region.
TypeName("reactingOneDim")
Runtime type information.
scalar minimumDelta_
Minimum delta for combustion.
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential)
autoPtr< solidReactionThermo > solidThermo_
Reference to solid thermo.
PtrList< volScalarField > & Ys_
List of solid components.
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
virtual const volScalarField & T() const
Return const temperature [K].
dimensionedScalar addedGasMass_
Cumulative mass generation of the gas phase [kg].
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
volScalarField qr_
Coupled region radiative heat flux [W/m2].
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
virtual scalar addMassSources(const label patchi, const label facei)
External hook to add mass to the primary region.
void updatePhiGas()
Update enthalpy flux for pyrolysis gases.
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
label nNonOrthCorr() const
Return the number of non-orthogonal correctors.
void updateMesh(const scalarField &mass0)
Update/move mesh based on change in mass.
virtual void preEvolveRegion()
Pre-evolve region.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m2/s].
virtual void info()
Provide some feedback.
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
bool gasHSource_
Add gas enthalpy source term.
label nNonOrthCorr_
Number of non-orthogonal correctors.
bool read()
Read control parameters from dictionary.
volScalarField chemistryQdot_
Heat release rate [J/s/m3].
virtual void evolveRegion()
Evolve the pyrolysis equations.
scalar totalGasMassFlux_
Total mass gas flux at the pyrolysing walls [kg/s].
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73