phaseChangeTwoPhaseMixture.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) 2011-2018 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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(phaseChangeTwoPhaseMixture, 0);
33  defineRunTimeSelectionTable(phaseChangeTwoPhaseMixture, components);
34 }
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const word& type,
41  const volVectorField& U,
42  const surfaceScalarField& phi
43 )
44 :
45  incompressibleTwoPhaseMixture(U, phi),
46  phaseChangeTwoPhaseMixtureCoeffs_(optionalSubDict(type + "Coeffs")),
47  pSat_("pSat", dimPressure, lookup("pSat"))
48 {}
49 
50 
51 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
52 
55 {
56  volScalarField alphalCoeff(1.0/rho1() - alpha1_*(1.0/rho1() - 1.0/rho2()));
57  Pair<tmp<volScalarField>> mDotAlphal = this->mDotAlphal();
58 
59  return Pair<tmp<volScalarField>>
60  (
61  alphalCoeff*mDotAlphal[0],
62  alphalCoeff*mDotAlphal[1]
63  );
64 }
65 
68 {
69  dimensionedScalar pCoeff(1.0/rho1() - 1.0/rho2());
70  Pair<tmp<volScalarField>> mDotP = this->mDotP();
71 
72  return Pair<tmp<volScalarField>>(pCoeff*mDotP[0], pCoeff*mDotP[1]);
73 }
74 
75 
77 {
79  {
81  lookup("pSat") >> pSat_;
82 
83  return true;
84  }
85  else
86  {
87  return false;
88  }
89 }
90 
91 
92 // ************************************************************************* //
const dimensionedScalar & rho2() const
Return const-access to phase2 density.
virtual bool read()=0
Read the transportProperties dictionary and update.
surfaceScalarField & phi
virtual Pair< tmp< volScalarField > > mDotP() const
Return the mass condensation and vaporisation rates as coefficients.
Pair< tmp< volScalarField > > vDotAlphal() const
Return the volumetric condensation and vaporisation rates as a.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
phaseChangeTwoPhaseMixture(const phaseChangeTwoPhaseMixture &)
Disallow copy construct.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
volScalarField alpha1_
const dimensionSet dimPressure
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Pair< tmp< volScalarField > > vDotP() const
Return the volumetric condensation and vaporisation rates as.
dimensionedScalar pSat_
Saturation vapour pressure.
virtual bool read()
Read base transportProperties dictionary.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const dimensionedScalar & rho1() const
Return const-access to phase1 density.