MulticomponentPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2024 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 \*---------------------------------------------------------------------------*/
25 
27 
28 #include "phaseSystem.H"
29 
30 #include "fvmDdt.H"
31 #include "fvmDiv.H"
32 #include "fvmSup.H"
33 #include "fvmLaplacian.H"
34 #include "fvcDdt.H"
35 #include "fvcDiv.H"
36 #include "fvMatrix.H"
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasePhaseModel>
42 (
43  const phaseSystem& fluid,
44  const word& phaseName,
45  const bool referencePhase,
46  const label index
47 )
48 :
49  BasePhaseModel(fluid, phaseName, referencePhase, index)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54 
55 template<class BasePhaseModel>
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 template<class BasePhaseModel>
64 {
65  this->thermo_->normaliseY();
66  BasePhaseModel::correctSpecies();
67 }
68 
69 
70 template<class BasePhaseModel>
72 {
73  return false;
74 }
75 
76 
77 template<class BasePhaseModel>
80 {
81  const volScalarField& alpha = *this;
82  const volScalarField& rho = this->rho();
83 
84  const tmp<surfaceScalarField> talphaRhoPhi(this->alphaRhoPhi());
85  const surfaceScalarField& alphaRhoPhi(talphaRhoPhi());
86 
87  return
88  (
89  fvm::ddt(alpha, rho, Yi)
90  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
91  + this->divj(Yi)
92  ==
93  alpha*this->R(Yi)
94 
95  - correction
96  (
97  fvm::Sp
98  (
99  max(this->residualAlpha() - alpha, scalar(0))*rho
100  /this->mesh().time().deltaT(),
101  Yi
102  )
103  )
104  );
105 }
106 
107 
108 template<class BasePhaseModel>
111 {
112  return this->thermo_->Y();
113 }
114 
115 
116 template<class BasePhaseModel>
119 {
120  return this->thermo_->Y(name);
121 }
122 
123 
124 template<class BasePhaseModel>
127 {
128  return this->thermo_->Y();
129 }
130 
131 
132 template<class BasePhaseModel>
134 (
135  const label speciei
136 ) const
137 {
138  return this->thermo_->solveSpecie(speciei);
139 }
140 
141 
142 // ************************************************************************* //
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
virtual void correctSpecies()
Correct the species fractions.
virtual bool solveSpecie(const label speciei) const
Should the given specie be solved for? I.e., is it active and.
MulticomponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
virtual ~MulticomponentPhaseModel()
Destructor.
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.