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 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 (
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 }
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{}
185
186// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
187
188template<class Thermo, class OtherThermo>
191::Kexp(const volScalarField& T)
192{
193
194 const fvMesh& mesh = this->mesh_;
195
196 updateInterface(T);
197
198 auto tdeltaT = volScalarField::New
199 (
200 "tdeltaT",
202 mesh,
204 );
205 auto& deltaT = tdeltaT.ref();
206
208
209 if (sign(R_.value()) > 0)
210 {
211 deltaT = max(T - Tactivate_, T0);
212 }
213 else
214 {
215 deltaT = max(Tactivate_ - T, T0);
216 }
217
218 word fullSpeciesName = this->transferSpecie();
219 auto tempOpen = fullSpeciesName.find('.');
220 const word speciesName(fullSpeciesName.substr(0, tempOpen));
221
222 tmp<volScalarField> L = mag(this->L(speciesName, T));
223
224 htc_ = R_/L();
225
226 const volScalarField& to = this->pair().to();
227 const volScalarField& from = this->pair().from();
228
230 (
231 "D",
232 dimArea,
233 spread_/sqr(gAverage(this->mesh_.nonOrthDeltaCoeffs()))
234 );
235
236 const dimensionedScalar MdotMin("MdotMin", mDotc_.dimensions(), 1e-3);
237
238 if (max(mDotc_) > MdotMin)
239 {
241 (
242 mDotcSpread_,
243 mDotc_,
244 from,
245 to,
246 D,
247 1e-3
248 );
249 }
250
251 mDotc_ = interfaceArea_*htc_*deltaT;
252
254}
255
256
257template<class Thermo, class OtherThermo>
260::KSp
261(
262 label variable,
263 const volScalarField& refValue
264)
265{
266 if (this->modelVariable_ == variable)
267 {
268 const volScalarField coeff(htc_*interfaceArea_);
269
270 if (sign(R_.value()) > 0)
271 {
272 return(coeff*pos(refValue - Tactivate_));
273 }
274 else
275 {
276 return(coeff*pos(Tactivate_ - refValue));
277 }
278 }
279
280 return nullptr;
281}
282
283
284template<class Thermo, class OtherThermo>
287::KSu
288(
289 label variable,
290 const volScalarField& refValue
291)
292{
293 if (this->modelVariable_ == variable)
294 {
295 const volScalarField coeff(htc_*interfaceArea_*Tactivate_);
296
297 if (sign(R_.value()) > 0)
298 {
299 return(-coeff*pos(refValue - Tactivate_));
300 }
301 else
302 {
303 return(coeff*pos(Tactivate_ - refValue));
304 }
305 }
306 else if (interfaceCompositionModel::P == variable)
307 {
308 return tmp<volScalarField>::New(mDotcSpread_);
309 }
310
311 return nullptr;
312}
313
314
315// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Base class for interface composition models, templated on the two thermodynamic models either side of...
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo).
static FOAM_NO_DANGLING_REFERENCE const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Class for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with an isosurface defined ...
Definition cutCellIso.H:75
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
Definition cutCell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
const word transferSpecie() const
Return the transferring species name.
const word & variable() const
Returns the variable on which the model is based.
modelVariable modelVariable_
Enumeration for the model variable.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
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:318
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
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:40
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)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & alpha
dictionary dict
const dimensionedScalar & D
scalar T0
volScalarField & e
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const vector L(dict.get< vector >("L"))