reactingOneDim.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::regionModels::pyrolysisModels::reactingOneDim
26 
27 Description
28  Reacting, 1-D pyrolysis model
29 
30 SourceFiles
31  reactingOneDim.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef reactingOneDim_H
36 #define reactingOneDim_H
37 
38 #include "pyrolysisModel.H"
40 #include "radiationModel.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace regionModels
47 {
48 namespace pyrolysisModels
49 {
50 
51 
52 /*---------------------------------------------------------------------------*\
53  Class reactingOneDim Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class reactingOneDim
57 :
58  public pyrolysisModel
59 {
60  // Private Member Functions
61 
62  //- Read model controls
63  void readReactingOneDimControls();
64 
65 
66 protected:
67 
68  // Protected data
69 
70  //- Reference to solid thermo
72 
73  //- Reference to the solid chemistry model
75 
76  //- Pointer to radiation model
78 
79 
80  // Reference to solid thermo properties
81 
82  //- Density [kg/m^3]
84 
85  //- List of solid components
87 
88  // Non-const access to enthalpy
90 
91 
92  // Solution parameters
93 
94  //- Number of non-orthogonal correctors
96 
97  //- Maximum diffussivity
98  scalar maxDiff_;
99 
100  //- Minimum delta for combustion
101  scalar minimumDelta_;
102 
103 
104  // Fields
105 
106  //- Total gas mass flux to the primary region [kg/m^2/s]
108 
109  //- Sensible enthalpy gas flux [J/m2/s]
111 
112  //- Heat release rate [J/s/m^3]
114 
115 
116  // Source term fields
117 
118  //- Coupled region radiative heat flux [W/m^2]
119  // Requires user to input mapping info for coupled patches
120  // volScalarField qrCoupled_;
121 
122  //- In depth radiative heat flux [W/m^2]
124 
125 
126  // Checks
127 
128  //- Cumulative lost mass of the condensed phase [kg]
130 
131  //- Cumulative mass generation of the gas phase [kg]
133 
134  //- Total mass gas flux at the pyrolysing walls [kg/s]
135  scalar totalGasMassFlux_;
136 
137  //- Total heat release rate [J/s]
139 
140 
141  // Options
142 
143  //- Add gas enthalpy source term
144  bool gasHSource_;
145 
146  //- Add in depth radiation source term
147  bool qrHSource_;
148 
149  //- Use chemistry solvers (ode or sequential)
151 
152 
153  // Protected member functions
154 
155  //- Read control parameters from dictionary
156  bool read();
157 
158  //- Read control parameters from dict
159  bool read(const dictionary& dict);
160 
161  //- Update submodels
162  void updateFields();
163 
164  //- Update/move mesh based on change in mass
165  void updateMesh(const scalarField& mass0);
166 
167  //- Update radiative flux in pyrolysis region
168  void updateqr();
169 
170  //- Update enthalpy flux for pyrolysis gases
171  void updatePhiGas();
172 
173  //- Mass check
174  void calculateMassTransfer();
175 
176 
177  // Equations
178 
179  //- Solve continuity equation
180  void solveContinuity();
181 
182  //- Solve energy
183  void solveEnergy();
184 
185  //- Solve solid species mass conservation
186  void solveSpeciesMass();
187 
188 
189 public:
190 
191  //- Runtime type information
192  TypeName("reactingOneDim");
193 
194 
195  // Constructors
196 
197  //- Construct from type name and mesh
199  (
200  const word& modelType,
201  const fvMesh& mesh,
202  const word& regionType
203  );
204 
205  //- Construct from type name, mesh and dictionary
207  (
208  const word& modelType,
209  const fvMesh& mesh,
210  const dictionary& dict,
211  const word& regionType
212  );
213 
214  //- Disallow default bitwise copy construction
215  reactingOneDim(const reactingOneDim&) = delete;
216 
217 
218  //- Destructor
219  virtual ~reactingOneDim();
220 
221 
222  // Member Functions
223 
224  // Access
225 
226  //- Fields
227 
228  //- Return const density [Kg/m^3]
229  const volScalarField& rho() const;
230 
231  //- Return const temperature [K]
232  virtual const volScalarField& T() const;
233 
234  //- Return specific heat capacity [J/kg/K]
235  virtual const tmp<volScalarField> Cp() const;
236 
237  //- Return the region absorptivity [1/m]
238  virtual tmp<volScalarField> kappaRad() const;
239 
240  //- Return the region thermal conductivity [W/m/k]
241  virtual tmp<volScalarField> kappa() const;
242 
243  //- Return the total gas mass flux to primary region [kg/m^2/s]
244  virtual const surfaceScalarField& phiGas() const;
245 
246 
247  // Solution parameters
248 
249  //- Return the number of non-orthogonal correctors
250  inline label nNonOrthCorr() const;
251 
252  //- Return max diffusivity allowed in the solid
253  virtual scalar maxDiff() const;
254 
255 
256  // Helper functions
257 
258  //- External hook to add mass to the primary region
259  virtual scalar addMassSources
260  (
261  const label patchi, // patchi on primary region
262  const label facei // facei of patchi
263  );
264 
265  //- Mean diffusion number of the solid region
266  virtual scalar solidRegionDiffNo() const;
267 
268 
269  // Evolution
270 
271  //- Pre-evolve region
272  virtual void preEvolveRegion();
273 
274  //- Evolve the pyrolysis equations
275  virtual void evolveRegion();
276 
277 
278  // I-O
279 
280  //- Provide some feedback
281  virtual void info();
282 
283 
284  // Member Operators
285 
286  //- Disallow default bitwise assignment
287  void operator=(const reactingOneDim&) = delete;
288 };
289 
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace pyrolysisModels
294 } // End namespace regionModels
295 } // End namespace Foam
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
299 #include "reactingOneDimI.H"
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 #endif
304 
305 // ************************************************************************* //
virtual const volScalarField & T() const
Return const temperature [K].
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential)
dictionary dict
volScalarField chemistryQdot_
Heat release rate [J/s/m^3].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedScalar addedGasMass_
Cumulative mass generation of the gas phase [kg].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
scalar totalGasMassFlux_
Total mass gas flux at the pyrolysing walls [kg/s].
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m^2/s].
autoPtr< solidReactionThermo > solidThermo_
Reference to solid thermo.
Generic GeometricField class.
TypeName("reactingOneDim")
Runtime type information.
scalar minimumDelta_
Minimum delta for combustion.
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
void solveSpeciesMass()
Solve solid species mass conservation.
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
const volScalarField & rho() const
Fields.
label nNonOrthCorr() const
Return the number of non-orthogonal correctors.
dynamicFvMesh & mesh
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m^2/s].
reactingOneDim(const word &modelType, const fvMesh &mesh, const word &regionType)
Construct from type name and mesh.
A class for handling words, derived from string.
Definition: word.H:59
virtual void info()
Provide some feedback.
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
void updateqr()
Update radiative flux in pyrolysis region.
bool qrHSource_
Add in depth radiation source term.
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
void updateMesh(const scalarField &mass0)
Update/move mesh based on change in mass.
dimensionedScalar lostSolidMass_
Cumulative lost mass of the condensed phase [kg].
virtual scalar addMassSources(const label patchi, const label facei)
External hook to add mass to the primary region.
autoPtr< radiationModel > radiation_
Pointer to radiation model.
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void operator=(const reactingOneDim &)=delete
Disallow default bitwise assignment.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void evolveRegion()
Evolve the pyrolysis equations.
void updatePhiGas()
Update enthalpy flux for pyrolysis gases.
label nNonOrthCorr_
Number of non-orthogonal correctors.
PtrList< volScalarField > & Ys_
List of solid components.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
void solveContinuity()
Solve continuity equation.
A class for managing temporary objects.
Definition: PtrList.H:53
bool gasHSource_
Add gas enthalpy source term.
virtual void preEvolveRegion()
Pre-evolve region.
volScalarField qr_
Coupled region radiative heat flux [W/m^2].
bool read()
Read control parameters from dictionary.
Namespace for OpenFOAM.