thermo.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) 2016-2019 OpenCFD Ltd
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "thermo.H"
30#include "volFields.H"
32#include "fvm.H"
33#include "fvcLaplacian.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace regionModels
40{
41namespace pyrolysisModels
42{
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
49
50// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51
53{
54 const dictionary& solution = this->solution().subDict("SIMPLE");
55 solution.readEntry("nNonOrthCorr", nNonOrthCorr_);
56 time().controlDict().readEntry("maxDi", maxDiff_);
57}
58
59
61{
63 {
65 return true;
66 }
67
68 return false;
69}
70
71
73{
75 {
77 return true;
78 }
79
80 return false;
81}
82
83
84// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85
87(
88 const word& modelType,
89 const fvMesh& mesh,
90 const word& regionType
91)
92:
93 pyrolysisModel(modelType, mesh, regionType),
94 solidThermo_(solidThermo::New(regionMesh())),
95 radiation_(radiation::radiationModel::New(solidThermo_->T())),
96 nNonOrthCorr_(-1),
97 maxDiff_(10)
98{
99 if (active())
100 {
101 readControls();
102 }
103}
104
105
107(
108 const word& modelType,
109 const fvMesh& mesh,
110 const dictionary& dict,
111 const word& regionType
112)
113:
114 pyrolysisModel(modelType, mesh, dict, regionType),
115 solidThermo_(solidThermo::New(regionMesh())),
116 radiation_(radiation::radiationModel::New(solidThermo_->T())),
117 nNonOrthCorr_(-1),
118 maxDiff_(10)
119{
120 if (active())
121 {
122 readControls();
123 }
124}
125
126
127// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
128
130{}
131
132
133// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
134
136{
137 if (active())
138 {
140 }
141}
142
143
145{
146 if (active())
147 {
148 Info<< "\nEvolving energy in region: " << regionMesh().name() << endl;
149
151
152 for (int nonOrth=0; nonOrth<=nNonOrthCorr_; nonOrth++)
153 {
155 fvScalarMatrix hEqn
156 (
157 fvm::ddt(rho(), h)
160 - fvc::laplacian(solidThermo_->kappa(), T())
161 );
162
163 hEqn.relax();
164 hEqn.solve();
165 }
166
167 solidThermo_->correct();
168 }
169}
170
171
173{
174 return solidThermo_->rho();
175}
176
177
179{
180 return solidThermo_->T();
181}
182
183
185{
186 return solidThermo_->Cp();
187}
188
189
191{
192 return radiation_->absorptionEmission().a();
193}
194
195
197{
198 return solidThermo_->kappa();
199}
200
201
203{
205 << "phiGas field not available for " << type() << abort(FatalError);
207}
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212} // End namespace surfaceFilmModels
213} // End namespace regionModels
214} // End namespace Foam
215
216// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static const GeometricField< scalar, fvsPatchField, surfaceMesh > & null()
Return a null geometric field.
const dictionary & controlDict() const
Return read access to the controlDict dictionary.
Definition: Time.H:364
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const word & name() const
Return reference to name.
Definition: fvMesh.H:310
Top level model for radiation modelling.
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
virtual bool read()
Read control parameters.
Pyrolysis model which solves only the energy equation in the region.
Definition: thermo.H:61
scalar maxDiff_
Maximum diffusivity.
Definition: thermo.H:102
virtual const volScalarField & rho() const
Return density [kg/m3].
Definition: thermo.C:172
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
Definition: thermo.C:196
virtual const volScalarField & T() const
Return const temperature [K].
Definition: thermo.C:178
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
Definition: thermo.C:190
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: thermo.C:184
virtual void preEvolveRegion()
Pre-evolve region.
Definition: thermo.C:135
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
Definition: thermo.H:93
void readControls()
Read control options.
Definition: thermo.C:52
autoPtr< solidThermo > solidThermo_
Pointer to the solid chemistry model.
Definition: thermo.H:90
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
Definition: thermo.C:202
label nNonOrthCorr_
Number of non-orthogonal correctors.
Definition: thermo.H:99
virtual bool read()
Read control parameters from dictionary.
Definition: thermo.C:60
virtual void evolveRegion()
Evolve the pyrolysis equations.
Definition: thermo.C:144
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:99
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Switch active() const
Return the active flag.
Definition: regionModelI.H:46
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:55
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:66
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Calculate the laplacian of the given field.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
volScalarField & alpha
dictionary dict
volScalarField & h