BlendedInterfacialModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2018 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(pair_.phase1().phi().boundaryField(), patchi)
43  {
44  if
45  (
46  isA<fixedValueFvsPatchScalarField>
47  (
48  pair_.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 phasePair::dictTable& modelTable,
64  const blendingMethod& blending,
65  const phasePair& pair,
66  const orderedPhasePair& pair1In2,
67  const orderedPhasePair& pair2In1,
68  const bool correctFixedFluxBCs
69 )
70 :
71  pair_(pair),
72  pair1In2_(pair1In2),
73  pair2In1_(pair2In1),
74  blending_(blending),
75  correctFixedFluxBCs_(correctFixedFluxBCs)
76 {
77  if (modelTable.found(pair_))
78  {
79  model_.set
80  (
82  (
83  modelTable[pair_],
84  pair_
85  ).ptr()
86  );
87  }
88 
89  if (modelTable.found(pair1In2_))
90  {
91  model1In2_.set
92  (
94  (
95  modelTable[pair1In2_],
96  pair1In2_
97  ).ptr()
98  );
99  }
100 
101  if (modelTable.found(pair2In1_))
102  {
103  model2In1_.set
104  (
106  (
107  modelTable[pair2In1_],
108  pair2In1_
109  ).ptr()
110  );
111  }
112 }
113 
114 
115 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
116 
117 template<class modelType>
119 {}
120 
121 
122 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
123 
124 template<class modelType>
127 {
128  tmp<volScalarField> f1, f2;
129 
130  if (model_.valid() || model1In2_.valid())
131  {
132  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
133  }
134 
135  if (model_.valid() || model2In1_.valid())
136  {
137  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
138  }
139 
140  tmp<volScalarField> x
141  (
143  (
144  modelType::typeName + ":K",
145  pair_.phase1().mesh(),
146  dimensionedScalar(modelType::dimK, 0)
147  )
148  );
149 
150  if (model_.valid())
151  {
152  x.ref() += model_->K()*(f1() - f2());
153  }
154 
155  if (model1In2_.valid())
156  {
157  x.ref() += model1In2_->K()*(1 - f1);
158  }
159 
160  if (model2In1_.valid())
161  {
162  x.ref() += model2In1_->K()*f2;
163  }
164 
165  if
166  (
167  correctFixedFluxBCs_
168  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
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_.valid() || model1In2_.valid())
185  {
186  f1 = fvc::interpolate
187  (
188  blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
189  );
190  }
191 
192  if (model_.valid() || model2In1_.valid())
193  {
194  f2 = fvc::interpolate
195  (
196  blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
197  );
198  }
199 
200  tmp<surfaceScalarField> x
201  (
203  (
204  modelType::typeName + ":Kf",
205  pair_.phase1().mesh(),
206  dimensionedScalar(modelType::dimK, 0)
207  )
208  );
209 
210  if (model_.valid())
211  {
212  x.ref() += model_->Kf()*(f1() - f2());
213  }
214 
215  if (model1In2_.valid())
216  {
217  x.ref() += model1In2_->Kf()*(1 - f1);
218  }
219 
220  if (model2In1_.valid())
221  {
222  x.ref() += model2In1_->Kf()*f2;
223  }
224 
225  if
226  (
227  correctFixedFluxBCs_
228  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
229  )
230  {
231  correctFixedFluxBCs(x.ref());
232  }
233 
234  return x;
235 }
236 
237 
238 template<class modelType>
239 template<class Type>
242 {
243  tmp<volScalarField> f1, f2;
244 
245  if (model_.valid() || model1In2_.valid())
246  {
247  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
248  }
249 
250  if (model_.valid() || model2In1_.valid())
251  {
252  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
253  }
254 
255  tmp<GeometricField<Type, fvPatchField, volMesh>> x
256  (
257  new GeometricField<Type, fvPatchField, volMesh>
258  (
259  IOobject
260  (
261  modelType::typeName + ":F",
262  pair_.phase1().mesh().time().timeName(),
263  pair_.phase1().mesh(),
266  false
267  ),
268  pair_.phase1().mesh(),
269  dimensioned<Type>("zero", modelType::dimF, Zero)
270  )
271  );
272 
273  if (model_.valid())
274  {
275  x.ref() += model_->F()*(f1() - f2());
276  }
277 
278  if (model1In2_.valid())
279  {
280  x.ref() += model1In2_->F()*(1 - f1);
281  }
282 
283  if (model2In1_.valid())
284  {
285  x.ref() -= model2In1_->F()*f2; // note : subtraction
286  }
287 
288  if
289  (
290  correctFixedFluxBCs_
291  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
292  )
293  {
294  correctFixedFluxBCs(x.ref());
295  }
296 
297  return x;
298 }
299 
300 
301 template<class modelType>
304 {
305  tmp<surfaceScalarField> f1, f2;
306 
307  if (model_.valid() || model1In2_.valid())
308  {
309  f1 = fvc::interpolate
310  (
311  blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed())
312  );
313  }
314 
315  if (model_.valid() || model2In1_.valid())
316  {
317  f2 = fvc::interpolate
318  (
319  blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed())
320  );
321  }
322 
323  tmp<surfaceScalarField> x
324  (
326  (
327  modelType::typeName + ":Ff",
328  pair_.phase1().mesh(),
329  dimensionedScalar(modelType::dimF*dimArea, 0)
330  )
331  );
332 
333  if (model_.valid())
334  {
335  x.ref() += model_->Ff()*(f1() - f2());
336  }
337 
338  if (model1In2_.valid())
339  {
340  x.ref() += model1In2_->Ff()*(1 - f1);
341  }
342 
343  if (model2In1_.valid())
344  {
345  x.ref() -= model2In1_->Ff()*f2; // note : subtraction
346  }
347 
348  if
349  (
350  correctFixedFluxBCs_
351  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
352  )
353  {
354  correctFixedFluxBCs(x.ref());
355  }
356 
357  return x;
358 }
359 
360 
361 template<class modelType>
364 {
365  tmp<volScalarField> f1, f2;
366 
367  if (model_.valid() || model1In2_.valid())
368  {
369  f1 = blending_.f1(pair1In2_.dispersed(), pair2In1_.dispersed());
370  }
371 
372  if (model_.valid() || model2In1_.valid())
373  {
374  f2 = blending_.f2(pair1In2_.dispersed(), pair2In1_.dispersed());
375  }
376 
377  tmp<volScalarField> x
378  (
380  (
381  modelType::typeName + ":D",
382  pair_.phase1().mesh(),
383  dimensionedScalar(modelType::dimD, 0)
384  )
385  );
386 
387  if (model_.valid())
388  {
389  x.ref() += model_->D()*(f1() - f2());
390  }
391 
392  if (model1In2_.valid())
393  {
394  x.ref() += model1In2_->D()*(1 - f1);
395  }
396 
397  if (model2In1_.valid())
398  {
399  x.ref() += model2In1_->D()*f2;
400  }
401 
402  if
403  (
404  correctFixedFluxBCs_
405  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
406  )
407  {
408  correctFixedFluxBCs(x.ref());
409  }
410 
411  return x;
412 }
413 
414 
415 template<class modelType>
417 (
418  const class phaseModel& phase
419 ) const
420 {
421  return
422  &phase == &(pair_.phase1())
423  ? model1In2_.valid()
424  : model2In1_.valid();
425 }
426 
427 
428 template<class modelType>
430 (
431  const class phaseModel& phase
432 ) const
433 {
434  return &phase == &(pair_.phase1()) ? model1In2_ : model2In1_;
435 }
436 
437 
438 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
const modelType & phaseModel(const phaseModel &phase) const
Return the model for the supplied phase.
~BlendedInterfacialModel()
Destructor.
tmp< volScalarField > K() const
Return the blended force coefficient.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
scalar f1
Definition: createFields.H:28
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
BlendedInterfacialModel(const phaseModel &phase1, const phaseModel &phase2, const blendingMethod &blending, autoPtr< ModelType > model, autoPtr< ModelType > model1In2, autoPtr< ModelType > model2In1, const bool correctFixedFluxBCs=true)
Construct from two phases, blending method and three models.
static const zero Zero
Definition: zero.H:97
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.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< volScalarField > D() const
Return the blended diffusivity.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57