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-2020 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 
28 #include "surfaceInterpolate.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace blendedInterfacialModel
35 {
36 
37 template<class GeoField>
38 inline tmp<GeoField> interpolate(tmp<volScalarField> f);
39 
40 template<>
41 inline tmp<Foam::volScalarField> interpolate(tmp<volScalarField> f)
42 {
43  return f;
44 }
45 
46 template<>
47 inline tmp<Foam::surfaceScalarField> interpolate(tmp<volScalarField> f)
48 {
49  return fvc::interpolate(f);
50 }
51 
52 } // End namespace blendedInterfacialModel
53 } // End namespace Foam
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 template<class ModelType>
58 template<template<class> class PatchField, class GeoMesh>
60 (
61  tmp<GeometricField<scalar, PatchField, GeoMesh>>& f1,
62  tmp<GeometricField<scalar, PatchField, GeoMesh>>& f2,
63  const bool subtract
64 ) const
65 {
66  typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
67 
68  if (model_.valid() && subtract)
69  {
71  << "Cannot treat an interfacial model with no distinction between "
72  << "continuous and dispersed phases as signed"
73  << exit(FatalError);
74  }
75 
76  if (model_.valid() || model1In2_.valid())
77  {
78  f1 =
79  blendedInterfacialModel::interpolate<scalarGeoField>
80  (
81  blending_.f1(phase1_, phase2_)
82  );
83  }
84 
85  if (model_.valid() || model2In1_.valid())
86  {
87  f2 =
88  (subtract ? -1 : +1)
89  *blendedInterfacialModel::interpolate<scalarGeoField>
90  (
91  blending_.f2(phase1_, phase2_)
92  );
93  }
94 }
95 
96 
97 template<class ModelType>
98 template<class Type, template<class> class PatchField, class GeoMesh>
100 (
101  GeometricField<Type, PatchField, GeoMesh>& field
102 ) const
103 {
104  typedef GeometricField<Type, PatchField, GeoMesh> typeGeoField;
105 
106  typename typeGeoField::Boundary& fieldBf = field.boundaryFieldRef();
107 
108  forAll(fieldBf, patchi)
109  {
110  if
111  (
112  (
113  !phase1_.stationary()
114  && isA<fixedValueFvsPatchScalarField>
115  (
116  phase1_.phi()().boundaryField()[patchi]
117  )
118  )
119  || (
120  !phase2_.stationary()
121  && isA<fixedValueFvsPatchScalarField>
122  (
123  phase2_.phi()().boundaryField()[patchi]
124  )
125  )
126  )
127  {
128  fieldBf[patchi] = Zero;
129  }
130  }
131 }
132 
133 
134 template<class ModelType>
135 template
136 <
137  class Type,
138  template<class> class PatchField,
139  class GeoMesh,
140  class ... Args
141 >
144 (
145  tmp<GeometricField<Type, PatchField, GeoMesh>>
146  (ModelType::*method)(Args ...) const,
147  const word& name,
148  const dimensionSet& dims,
149  const bool subtract,
150  Args ... args
151 ) const
152 {
153  typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
154  typedef GeometricField<Type, PatchField, GeoMesh> typeGeoField;
155 
156  tmp<scalarGeoField> f1, f2;
157  calculateBlendingCoeffs(f1, f2, subtract);
158 
159  tmp<typeGeoField> x =
161  (
162  ModelType::typeName + ":"
163  + IOobject::groupName(name, phasePair(phase1_, phase2_).name()),
164  phase1_.mesh(),
165  dimensioned<Type>(dims, Zero)
166  );
167 
168  if (model_.valid())
169  {
170  x.ref() += (scalar(1) - f1() - f2())*(model_().*method)(args ...);
171  }
172 
173  if (model1In2_.valid())
174  {
175  x.ref() += f1*(model1In2_().*method)(args ...);
176  }
177 
178  if (model2In1_.valid())
179  {
180  x.ref() += f2*(model2In1_().*method)(args ...);
181  }
182 
183  if
184  (
185  correctFixedFluxBCs_
186  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
187  )
188  {
189  correctFixedFluxBCs(x.ref());
190  }
191 
192  return x;
193 }
194 
195 
196 template<class ModelType>
197 template
198 <
199  class Type,
200  template<class> class PatchField,
201  class GeoMesh,
202  class ... Args
203 >
206 (
207  HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>
208  (ModelType::*method)(Args ...) const,
209  const word& name,
210  const dimensionSet& dims,
211  const bool subtract,
212  Args ... args
213 ) const
214 {
215  typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
216  typedef GeometricField<Type, PatchField, GeoMesh> typeGeoField;
217 
218  tmp<scalarGeoField> f1, f2;
219  calculateBlendingCoeffs(f1, f2, subtract);
220 
221  HashPtrTable<typeGeoField> xs;
222 
223  auto addToXs = [&]
224  (
225  const scalarGeoField& f,
226  const HashPtrTable<typeGeoField>& dxs
227  )
228  {
229  forAllConstIter(typename HashPtrTable<typeGeoField>, dxs, dxIter)
230  {
231  if (xs.found(dxIter.key()))
232  {
233  *xs[dxIter.key()] += f**dxIter();
234  }
235  else
236  {
237  xs.insert
238  (
239  dxIter.key(),
241  (
242  ModelType::typeName + ':'
244  (
245  IOobject::groupName(name, dxIter.key()),
246  phasePair(phase1_, phase2_).name()
247  ),
248  f**dxIter()
249  ).ptr()
250  );
251  }
252  }
253  };
254 
255  if (model_.valid())
256  {
257  addToXs(scalar(1) - f1() - f2(), (model_().*method)(args ...));
258  }
259 
260  if (model1In2_.valid())
261  {
262  addToXs(f1, (model1In2_().*method)(args ...));
263  }
264 
265  if (model2In1_.valid())
266  {
267  addToXs(f2, (model2In1_().*method)(args ...));
268  }
269 
270  if
271  (
272  correctFixedFluxBCs_
273  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
274  )
275  {
276  forAllIter(typename HashPtrTable<typeGeoField>, xs, xIter)
277  {
278  correctFixedFluxBCs(*xIter());
279  }
280  }
281 
282  return xs;
283 }
284 
285 
286 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
287 
288 template<class ModelType>
290 (
291  const phaseModel& phase1,
292  const phaseModel& phase2,
293  const blendingMethod& blending,
294  autoPtr<ModelType> model,
295  autoPtr<ModelType> model1In2,
296  autoPtr<ModelType> model2In1,
297  const bool correctFixedFluxBCs
298 )
299 :
300  regIOobject
301  (
302  IOobject
303  (
304  IOobject::groupName(typeName, phasePair(phase1, phase2).name()),
305  phase1.mesh().time().timeName(),
306  phase1.mesh()
307  )
308  ),
309  phase1_(phase1),
310  phase2_(phase2),
311  blending_(blending),
312  model_(model),
313  model1In2_(model1In2),
314  model2In1_(model2In1),
315  correctFixedFluxBCs_(correctFixedFluxBCs)
316 {}
317 
318 
319 template<class ModelType>
321 (
322  const phasePair::dictTable& modelTable,
323  const blendingMethod& blending,
324  const phasePair& pair,
325  const orderedPhasePair& pair1In2,
326  const orderedPhasePair& pair2In1,
327  const bool correctFixedFluxBCs
328 )
329 :
330  regIOobject
331  (
332  IOobject
333  (
334  IOobject::groupName(typeName, pair.name()),
335  pair.phase1().mesh().time().timeName(),
336  pair.phase1().mesh()
337  )
338  ),
339  phase1_(pair.phase1()),
340  phase2_(pair.phase2()),
341  blending_(blending),
342  correctFixedFluxBCs_(correctFixedFluxBCs)
343 {
344  if (modelTable.found(pair))
345  {
346  model_.set
347  (
349  (
350  modelTable[pair],
351  pair
352  ).ptr()
353  );
354  }
355 
356  if (modelTable.found(pair1In2))
357  {
358  model1In2_.set
359  (
361  (
362  modelTable[pair1In2],
363  pair1In2
364  ).ptr()
365  );
366  }
367 
368  if (modelTable.found(pair2In1))
369  {
370  model2In1_.set
371  (
373  (
374  modelTable[pair2In1],
375  pair2In1
376  ).ptr()
377  );
378  }
379 }
380 
381 
382 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
383 
384 template<class ModelType>
386 {}
387 
388 
389 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
390 
391 template<class ModelType>
394 {
395  tmp<volScalarField> (ModelType::*k)() const = &ModelType::K;
396 
397  return evaluate(k, "K", ModelType::dimK, false);
398 }
399 
400 
401 template<class ModelType>
403 Foam::BlendedInterfacialModel<ModelType>::K(const scalar residualAlpha) const
404 {
405  tmp<volScalarField> (ModelType::*k)(const scalar) const = &ModelType::K;
406 
407  return evaluate(k, "K", ModelType::dimK, false, residualAlpha);
408 }
409 
410 
411 template<class ModelType>
414 {
415  return evaluate(&ModelType::Kf, "Kf", ModelType::dimK, false);
416 }
417 
418 
419 template<class ModelType>
420 template<class Type>
423 {
424  return evaluate(&ModelType::F, "F", ModelType::dimF, true);
425 }
426 
427 
428 template<class ModelType>
431 {
432  return evaluate(&ModelType::Ff, "Ff", ModelType::dimF*dimArea, true);
433 }
434 
435 
436 template<class ModelType>
439 {
440  return evaluate(&ModelType::D, "D", ModelType::dimD, false);
441 }
442 
443 
444 template<class ModelType>
446 {
447  return
448  (model1In2_.valid() && model1In2_->mixture())
449  || (model2In1_.valid() && model2In1_->mixture())
450  || (model_.valid() && model_->mixture());
451 }
452 
453 
454 template<class ModelType>
457 {
458  return evaluate(&ModelType::dmdtf, "dmdtf", ModelType::dimDmdt, true);
459 }
460 
461 
462 template<class ModelType>
464 {
465  wordList species;
466 
467  if (model1In2_.valid())
468  {
469  species.append(model1In2_->species());
470  }
471  if (model2In1_.valid())
472  {
473  species.append(model2In1_->species());
474  }
475  if (model_.valid())
476  {
477  species.append(model_->species());
478  }
479 
480  return hashedWordList(move(species));
481 }
482 
483 
484 template<class ModelType>
487 {
488  return evaluate(&ModelType::dmidtf, "dmidtf", ModelType::dimDmdt, true);
489 }
490 
491 
492 template<class ModelType>
494 {
495  return os.good();
496 }
497 
498 
499 // ************************************************************************* //
tmp< volScalarField > dmdtf() const
Return the blended mass transfer rate.
const dimensionedScalar & F
Faraday constant: default SI units: [C/mol].
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool writeData(Ostream &os) const
Dummy write for regIOobject.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
~BlendedInterfacialModel()
Destructor.
bool mixture() const
Return the list of individual species that are transferred.
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.
label k
Boltzmann constant.
A HashTable specialization for hashing pointers.
Definition: HashPtrTable.H:50
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
CGAL::Exact_predicates_exact_constructions_kernel K
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
scalar f1
Definition: createFields.H:28
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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 word groupName(Name name, const word &group)
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
hashedWordList species() const
Return the list of individual species that are transferred.
labelList f(nPoints)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
A wordList with hashed indices for faster lookup by name.
tmp< volScalarField > D() const
Return the blended diffusivity.
HashPtrTable< volScalarField > dmidtf() const
Return the blended mass transfer rates for individual species.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.