coefficientPhaseChange.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) 2021-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 
26 #include "coefficientPhaseChange.H"
27 #include "fvcGrad.H"
28 #include "multicomponentThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
39 }
40 }
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::fv::coefficientPhaseChange::readCoeffs()
45 {
46  C_.read(coeffs());
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& name,
55  const word& modelType,
56  const fvMesh& mesh,
57  const dictionary& dict
58 )
59 :
61  (
62  name,
63  modelType,
64  mesh,
65  dict,
66  {false, false},
67  {false, false}
68  ),
69  C_("C", dimMass/dimArea/dimTime, NaN)
70 {
71  readCoeffs();
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
79 {
80  const volScalarField& alpha1 =
81  mesh().lookupObject<volScalarField>(alphaNames().first());
82 
84  C_*alpha1()*mag(fvc::grad(alpha1))()();
85 
86  if (specieis().first() != -1)
87  {
88  tmDot.ref() *= specieThermos().first().Y()[specieis().first()];
89  }
90 
91  return tmDot;
92 }
93 
94 
96 (
97  const volScalarField& alpha,
98  const volScalarField& rho,
99  fvMatrix<scalar>& eqn
100 ) const
101 {
102  const label i = index(alphaNames(), eqn.psi().name());
103 
104  if (i != -1)
105  {
106  const volScalarField& alpha1 =
107  mesh().lookupObject<volScalarField>(alphaNames().first());
108 
109  volScalarField::Internal mDotByAlpha1(C_*mag(fvc::grad(alpha1)));
110 
111  if (specieis().first() != -1)
112  {
113  mDotByAlpha1 *= specieThermos().first().Y()[specieis().first()];
114  }
115 
116  if (i == 0)
117  {
118  eqn -= fvm::Sp(mDotByAlpha1, eqn.psi());
119  }
120  else
121  {
122  eqn +=
123  mDotByAlpha1*alpha1
124  - correction(fvm::Sp(mDotByAlpha1, eqn.psi()));
125  }
126  }
127  else
128  {
130  }
131 }
132 
133 
135 (
136  const volScalarField& alpha,
137  const volScalarField& rho,
138  const volScalarField& Yi,
139  fvMatrix<scalar>& eqn
140 ) const
141 {
142  const label i = index(alphaNames(), eqn.psi().name());
143 
144  if (i == 0 && specieis().first() != -1 && Yi.member() == specie())
145  {
146  eqn -= fvm::Sp(C_*alpha()*mag(fvc::grad(alpha))()(), Yi);
147  }
148  else
149  {
151  }
152 }
153 
154 
156 {
158  {
159  readCoeffs();
160  return true;
161  }
162  else
163  {
164  return false;
165  }
166 }
167 
168 
169 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:149
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
void read(const dictionary &)
Update the value of dimensioned<Type>
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
This simple model generates a phase change between two phases calculated from the following expressio...
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the phase change rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
coefficientPhaseChange(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Override the compressible continuity equation to add.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:218
Base class for phase change models in which only a single component changes phase....
virtual bool read(const dictionary &dict)
Read source dictionary.
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Override the energy equation to add the phase change heat, or.
Base class of the thermophysical property types.
Definition: specie.H:66
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the gradient of the given field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
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
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
const dimensionSet dimMass
const dimensionSet dimArea
labelList fv(nPoints)
dictionary dict