Loading...
Searching...
No Matches
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 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
29#include "solidChemistryModel.H"
30#include "reactingMixture.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class CompType, class SolidThermo>
35Foam::solidChemistryModel<CompType, SolidThermo>::solidChemistryModel
36(
37 typename CompType::reactionThermo& thermo
38)
39:
40 CompType(thermo),
41 ODESystem(),
42 Ys_(this->solidThermo().composition().Y()),
44 (
45 dynamic_cast<const reactingMixture<SolidThermo>&>
46 (
47 this->solidThermo()
48 )
49 ),
51 (
52 dynamic_cast<const reactingMixture<SolidThermo>&>
53 (
54 this->solidThermo()
55 ).speciesData()
56 ),
57 nSolids_(Ys_.size()),
58 nReaction_(reactions_.size()),
60 reactingCells_(this->mesh().nCells(), true)
61{
62 // create the fields for the chemistry sources
63 forAll(RRs_, fieldi)
64 {
65 RRs_.set
66 (
67 fieldi,
69 (
71 (
72 "RRs." + Ys_[fieldi].name(),
73 this->mesh().time().timeName(),
74 this->mesh(),
77 ),
78 this->mesh(),
80 )
81 );
82 }
83}
84
85
86// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87
88template<class CompType, class SolidThermo>
91{}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
96template<class CompType, class SolidThermo>
98(
99 const scalarField& deltaT
100)
101{
103 return 0;
104}
105
106
107template<class CompType, class SolidThermo>
110{
112 return volScalarField::null();
113}
114
115
116template<class CompType, class SolidThermo>
119{
120 auto tQdot = volScalarField::New
121 (
122 "Qdot",
124 this->mesh_,
126 );
127 scalarField& Qdot = tQdot.ref();
128
129 if (this->chemistry_)
130 {
131 forAll(Ys_, i)
132 {
133 forAll(Qdot, celli)
134 {
135 scalar hf = solidThermo_[i].Hc();
136 Qdot[celli] -= hf*RRs_[i][celli];
137 }
138 }
140
141 return tQdot;
142}
143
144
145template<class CompType, class SolidThermo>
147(
148 const label celli,
149 const bool active
150)
151{
152 reactingCells_[celli] = active;
153}
154
155// ************************************************************************* //
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
ODESystem()
Construct null.
Definition ODESystem.H:53
Foam::reactingMixture.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
virtual ~solidChemistryModel()
Destructor.
PtrList< volScalarField > & Ys_
Reference to solid mass fractions.
const PtrList< Reaction< SolidThermo > > & reactions_
Reactions.
label nSolids_
Number of solid components.
virtual scalar solve(const scalar deltaT)=0
Solve the reaction system for the given time step.
PtrList< volScalarField::Internal > RRs_
List of reaction rate per solid [kg/m3/s].
const PtrList< SolidThermo > & solidThermo_
Thermodynamic data of solids.
List< bool > reactingCells_
List of active reacting cells.
label nReaction_
Number of solid reactions.
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:51
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:75
basicSpecieMixture & composition
PtrList< volScalarField > & Y
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
auto & name
word timeName
Definition getTimeIndex.H:3
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
scalar Qdot
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299