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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2018 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "surfaceInterpolate.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39 
40 template<>
43 {
44  return f;
45 }
46 
47 
48 template<>
51 {
52  return fvc::interpolate(f);
53 }
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 } // End namespace Foam
58 
59 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60 
61 template<class ModelType>
62 template<class GeoField>
64 (
65  GeoField& field
66 ) const
67 {
68  typename GeoField::Boundary& fieldBf = field.boundaryFieldRef();
69 
70  forAll(phase1_.phi()().boundaryField(), patchi)
71  {
72  if
73  (
74  isA<fixedValueFvsPatchScalarField>
75  (
76  phase1_.phi()().boundaryField()[patchi]
77  )
78  )
79  {
80  fieldBf[patchi] = Zero;
81  }
82  }
83 }
84 
85 
86 template<class ModelType>
87 template
88 <
89  class Type,
90  template<class> class PatchField,
91  class GeoMesh,
92  class ... Args
93 >
96 (
98  (ModelType::*method)(Args ...) const,
99  const word& name,
100  const dimensionSet& dims,
101  const bool subtract,
102  Args ... args
103 ) const
104 {
105  typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
106  typedef GeometricField<Type, PatchField, GeoMesh> typeGeoField;
107 
108  tmp<scalarGeoField> f1, f2;
109 
110  if (model_ || model1In2_)
111  {
112  f1 =
113  blendedInterfacialModel::interpolate<scalarGeoField>
114  (
115  blending_.f1(phase1_, phase2_)
116  );
117  }
118 
119  if (model_ || model2In1_)
120  {
121  f2 =
122  blendedInterfacialModel::interpolate<scalarGeoField>
123  (
124  blending_.f2(phase1_, phase2_)
125  );
126  }
127 
129  (
130  new typeGeoField
131  (
132  IOobject
133  (
134  ModelType::typeName + ":" + name,
135  phase1_.mesh().time().timeName(),
136  phase1_.mesh(),
137  IOobject::NO_READ,
138  IOobject::NO_WRITE,
139  IOobject::NO_REGISTER
140  ),
141  phase1_.mesh(),
142  dimensioned<Type>(dims, Zero)
143  )
144  );
145 
146  if (model_)
147  {
148  if (subtract)
149  {
151  << "Cannot treat an interfacial model with no distinction "
152  << "between continuous and dispersed phases as signed"
153  << exit(FatalError);
154  }
155 
156  x.ref() += (model_().*method)(args ...)*(scalar(1) - f1() - f2());
157  }
158 
159  if (model1In2_)
160  {
161  x.ref() += (model1In2_().*method)(args ...)*f1;
162  }
163 
164  if (model2In1_)
165  {
166  tmp<typeGeoField> dx = (model2In1_().*method)(args ...)*f2;
167 
168  if (subtract)
169  {
170  x.ref() -= dx;
171  }
172  else
173  {
174  x.ref() += dx;
175  }
176  }
177 
178  if
179  (
180  correctFixedFluxBCs_
181  && (model_ || model1In2_ || model2In1_)
182  )
183  {
184  correctFixedFluxBCs(x.ref());
185  }
186 
187  return x;
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
193 template<class ModelType>
195 (
196  const phaseModel& phase1,
197  const phaseModel& phase2,
198  const blendingMethod& blending,
199  autoPtr<ModelType> model,
200  autoPtr<ModelType> model1In2,
201  autoPtr<ModelType> model2In1,
202  const bool correctFixedFluxBCs
203 )
204 :
206  (
207  IOobject
208  (
209  IOobject::groupName(typeName, phasePair(phase1, phase2).name()),
210  phase1.mesh().time().timeName(),
211  phase1.mesh()
212  )
213  ),
214  phase1_(phase1),
215  phase2_(phase2),
216  blending_(blending),
217  model_(model),
218  model1In2_(model1In2),
219  model2In1_(model2In1),
220  correctFixedFluxBCs_(correctFixedFluxBCs)
221 {}
222 
223 
224 template<class ModelType>
226 (
227  const phasePair::dictTable& modelTable,
228  const blendingMethod& blending,
229  const phasePair& pair,
230  const orderedPhasePair& pair1In2,
231  const orderedPhasePair& pair2In1,
232  const bool correctFixedFluxBCs
233 )
234 :
236  (
237  IOobject
238  (
239  IOobject::groupName(typeName, pair.name()),
240  pair.phase1().mesh().time().timeName(),
241  pair.phase1().mesh()
242  )
243  ),
244  phase1_(pair.phase1()),
245  phase2_(pair.phase2()),
246  blending_(blending),
247  correctFixedFluxBCs_(correctFixedFluxBCs)
248 {
249  if (modelTable.found(pair))
250  {
251  model_.set
252  (
254  (
255  modelTable[pair],
256  pair
257  ).ptr()
258  );
259  }
260 
261  if (modelTable.found(pair1In2))
262  {
263  model1In2_.set
264  (
266  (
267  modelTable[pair1In2],
268  pair1In2
269  ).ptr()
270  );
271  }
272 
273  if (modelTable.found(pair2In1))
274  {
275  model2In1_.set
276  (
278  (
279  modelTable[pair2In1],
280  pair2In1
281  ).ptr()
282  );
283  }
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
288 
289 template<class ModelType>
291 {}
292 
293 
294 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
295 
296 template<class ModelType>
298 (
299  const class phaseModel& phase
300 ) const
301 {
302  return
303  &phase == &(phase1_)
304  ? bool(model1In2_)
305  : bool(model2In1_);
306 }
307 
308 
309 template<class ModelType>
311 (
312  const class phaseModel& phase
313 ) const
314 {
315  return &phase == &(phase1_) ? model1In2_ : model2In1_;
316 }
317 
318 
319 template<class ModelType>
322 {
323  tmp<volScalarField> (ModelType::*k)() const = &ModelType::K;
325  return evaluate(k, "K", ModelType::dimK, false);
326 }
327 
328 
329 template<class ModelType>
331 Foam::BlendedInterfacialModel<ModelType>::K(const scalar residualAlpha) const
332 {
333  tmp<volScalarField> (ModelType::*k)(const scalar) const = &ModelType::K;
335  return evaluate(k, "K", ModelType::dimK, false, residualAlpha);
336 }
337 
338 
339 template<class ModelType>
342 {
343  return evaluate(&ModelType::Kf, "Kf", ModelType::dimK, false);
344 }
345 
346 
347 template<class ModelType>
348 template<class Type>
351 {
352  return evaluate(&ModelType::F, "F", ModelType::dimF, true);
353 }
354 
355 
356 template<class ModelType>
359 {
360  return evaluate(&ModelType::Ff, "Ff", ModelType::dimF*dimArea, true);
361 }
362 
363 
364 template<class ModelType>
367 {
368  return evaluate(&ModelType::D, "D", ModelType::dimD, false);
369 }
370 
371 
372 template<class ModelType>
375 {
376  return evaluate(&ModelType::dmdt, "dmdt", ModelType::dimDmdt, false);
377 }
378 
379 
380 template<class ModelType>
382 {
383  return os.good();
384 }
385 
386 
387 // ************************************************************************* //
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
rDeltaTY field()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1354
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
tmp< volScalarField > K() const
Return the blended force coefficient.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
label k
Boltzmann constant.
Generic dimensioned Type class.
tmp< surfaceScalarField > Ff() const
Return the face blended force.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
phaseModel & phase2
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)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
static tmp< GeoField > interpolate(tmp< volScalarField > f)
Convenience function to interpolate blending values. Needs to be.
volVectorField F(fluid.F())
A HashTable similar to std::unordered_map.
Definition: HashTable.H:108
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
phaseModel & phase1
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
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< surfaceScalarField > Kf() const
Return the face blended force coefficient.
tmp< volScalarField > D() const
Return the blended diffusivity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
const dimensionedScalar & D
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
A class for managing temporary objects.
Definition: HashPtrTable.H:50
surfaceScalarField Ff(fluid.Ff())
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
tmp< volScalarField > dmdt() const
Return the blended mass transfer rate.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127