solidChemistryModel.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) 2011-2017 OpenFOAM Foundation
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 "solidChemistryModel.H"
29#include "reactingMixture.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CompType, class SolidThermo>
35(
37)
38:
39 CompType(thermo),
40 ODESystem(),
41 Ys_(this->solidThermo().composition().Y()),
42 reactions_
43 (
44 dynamic_cast<const reactingMixture<SolidThermo>&>
45 (
46 this->solidThermo()
47 )
48 ),
49 solidThermo_
50 (
51 dynamic_cast<const reactingMixture<SolidThermo>&>
52 (
53 this->solidThermo()
54 ).speciesData()
55 ),
56 nSolids_(Ys_.size()),
57 nReaction_(reactions_.size()),
58 RRs_(nSolids_),
59 reactingCells_(this->mesh().nCells(), true)
60{
61 // create the fields for the chemistry sources
62 forAll(RRs_, fieldi)
63 {
64 RRs_.set
65 (
66 fieldi,
68 (
70 (
71 "RRs." + Ys_[fieldi].name(),
72 this->mesh().time().timeName(),
73 this->mesh(),
76 ),
77 this->mesh(),
79 )
80 );
81 }
82}
83
84
85// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
86
87template<class CompType, class SolidThermo>
90{}
91
92
93// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94
95template<class CompType, class SolidThermo>
97(
98 const scalarField& deltaT
99)
100{
102 return 0;
103}
104
105
106template<class CompType, class SolidThermo>
109{
111 return volScalarField::null();
112}
113
114
115template<class CompType, class SolidThermo>
118{
120 (
122 (
124 (
125 "Qdot",
126 this->mesh_.time().timeName(),
127 this->mesh_,
130 false
131 ),
132 this->mesh_,
134 )
135 );
136
137 if (this->chemistry_)
138 {
139 scalarField& Qdot = tQdot.ref();
140
141 forAll(Ys_, i)
142 {
143 forAll(Qdot, celli)
144 {
145 scalar hf = solidThermo_[i].Hc();
146 Qdot[celli] -= hf*RRs_[i][celli];
147 }
148 }
149 }
150
151 return tQdot;
152}
153
154
155template<class CompType, class SolidThermo>
157(
158 const label celli,
159 const bool active
160)
161{
162 reactingCells_[celli] = active;
163}
164
165// ************************************************************************* //
ReactionThermo reactionThermo
Thermo type.
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:50
Foam::reactingMixture.
Extends base solid chemistry model by adding a thermo package, and ODE functions.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
virtual ~solidChemistryModel()
Destructor.
PtrList< volScalarField > & Ys_
Reference to solid mass fractions.
PtrList< volScalarField::Internal > RRs_
List of reaction rate per solid [kg/m3/s].
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s3].
void setCellReacting(const label celli, const bool active)
Set reacting status of cell, celli.
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:55
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
basicSpecieMixture & composition
PtrList< volScalarField > & Y
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
word timeName
Definition: getTimeIndex.H:3
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
CEqn solve()
scalar Qdot
Definition: solveChemistry.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333