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-2023 OpenCFD Ltd.
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 
31 #include "surfaceInterpolate.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class modelType>
36 template<class GeometricField>
38 (
39  GeometricField& field
40 ) const
41 {
42  auto& fieldBf = field.boundaryFieldRef();
43 
44  forAll(pair_.phase1().phi().boundaryField(), patchi)
45  {
46  if
47  (
48  isA<fixedValueFvsPatchScalarField>
49  (
50  pair_.phase1().phi().boundaryField()[patchi]
51  )
52  )
53  {
54  fieldBf[patchi] = Zero;
55  }
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
62 template<class modelType>
64 (
65  const phasePair::dictTable& modelTable,
66  const blendingMethod& blending,
67  const phasePair& pair,
68  const orderedPhasePair& pair1In2,
69  const orderedPhasePair& pair2In1,
70  const bool correctFixedFluxBCs
71 )
72 :
73  pair_(pair),
74  pair1In2_(pair1In2),
75  pair2In1_(pair2In1),
76  blending_(blending),
77  correctFixedFluxBCs_(correctFixedFluxBCs)
78 {
79  if (modelTable.found(pair_))
80  {
81  model_.reset
82  (
84  (
85  modelTable[pair_],
86  pair_
87  ).ptr()
88  );
89  }
90 
91  if (modelTable.found(pair1In2_))
92  {
93  model1In2_.reset
94  (
96  (
97  modelTable[pair1In2_],
98  pair1In2_
99  ).ptr()
100  );
101  }
102 
103  if (modelTable.found(pair2In1_))
104  {
105  model2In1_.reset
106  (
108  (
109  modelTable[pair2In1_],
110  pair2In1_
111  ).ptr()
112  );
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
118 
119 template<class modelType>
121 {}
122 
123 
124 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125 
126 template<class modelType>
129 {
130  tmp<volScalarField> f1, f2;
131 
132  if (model_ || model1In2_)
133  {
134  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
135  }
136 
137  if (model_ || model2In1_)
138  {
139  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
140  }
141 
142  auto x = volScalarField::New
143  (
144  IOobject::scopedName(modelType::typeName, "K"),
145  IOobject::NO_REGISTER,
146  pair_.phase1().mesh(),
147  dimensionedScalar(modelType::dimK, Zero)
148  );
149 
150  if (model_)
151  {
152  x.ref() += model_->K()*(f1() - f2());
153  }
154 
155  if (model1In2_)
156  {
157  x.ref() += model1In2_->K()*(1 - f1);
158  }
159 
160  if (model2In1_)
161  {
162  x.ref() += model2In1_->K()*f2;
163  }
164 
165  if
166  (
167  correctFixedFluxBCs_
168  && (model_ || model1In2_ || model2In1_)
169  )
170  {
171  correctFixedFluxBCs(x.ref());
172  }
173 
174  return x;
175 }
176 
177 
178 template<class modelType>
181 {
182  tmp<surfaceScalarField> f1, f2;
183 
184  if (model_ || model1In2_)
185  {
186  f1 = fvc::interpolate
187  (
188  blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
189  );
190  }
191 
192  if (model_ || model2In1_)
193  {
194  f2 = fvc::interpolate
195  (
196  blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
197  );
198  }
199 
201  (
202  IOobject::scopedName(modelType::typeName, "Kf"),
203  pair_.phase1().mesh(),
204  dimensionedScalar(modelType::dimK, Zero)
205  );
206 
207  if (model_)
208  {
209  x.ref() += model_->Kf()*(f1() - f2());
210  }
211 
212  if (model1In2_)
213  {
214  x.ref() += model1In2_->Kf()*(1 - f1);
215  }
216 
217  if (model2In1_)
218  {
219  x.ref() += model2In1_->Kf()*f2;
220  }
221 
222  if
223  (
224  correctFixedFluxBCs_
225  && (model_ || model1In2_ || model2In1_)
226  )
227  {
228  correctFixedFluxBCs(x.ref());
229  }
230 
231  return x;
232 }
233 
234 
235 template<class modelType>
236 template<class Type>
239 {
240  tmp<volScalarField> f1, f2;
241 
242  if (model_ || model1In2_)
243  {
244  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
245  }
246 
247  if (model_ || model2In1_)
248  {
249  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
250  }
251 
253  (
254  IOobject::scopedName(modelType::typeName, "F"),
255  IOobject::NO_REGISTER,
256  pair_.phase1().mesh(),
257  dimensioned<Type>(modelType::dimF, Zero)
258  );
259 
260  if (model_)
261  {
262  x.ref() += model_->F()*(f1() - f2());
263  }
264 
265  if (model1In2_)
266  {
267  x.ref() += model1In2_->F()*(1 - f1);
268  }
269 
270  if (model2In1_)
271  {
272  x.ref() -= model2In1_->F()*f2; // note : subtraction
273  }
274 
275  if
276  (
277  correctFixedFluxBCs_
278  && (model_ || model1In2_ || model2In1_)
279  )
280  {
281  correctFixedFluxBCs(x.ref());
282  }
283 
284  return x;
285 }
286 
287 
288 template<class modelType>
291 {
292  tmp<surfaceScalarField> f1, f2;
293 
294  if (model_ || model1In2_)
295  {
296  f1 = fvc::interpolate
297  (
298  blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
299  );
300  }
301 
302  if (model_ || model2In1_)
303  {
304  f2 = fvc::interpolate
305  (
306  blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
307  );
308  }
309 
311  (
312  IOobject::scopedName(modelType::typeName, "Ff"),
313  IOobject::NO_REGISTER,
314  pair_.phase1().mesh(),
315  dimensionedScalar(modelType::dimF*dimArea, Zero)
316  );
317 
318  x.ref().setOriented();
319 
320  if (model_)
321  {
322  x.ref() += model_->Ff()*(f1() - f2());
323  }
324 
325  if (model1In2_)
326  {
327  x.ref() += model1In2_->Ff()*(1 - f1);
328  }
329 
330  if (model2In1_)
331  {
332  x.ref() -= model2In1_->Ff()*f2; // note : subtraction
333  }
334 
335  if
336  (
337  correctFixedFluxBCs_
338  && (model_ || model1In2_ || model2In1_)
339  )
340  {
341  correctFixedFluxBCs(x.ref());
342  }
343 
344  return x;
345 }
346 
347 
348 template<class modelType>
351 {
352  tmp<volScalarField> f1, f2;
353 
354  if (model_ || model1In2_)
355  {
356  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
357  }
358 
359  if (model_ || model2In1_)
360  {
361  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
362  }
363 
364  auto x = volScalarField::New
365  (
366  IOobject::scopedName(modelType::typeName, "D"),
367  IOobject::NO_REGISTER,
368  pair_.phase1().mesh(),
369  dimensionedScalar(modelType::dimD, Zero)
370  );
371 
372  if (model_)
373  {
374  x.ref() += model_->D()*(f1() - f2());
375  }
376 
377  if (model1In2_)
378  {
379  x.ref() += model1In2_->D()*(1 - f1);
380  }
381 
382  if (model2In1_)
383  {
384  x.ref() += model2In1_->D()*f2;
385  }
386 
387  if
388  (
389  correctFixedFluxBCs_
390  && (model_ || model1In2_ || model2In1_)
391  )
392  {
393  correctFixedFluxBCs(x.ref());
394  }
395 
396  return x;
397 }
398 
399 
400 template<class modelType>
402 (
403  const class phaseModel& phase
404 ) const
405 {
406  return
407  (
408  &phase == &(pair_.phase1())
409  ? bool(model1In2_)
410  : bool(model2In1_)
411  );
412 }
413 
414 
415 template<class modelType>
417 (
418  const class phaseModel& phase
419 ) const
420 {
421  return &phase == &(pair_.phase1()) ? *model1In2_ : *model2In1_;
422 }
423 
424 
425 // ************************************************************************* //
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
rDeltaTY field()
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
tmp< volScalarField > K() const
Return the blended force coefficient.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< volScalarField > D() const
Return the blended diffusivity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127