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(pair_.phase1().phi().boundaryField(), patchI)
40  {
41  if
42  (
43  isA<fixedValueFvsPatchScalarField>
44  (
45  pair_.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 phasePair::dictTable& modelTable,
62  const blendingMethod& blending,
63  const phasePair& pair,
64  const orderedPhasePair& pair1In2,
65  const orderedPhasePair& pair2In1,
66  const bool correctFixedFluxBCs
67 )
68 :
69  pair_(pair),
70  pair1In2_(pair1In2),
71  pair2In1_(pair2In1),
72  blending_(blending),
73  correctFixedFluxBCs_(correctFixedFluxBCs)
74 {
75  if (modelTable.found(pair_))
76  {
77  model_.set
78  (
80  (
81  modelTable[pair_],
82  pair_
83  ).ptr()
84  );
85  }
86 
87  if (modelTable.found(pair1In2_))
88  {
89  model1In2_.set
90  (
92  (
93  modelTable[pair1In2_],
94  pair1In2_
95  ).ptr()
96  );
97  }
98 
99  if (modelTable.found(pair2In1_))
100  {
101  model2In1_.set
102  (
104  (
105  modelTable[pair2In1_],
106  pair2In1_
107  ).ptr()
108  );
109  }
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
114 
115 template<class modelType>
117 {}
118 
119 
120 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
121 
122 template<class modelType>
125 {
126  tmp<volScalarField> f1, f2;
127 
128  if (model_.valid() || model1In2_.valid())
129  {
130  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
131  }
132 
133  if (model_.valid() || model2In1_.valid())
134  {
135  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
136  }
137 
138  tmp<volScalarField> x
139  (
140  new volScalarField
141  (
142  IOobject
143  (
144  modelType::typeName + ":K",
145  pair_.phase1().mesh().time().timeName(),
146  pair_.phase1().mesh(),
149  false
150  ),
151  pair_.phase1().mesh(),
152  dimensionedScalar("zero", modelType::dimK, 0)
153  )
154  );
155 
156  if (model_.valid())
157  {
158  x() += model_->K()*(f1() - f2());
159  }
160 
161  if (model1In2_.valid())
162  {
163  x() += model1In2_->K()*(1 - f1);
164  }
165 
166  if (model2In1_.valid())
167  {
168  x() += model2In1_->K()*f2;
169  }
170 
171  if
172  (
173  correctFixedFluxBCs_
174  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
175  )
176  {
177  correctFixedFluxBCs(x());
178  }
179 
180  return x;
181 }
182 
183 
184 template<class modelType>
187 {
188  tmp<surfaceScalarField> f1, f2;
189 
190  if (model_.valid() || model1In2_.valid())
191  {
193  (
194  blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
195  );
196  }
197 
198  if (model_.valid() || model2In1_.valid())
199  {
200  f2 = fvc::interpolate
201  (
202  blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
203  );
204  }
205 
206  tmp<surfaceScalarField> x
207  (
209  (
210  IOobject
211  (
212  modelType::typeName + ":Kf",
213  pair_.phase1().mesh().time().timeName(),
214  pair_.phase1().mesh(),
217  false
218  ),
219  pair_.phase1().mesh(),
220  dimensionedScalar("zero", modelType::dimK, 0)
221  )
222  );
223 
224  if (model_.valid())
225  {
226  x() += model_->Kf()*(f1() - f2());
227  }
228 
229  if (model1In2_.valid())
230  {
231  x() += model1In2_->Kf()*(1 - f1);
232  }
233 
234  if (model2In1_.valid())
235  {
236  x() += model2In1_->Kf()*f2;
237  }
238 
239  if
240  (
241  correctFixedFluxBCs_
242  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
243  )
244  {
245  correctFixedFluxBCs(x());
246  }
247 
248  return x;
249 }
250 
251 
252 template<class modelType>
253 template<class Type>
256 {
257  tmp<volScalarField> f1, f2;
258 
259  if (model_.valid() || model1In2_.valid())
260  {
261  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
262  }
263 
264  if (model_.valid() || model2In1_.valid())
265  {
266  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
267  }
268 
269  tmp<GeometricField<Type, fvPatchField, volMesh> > x
270  (
271  new GeometricField<Type, fvPatchField, volMesh>
272  (
273  IOobject
274  (
275  modelType::typeName + ":F",
276  pair_.phase1().mesh().time().timeName(),
277  pair_.phase1().mesh(),
280  false
281  ),
282  pair_.phase1().mesh(),
283  dimensioned<Type>("zero", modelType::dimF, pTraits<Type>::zero)
284  )
285  );
286 
287  if (model_.valid())
288  {
289  x() += model_->F()*(f1() - f2());
290  }
291 
292  if (model1In2_.valid())
293  {
294  x() += model1In2_->F()*(1 - f1);
295  }
296 
297  if (model2In1_.valid())
298  {
299  x() -= model2In1_->F()*f2; // note : subtraction
300  }
301 
302  if
303  (
304  correctFixedFluxBCs_
305  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
306  )
307  {
308  correctFixedFluxBCs(x());
309  }
310 
311  return x;
312 }
313 
314 
315 template<class modelType>
318 {
319  tmp<surfaceScalarField> f1, f2;
320 
321  if (model_.valid() || model1In2_.valid())
322  {
324  (
325  blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
326  );
327  }
328 
329  if (model_.valid() || model2In1_.valid())
330  {
331  f2 = fvc::interpolate
332  (
333  blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
334  );
335  }
336 
337  tmp<surfaceScalarField> x
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("zero", modelType::dimF*dimArea, 0)
352  )
353  );
354 
355  if (model_.valid())
356  {
357  x() += model_->Ff()*(f1() - f2());
358  }
359 
360  if (model1In2_.valid())
361  {
362  x() += model1In2_->Ff()*(1 - f1);
363  }
364 
365  if (model2In1_.valid())
366  {
367  x() -= model2In1_->Ff()*f2; // note : subtraction
368  }
369 
370  if
371  (
372  correctFixedFluxBCs_
373  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
374  )
375  {
376  correctFixedFluxBCs(x());
377  }
378 
379  return x;
380 }
381 
382 
383 template<class modelType>
386 {
387  tmp<volScalarField> f1, f2;
388 
389  if (model_.valid() || model1In2_.valid())
390  {
391  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
392  }
393 
394  if (model_.valid() || model2In1_.valid())
395  {
396  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
397  }
398 
399  tmp<volScalarField> x
400  (
401  new volScalarField
402  (
403  IOobject
404  (
405  modelType::typeName + ":D",
406  pair_.phase1().mesh().time().timeName(),
407  pair_.phase1().mesh(),
410  false
411  ),
412  pair_.phase1().mesh(),
413  dimensionedScalar("zero", modelType::dimD, 0)
414  )
415  );
416 
417  if (model_.valid())
418  {
419  x() += model_->D()*(f1() - f2());
420  }
421 
422  if (model1In2_.valid())
423  {
424  x() += model1In2_->D()*(1 - f1);
425  }
426 
427  if (model2In1_.valid())
428  {
429  x() += model2In1_->D()*f2;
430  }
431 
432  if
433  (
434  correctFixedFluxBCs_
435  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
436  )
437  {
438  correctFixedFluxBCs(x());
439  }
440 
441  return x;
442 }
443 
444 
445 template<class modelType>
447 (
448  const class phaseModel& phase
449 ) const
450 {
451  return
452  &phase == &(pair_.phase1())
453  ? model1In2_.valid()
454  : model2In1_.valid();
455 }
456 
457 
458 template<class modelType>
460 (
461  const class phaseModel& phase
462 ) const
463 {
464  return &phase == &(pair_.phase1()) ? model1In2_ : model2In1_;
465 }
466 
467 
468 // ************************************************************************* //
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
Foam::BlendedInterfacialModel::phaseModel
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
fixedValueFvsPatchFields.H
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::BlendedInterfacialModel::Kf
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::BlendedInterfacialModel::K
tmp< volScalarField > K() const
Return the blended force coefficient.
Foam::wallLubricationModel::pair_
const phasePair & pair_
Phase pair.
Definition: wallLubricationModel.H:62
Foam::BlendedInterfacialModel::correctFixedFluxBCs
void correctFixedFluxBCs(GeometricField &field) const
Correct coeff/value on fixed flux boundary conditions.