phaseSystemTemplates.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) 2015-2019 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 
28 #include "BlendedInterfacialModel.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class modelType>
34 (
35  const dictTable& modelDicts,
36  HashTable
37  <
41  >& models
42 )
43 {
44  forAllConstIter(dictTable, modelDicts, iter)
45  {
46  const phasePairKey& key = iter.key();
47 
48  models.insert
49  (
50  key,
52  (
53  *iter,
54  phasePairs_[key]
55  )
56  );
57  }
58 }
59 
60 
61 template<class modelType>
63 (
64  const word& modelName,
65  HashTable
66  <
70  >& models
71 )
72 {
73  dictTable modelDicts(lookup(modelName));
74 
75  generatePairs(modelDicts);
76 
77  createSubModels(modelDicts, models);
78 }
79 
80 
81 template<class modelType>
83 (
84  const word& modelName,
85  HashTable
86  <
90  >& models,
91  const bool correctFixedFluxBCs
92 )
93 {
94  typedef
96  modelTypeTable;
97 
98  modelTypeTable tempModels;
99  generatePairsAndSubModels(modelName, tempModels);
100 
101  const blendingMethod& blending
102  (
103  blendingMethods_.found(modelName)
104  ? blendingMethods_[modelName]
105  : blendingMethods_.found(member(modelName))
106  ? blendingMethods_[member(modelName)]
107  : blendingMethods_["default"]
108  );
109 
110  autoPtr<modelType> noModel(nullptr);
111 
112  forAllConstIter(typename modelTypeTable, tempModels, iter)
113  {
114  if (!iter().valid())
115  {
116  continue;
117  }
118 
119  const phasePairKey key(iter.key().first(), iter.key().second());
120  const phasePairKey key1In2(key.first(), key.second(), true);
121  const phasePairKey key2In1(key.second(), key.first(), true);
122 
123  models.insert
124  (
125  key,
127  (
129  (
130  phaseModels_[key.first()],
131  phaseModels_[key.second()],
132  blending,
133  tempModels.found(key ) ? tempModels[key ] : noModel,
134  tempModels.found(key1In2) ? tempModels[key1In2] : noModel,
135  tempModels.found(key2In1) ? tempModels[key2In1] : noModel,
136  correctFixedFluxBCs
137  )
138  )
139  );
140 
141  if (!phasePairs_.found(key))
142  {
143  phasePairs_.insert
144  (
145  key,
147  (
148  new phasePair
149  (
150  phaseModels_[key.first()],
151  phaseModels_[key.second()]
152  )
153  )
154  );
155  }
156  }
157 }
158 
159 
160 template<class modelType>
162 (
163  const word& modelName,
164  HashTable
165  <
167  phasePairKey,
169  >& models,
170  const bool correctFixedFluxBCs
171 )
172 {
173  typedef
175  modelTypeTable;
176 
177  forAll(phaseModels_, phasei)
178  {
179  const phaseModel& phase = phaseModels_[phasei];
180 
181  modelTypeTable tempModels;
182  generatePairsAndSubModels
183  (
184  IOobject::groupName(modelName, phase.name()),
185  tempModels,
186  correctFixedFluxBCs
187  );
188 
189  forAllIter(typename modelTypeTable, tempModels, tempModelIter)
190  {
191  const phasePairKey& key(tempModelIter.key());
192 
193  if (!models.found(key))
194  {
195  models.insert
196  (
197  key,
199  );
200  }
201 
202  const phasePair& pair = phasePairs_[key];
203 
204  if (!pair.contains(phase))
205  {
207  << "A two-sided " << modelType::typeName << " was "
208  << "specified for the " << phase.name() << " side of the "
209  << pair << " pair, but that phase is not part of that pair."
210  << exit(FatalError);
211  }
212 
213  models[key][pair.index(phase)].set(tempModelIter().ptr());
214  }
215  }
216 }
217 
218 
219 template<class GeoField>
221 (
222  const phaseModel& phase,
223  const word& fieldName,
225  PtrList<GeoField>& fieldList
226 ) const
227 {
228  if (fieldList.set(phase.index()))
229  {
230  fieldList[phase.index()] += field;
231  }
232  else
233  {
234  fieldList.set
235  (
236  phase.index(),
237  new GeoField
238  (
239  IOobject::groupName(fieldName, phase.name()),
240  field
241  )
242  );
243  }
244 }
245 
246 
247 template<class GeoField>
249 (
250  const phaseModel& phase,
251  const word& fieldName,
252  const GeoField& field,
253  PtrList<GeoField>& fieldList
254 ) const
255 {
256  addField(phase, fieldName, tmp<GeoField>(field), fieldList);
257 }
258 
259 
260 template<class GeoField>
262 (
263  const phaseModel& phase,
264  const word& fieldName,
266  HashPtrTable<GeoField>& fieldTable
267 ) const
268 {
269  if (fieldTable.found(phase.name()))
270  {
271  *fieldTable[phase.name()] += field;
272  }
273  else
274  {
275  fieldTable.set
276  (
277  phase.name(),
278  new GeoField
279  (
280  IOobject::groupName(fieldName, phase.name()),
281  field
282  )
283  );
284  }
285 }
286 
287 
288 template<class GeoField>
290 (
291  const phaseModel& phase,
292  const word& fieldName,
293  const GeoField& field,
294  HashPtrTable<GeoField>& fieldTable
295 ) const
296 {
297  addField(phase, fieldName, tmp<GeoField>(field), fieldTable);
298 }
299 
300 
301 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
302 
303 template<class Type, template<class> class PatchField, class GeoMesh>
305 (
306  const word& name,
307  const dimensionSet& dims,
309 ) const
310 {
311  forAll(this->phaseModels_, phasei)
312  {
313  if (fieldList.set(phasei))
314  {
315  continue;
316  }
317 
318  const phaseModel& phase = this->phaseModels_[phasei];
319 
320  fieldList.set
321  (
322  phasei,
324  (
325  IOobject
326  (
327  IOobject::groupName(name, phase.name()),
328  this->mesh_.time().timeName(),
329  this->mesh_
330  ),
331  this->mesh_,
333  )
334  );
335  }
336 }
337 
338 
339 template<class Type, template<class> class PatchField, class GeoMesh>
341 (
342  const word& name,
343  const dimensionSet& dims,
345 ) const
346 {
347  forAll(this->phaseModels_, phasei)
348  {
349  const phaseModel& phase = this->phaseModels_[phasei];
350 
351  if (fieldTable.set(phase.name()))
352  {
353  continue;
354  }
355 
356  fieldTable.set
357  (
358  phase.name(),
360  (
361  IOobject
362  (
363  IOobject::groupName(name, phase.name()),
364  this->mesh_.time().timeName(),
365  this->mesh_
366  ),
367  this->mesh_,
369  )
370  );
371  }
372 }
373 
374 
375 template<class modelType>
377 {
378  const word name(IOobject::groupName(modelType::typeName, key.name()));
379 
380  if (key.ordered())
381  {
382  if (mesh().foundObject<modelType>(name))
383  {
384  return true;
385  }
386  else
387  {
388  return false;
389  }
390  }
391  else
392  {
393  if
394  (
395  mesh().foundObject<modelType>(name)
396  ||
397  mesh().foundObject<modelType>
398  (
399  IOobject::groupName(modelType::typeName, key.otherName())
400  )
401  )
402  {
403  return true;
404  }
405  else
406  {
407  return false;
408  }
409  }
410 }
411 
412 
413 template<class modelType>
414 const modelType& Foam::phaseSystem::lookupSubModel(const phasePair& key) const
415 {
416  const word name(IOobject::groupName(modelType::typeName, key.name()));
417 
418  if (key.ordered() || mesh().foundObject<modelType>(name))
419  {
420  return mesh().lookupObject<modelType>(name);
421  }
422  else
423  {
424  return
425  mesh().lookupObject<modelType>
426  (
427  IOobject::groupName(modelType::typeName, key.otherName())
428  );
429  }
430 }
431 
432 
433 template<class modelType>
435 (
436  const phaseModel& dispersed,
437  const phaseModel& continuous
438 ) const
439 {
440  return foundSubModel<modelType>(orderedPhasePair(dispersed, continuous));
441 }
442 
443 
444 template<class modelType>
445 const modelType& Foam::phaseSystem::lookupSubModel
446 (
447  const phaseModel& dispersed,
448  const phaseModel& continuous
449 ) const
450 {
451  return lookupSubModel<modelType>(orderedPhasePair(dispersed, continuous));
452 }
453 
454 
455 template<class modelType>
457 {
458  if
459  (
461  (
463  (
465  key.name()
466  )
467  )
469  (
471  (
473  key.otherName()
474  )
475  )
476  )
477  {
478  return true;
479  }
480  else
481  {
482  return false;
483  }
484 }
485 
486 
487 template<class modelType>
490 {
491  const word name
492  (
494  (
496  key.name()
497  )
498  );
499 
500  if (mesh().foundObject<BlendedInterfacialModel<modelType>>(name))
501  {
503  }
504  else
505  {
506  return
508  (
510  (
512  key.otherName()
513  )
514  );
515  }
516 }
517 
518 
519 // ************************************************************************* //
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::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:51
Foam::phasePair::otherName
virtual word otherName() const
Other pair name.
Definition: phasePair.C:100
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:319
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
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::phasePairKey::hash
Ordered or unordered hashing of word pair.
Definition: phasePairKey.H:65
Foam::phaseSystem::fillFields
void fillFields(const word &name, const dimensionSet &dims, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fieldList) const
Fill up gaps in a phase-indexed list of fields with zeros.
Definition: phaseSystemTemplates.C:305
Foam::phaseSystem::createSubModels
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
Definition: phaseSystemTemplates.C:34
Foam::phaseSystem::addField
void addField(const phaseModel &phase, const word &fieldName, tmp< GeoField > field, PtrList< GeoField > &fieldList) const
Add the field to a phase-indexed list, with the given name,.
Definition: phaseSystemTemplates.C:221
Foam::phasePair::index
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
Definition: phasePairI.H:72
phasei
label phasei
Definition: pEqn.H:27
Foam::phaseSystem::lookupSubModel
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
Definition: phaseSystemTemplates.C:414
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::baseIOdictionary::name
const word & name() const
Name function is needed to disambiguate those inherited.
Definition: baseIOdictionary.C:88
Foam::phaseSystem::foundSubModel
bool foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
Definition: phaseSystemTemplates.C:376
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
forAllConstIter
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:338
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::HashPtrTable::set
bool set(const Key &key, T *ptr)
Assign a new entry, overwriting existing entries.
Definition: HashPtrTableI.H:84
Foam::phaseSystem::generatePairsAndSubModels
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
Definition: phaseSystemTemplates.C:63
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::cellModeller::lookup
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
Definition: cellModeller.H:46
Foam::orderedPhasePair
Definition: orderedPhasePair.H:49
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::phasePairKey
Definition: phasePairKey.H:59
Foam::phaseSystem::lookupBlendedSubModel
const BlendedInterfacialModel< modelType > & lookupBlendedSubModel(const phasePair &key) const
Return a blended sub model between a phase pair.
Foam::FatalError
error FatalError
Foam::phasePair::contains
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Definition: phasePairI.H:42
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:48
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: PtrListI.H:143
Foam::phasePairKey::ordered
bool ordered() const
Return the ordered flag.
Definition: phasePairKey.C:60
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::phaseSystem::mesh
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:30
Foam::phase::name
const word & name() const
Definition: phase.H:111
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 >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::HashPtrTable
A HashTable of pointers to objects of type <T>.
Definition: HashPtrTable.H:54
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
Foam::BlendedInterfacialModel
Definition: BlendedInterfacialModel.H:68
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
Foam::phasePair::name
virtual word name() const
Pair name.
Definition: phasePair.C:92
Foam::phaseSystem::foundBlendedSubModel
bool foundBlendedSubModel(const phasePair &key) const
Check availability of a blended sub model for a given phase pair.
Definition: phaseSystemTemplates.C:456
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53