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-2021 OpenCFD Ltd.
9  Copyright (C) 2020 Henning Scheufler
10 -------------------------------------------------------------------------------
11 License
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"
31 #include "twoPhaseMixtureEThermo.H"
32 #include "cutCellIso.H"
33 #include "volPointInterpolation.H"
35 #include "wallPolyPatch.H"
36 #include "fvcSmooth.H"
37 #include "fvmSup.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace temperaturePhaseChangeTwoPhaseMixtures
44 {
47  (
48  temperaturePhaseChangeTwoPhaseMixture,
50  components
51  );
52 }
53 }
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
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 
154 vDotAlphal() 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 
168 mDotAlphal() 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 
190 mDot() const
191 {
192  return Pair<tmp<volScalarField>>
193  (
194  tmp<volScalarField>(mDotc_),
195  tmp<volScalarField>(mDote_)
196  );
197 }
198 
199 
202 mDotDeltaT() 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 
226 TSource() const
227 {
228  const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
229 
230  auto tTSource = tmp<fvScalarMatrix>::New(T, dimEnergy/dimTime);
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 
251 correct()
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  dimensionedScalar L = 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  // Limiting max condensation
275  forAll(mDotc_, celli)
276  {
277  scalar rhobyDt = mixture_.rho1().value()/mesh_.time().deltaTValue();
278  scalar maxEvap = mixture_.alpha1()[celli]*rhobyDt; // positive
279  scalar maxCond = -mixture_.alpha2()[celli]*rhobyDt; // negative
280  mDotc_[celli] = min(max(mDotc_[celli], maxCond), maxEvap);
281  }
282 
283  // Calculate the spread sources
284 
286  (
287  "D",
288  dimArea,
289  spread_/sqr(gAverage(mesh_.nonOrthDeltaCoeffs()))
290  );
291 
292 
293  const volScalarField& alpha1 = mixture_.alpha1();
294  const volScalarField& alpha2 = mixture_.alpha2();
295 
296  const dimensionedScalar MDotMin("MdotMin", mDotc_.dimensions(), 1e-3);
297 
298  if (max(mDotc_) > MDotMin)
299  {
301  (
302  mDotcSpread_,
303  mDotc_,
304  alpha1,
305  alpha2,
306  D,
307  1e-3
308  );
309  }
310 
311  if (max(mDote_) > MDotMin)
312  {
314  (
315  mDoteSpread_,
316  mDote_,
317  alpha1,
318  alpha2,
319  D,
320  1e-3
321  );
322  }
323 }
324 
325 
326 void Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::
327 updateInterface()
328 {
329 
330  // interface heat resistance
331  // Interpolating alpha1 cell centre values to mesh points (vertices)
332  scalarField ap
333  (
334  volPointInterpolation::New(mesh_).interpolate(mixture_.alpha1())
335  );
336 
337  cutCellIso cutCell(mesh_, ap);
338 
339  forAll(interfaceArea_, celli)
340  {
341  label status = cutCell.calcSubCell(celli, 0.5);
342  interfaceArea_[celli] = 0;
343  if (status == 0) // cell is cut
344  {
345  interfaceArea_[celli] =
346  mag(cutCell.faceArea())/mesh_.V()[celli];
347  }
348  }
349 }
350 
351 
354 vDot() const
355 {
356 
357  dimensionedScalar pCoeff(1.0/mixture_.rho1() - 1.0/mixture_.rho2());
358 
359  return Pair<tmp<volScalarField>>
360  (
361  pCoeff*mDotcSpread_,
362  -pCoeff*mDoteSpread_
363  );
364 }
365 
366 
368 read()
369 {
371  {
372  optionalSubDict(type() + "Coeffs").readEntry("R", R_);
373  optionalSubDict(type() + "Coeffs").readEntry("spread", spread_);
374 
375  return true;
376  }
377 
378  return false;
379 }
380 
381 
382 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::vDot
virtual Pair< tmp< volScalarField > > vDot() const
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
L
const vector L(dict.get< vector >("L"))
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::interfaceHeatResistance
interfaceHeatResistance(const thermoIncompressibleTwoPhaseMixture &mixture, const fvMesh &mesh)
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::mDot
virtual Pair< tmp< volScalarField > > mDot() const
twoPhaseMixtureEThermo.H
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:597
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
Foam::basicThermo::dictName
static const word dictName
Definition: basicThermo.H:252
calculatedFvPatchFields.H
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Definition: MeshObject.C:41
wallPolyPatch.H
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::vDotAlphal
virtual Pair< tmp< volScalarField > > vDotAlphal() const
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::TSource
virtual tmp< fvScalarMatrix > TSource() const
Foam::fvc::spreadSource
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
Definition: fvcSmooth.C:318
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
cutCellIso.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
interfaceHeatResistance
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::read
virtual bool read()
correct
fvOptions correct(rho)
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::dimPower
const dimensionSet dimPower
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::mDotAlphal
virtual Pair< tmp< volScalarField > > mDotAlphal() const
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
fvmSup.H
Calculate the matrix for implicit and explicit sources.
volPointInterpolation.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:50
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::mDotDeltaT
virtual Pair< tmp< volScalarField > > mDotDeltaT() const
Foam::tmp::New
static tmp< T > New(Args &&... args)
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::roots::type
type
Definition: Roots.H:52
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
T0
scalar T0
Definition: createFields.H:22
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::correct
virtual void correct()
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Foam::dimless
const dimensionSet dimless
Definition: dimensionSets.C:182
Foam::temperaturePhaseChangeTwoPhaseMixture::read
virtual bool read()
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:170