BlendedInterfacialModel.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) 2014-2016 OpenFOAM Foundation
9 Copyright (C) 2020 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
31#include "surfaceInterpolate.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35template<class modelType>
36template<class GeometricField>
38(
39 GeometricField& field
40) const
41{
42 typename GeometricField::Boundary& fieldBf =
43 field.boundaryFieldRef();
44
45 forAll(pair_.phase1().phi().boundaryField(), patchi)
46 {
47 if
48 (
49 isA<fixedValueFvsPatchScalarField>
50 (
51 pair_.phase1().phi().boundaryField()[patchi]
52 )
53 )
54 {
55 fieldBf[patchi] = Zero;
56 }
57 }
58}
59
60
61// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62
63template<class modelType>
65(
66 const phasePair::dictTable& modelTable,
67 const blendingMethod& blending,
68 const phasePair& pair,
69 const orderedPhasePair& pair1In2,
70 const orderedPhasePair& pair2In1,
71 const bool correctFixedFluxBCs
72)
73:
74 pair_(pair),
75 pair1In2_(pair1In2),
76 pair2In1_(pair2In1),
77 blending_(blending),
78 correctFixedFluxBCs_(correctFixedFluxBCs)
79{
80 if (modelTable.found(pair_))
81 {
82 model_.reset
83 (
84 modelType::New
85 (
86 modelTable[pair_],
87 pair_
88 ).ptr()
89 );
90 }
91
92 if (modelTable.found(pair1In2_))
93 {
94 model1In2_.reset
95 (
96 modelType::New
97 (
98 modelTable[pair1In2_],
99 pair1In2_
100 ).ptr()
101 );
102 }
103
104 if (modelTable.found(pair2In1_))
105 {
106 model2In1_.reset
107 (
108 modelType::New
109 (
110 modelTable[pair2In1_],
111 pair2In1_
112 ).ptr()
113 );
114 }
115}
116
117
118// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
119
120template<class modelType>
122{}
123
124
125// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
126
127template<class modelType>
130{
131 tmp<volScalarField> f1, f2;
132
133 if (model_ || model1In2_)
134 {
135 f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
136 }
137
138 if (model_ || model2In1_)
139 {
140 f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
141 }
142
144 (
146 (
148 (
149 modelType::typeName + ":K",
150 pair_.phase1().mesh().time().timeName(),
151 pair_.phase1().mesh(),
154 false
155 ),
156 pair_.phase1().mesh(),
157 dimensionedScalar(modelType::dimK, Zero)
158 )
159 );
160
161 if (model_)
162 {
163 x.ref() += model_->K()*(f1() - f2());
164 }
165
166 if (model1In2_)
167 {
168 x.ref() += model1In2_->K()*(1 - f1);
169 }
170
171 if (model2In1_)
172 {
173 x.ref() += model2In1_->K()*f2;
174 }
175
176 if
177 (
178 correctFixedFluxBCs_
179 && (model_ || model1In2_ || model2In1_)
180 )
181 {
182 correctFixedFluxBCs(x.ref());
183 }
184
185 return x;
186}
187
188
189template<class modelType>
192{
194
195 if (model_ || model1In2_)
196 {
198 (
199 blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
200 );
201 }
202
203 if (model_ || model2In1_)
204 {
206 (
207 blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
208 );
209 }
210
211 tmp<surfaceScalarField> x
212 (
214 (
215 IOobject
216 (
217 modelType::typeName + ":Kf",
218 pair_.phase1().mesh().time().timeName(),
219 pair_.phase1().mesh(),
222 false
223 ),
224 pair_.phase1().mesh(),
225 dimensionedScalar(modelType::dimK, Zero)
226 )
227 );
228
229 if (model_)
230 {
231 x.ref() += model_->Kf()*(f1() - f2());
232 }
233
234 if (model1In2_)
235 {
236 x.ref() += model1In2_->Kf()*(1 - f1);
237 }
238
239 if (model2In1_)
240 {
241 x.ref() += model2In1_->Kf()*f2;
242 }
243
244 if
245 (
246 correctFixedFluxBCs_
247 && (model_ || model1In2_ || model2In1_)
248 )
249 {
250 correctFixedFluxBCs(x.ref());
251 }
252
253 return x;
254}
255
256
257template<class modelType>
258template<class Type>
261{
262 tmp<volScalarField> f1, f2;
263
264 if (model_ || model1In2_)
265 {
266 f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
267 }
268
269 if (model_ || model2In1_)
270 {
271 f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
272 }
273
274 auto x = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
275 (
276 IOobject
277 (
278 modelType::typeName + ":F",
279 pair_.phase1().mesh().time().timeName(),
280 pair_.phase1().mesh(),
283 false
284 ),
285 pair_.phase1().mesh(),
286 dimensioned<Type>(modelType::dimF, Zero)
287 );
288
289 if (model_)
290 {
291 x.ref() += model_->F()*(f1() - f2());
292 }
293
294 if (model1In2_)
295 {
296 x.ref() += model1In2_->F()*(1 - f1);
297 }
298
299 if (model2In1_)
300 {
301 x.ref() -= model2In1_->F()*f2; // note : subtraction
302 }
303
304 if
305 (
306 correctFixedFluxBCs_
307 && (model_ || model1In2_ || model2In1_)
308 )
309 {
310 correctFixedFluxBCs(x.ref());
311 }
312
313 return x;
314}
315
316
317template<class modelType>
320{
321 tmp<surfaceScalarField> f1, f2;
322
323 if (model_ || model1In2_)
324 {
326 (
327 blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
328 );
329 }
330
331 if (model_ || model2In1_)
332 {
334 (
335 blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
336 );
337 }
338
340 (
341 IOobject
342 (
343 modelType::typeName + ":Ff",
344 pair_.phase1().mesh().time().timeName(),
345 pair_.phase1().mesh(),
348 false
349 ),
350 pair_.phase1().mesh(),
351 dimensionedScalar(modelType::dimF*dimArea, Zero)
352 );
353
354 x.ref().setOriented();
355
356 if (model_)
357 {
358 x.ref() += model_->Ff()*(f1() - f2());
359 }
360
361 if (model1In2_)
362 {
363 x.ref() += model1In2_->Ff()*(1 - f1);
364 }
365
366 if (model2In1_)
367 {
368 x.ref() -= model2In1_->Ff()*f2; // note : subtraction
369 }
370
371 if
372 (
373 correctFixedFluxBCs_
374 && (model_ || model1In2_ || model2In1_)
375 )
376 {
377 correctFixedFluxBCs(x.ref());
378 }
379
380 return x;
381}
382
383
384template<class modelType>
387{
388 tmp<volScalarField> f1, f2;
389
390 if (model_ || model1In2_)
391 {
392 f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
393 }
394
395 if (model_ || model2In1_)
396 {
397 f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
398 }
399
400 tmp<volScalarField> x
401 (
403 (
404 IOobject
405 (
406 modelType::typeName + ":D",
407 pair_.phase1().mesh().time().timeName(),
408 pair_.phase1().mesh(),
411 false
412 ),
413 pair_.phase1().mesh(),
414 dimensionedScalar(modelType::dimD, Zero)
415 )
416 );
417
418 if (model_)
419 {
420 x.ref() += model_->D()*(f1() - f2());
421 }
422
423 if (model1In2_)
424 {
425 x.ref() += model1In2_->D()*(1 - f1);
426 }
427
428 if (model2In1_)
429 {
430 x.ref() += model2In1_->D()*f2;
431 }
432
433 if
434 (
435 correctFixedFluxBCs_
436 && (model_ || model1In2_ || model2In1_)
437 )
438 {
439 correctFixedFluxBCs(x.ref());
440 }
441
442 return x;
443}
444
445
446template<class modelType>
448(
449 const class phaseModel& phase
450) const
451{
452 return
453 (
454 &phase == &(pair_.phase1())
455 ? bool(model1In2_)
456 : bool(model2In1_)
457 );
458}
459
460
461template<class modelType>
463(
464 const class phaseModel& phase
465) const
466{
467 return &phase == &(pair_.phase1()) ? *model1In2_ : *model2In1_;
468}
469
470
471// ************************************************************************* //
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< volScalarField > D() const
Return the blended diffusivity.
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
tmp< volScalarField > K() const
Return the blended force coefficient.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
rDeltaTY field()
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333