Loading...
Searching...
No Matches
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
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{
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
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
166Foam::Pair<Foam::tmp<Foam::volScalarField>>
168mDotAlphal() const
169{
170 return Pair<tmp<volScalarField>>
171 (
172 (mDotc_/clamp(mixture_.alpha2(), scalarMinMax(SMALL, 1))),
173 -(mDote_/clamp(mixture_.alpha1(), scalarMinMax(SMALL, 1)))
174 );
175}
176
177
178Foam::Pair<Foam::tmp<Foam::volScalarField>>
180mDot() const
181{
182 return Pair<tmp<volScalarField>>
183 (
184 tmp<volScalarField>(mDotc_),
185 tmp<volScalarField>(mDote_)
186 );
187}
188
189
190Foam::Pair<Foam::tmp<Foam::volScalarField>>
192mDotDeltaT() const
193{
194 const twoPhaseMixtureEThermo& thermo =
196 (
197 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
198 );
199
200 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
201
202 const dimensionedScalar& TSat = thermo.TSat();
203
204 Pair<tmp<volScalarField>> mDotce(mDot());
205
206 return Pair<tmp<volScalarField>>
207 (
208 mDotc_*pos(TSat - T.oldTime())/(TSat - T.oldTime()),
209 -mDote_*pos(T.oldTime() - TSat)/(T.oldTime() - TSat)
210 );
211}
212
213
214Foam::tmp<Foam::fvScalarMatrix>
216TSource() const
217{
218 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
219
221 auto& TSource = tTSource.ref();
222
223 const twoPhaseMixtureEThermo& thermo =
225 (
226 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
227 );
228
229 const dimensionedScalar& TSat = thermo.TSat();
230
231 // interface heat resistance
232 volScalarField IHRcoeff(interfaceArea_*R_);
233
234 TSource = fvm::Sp(IHRcoeff, T) - IHRcoeff*TSat;
235
236 return tTSource;
237}
238
239
241correct()
242{
243 // Update Interface
244 updateInterface();
245
246 // Update mDotc_ and mDote_
247 const volScalarField& T = mesh_.lookupObject<volScalarField>("T");
248
249 const twoPhaseMixtureEThermo& thermo =
251 (
252 mesh_.lookupObject<basicThermo>(basicThermo::dictName)
253 );
254
255 const dimensionedScalar& TSat = thermo.TSat();
257
258 const dimensionedScalar L(mag(mixture_.Hf2() - mixture_.Hf1()));
259
260 // interface heat resistance
261 mDotc_ = interfaceArea_*R_*max(TSat - T, T0)/L;
262 mDote_ = interfaceArea_*R_*max(T - TSat, T0)/L;
263
264 // Calculate the spread sources
266 (
267 "D",
268 dimArea,
269 spread_/sqr(gAverage(mesh_.nonOrthDeltaCoeffs()))
270 );
271
272
273 const volScalarField& alpha1 = mixture_.alpha1();
274 const volScalarField& alpha2 = mixture_.alpha2();
275
276 const dimensionedScalar MDotMin("MdotMin", mDotc_.dimensions(), 1e-3);
277
278 if (max(mDotc_) > MDotMin)
279 {
281 (
282 mDotcSpread_,
283 mDotc_,
284 alpha1,
285 alpha2,
286 D,
287 1e-3
288 );
289 }
290
291 if (max(mDote_) > MDotMin)
292 {
294 (
295 mDoteSpread_,
296 mDote_,
297 alpha1,
298 alpha2,
299 D,
300 1e-3
301 );
302 }
303}
304
305
306void Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::
307updateInterface()
308{
309
310 // interface heat resistance
311 // Interpolating alpha1 cell centre values to mesh points (vertices)
312 scalarField ap
313 (
314 volPointInterpolation::New(mesh_).interpolate(mixture_.alpha1())
315 );
316
317 cutCellIso cutCell(mesh_, ap);
318
319 forAll(interfaceArea_, celli)
320 {
321 label status = cutCell.calcSubCell(celli, 0.5);
322 interfaceArea_[celli] = 0;
323 if (status == 0) // cell is cut
324 {
325 interfaceArea_[celli] =
326 mag(cutCell.faceArea())/mesh_.V()[celli];
327 }
328 }
329}
330
331
332Foam::Pair<Foam::tmp<Foam::volScalarField>>
334vDot() const
335{
336
337 dimensionedScalar pCoeff(1.0/mixture_.rho1() - 1.0/mixture_.rho2());
338
339 return Pair<tmp<volScalarField>>
340 (
341 pCoeff*mDotcSpread_,
342 -pCoeff*mDoteSpread_
343 );
344}
345
346
348read()
349{
351 {
352 optionalSubDict(type() + "Coeffs").readEntry("R", R_);
353 optionalSubDict(type() + "Coeffs").readEntry("spread", spread_);
354
355 return true;
356 }
357
358 return false;
359}
360
361
362// ************************************************************************* //
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
static FOAM_NO_DANGLING_REFERENCE const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
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.
interfaceHeatResistance(const thermoIncompressibleTwoPhaseMixture &mixture, const fvMesh &mesh)
Construct from components.
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.
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
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)
List of values generated by applying the access operation to each list item.
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
Definition fvcSmooth.C:318
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.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
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 dimPower
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
const dimensionedScalar & D
psiReactionThermo & thermo
scalar T0
volScalarField & e
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const vector L(dict.get< vector >("L"))