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-2021 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  residualAlpha_
51  (
52  "residualAlpha",
53  dimless,
54  fluid.mesh().solverDict("Yi")
55  )
56 {
57  PtrList<volScalarField>& Y = this->thermo_->composition().Y();
58 
59  forAll(Y, i)
60  {
61  if (this->thermo_->composition().solve(i))
62  {
63  const label j = YActive_.size();
64  YActive_.resize(j + 1);
65  YActive_.set(j, &Y[i]);
66  }
67  }
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
73 template<class BasePhaseModel>
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class BasePhaseModel>
82 {
83  this->thermo_->composition().normalise();
85 }
86 
87 
88 template<class BasePhaseModel>
90 {
91  return false;
92 }
93 
94 
95 template<class BasePhaseModel>
98 {
99  const volScalarField& alpha = *this;
100  const volScalarField& rho = this->thermo().rho();
101 
102  const tmp<surfaceScalarField> talphaRhoPhi(this->alphaRhoPhi());
103  const surfaceScalarField& alphaRhoPhi(talphaRhoPhi());
104 
105  return
106  (
107  fvm::ddt(alpha, rho, Yi)
108  + fvm::div(alphaRhoPhi, Yi, "div(" + alphaRhoPhi.name() + ",Yi)")
109  + this->divj(Yi)
110  ==
111  alpha*this->R(Yi)
112 
113  + fvc::ddt(residualAlpha_*rho, Yi)
114  - fvm::ddt(residualAlpha_*rho, Yi)
115  );
116 }
117 
118 
119 template<class BasePhaseModel>
122 {
123  return this->thermo_->composition().Y();
124 }
125 
126 
127 template<class BasePhaseModel>
130 {
131  return this->thermo_->composition().Y(name);
132 }
133 
134 
135 template<class BasePhaseModel>
138 {
139  return this->thermo_->composition().Y();
140 }
141 
142 
143 template<class BasePhaseModel>
146 {
147  return YActive_;
148 }
149 
150 
151 template<class BasePhaseModel>
154 {
155  return YActive_;
156 }
157 
158 
159 // ************************************************************************* //
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
virtual void correctSpecies()
Correct the species fractions.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
Calculate the matrix for the laplacian of the field.
const dimensionSet dimless
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
MultiComponentPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
fluid correctSpecies()
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Calculate the first temporal derivative.
dynamicFvMesh & mesh
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
virtual ~MultiComponentPhaseModel()
Destructor.
Calculate the matrix for the divergence of the given field and flux.
#define R(A, B, C, D, E, F, K, M)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
Calculate the matrix for implicit and explicit sources.