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 | Copyright (C) 2011-2016 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::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:
61 
62  // Private member functions
63 
64  //- Disallow default bitwise copy construct
66 
67  //- Disallow default bitwise assignment
68  void operator=(const reactingOneDim&);
69 
70  //- Read model controls
71  void readReactingOneDimControls();
72 
73 
74 protected:
75 
76  // Protected data
77 
78  //- Reference to the solid chemistry model
80 
81  //- Reference to solid thermo
83 
84  //- Pointer to radiation model
86 
87 
88  // Reference to solid thermo properties
89 
90  //- Density [kg/m3]
92 
93  //- List of solid components
95 
96  // Non-const access to enthalpy
98 
99 
100  // Solution parameters
101 
102  //- Number of non-orthogonal correctors
104 
105  //- Maximum diffussivity
106  scalar maxDiff_;
107 
108  //- Minimum delta for combustion
109  scalar minimumDelta_;
110 
111 
112  // Fields
113 
114  //- Total gas mass flux to the primary region [kg/m2/s]
116 
117  //- Sensible enthalpy gas flux [J/m2/s]
119 
120  //- Heat release [J/s/m3]
122 
123 
124  // Source term fields
125 
126  //- Coupled region radiative heat flux [W/m2]
127  // Requires user to input mapping info for coupled patches
128  //volScalarField QrCoupled_;
129 
130  //- In depth radiative heat flux [W/m2]
132 
133 
134  // Checks
135 
136  //- Cumulative lost mass of the condensed phase [kg]
138 
139  //- Cumulative mass generation of the gas phase [kg]
141 
142  //- Total mass gas flux at the pyrolysing walls [kg/s]
143  scalar totalGasMassFlux_;
144 
145  //- Total heat release rate [J/s]
147 
148 
149  // Options
150 
151  //- Add gas enthalpy source term
152  bool gasHSource_;
153 
154  //- Add in depth radiation source term
155  bool QrHSource_;
156 
157  //- Use chemistry solvers (ode or sequential)
159 
160 
161  // Protected member functions
162 
163  //- Read control parameters from dictionary
164  bool read();
165 
166  //- Read control parameters from dict
167  bool read(const dictionary& dict);
168 
169  //- Update submodels
170  void updateFields();
171 
172  //- Update/move mesh based on change in mass
173  void updateMesh(const scalarField& mass0);
174 
175  //- Update radiative flux in pyrolysis region
176  void updateQr();
177 
178  //- Update enthalpy flux for pyrolysis gases
179  void updatePhiGas();
180 
181  //- Mass check
182  void calculateMassTransfer();
183 
184 
185  // Equations
186 
187  //- Solve continuity equation
188  void solveContinuity();
189 
190  //- Solve energy
191  void solveEnergy();
192 
193  //- Solve solid species mass conservation
194  void solveSpeciesMass();
195 
196 
197 public:
198 
199  //- Runtime type information
200  TypeName("reactingOneDim");
201 
202 
203  // Constructors
204 
205  //- Construct from type name and mesh
207  (
208  const word& modelType,
209  const fvMesh& mesh,
210  const word& regionType
211  );
212 
213  //- Construct from type name, mesh and dictionary
215  (
216  const word& modelType,
217  const fvMesh& mesh,
218  const dictionary& dict,
219  const word& regionType
220  );
221 
222 
223 
224 
225  //- Destructor
226  virtual ~reactingOneDim();
227 
228 
229  // Member Functions
230 
231  // Access
232 
233  //- Fields
234 
235  //- Return const density [Kg/m3]
236  const volScalarField& rho() const;
237 
238  //- Return const temperature [K]
239  virtual const volScalarField& T() const;
240 
241  //- Return specific heat capacity [J/kg/K]
242  virtual const tmp<volScalarField> Cp() const;
243 
244  //- Return the region absorptivity [1/m]
245  virtual tmp<volScalarField> kappaRad() const;
246 
247  //- Return the region thermal conductivity [W/m/k]
248  virtual tmp<volScalarField> kappa() const;
249 
250  //- Return the total gas mass flux to primary region [kg/m2/s]
251  virtual const surfaceScalarField& phiGas() const;
252 
253 
254  // Solution parameters
255 
256  //- Return the number of non-orthogonal correctors
257  inline label nNonOrthCorr() const;
258 
259  //- Return max diffusivity allowed in the solid
260  virtual scalar maxDiff() const;
261 
262 
263  // Helper functions
264 
265  //- External hook to add mass to the primary region
266  virtual scalar addMassSources
267  (
268  const label patchi, // patchi on primary region
269  const label facei // facei of patchi
270  );
271 
272  //- Mean diffusion number of the solid region
273  virtual scalar solidRegionDiffNo() const;
274 
275 
276  // Evolution
277 
278  //- Pre-evolve region
279  virtual void preEvolveRegion();
280 
281  //- Evolve the pyrolysis equations
282  virtual void evolveRegion();
283 
284 
285  // I-O
286 
287  //- Provide some feedback
288  virtual void info();
289 };
290 
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 } // End namespace pyrolysisModels
295 } // End namespace regionModels
296 } // End namespace Foam
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 #include "reactingOneDimI.H"
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #endif
305 
306 // ************************************************************************* //
virtual scalar solidRegionDiffNo() const
Mean diffusion number of the solid region.
bool useChemistrySolvers_
Use chemistry solvers (ode or sequential)
solidReactionThermo & solidThermo_
Reference to solid thermo.
dictionary dict
volScalarField chemistrySh_
Heat release [J/s/m3].
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:137
Foam::solidReactionThermo.
scalar totalGasMassFlux_
Total mass gas flux at the pyrolysing walls [kg/s].
Generic GeometricField class.
virtual scalar maxDiff() const
Return max diffusivity allowed in the solid.
TypeName("reactingOneDim")
Runtime type information.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
scalar minimumDelta_
Minimum delta for combustion.
dimensionedScalar totalHeatRR_
Total heat release rate [J/s].
void solveSpeciesMass()
Solve solid species mass conservation.
virtual const volScalarField & T() const
Return const temperature [K].
autoPtr< basicSolidChemistryModel > solidChemistry_
Reference to the solid chemistry model.
dynamicFvMesh & mesh
surfaceScalarField phiGas_
Total gas mass flux to the primary region [kg/m2/s].
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
A class for handling words, derived from string.
Definition: word.H:59
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
virtual void info()
Provide some feedback.
volScalarField phiHsGas_
Sensible enthalpy gas flux [J/m2/s].
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
volScalarField Qr_
Coupled region radiative heat flux [W/m2].
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.
void updateQr()
Update radiative flux in pyrolysis region.
label patchi
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
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.
const volScalarField & rho() const
Fields.
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:53
label nNonOrthCorr() const
Return the number of non-orthogonal correctors.
void solveContinuity()
Solve continuity equation.
A class for managing temporary objects.
Definition: PtrList.H:54
bool gasHSource_
Add gas enthalpy source term.
virtual void preEvolveRegion()
Pre-evolve region.
bool read()
Read control parameters from dictionary.
bool QrHSource_
Add in depth radiation source term.
Namespace for OpenFOAM.