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