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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36 
37 template<>
38 inline tmp<Foam::volScalarField>
39 blendedInterfacialModel::interpolate(tmp<volScalarField> f)
40 {
41  return f;
42 }
43 
44 
45 template<>
46 inline tmp<Foam::surfaceScalarField>
47 blendedInterfacialModel::interpolate(tmp<volScalarField> f)
48 {
49  return fvc::interpolate(f);
50 }
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 } // End namespace Foam
55 
56 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
57 
58 template<class ModelType>
59 template<class GeoField>
61 (
62  GeoField& field
63 ) const
64 {
65  typename GeoField::Boundary& fieldBf = field.boundaryFieldRef();
66 
67  forAll(phase1_.phi()().boundaryField(), patchi)
68  {
69  if
70  (
71  isA<fixedValueFvsPatchScalarField>
72  (
73  phase1_.phi()().boundaryField()[patchi]
74  )
75  )
76  {
77  fieldBf[patchi] = Zero;
78  }
79  }
80 }
81 
82 
83 template<class ModelType>
84 template
85 <
86  class Type,
87  template<class> class PatchField,
88  class GeoMesh,
89  class ... Args
90 >
93 (
94  tmp<GeometricField<Type, PatchField, GeoMesh>>
95  (ModelType::*method)(Args ...) const,
96  const word& name,
97  const dimensionSet& dims,
98  const bool subtract,
99  Args ... args
100 ) const
101 {
102  typedef GeometricField<scalar, PatchField, GeoMesh> scalarGeoField;
103  typedef GeometricField<Type, PatchField, GeoMesh> typeGeoField;
104 
105  tmp<scalarGeoField> f1, f2;
106 
107  if (model_.valid() || model1In2_.valid())
108  {
109  f1 =
110  blendedInterfacialModel::interpolate<scalarGeoField>
111  (
112  blending_.f1(phase1_, phase2_)
113  );
114  }
115 
116  if (model_.valid() || model2In1_.valid())
117  {
118  f2 =
119  blendedInterfacialModel::interpolate<scalarGeoField>
120  (
121  blending_.f2(phase1_, phase2_)
122  );
123  }
124 
125  tmp<typeGeoField> x
126  (
127  new typeGeoField
128  (
129  IOobject
130  (
131  ModelType::typeName + ":" + name,
132  phase1_.mesh().time().timeName(),
133  phase1_.mesh(),
136  false
137  ),
138  phase1_.mesh(),
139  dimensioned<Type>("zero", dims, Zero)
140  )
141  );
142 
143  if (model_.valid())
144  {
145  if (subtract)
146  {
148  << "Cannot treat an interfacial model with no distinction "
149  << "between continuous and dispersed phases as signed"
150  << exit(FatalError);
151  }
152 
153  x.ref() += (model_().*method)(args ...)*(scalar(1) - f1() - f2());
154  }
155 
156  if (model1In2_.valid())
157  {
158  x.ref() += (model1In2_().*method)(args ...)*f1;
159  }
160 
161  if (model2In1_.valid())
162  {
163  tmp<typeGeoField> dx = (model2In1_().*method)(args ...)*f2;
164 
165  if (subtract)
166  {
167  x.ref() -= dx;
168  }
169  else
170  {
171  x.ref() += dx;
172  }
173  }
174 
175  if
176  (
177  correctFixedFluxBCs_
178  && (model_.valid() || model1In2_.valid() || model2In1_.valid())
179  )
180  {
181  correctFixedFluxBCs(x.ref());
182  }
183 
184  return x;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189 
190 template<class ModelType>
192 (
193  const phaseModel& phase1,
194  const phaseModel& phase2,
195  const blendingMethod& blending,
196  autoPtr<ModelType> model,
197  autoPtr<ModelType> model1In2,
198  autoPtr<ModelType> model2In1,
199  const bool correctFixedFluxBCs
200 )
201 :
202  regIOobject
203  (
204  IOobject
205  (
206  IOobject::groupName(typeName, phasePair(phase1, phase2).name()),
207  phase1.mesh().time().timeName(),
208  phase1.mesh()
209  )
210  ),
211  phase1_(phase1),
212  phase2_(phase2),
213  blending_(blending),
214  model_(model),
215  model1In2_(model1In2),
216  model2In1_(model2In1),
217  correctFixedFluxBCs_(correctFixedFluxBCs)
218 {}
219 
220 
221 template<class ModelType>
223 (
224  const phasePair::dictTable& modelTable,
225  const blendingMethod& blending,
226  const phasePair& pair,
227  const orderedPhasePair& pair1In2,
228  const orderedPhasePair& pair2In1,
229  const bool correctFixedFluxBCs
230 )
231 :
232  regIOobject
233  (
234  IOobject
235  (
236  IOobject::groupName(typeName, pair.name()),
237  pair.phase1().mesh().time().timeName(),
238  pair.phase1().mesh()
239  )
240  ),
241  phase1_(pair.phase1()),
242  phase2_(pair.phase2()),
243  blending_(blending),
244  correctFixedFluxBCs_(correctFixedFluxBCs)
245 {
246  if (modelTable.found(pair))
247  {
248  model_.set
249  (
251  (
252  modelTable[pair],
253  pair
254  ).ptr()
255  );
256  }
257 
258  if (modelTable.found(pair1In2))
259  {
260  model1In2_.set
261  (
263  (
264  modelTable[pair1In2],
265  pair1In2
266  ).ptr()
267  );
268  }
269 
270  if (modelTable.found(pair2In1))
271  {
272  model2In1_.set
273  (
275  (
276  modelTable[pair2In1],
277  pair2In1
278  ).ptr()
279  );
280  }
281 }
282 
283 
284 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
285 
286 template<class ModelType>
288 {}
289 
290 
291 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
292 
293 template<class ModelType>
295 (
296  const class phaseModel& phase
297 ) const
298 {
299  return
300  &phase == &(phase1_)
301  ? model1In2_.valid()
302  : model2In1_.valid();
303 }
304 
305 
306 template<class ModelType>
308 (
309  const class phaseModel& phase
310 ) const
311 {
312  return &phase == &(phase1_) ? model1In2_ : model2In1_;
313 }
314 
315 
316 template<class ModelType>
319 {
320  tmp<volScalarField> (ModelType::*k)() const = &ModelType::K;
321 
322  return evaluate(k, "K", ModelType::dimK, false);
323 }
324 
325 
326 template<class ModelType>
328 Foam::BlendedInterfacialModel<ModelType>::K(const scalar residualAlpha) const
329 {
330  tmp<volScalarField> (ModelType::*k)(const scalar) const = &ModelType::K;
331 
332  return evaluate(k, "K", ModelType::dimK, false, residualAlpha);
333 }
334 
335 
336 template<class ModelType>
339 {
340  return evaluate(&ModelType::Kf, "Kf", ModelType::dimK, false);
341 }
342 
343 
344 template<class ModelType>
345 template<class Type>
348 {
349  return evaluate(&ModelType::F, "F", ModelType::dimF, true);
350 }
351 
352 
353 template<class ModelType>
356 {
357  return evaluate(&ModelType::Ff, "Ff", ModelType::dimF*dimArea, true);
358 }
359 
360 
361 template<class ModelType>
364 {
365  return evaluate(&ModelType::D, "D", ModelType::dimD, false);
366 }
367 
368 
369 template<class ModelType>
372 {
373  return evaluate(&ModelType::dmdt, "dmdt", ModelType::dimDmdt, false);
374 }
375 
376 
377 template<class ModelType>
379 {
380  return os.good();
381 }
382 
383 
384 // ************************************************************************* //
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool hasModel(const phaseModel &phase) const
Return true if a model is specified for the supplied phase.
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
~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< volScalarField > dmdt() const
Return the blended mass transfer rate.
label k
Boltzmann constant.
tmp< GeometricField< Type, fvPatchField, volMesh > > F() const
Return the blended force.
CGAL::Exact_predicates_exact_constructions_kernel K
scalar f1
Definition: createFields.H:28
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const ModelType & model(const phaseModel &phase) const
Return the model for the supplied phase.
dynamicFvMesh & mesh
static tmp< GeoField > interpolate(tmp< volScalarField > f)
Convenience function to interpolate blending values. Needs to be.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
volVectorField F(fluid.F())
phaseModel & phase1
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.
label patchi
tmp< volScalarField > D() const
Return the blended diffusivity.
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:53
surfaceScalarField Ff(fluid.Ff())
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.