constant.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) 2016-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "constant.H"
30 #include "fvcGrad.H"
31 #include "twoPhaseMixtureEThermo.H"
32 #include "fvmSup.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace temperaturePhaseChangeTwoPhaseMixtures
39 {
40  defineTypeNameAndDebug(constant, 0);
42  (
43  temperaturePhaseChangeTwoPhaseMixture,
44  constant,
45  components
46  );
47 }
48 }
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const thermoIncompressibleTwoPhaseMixture& mixture,
55  const fvMesh& mesh
56 )
57 :
58  temperaturePhaseChangeTwoPhaseMixture(mixture, mesh),
59  coeffC_
60  (
61  "coeffC",
63  optionalSubDict(type() + "Coeffs")
64  ),
65  coeffE_
66  (
67  "coeffE",
69  optionalSubDict(type() + "Coeffs")
70  )
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
78 {
79  const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
80 
81  const twoPhaseMixtureEThermo& thermo =
82  refCast<const twoPhaseMixtureEThermo>
83  (
84  mesh_.lookupObject<basicThermo>(basicThermo::dictName)
85  );
86 
87  const dimensionedScalar& TSat = thermo.TSat();
88 
90 
91  return Pair<tmp<volScalarField>>
92  (
93  coeffC_*mixture_.rho2()*max(TSat - T, T0),
94  -coeffE_*mixture_.rho1()*max(T - TSat, T0)
95  );
96 }
97 
98 
101 {
102  const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
103 
104  const twoPhaseMixtureEThermo& thermo =
105  refCast<const twoPhaseMixtureEThermo>
106  (
107  mesh_.lookupObject<basicThermo>(basicThermo::dictName)
108  );
109 
110  const dimensionedScalar& TSat = thermo.TSat();
111 
113 
114  volScalarField mDotE
115  (
116  "mDotE",
117  coeffE_*mixture_.rho1()*clamp(mixture_.alpha1(), zero_one{})
118  * max(T - TSat, T0)
119  );
120  volScalarField mDotC
121  (
122  "mDotC",
123  coeffC_*mixture_.rho2()*clamp(mixture_.alpha2(), zero_one{})
124  * max(TSat - T, T0)
125  );
126 
127  if (mesh_.time().outputTime())
128  {
129  mDotC.write();
130  mDotE.write();
131  }
132 
133  return Pair<tmp<volScalarField>>
134  (
135  tmp<volScalarField>(new volScalarField(mDotC)),
136  tmp<volScalarField>(new volScalarField(-mDotE))
137  );
138 }
139 
140 
143 {
144  const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
145 
146  const twoPhaseMixtureEThermo& thermo =
147  refCast<const twoPhaseMixtureEThermo>
148  (
149  mesh_.lookupObject<basicThermo>(basicThermo::dictName)
150  );
151 
152  const dimensionedScalar& TSat = thermo.TSat();
153 
154  return Pair<tmp<volScalarField>>
155  (
156  (
157  coeffC_*mixture_.rho2()*clamp(mixture_.alpha2(), zero_one{})
158  * pos(TSat - T)
159  ),
160  (
161  coeffE_*mixture_.rho1()*clamp(mixture_.alpha1(), zero_one{})
162  * pos(T - TSat)
163  )
164  );
165 }
166 
167 
170 {
171 
172  const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
173 
174  tmp<fvScalarMatrix> tTSource
175  (
176  new fvScalarMatrix
177  (
178  T,
180  )
181  );
182 
183  fvScalarMatrix& TSource = tTSource.ref();
184 
185  const twoPhaseMixtureEThermo& thermo =
186  refCast<const twoPhaseMixtureEThermo>
187  (
188  mesh_.lookupObject<basicThermo>(basicThermo::dictName)
189  );
190 
191  const dimensionedScalar& TSat = thermo.TSat();
192 
193  const dimensionedScalar L = mixture_.Hf2() - mixture_.Hf1();
194 
195  const volScalarField Vcoeff
196  (
197  coeffE_*mixture_.rho1()*clamp(mixture_.alpha1(), zero_one{})
198  * L*pos(T - TSat)
199  );
200  const volScalarField Ccoeff
201  (
202  coeffC_*mixture_.rho2()*clamp(mixture_.alpha2(), zero_one{})
203  * L*pos(TSat - T)
204  );
205 
206  TSource =
207  fvm::Sp(Vcoeff, T) - Vcoeff*TSat
208  + fvm::Sp(Ccoeff, T) - Ccoeff*TSat;
209 
210  return tTSource;
211 }
212 
213 
215 {
216 }
217 
218 
220 {
222  {
223  subDict(type() + "Coeffs").readEntry("coeffC", coeffC_);
224  subDict(type() + "Coeffs").readEntry("coeffE", coeffE_);
225 
226  return true;
227  }
228 
229  return false;
230 }
231 
232 
233 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
virtual bool read()
Read the transportProperties dictionary and update.
virtual void correct()
Correct the constant phaseChange model.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
type
Types of root.
Definition: Roots.H:52
const vector L(dict.get< vector >("L"))
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
virtual Pair< tmp< volScalarField > > mDot() const
Return the mass condensation and vaporisation rates as coefficients.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
dynamicFvMesh & mesh
virtual Pair< tmp< volScalarField > > mDotDeltaT() const
Return the mass condensation and vaporisation rates as a.
Calculate the gradient of the given field.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual bool read()
Read the transportProperties dictionary and update.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual tmp< fvScalarMatrix > TSource() const
Source for T equarion.
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Calculate the finiteVolume matrix for implicit and explicit sources.
constant(const thermoIncompressibleTwoPhaseMixture &mixture, const fvMesh &mesh)
Construct from components.
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
scalar T0
Definition: createFields.H:22
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133