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 Henning Scheufler
9 Copyright (C) 2020-2022 OpenCFD Ltd.
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
30#include "constants.H"
31#include "cutCellIso.H"
33#include "wallPolyPatch.H"
34#include "fvcSmooth.H"
35
36using namespace Foam::constant;
37
38
39// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
40
41template<class Thermo, class OtherThermo>
42void Foam::meltingEvaporationModels::
43interfaceHeatResistance<Thermo, OtherThermo>
44::updateInterface(const volScalarField& T)
45{
46 const fvMesh& mesh = this->mesh_;
47
48 const volScalarField& alpha = this->pair().from();
49
51 (
52 volPointInterpolation::New(mesh).interpolate(alpha)
53 );
54
56
57 forAll(interfaceArea_, celli)
58 {
59 label status = cutCell.calcSubCell(celli, isoAlpha_);
60 interfaceArea_[celli] = 0;
61 if (status == 0) // cell is cut
62 {
63 interfaceArea_[celli] =
64 mag(cutCell.faceArea())/mesh.V()[celli];
65 }
66 }
67
68 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
69
70 forAll(pbm, patchi)
71 {
72 if (isA<wallPolyPatch>(pbm[patchi]))
73 {
74 const polyPatch& pp = pbm[patchi];
75 forAll(pp.faceCells(), faceI)
76 {
77 const label pCelli = pp.faceCells()[faceI];
78 bool interface(false);
79 if
80 (
81 sign(R_.value()) > 0
82 && (T[pCelli] - Tactivate_.value()) > 0
83 )
84 {
85 interface = true;
86 }
87
88 if
89 (
90 sign(R_.value()) < 0
91 && (T[pCelli] - Tactivate_.value()) < 0
92 )
93 {
94 interface = true;
95 }
96
97 if
98 (
100 &&
101 (
102 alpha[pCelli] < 2*isoAlpha_
103 && alpha[pCelli] > 0.5*isoAlpha_
104 )
105 )
106 {
107 interfaceArea_[pCelli] =
108 mag(pp.faceAreas()[faceI])/mesh.V()[pCelli];
109 }
110 }
111 }
112 }
113}
114
115
116// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117
118template<class Thermo, class OtherThermo>
121(
122 const dictionary& dict,
123 const phasePair& pair
124)
125:
126 InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
128 Tactivate_("Tactivate", dimTemperature, dict),
129 interfaceArea_
130 (
132 (
133 "interfaceArea",
134 this->mesh_.time().timeName(),
135 this->mesh_,
136 IOobject::NO_READ,
137 IOobject::NO_WRITE
138 ),
139 this->mesh_,
141 ),
142 mDotc_
143 (
145 (
146 "mDotc",
147 this->mesh_.time().timeName(),
148 this->mesh_,
149 IOobject::NO_READ,
150 IOobject::AUTO_WRITE
151 ),
152 this->mesh_,
154 ),
155 mDotcSpread_
156 (
158 (
159 "mDotcSpread",
160 this->mesh_.time().timeName(),
161 this->mesh_,
162 IOobject::NO_READ,
163 IOobject::NO_WRITE
164 ),
165 this->mesh_,
167 ),
168 htc_
169 (
171 (
172 "htc",
173 this->mesh_.time().timeName(),
174 this->mesh_,
175 IOobject::NO_READ,
176 IOobject::NO_WRITE
177 ),
178 this->mesh_,
180 ),
181 isoAlpha_(dict.getOrDefault<scalar>("isoAlpha", 0.5)),
182 spread_(dict.getOrDefault<scalar>("spread", 3))
183{}
184
185
186// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
187
188template<class Thermo, class OtherThermo>
192{
193
194 const fvMesh& mesh = this->mesh_;
195
196 updateInterface(T);
197
198 auto tdeltaT = tmp<volScalarField>::New
199 (
201 (
202 "tdeltaT",
203 mesh.time().timeName(),
204 mesh
205 ),
206 mesh,
208 );
209 auto& deltaT = tdeltaT.ref();
210
212
213 if (sign(R_.value()) > 0)
214 {
215 deltaT = max(T - Tactivate_, T0);
216 }
217 else
218 {
219 deltaT = max(Tactivate_ - T, T0);
220 }
221
222 word fullSpeciesName = this->transferSpecie();
223 auto tempOpen = fullSpeciesName.find('.');
224 const word speciesName(fullSpeciesName.substr(0, tempOpen));
225
226 tmp<volScalarField> L = mag(this->L(speciesName, T));
227
228 htc_ = R_/L();
229
230 const volScalarField& to = this->pair().to();
231 const volScalarField& from = this->pair().from();
232
234 (
235 "D",
236 dimArea,
237 spread_/sqr(gAverage(this->mesh_.nonOrthDeltaCoeffs()))
238 );
239
240 const dimensionedScalar MdotMin("MdotMin", mDotc_.dimensions(), 1e-3);
241
242 if (max(mDotc_) > MdotMin)
243 {
245 (
246 mDotcSpread_,
247 mDotc_,
248 from,
249 to,
250 D,
251 1e-3
252 );
253 }
254
255 mDotc_ = interfaceArea_*htc_*deltaT;
256
257 return tmp<volScalarField>::New(mDotc_);
258}
259
260
261template<class Thermo, class OtherThermo>
265(
266 label variable,
267 const volScalarField& refValue
268)
269{
270 if (this->modelVariable_ == variable)
271 {
272 const volScalarField coeff(htc_*interfaceArea_);
273
274 if (sign(R_.value()) > 0)
275 {
276 return(coeff*pos(refValue - Tactivate_));
277 }
278 else
279 {
280 return(coeff*pos(Tactivate_ - refValue));
281 }
282 }
283
284 return nullptr;
285}
286
287
288template<class Thermo, class OtherThermo>
292(
293 label variable,
294 const volScalarField& refValue
295)
296{
297 if (this->modelVariable_ == variable)
298 {
299 const volScalarField coeff(htc_*interfaceArea_*Tactivate_);
300
301 if (sign(R_.value()) > 0)
302 {
303 return(-coeff*pos(refValue - Tactivate_));
304 }
305 else
306 {
307 return(coeff*pos(Tactivate_ - refValue));
308 }
309 }
310 else if (interfaceCompositionModel::P == variable)
311 {
312 return tmp<volScalarField>::New(mDotcSpread_);
313 }
314
315 return nullptr;
316}
317
318
319// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Base class for interface composition models, templated on the two thermodynamic models either side of...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Definition: cutCellIso.H:79
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition: cutCell.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:327
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & T
mesh interpolate(rAU)
dynamicFvMesh & mesh
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
word timeName
Definition: getTimeIndex.H:3
Different types of constants.
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
Definition: fvcSmooth.C:325
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.
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & alpha
dictionary dict
const dimensionedScalar & D
scalar T0
Definition: createFields.H:22
volScalarField & e
Definition: createFields.H:11
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const vector L(dict.get< vector >("L"))