interfaceHeatResistance.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) 2020-2022 OpenCFD Ltd.
9 Copyright (C) 2020 Henning Scheufler
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "interfaceHeatResistance.H"
32#include "cutCellIso.H"
35#include "wallPolyPatch.H"
36#include "fvcSmooth.H"
37#include "fvmSup.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace temperaturePhaseChangeTwoPhaseMixtures
44{
47 (
48 temperaturePhaseChangeTwoPhaseMixture,
50 components
51 );
52}
53}
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
57Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::
58interfaceHeatResistance
59(
60 const thermoIncompressibleTwoPhaseMixture& mixture,
61 const fvMesh& mesh
62)
63:
64 temperaturePhaseChangeTwoPhaseMixture(mixture, mesh),
65 R_
66 (
67 "R",
68 dimPower/dimArea/dimTemperature, optionalSubDict(type() + "Coeffs")
69 ),
70
71 interfaceArea_
72 (
73 IOobject
74 (
75 "interfaceArea",
76 mesh_.time().timeName(),
77 mesh_,
78 IOobject::NO_READ,
79 IOobject::AUTO_WRITE
80 ),
81 mesh_,
83 ),
84
85 mDotc_
86 (
87 IOobject
88 (
89 "mDotc",
90 mesh_.time().timeName(),
91 mesh_,
92 IOobject::NO_READ,
93 IOobject::AUTO_WRITE
94 ),
95 mesh_,
97 ),
98
99 mDote_
100 (
101 IOobject
102 (
103 "mDote",
104 mesh_.time().timeName(),
105 mesh_,
106 IOobject::NO_READ,
107 IOobject::AUTO_WRITE
108 ),
109 mesh_,
111 ),
112
113 mDotcSpread_
114 (
115 IOobject
116 (
117 "mDotcSpread",
118 mesh_.time().timeName(),
119 mesh_,
120 IOobject::NO_READ,
121 IOobject::NO_WRITE
122 ),
123 mesh_,
125 ),
126
127 mDoteSpread_
128 (
129 IOobject
130 (
131 "mDoteSpread",
132 mesh_.time().timeName(),
133 mesh_,
134 IOobject::NO_READ,
135 IOobject::NO_WRITE
136 ),
137 mesh_,
139 ),
140
141 spread_
142 (
143 optionalSubDict(type() + "Coeffs").get<scalar>("spread")
144 )
145{
146 correct();
147}
148
149
150// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
151
154vDotAlphal() const
155{
156 dimensionedScalar alphalCoeff(1.0/mixture_.rho1());
157
158 return Pair<tmp<volScalarField>>
159 (
160 (alphalCoeff*mDotc_)/(mixture_.alpha2() + SMALL),
161 -(alphalCoeff*mDote_)/(mixture_.alpha1() + SMALL)
162 );
163}
164
165
168mDotAlphal() const
169{
170 volScalarField limitedAlpha1
171 (
172 min(max(mixture_.alpha1(), scalar(0)), scalar(1))
173 );
174
175 volScalarField limitedAlpha2
176 (
177 min(max(mixture_.alpha2(), scalar(0)), scalar(1))
178 );
179
180 return Pair<tmp<volScalarField>>
181 (
182 (mDotc_/(limitedAlpha2 + SMALL)),
183 -(mDote_/(limitedAlpha1 + SMALL))
184 );
185}
186
187
190mDot() const
191{
192 return Pair<tmp<volScalarField>>
193 (
194 tmp<volScalarField>(mDotc_),
195 tmp<volScalarField>(mDote_)
196 );
197}
198
199
202mDotDeltaT() const
203{
204 const twoPhaseMixtureEThermo& thermo =
205 refCast<const twoPhaseMixtureEThermo>
206 (
207 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
208 );
209
210 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
211
212 const dimensionedScalar& TSat = thermo.TSat();
213
214 Pair<tmp<volScalarField>> mDotce(mDot());
215
216 return Pair<tmp<volScalarField>>
217 (
218 mDotc_*pos(TSat - T.oldTime())/(TSat - T.oldTime()),
219 -mDote_*pos(T.oldTime() - TSat)/(T.oldTime() - TSat)
220 );
221}
222
223
226TSource() const
227{
228 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
229
231 auto& TSource = tTSource.ref();
232
233 const twoPhaseMixtureEThermo& thermo =
234 refCast<const twoPhaseMixtureEThermo>
235 (
236 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
237 );
238
239 const dimensionedScalar& TSat = thermo.TSat();
240
241 // interface heat resistance
242 volScalarField IHRcoeff(interfaceArea_*R_);
243
244 TSource = fvm::Sp(IHRcoeff, T) - IHRcoeff*TSat;
245
246 return tTSource;
247}
248
249
251correct()
252{
253 // Update Interface
254 updateInterface();
255
256 // Update mDotc_ and mDote_
257 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
258
259 const twoPhaseMixtureEThermo& thermo =
260 refCast<const twoPhaseMixtureEThermo>
261 (
262 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
263 );
264
265 const dimensionedScalar& TSat = thermo.TSat();
267
268 const dimensionedScalar L(mag(mixture_.Hf2() - mixture_.Hf1()));
269
270 // interface heat resistance
271 mDotc_ = interfaceArea_*R_*max(TSat - T, T0)/L;
272 mDote_ = interfaceArea_*R_*max(T - TSat, T0)/L;
273
274 // Calculate the spread sources
276 (
277 "D",
278 dimArea,
279 spread_/sqr(gAverage(mesh_.nonOrthDeltaCoeffs()))
280 );
281
282
283 const volScalarField& alpha1 = mixture_.alpha1();
284 const volScalarField& alpha2 = mixture_.alpha2();
285
286 const dimensionedScalar MDotMin("MdotMin", mDotc_.dimensions(), 1e-3);
287
288 if (max(mDotc_) > MDotMin)
289 {
291 (
292 mDotcSpread_,
293 mDotc_,
294 alpha1,
295 alpha2,
296 D,
297 1e-3
298 );
299 }
300
301 if (max(mDote_) > MDotMin)
302 {
304 (
305 mDoteSpread_,
306 mDote_,
307 alpha1,
308 alpha2,
309 D,
310 1e-3
311 );
312 }
313}
314
315
316void Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::
317updateInterface()
318{
319
320 // interface heat resistance
321 // Interpolating alpha1 cell centre values to mesh points (vertices)
322 scalarField ap
323 (
324 volPointInterpolation::New(mesh_).interpolate(mixture_.alpha1())
325 );
326
327 cutCellIso cutCell(mesh_, ap);
328
329 forAll(interfaceArea_, celli)
330 {
331 label status = cutCell.calcSubCell(celli, 0.5);
332 interfaceArea_[celli] = 0;
333 if (status == 0) // cell is cut
334 {
335 interfaceArea_[celli] =
336 mag(cutCell.faceArea())/mesh_.V()[celli];
337 }
338 }
339}
340
341
344vDot() const
345{
346
347 dimensionedScalar pCoeff(1.0/mixture_.rho1() - 1.0/mixture_.rho2());
348
349 return Pair<tmp<volScalarField>>
350 (
351 pCoeff*mDotcSpread_,
352 -pCoeff*mDoteSpread_
353 );
354}
355
356
358read()
359{
361 {
362 optionalSubDict(type() + "Coeffs").readEntry("R", R_);
363 optionalSubDict(type() + "Coeffs").readEntry("spread", spread_);
364
365 return true;
366 }
367
368 return false;
369}
370
371
372// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
const volScalarField & alpha2
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static const word dictName
Definition: basicThermo.H:256
const dimensionedScalar & rho1() const
Return const-access to phase1 density.
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
const thermoIncompressibleTwoPhaseMixture & mixture_
Reference to the thermoIncompressibleTwoPhaseMixture.
virtual bool read()
Read the transportProperties dictionary and update.
virtual Pair< tmp< volScalarField > > vDot() const
Return the volumetric condensation and vaporisation rates as.
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 interfaceHeatResistance phaseChange model.
virtual Pair< tmp< volScalarField > > vDotAlphal() const
Volumetric source for alpha (used by alphaEq)
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
const volScalarField & alpha2() const
Return the phase-fraction of phase 2.
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
thermo correct()
const volScalarField & T
dynamicFvMesh & mesh
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Calculate the finiteVolume matrix for implicit and explicit sources.
word timeName
Definition: getTimeIndex.H:3
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
Definition: fvcSmooth.C:325
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 dimPower
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
const dimensionSet dimDensity
const dimensionedScalar & D
scalar T0
Definition: createFields.H:22
volScalarField & e
Definition: createFields.H:11
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const vector L(dict.get< vector >("L"))