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 -------------------------------------------------------------------------------
10 License
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 
37 namespace Foam
38 {
39 namespace regionModels
40 {
41 namespace pyrolysisModels
42 {
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
47 addToRunTimeSelectionTable(pyrolysisModel, thermo, mesh);
48 addToRunTimeSelectionTable(pyrolysisModel, thermo, dictionary);
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  {
64  readControls();
65  return true;
66  }
67 
68  return false;
69 }
70 
71 
73 {
75  {
76  readControls();
77  return true;
78  }
79 
80  return false;
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
86 thermo::thermo
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 
106 thermo::thermo
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 
150  volScalarField& h = solidThermo_->he();
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 
178 const volScalarField& thermo::T() const
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);
206  return surfaceScalarField::null();
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 } // End namespace surfaceFilmModels
213 } // End namespace regionModels
214 } // End namespace Foam
215 
216 // ************************************************************************* //
Foam::regionModels::pyrolysisModels::pyrolysisModel
Base class for pyrolysis models.
Definition: pyrolysisModel.H:61
volFields.H
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:55
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::regionModels::regionModel::active
Switch active() const
Return the active flag.
Definition: regionModelI.H:46
Foam::regionModels::pyrolysisModels::thermo::~thermo
virtual ~thermo()
Destructor.
Definition: thermo.C:129
Foam::solidThermo::New
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Return a pointer to a new solidThermo created from.
Definition: solidThermo.C:82
Foam::regionModels::pyrolysisModels::thermo::phiGas
virtual const surfaceScalarField & phiGas() const
Return the total gas mass flux to primary region [kg/m2/s].
Definition: thermo.C:202
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::regionModels::pyrolysisModels::thermo::Cp
virtual const tmp< volScalarField > Cp() const
Return specific heat capacity [J/kg/K].
Definition: thermo.C:184
Foam::regionModels::pyrolysisModels::thermo::radiation_
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
Definition: thermo.H:93
Foam::regionModels::pyrolysisModels::thermo::evolveRegion
virtual void evolveRegion()
Evolve the pyrolysis equations.
Definition: thermo.C:144
Foam::regionModels::pyrolysisModels::thermo::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: thermo.C:135
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
Foam::regionModels::pyrolysisModels::defineTypeNameAndDebug
defineTypeNameAndDebug(noPyrolysis, 0)
Foam::regionModels::pyrolysisModels::thermo::read
virtual bool read()
Read control parameters from dictionary.
Definition: thermo.C:60
Foam::regionModels::pyrolysisModels::thermo::readControls
void readControls()
Read control options.
Definition: thermo.C:52
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::radiation::radiationModel::New
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
Definition: radiationModelNew.C:36
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::regionModels::pyrolysisModels::thermo::solidThermo_
autoPtr< solidThermo > solidThermo_
Pointer to the solid chemistry model.
Definition: thermo.H:90
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Foam::regionModels::pyrolysisModels::thermo::kappaRad
virtual tmp< volScalarField > kappaRad() const
Return the region absorptivity [1/m].
Definition: thermo.C:190
Foam::regionModels::regionModel::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:500
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::regionModels::pyrolysisModels::pyrolysisModel::read
virtual bool read()
Read control parameters.
Definition: pyrolysisModel.C:52
Foam::regionModels::regionModel::solution
const dictionary & solution() const
Return the solution dictionary.
Definition: regionModelI.H:99
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::regionModels::pyrolysisModels::thermo::nNonOrthCorr_
label nNonOrthCorr_
Number of non-orthogonal correctors.
Definition: thermo.H:99
fvm.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::Time::controlDict
const dictionary & controlDict() const
Return read access to the controlDict dictionary.
Definition: Time.H:364
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
thermo.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
fvcLaplacian.H
Calculate the laplacian of the given field.
absorptionEmissionModel.H
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >::null
static const GeometricField< scalar, fvsPatchField, surfaceMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:32
Foam::regionModels::pyrolysisModels::thermo::rho
virtual const volScalarField & rho() const
Return density [kg/m3].
Definition: thermo.C:172
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::regionModels::pyrolysisModels::thermo::maxDiff_
scalar maxDiff_
Maximum diffusivity.
Definition: thermo.H:102
Foam::regionModels::pyrolysisModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(pyrolysisModel, noPyrolysis, mesh)
Foam::regionModels::pyrolysisModels::thermo::T
virtual const volScalarField & T() const
Return const temperature [K].
Definition: thermo.C:178
Foam::regionModels::pyrolysisModels::thermo::kappa
virtual tmp< volScalarField > kappa() const
Return the region thermal conductivity [W/m/k].
Definition: thermo.C:196
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:300