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 | Copyright (C) 2014-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "BlendedInterfacialModel.H"
28 #include "surfaceInterpolate.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class ModelType>
33 template<class GeometricField>
35 (
36  GeometricField& field
37 ) const
38 {
39  typename GeometricField::Boundary& fieldBf =
40  field.boundaryFieldRef();
41 
42  forAll(phase1_.phi()().boundaryField(), patchi)
43  {
44  if
45  (
46  isA<fixedValueFvsPatchScalarField>
47  (
48  phase1_.phi()().boundaryField()[patchi]
49  )
50  )
51  {
52  fieldBf[patchi] = Zero;
53  }
54  }
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 template<class ModelType>
62 (
63  const phaseModel& phase1,
64  const phaseModel& phase2,
65  const blendingMethod& blending,
66  autoPtr<ModelType> model,
67  autoPtr<ModelType> model1In2,
68  autoPtr<ModelType> model2In1,
69  const bool correctFixedFluxBCs
70 )
71 :
72  phase1_(phase1),
73  phase2_(phase2),
74  blending_(blending),
75  model_(model),
76  model1In2_(model1In2),
77  model2In1_(model2In1),
78  correctFixedFluxBCs_(correctFixedFluxBCs)
79 {}
80 
81 
82 template<class ModelType>
84 (
85  const phasePair::dictTable& modelTable,
86  const blendingMethod& blending,
87  const phasePair& pair,
88  const orderedPhasePair& pair1In2,
89  const orderedPhasePair& pair2In1,
90  const bool correctFixedFluxBCs
91 )
92 :
93  phase1_(pair.phase1()),
94  phase2_(pair.phase2()),
95  blending_(blending),
96  correctFixedFluxBCs_(correctFixedFluxBCs)
97 {
98  if (modelTable.found(pair))
99  {
100  model_.set
101  (
103  (
104  modelTable[pair],
105  pair
106  ).ptr()
107  );
108  }
109 
110  if (modelTable.found(pair1In2))
111  {
112  model1In2_.set
113  (
115  (
116  modelTable[pair1In2],
117  pair1In2
118  ).ptr()
119  );
120  }
121 
122  if (modelTable.found(pair2In1))
123  {
124  model2In1_.set
125  (
127  (
128  modelTable[pair2In1],
129  pair2In1
130  ).ptr()
131  );
132  }
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
138 template<class ModelType>
140 {}
141 
142 
143 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
144 
145 template<class ModelType>
147 (
148  const class phaseModel& phase
149 ) const
150 {
151  return
152  &phase == &(phase1_)
153  ? model1In2_.valid()
154  : model2In1_.valid();
155 }
156 
157 
158 template<class ModelType>
160 (
161  const class phaseModel& phase
162 ) const
163 {
164  return &phase == &(phase1_) ? model1In2_ : model2In1_;
165 }
166 
167 
168 template<class ModelType>
171 {
172  tmp<volScalarField> f1, f2;
173 
174  if (model_.valid() || model1In2_.valid())
175  {
176  f1 = blending_.f1(phase1_, phase2_);
177  }
178 
179  if (model_.valid() || model2In1_.valid())
180  {
181  f2 = blending_.f2(phase1_, phase2_);
182  }
183 
184  tmp<volScalarField> x
185  (
186  new volScalarField
187  (
188  IOobject
189  (
190  ModelType::typeName + ":K",
191  phase1_.mesh().time().timeName(),
192  phase1_.mesh(),
195  false
196  ),
197  phase1_.mesh(),
198  dimensionedScalar("zero", ModelType::dimK, 0)
199  )
200  );
201 
202  if (model_.valid())
203  {
204  x.ref() += model_->K()*(scalar(1) - f1() - f2());
205  }
206  if (model1In2_.valid())
207  {
208  x.ref() += model1In2_->K()*f1;
209  }
210  if (model2In1_.valid())
211  {
212  x.ref() += model2In1_->K()*f2;
213  }
214 
215  if
216  (
217  correctFixedFluxBCs_
218  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
219  )
220  {
221  correctFixedFluxBCs(x.ref());
222  }
223 
224  return x;
225 }
226 
227 
228 template<class ModelType>
230 Foam::BlendedInterfacialModel<ModelType>::K(const scalar residualAlpha) const
231 {
232  tmp<volScalarField> f1, f2;
233 
234  if (model_.valid() || model1In2_.valid())
235  {
236  f1 = blending_.f1(phase1_, phase2_);
237  }
238 
239  if (model_.valid() || model2In1_.valid())
240  {
241  f2 = blending_.f2(phase1_, phase2_);
242  }
243 
244  tmp<volScalarField> x
245  (
246  new volScalarField
247  (
248  IOobject
249  (
250  ModelType::typeName + ":K",
251  phase1_.mesh().time().timeName(),
252  phase1_.mesh(),
255  false
256  ),
257  phase1_.mesh(),
258  dimensionedScalar("zero", ModelType::dimK, 0)
259  )
260  );
261 
262  if (model_.valid())
263  {
264  x.ref() += model_->K(residualAlpha)*(scalar(1) - f1() - f2());
265  }
266  if (model1In2_.valid())
267  {
268  x.ref() += model1In2_->K(residualAlpha)*f1;
269  }
270  if (model2In1_.valid())
271  {
272  x.ref() += model2In1_->K(residualAlpha)*f2;
273  }
274 
275  if
276  (
277  correctFixedFluxBCs_
278  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
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_.valid() || model1In2_.valid())
295  {
296  f1 = fvc::interpolate
297  (
298  blending_.f1(phase1_, phase2_)
299  );
300  }
301 
302  if (model_.valid() || model2In1_.valid())
303  {
304  f2 = fvc::interpolate
305  (
306  blending_.f2(phase1_, phase2_)
307  );
308  }
309 
310  tmp<surfaceScalarField> x
311  (
313  (
314  IOobject
315  (
316  ModelType::typeName + ":Kf",
317  phase1_.mesh().time().timeName(),
318  phase1_.mesh(),
321  false
322  ),
323  phase1_.mesh(),
324  dimensionedScalar("zero", ModelType::dimK, 0)
325  )
326  );
327 
328  if (model_.valid())
329  {
330  x.ref() += model_->Kf()*(scalar(1) - f1() - f2());
331  }
332 
333  if (model1In2_.valid())
334  {
335  x.ref() += model1In2_->Kf()*f1;
336  }
337 
338  if (model2In1_.valid())
339  {
340  x.ref() += model2In1_->Kf()*f2;
341  }
342 
343  if
344  (
345  correctFixedFluxBCs_
346  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
347  )
348  {
349  correctFixedFluxBCs(x.ref());
350  }
351 
352  return x;
353 }
354 
355 
356 template<class ModelType>
357 template<class Type>
360 {
361  tmp<volScalarField> f1, f2;
362 
363  if (model_.valid() || model1In2_.valid())
364  {
365  f1 = blending_.f1(phase1_, phase2_);
366  }
367 
368  if (model_.valid() || model2In1_.valid())
369  {
370  f2 = blending_.f2(phase1_, phase2_);
371  }
372 
373  tmp<GeometricField<Type, fvPatchField, volMesh>> x
374  (
375  new GeometricField<Type, fvPatchField, volMesh>
376  (
377  IOobject
378  (
379  ModelType::typeName + ":F",
380  phase1_.mesh().time().timeName(),
381  phase1_.mesh(),
384  false
385  ),
386  phase1_.mesh(),
387  dimensioned<Type>("zero", ModelType::dimF, Zero)
388  )
389  );
390 
391  if (model_.valid())
392  {
393  x.ref() += model_->F()*(scalar(1) - f1() - f2());
394  }
395 
396  if (model1In2_.valid())
397  {
398  x.ref() += model1In2_->F()*f1;
399  }
400 
401  if (model2In1_.valid())
402  {
403  x.ref() -= model2In1_->F()*f2; // note : subtraction
404  }
405 
406  if
407  (
408  correctFixedFluxBCs_
409  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
410  )
411  {
412  correctFixedFluxBCs(x.ref());
413  }
414 
415  return x;
416 }
417 
418 
419 template<class ModelType>
422 {
423  tmp<surfaceScalarField> f1, f2;
424 
425  if (model_.valid() || model1In2_.valid())
426  {
427  f1 = fvc::interpolate
428  (
429  blending_.f1(phase1_, phase2_)
430  );
431  }
432 
433  if (model_.valid() || model2In1_.valid())
434  {
435  f2 = fvc::interpolate
436  (
437  blending_.f2(phase1_, phase2_)
438  );
439  }
440 
441  tmp<surfaceScalarField> x
442  (
444  (
445  IOobject
446  (
447  ModelType::typeName + ":Ff",
448  phase1_.mesh().time().timeName(),
449  phase1_.mesh(),
452  false
453  ),
454  phase1_.mesh(),
455  dimensionedScalar("zero", ModelType::dimF*dimArea, 0)
456  )
457  );
458 
459  if (model_.valid())
460  {
461  x.ref() += model_->Ff()*(scalar(1) - f1() - f2());
462  }
463 
464  if (model1In2_.valid())
465  {
466  x.ref() += model1In2_->Ff()*f1;
467  }
468 
469  if (model2In1_.valid())
470  {
471  x.ref() -= model2In1_->Ff()*f2; // note : subtraction
472  }
473 
474  if
475  (
476  correctFixedFluxBCs_
477  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
478  )
479  {
480  correctFixedFluxBCs(x.ref());
481  }
482 
483  return x;
484 }
485 
486 
487 template<class ModelType>
490 {
491  tmp<volScalarField> f1, f2;
492 
493  if (model_.valid() || model1In2_.valid())
494  {
495  f1 = blending_.f1(phase1_, phase2_);
496  }
497 
498  if (model_.valid() || model2In1_.valid())
499  {
500  f2 = blending_.f2(phase1_, phase2_);
501  }
502 
503  tmp<volScalarField> x
504  (
505  new volScalarField
506  (
507  IOobject
508  (
509  ModelType::typeName + ":D",
510  phase1_.mesh().time().timeName(),
511  phase1_.mesh(),
514  false
515  ),
516  phase1_.mesh(),
517  dimensionedScalar("zero", ModelType::dimD, 0)
518  )
519  );
520 
521  if (model_.valid())
522  {
523  x.ref() += model_->D()*(scalar(1) - f1() - f2());
524  }
525  if (model1In2_.valid())
526  {
527  x.ref() += model1In2_->D()*f1;
528  }
529  if (model2In1_.valid())
530  {
531  x.ref() += model2In1_->D()*f2;
532  }
533 
534  if
535  (
536  correctFixedFluxBCs_
537  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
538  )
539  {
540  correctFixedFluxBCs(x.ref());
541  }
542 
543  return x;
544 }
545 
546 
547 // ************************************************************************* //
tmp< volScalarField > D() const
Return the blended diffusivity.
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
~BlendedInterfacialModel()
Destructor.
scalar f1
Definition: createFields.H:28
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static const zero Zero
Definition: zero.H:91
phaseModel & phase1
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:54
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< volScalarField > K() const
Return the blended force coefficient.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.