Loading...
Searching...
No Matches
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{
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{
79 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
80
81 const twoPhaseMixtureEThermo& thermo =
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
99Foam::Pair<Foam::tmp<Foam::volScalarField>>
101{
102 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
103
104 const twoPhaseMixtureEThermo& thermo =
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().writeTime())
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
141Foam::Pair<Foam::tmp<Foam::volScalarField>>
143{
144 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
145
146 const twoPhaseMixtureEThermo& thermo =
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
168Foam::tmp<Foam::fvScalarMatrix>
170{
171
172 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
173
175 auto& TSource = tTSource.ref();
176
177 const twoPhaseMixtureEThermo& thermo =
179 (
180 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
181 );
182
183 const dimensionedScalar& TSat = thermo.TSat();
184
185 const dimensionedScalar L = mixture_.Hf2() - mixture_.Hf1();
186
187 const volScalarField Vcoeff
188 (
189 coeffE_*mixture_.rho1()*clamp(mixture_.alpha1(), zero_one{})
190 * L*pos(T - TSat)
191 );
192 const volScalarField Ccoeff
193 (
194 coeffC_*mixture_.rho2()*clamp(mixture_.alpha2(), zero_one{})
195 * L*pos(TSat - T)
196 );
197
198 TSource =
199 fvm::Sp(Vcoeff, T) - Vcoeff*TSat
200 + fvm::Sp(Ccoeff, T) - Ccoeff*TSat;
201
202 return tTSource;
203}
204
205
207{
208}
209
210
212{
214 {
215 subDict(type() + "Coeffs").readEntry("coeffC", coeffC_);
216 subDict(type() + "Coeffs").readEntry("coeffE", coeffE_);
217
218 return true;
219 }
220
221 return false;
222}
223
224
225// ************************************************************************* //
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:66
static const word dictName
The dictionary name ("thermophysicalProperties").
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.
constant(const thermoIncompressibleTwoPhaseMixture &mixture, const fvMesh &mesh)
Construct from components.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
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:53
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
psiReactionThermo & thermo
scalar T0
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
const vector L(dict.get< vector >("L"))