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-------------------------------------------------------------------------------
10License
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"
32#include "fvmSup.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace temperaturePhaseChangeTwoPhaseMixtures
39{
42 (
43 temperaturePhaseChangeTwoPhaseMixture,
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{
80
81 const twoPhaseMixtureEThermo& thermo =
82 refCast<const twoPhaseMixtureEThermo>
83 (
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
103 volScalarField limitedAlpha1
104 (
105 min(max(mixture_.alpha1(), scalar(0)), scalar(1))
106 );
107
108 volScalarField limitedAlpha2
109 (
110 min(max(mixture_.alpha2(), scalar(0)), scalar(1))
111 );
112
113 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
114
115 const twoPhaseMixtureEThermo& thermo =
116 refCast<const twoPhaseMixtureEThermo>
117 (
118 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
119 );
120
121 const dimensionedScalar& TSat = thermo.TSat();
122
124
125 volScalarField mDotE
126 (
127 "mDotE", coeffE_*mixture_.rho1()*limitedAlpha1*max(T - TSat, T0)
128 );
129 volScalarField mDotC
130 (
131 "mDotC", coeffC_*mixture_.rho2()*limitedAlpha2*max(TSat - T, T0)
132 );
133
134 if (mesh_.time().outputTime())
135 {
136 mDotC.write();
137 mDotE.write();
138 }
139
140 return Pair<tmp<volScalarField>>
141 (
142 tmp<volScalarField>(new volScalarField(mDotC)),
143 tmp<volScalarField>(new volScalarField(-mDotE))
144 );
145}
146
147
150{
151 volScalarField limitedAlpha1
152 (
153 min(max(mixture_.alpha1(), scalar(0)), scalar(1))
154 );
155
156 volScalarField limitedAlpha2
157 (
158 min(max(mixture_.alpha2(), scalar(0)), scalar(1))
159 );
160
161 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
162
163 const twoPhaseMixtureEThermo& thermo =
164 refCast<const twoPhaseMixtureEThermo>
165 (
166 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
167 );
168
169 const dimensionedScalar& TSat = thermo.TSat();
170
171 return Pair<tmp<volScalarField>>
172 (
173 coeffC_*mixture_.rho2()*limitedAlpha2*pos(TSat - T),
174 coeffE_*mixture_.rho1()*limitedAlpha1*pos(T - TSat)
175 );
176}
177
178
181{
182
183 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
184
185 tmp<fvScalarMatrix> tTSource
186 (
188 (
189 T,
191 )
192 );
193
194 fvScalarMatrix& TSource = tTSource.ref();
195
196 const twoPhaseMixtureEThermo& thermo =
197 refCast<const twoPhaseMixtureEThermo>
198 (
199 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
200 );
201
202 const dimensionedScalar& TSat = thermo.TSat();
203
204 dimensionedScalar L = mixture_.Hf2() - mixture_.Hf1();
205
206 volScalarField limitedAlpha1
207 (
208 min(max(mixture_.alpha1(), scalar(0)), scalar(1))
209 );
210
211 volScalarField limitedAlpha2
212 (
213 min(max(mixture_.alpha2(), scalar(0)), scalar(1))
214 );
215
216 const volScalarField Vcoeff
217 (
218 coeffE_*mixture_.rho1()*limitedAlpha1*L*pos(T - TSat)
219 );
220 const volScalarField Ccoeff
221 (
222 coeffC_*mixture_.rho2()*limitedAlpha2*L*pos(TSat - T)
223 );
224
225 TSource =
226 fvm::Sp(Vcoeff, T) - Vcoeff*TSat
227 + fvm::Sp(Ccoeff, T) - Ccoeff*TSat;
228
229 return tTSource;
230}
231
232
234{
235}
236
237
239{
241 {
242 subDict(type() + "Coeffs").readEntry("coeffC", coeffC_);
243 subDict(type() + "Coeffs").readEntry("coeffE", coeffE_);
244
245 return true;
246 }
247
248 return false;
249}
250
251
252// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
static const word dictName
Definition: basicThermo.H:256
const dimensionedScalar & rho1() const
Return const-access to phase1 density.
const dimensionedScalar & rho2() const
Return const-access to phase2 density.
Constant dispersed-phase particle diameter model.
const Type & lookupObject(const word &name, const bool recursive=false) const
constant condensation/saturation model.
const thermoIncompressibleTwoPhaseMixture & mixture_
Reference to the thermoIncompressibleTwoPhaseMixture.
virtual bool read()
Read the transportProperties dictionary and update.
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Return the mass condensation and vaporisation rates as a.
virtual Pair< tmp< volScalarField > > mDotDeltaT() const
Return the mass condensation and vaporisation rates as a.
virtual tmp< fvScalarMatrix > TSource() const
Source for T equarion.
virtual void correct()
Correct the constant phaseChange model.
virtual bool read()
Read the transportProperties dictionary and update.
virtual Pair< tmp< volScalarField > > mDot() const
Return the mass condensation and vaporisation rates as coefficients.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
Calculate the gradient of the given field.
Calculate the finiteVolume matrix for implicit and explicit sources.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
type
Types of root.
Definition: Roots.H:55
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:45
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
scalar T0
Definition: createFields.H:22
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
const vector L(dict.get< vector >("L"))