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-2015 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  forAll(phase1_.phi()->boundaryField(), patchI)
40  {
41  if
42  (
43  isA<fixedValueFvsPatchScalarField>
44  (
45  phase1_.phi()->boundaryField()[patchI]
46  )
47  )
48  {
49  field.boundaryField()[patchI]
50  = pTraits<typename GeometricField::value_type>::zero;
51  }
52  }
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
58 template<class ModelType>
60 (
61  const phaseModel& phase1,
62  const phaseModel& phase2,
63  const blendingMethod& blending,
64  autoPtr<ModelType> model,
65  autoPtr<ModelType> model1In2,
66  autoPtr<ModelType> model2In1,
67  const bool correctFixedFluxBCs
68 )
69 :
70  phase1_(phase1),
71  phase2_(phase2),
72  blending_(blending),
73  model_(model),
74  model1In2_(model1In2),
75  model2In1_(model2In1),
76  correctFixedFluxBCs_(correctFixedFluxBCs)
77 {}
78 
79 
80 template<class ModelType>
82 (
83  const phasePair::dictTable& modelTable,
84  const blendingMethod& blending,
85  const phasePair& pair,
86  const orderedPhasePair& pair1In2,
87  const orderedPhasePair& pair2In1,
88  const bool correctFixedFluxBCs
89 )
90 :
91  phase1_(pair.phase1()),
92  phase2_(pair.phase2()),
93  blending_(blending),
94  correctFixedFluxBCs_(correctFixedFluxBCs)
95 {
96  if (modelTable.found(pair))
97  {
98  model_.set
99  (
101  (
102  modelTable[pair],
103  pair
104  ).ptr()
105  );
106  }
107 
108  if (modelTable.found(pair1In2))
109  {
110  model1In2_.set
111  (
113  (
114  modelTable[pair1In2],
115  pair1In2
116  ).ptr()
117  );
118  }
119 
120  if (modelTable.found(pair2In1))
121  {
122  model2In1_.set
123  (
125  (
126  modelTable[pair2In1],
127  pair2In1
128  ).ptr()
129  );
130  }
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
136 template<class ModelType>
138 {}
139 
140 
141 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
142 
143 template<class ModelType>
145 (
146  const class phaseModel& phase
147 ) const
148 {
149  return
150  &phase == &(phase1_)
151  ? model1In2_.valid()
152  : model2In1_.valid();
153 }
154 
155 
156 template<class ModelType>
158 (
159  const class phaseModel& phase
160 ) const
161 {
162  return &phase == &(phase1_) ? model1In2_ : model2In1_;
163 }
164 
165 
166 template<class ModelType>
169 {
170  tmp<volScalarField> f1, f2;
171 
172  if (model_.valid() || model1In2_.valid())
173  {
174  f1 = blending_.f1(phase1_, phase2_);
175  }
176 
177  if (model_.valid() || model2In1_.valid())
178  {
179  f2 = blending_.f2(phase1_, phase2_);
180  }
181 
182  tmp<volScalarField> x
183  (
184  new volScalarField
185  (
186  IOobject
187  (
188  ModelType::typeName + ":K",
189  phase1_.mesh().time().timeName(),
190  phase1_.mesh(),
193  false
194  ),
195  phase1_.mesh(),
196  dimensionedScalar("zero", ModelType::dimK, 0)
197  )
198  );
199 
200  if (model_.valid())
201  {
202  x() += model_->K()*(scalar(1) - f1() - f2());
203  }
204  if (model1In2_.valid())
205  {
206  x() += model1In2_->K()*f1;
207  }
208  if (model2In1_.valid())
209  {
210  x() += model2In1_->K()*f2;
211  }
212 
213  if
214  (
215  correctFixedFluxBCs_
216  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
217  )
218  {
219  correctFixedFluxBCs(x());
220  }
221 
222  return x;
223 }
224 
225 
226 template<class ModelType>
228 Foam::BlendedInterfacialModel<ModelType>::K(const scalar residualAlpha) const
229 {
230  tmp<volScalarField> f1, f2;
231 
232  if (model_.valid() || model1In2_.valid())
233  {
234  f1 = blending_.f1(phase1_, phase2_);
235  }
236 
237  if (model_.valid() || model2In1_.valid())
238  {
239  f2 = blending_.f2(phase1_, phase2_);
240  }
241 
242  tmp<volScalarField> x
243  (
244  new volScalarField
245  (
246  IOobject
247  (
248  ModelType::typeName + ":K",
249  phase1_.mesh().time().timeName(),
250  phase1_.mesh(),
253  false
254  ),
255  phase1_.mesh(),
256  dimensionedScalar("zero", ModelType::dimK, 0)
257  )
258  );
259 
260  if (model_.valid())
261  {
262  x() += model_->K(residualAlpha)*(scalar(1) - f1() - f2());
263  }
264  if (model1In2_.valid())
265  {
266  x() += model1In2_->K(residualAlpha)*f1;
267  }
268  if (model2In1_.valid())
269  {
270  x() += model2In1_->K(residualAlpha)*f2;
271  }
272 
273  if
274  (
275  correctFixedFluxBCs_
276  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
277  )
278  {
279  correctFixedFluxBCs(x());
280  }
281 
282  return x;
283 }
284 
285 
286 template<class ModelType>
289 {
290  tmp<surfaceScalarField> f1, f2;
291 
292  if (model_.valid() || model1In2_.valid())
293  {
295  (
296  blending_.f1(phase1_, phase2_)
297  );
298  }
299 
300  if (model_.valid() || model2In1_.valid())
301  {
302  f2 = fvc::interpolate
303  (
304  blending_.f2(phase1_, phase2_)
305  );
306  }
307 
308  tmp<surfaceScalarField> x
309  (
311  (
312  IOobject
313  (
314  ModelType::typeName + ":Kf",
315  phase1_.mesh().time().timeName(),
316  phase1_.mesh(),
319  false
320  ),
321  phase1_.mesh(),
322  dimensionedScalar("zero", ModelType::dimK, 0)
323  )
324  );
325 
326  if (model_.valid())
327  {
328  x() += model_->Kf()*(scalar(1) - f1() - f2());
329  }
330 
331  if (model1In2_.valid())
332  {
333  x() += model1In2_->Kf()*f1;
334  }
335 
336  if (model2In1_.valid())
337  {
338  x() += model2In1_->Kf()*f2;
339  }
340 
341  if
342  (
343  correctFixedFluxBCs_
344  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
345  )
346  {
347  correctFixedFluxBCs(x());
348  }
349 
350  return x;
351 }
352 
353 
354 template<class ModelType>
355 template<class Type>
358 {
359  tmp<volScalarField> f1, f2;
360 
361  if (model_.valid() || model1In2_.valid())
362  {
363  f1 = blending_.f1(phase1_, phase2_);
364  }
365 
366  if (model_.valid() || model2In1_.valid())
367  {
368  f2 = blending_.f2(phase1_, phase2_);
369  }
370 
371  tmp<GeometricField<Type, fvPatchField, volMesh> > x
372  (
373  new GeometricField<Type, fvPatchField, volMesh>
374  (
375  IOobject
376  (
377  ModelType::typeName + ":F",
378  phase1_.mesh().time().timeName(),
379  phase1_.mesh(),
382  false
383  ),
384  phase1_.mesh(),
385  dimensioned<Type>("zero", ModelType::dimF, pTraits<Type>::zero)
386  )
387  );
388 
389  if (model_.valid())
390  {
391  x() += model_->F()*(scalar(1) - f1() - f2());
392  }
393 
394  if (model1In2_.valid())
395  {
396  x() += model1In2_->F()*f1;
397  }
398 
399  if (model2In1_.valid())
400  {
401  x() -= model2In1_->F()*f2; // note : subtraction
402  }
403 
404  if
405  (
406  correctFixedFluxBCs_
407  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
408  )
409  {
410  correctFixedFluxBCs(x());
411  }
412 
413  return x;
414 }
415 
416 
417 template<class ModelType>
420 {
421  tmp<surfaceScalarField> f1, f2;
422 
423  if (model_.valid() || model1In2_.valid())
424  {
426  (
427  blending_.f1(phase1_, phase2_)
428  );
429  }
430 
431  if (model_.valid() || model2In1_.valid())
432  {
433  f2 = fvc::interpolate
434  (
435  blending_.f2(phase1_, phase2_)
436  );
437  }
438 
439  tmp<surfaceScalarField> x
440  (
442  (
443  IOobject
444  (
445  ModelType::typeName + ":Ff",
446  phase1_.mesh().time().timeName(),
447  phase1_.mesh(),
450  false
451  ),
452  phase1_.mesh(),
453  dimensionedScalar("zero", ModelType::dimF*dimArea, 0)
454  )
455  );
456 
457  if (model_.valid())
458  {
459  x() += model_->Ff()*(scalar(1) - f1() - f2());
460  }
461 
462  if (model1In2_.valid())
463  {
464  x() += model1In2_->Ff()*f1;
465  }
466 
467  if (model2In1_.valid())
468  {
469  x() -= model2In1_->Ff()*f2; // note : subtraction
470  }
471 
472  if
473  (
474  correctFixedFluxBCs_
475  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
476  )
477  {
478  correctFixedFluxBCs(x());
479  }
480 
481  return x;
482 }
483 
484 
485 template<class ModelType>
488 {
489  tmp<volScalarField> f1, f2;
490 
491  if (model_.valid() || model1In2_.valid())
492  {
493  f1 = blending_.f1(phase1_, phase2_);
494  }
495 
496  if (model_.valid() || model2In1_.valid())
497  {
498  f2 = blending_.f2(phase1_, phase2_);
499  }
500 
501  tmp<volScalarField> x
502  (
503  new volScalarField
504  (
505  IOobject
506  (
507  ModelType::typeName + ":D",
508  phase1_.mesh().time().timeName(),
509  phase1_.mesh(),
512  false
513  ),
514  phase1_.mesh(),
515  dimensionedScalar("zero", ModelType::dimD, 0)
516  )
517  );
518 
519  if (model_.valid())
520  {
521  x() += model_->D()*(scalar(1) - f1() - f2());
522  }
523  if (model1In2_.valid())
524  {
525  x() += model1In2_->D()*f1;
526  }
527  if (model2In1_.valid())
528  {
529  x() += model2In1_->D()*f2;
530  }
531 
532  if
533  (
534  correctFixedFluxBCs_
535  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
536  )
537  {
538  correctFixedFluxBCs(x());
539  }
540 
541  return x;
542 }
543 
544 
545 // ************************************************************************* //
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::BlendedInterfacialModel::BlendedInterfacialModel
BlendedInterfacialModel(const BlendedInterfacialModel< ModelType > &)
Disallow default bitwise copy construct.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::BlendedInterfacialModel::Ff
tmp< surfaceScalarField > Ff() const
Return the face blended force.
Foam::BlendedInterfacialModel::~BlendedInterfacialModel
~BlendedInterfacialModel()
Destructor.
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::BlendedInterfacialModel::F
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
Foam::BlendedInterfacialModel::hasModel
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
f1
scalar f1
Definition: createFields.H:28
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
fixedValueFvsPatchFields.H
Foam::BlendedInterfacialModel::model
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
Foam::BlendedInterfacialModel::D
tmp< volScalarField > D() const
Return the blended diffusivity.
surfaceInterpolate.H
Surface Interpolation.
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::phasePair::dictTable
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:59
Foam::BlendedInterfacialModel::Kf
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
x
x
Definition: LISASMDCalcMethod2.H:52
phase2
phaseModel & phase2
Definition: createFields.H:13
Foam::BlendedInterfacialModel::K
tmp< volScalarField > K() const
Return the blended force coefficient.
phase1
phaseModel & phase1
Definition: createFields.H:12
Foam::BlendedInterfacialModel::correctFixedFluxBCs
void correctFixedFluxBCs(GeometricField &field) const
Correct coeff/value on fixed flux boundary conditions.