phaseSystemTemplates.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) 2015-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 
27 
28 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
29 
30 template<class modelType>
32 (
33  const dictTable& modelDicts,
34  HashTable
35  <
36  autoPtr<modelType>,
37  phasePairKey,
38  phasePairKey::hash
39  >& models
40 )
41 {
42  forAllConstIter(dictTable, modelDicts, iter)
43  {
44  const phasePairKey& key = iter.key();
45 
46  models.insert
47  (
48  key,
50  (
51  *iter,
52  phasePairs_[key]
53  )
54  );
55  }
56 }
57 
58 
59 template<class modelType>
61 (
62  const word& modelName,
63  HashTable
64  <
65  autoPtr<modelType>,
66  phasePairKey,
67  phasePairKey::hash
68  >& models,
69  const bool correctFixedFluxBCs
70 )
71 {
72  dictTable modelDicts(lookup(modelName));
73 
74  generatePairs(modelDicts);
75 
76  createSubModels(modelDicts, models);
77 }
78 
79 
80 template<class modelType>
82 (
83  const word& modelName,
84  HashTable
85  <
86  autoPtr<BlendedInterfacialModel<modelType>>,
87  phasePairKey,
88  phasePairKey::hash
89  >& models,
90  const bool correctFixedFluxBCs
91 )
92 {
93  typedef
94  HashTable<autoPtr<modelType>, phasePairKey, phasePairKey::hash>
95  modelTypeTable;
96 
97  modelTypeTable tempModels;
98  generatePairsAndSubModels(modelName, tempModels);
99 
100  const blendingMethod& blending
101  (
102  blendingMethods_.found(modelName)
104  : blendingMethods_.found(member(modelName))
105  ? blendingMethods_[member(modelName)]
106  : blendingMethods_["default"]
107  );
108 
109  autoPtr<modelType> noModel(nullptr);
110 
111  forAllConstIter(typename modelTypeTable, tempModels, iter)
112  {
113  if (!iter().valid())
114  {
115  continue;
116  }
117 
118  const phasePairKey key(iter.key().first(), iter.key().second());
119 
120  if (!phasePairs_.found(key))
121  {
123  (
124  key,
125  autoPtr<phasePair>
126  (
127  new phasePair
128  (
129  phaseModels_[key.first()],
130  phaseModels_[key.second()]
131  )
132  )
133  );
134  }
135 
136  const phasePair& pair = phasePairs_[key];
137  const phaseModel& phase1 = pair.phase1();
138  const phaseModel& phase2 = pair.phase2();
139  const phasePairKey key1In2(phase1.name(), phase2.name(), true);
140  const phasePairKey key2In1(phase2.name(), phase1.name(), true);
141 
142  models.insert
143  (
144  key,
145  autoPtr<BlendedInterfacialModel<modelType>>
146  (
147  new BlendedInterfacialModel<modelType>
148  (
149  phase1,
150  phase2,
151  blending,
152  tempModels.found(key) ? tempModels[key] : noModel,
153  tempModels.found(key1In2) ? tempModels[key1In2] : noModel,
154  tempModels.found(key2In1) ? tempModels[key2In1] : noModel,
155  correctFixedFluxBCs
156  )
157  )
158  );
159  }
160 }
161 
162 
163 template<class modelType>
165 (
166  const word& modelName,
167  HashTable
168  <
169  Pair<autoPtr<modelType>>,
170  phasePairKey,
171  phasePairKey::hash
172  >& models,
173  const bool correctFixedFluxBCs
174 )
175 {
176  typedef
177  HashTable<autoPtr<modelType>, phasePairKey, phasePairKey::hash>
178  modelTypeTable;
179 
181  {
182  const phaseModel& phase = phaseModels_[phasei];
183 
184  modelTypeTable tempModels;
186  (
187  IOobject::groupName(modelName, phase.name()),
188  tempModels,
189  correctFixedFluxBCs
190  );
191 
192  forAllIter(typename modelTypeTable, tempModels, tempModelIter)
193  {
194  const phasePairKey& key(tempModelIter.key());
195 
196  if (!models.found(key))
197  {
198  models.insert
199  (
200  key,
201  Pair<autoPtr<modelType>>()
202  );
203  }
204 
205  const phasePair& pair = phasePairs_[key];
206 
207  if (!pair.contains(phase))
208  {
210  << "A two-sided " << modelType::typeName << " was "
211  << "specified for the " << phase.name() << " side of the "
212  << pair << " pair, but that phase is not part of that pair."
213  << exit(FatalError);
214  }
215 
216  models[key][pair.index(phase)] = tempModelIter().ptr();
217  }
218  }
219 }
220 
221 
222 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
223 
224 template<class Type, template<class> class PatchField, class GeoMesh>
226 (
227  const word& name,
228  const dimensionSet& dims,
229  PtrList<GeometricField<Type, PatchField, GeoMesh>>& fieldList
230 ) const
231 {
232  forAll(this->phaseModels_, phasei)
233  {
234  if (fieldList.set(phasei))
235  {
236  continue;
237  }
238 
239  const phaseModel& phase = this->phaseModels_[phasei];
240 
241  fieldList.set
242  (
243  phasei,
244  new GeometricField<Type, PatchField, GeoMesh>
245  (
246  IOobject
247  (
248  IOobject::groupName(name, phase.name()),
249  this->mesh_.time().timeName(),
250  this->mesh_
251  ),
252  this->mesh_,
253  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
254  )
255  );
256  }
257 }
258 
259 
260 template<class Type, template<class> class PatchField, class GeoMesh>
262 (
263  const word& name,
264  const dimensionSet& dims,
265  HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>& fieldTable
266 ) const
267 {
268  forAll(this->phaseModels_, phasei)
269  {
270  const phaseModel& phase = this->phaseModels_[phasei];
271 
272  if (fieldTable.set(phase.name()))
273  {
274  continue;
275  }
276 
277  fieldTable.set
278  (
279  phase.name(),
280  new GeometricField<Type, PatchField, GeoMesh>
281  (
282  IOobject
283  (
284  IOobject::groupName(name, phase.name()),
285  this->mesh_.time().timeName(),
286  this->mesh_
287  ),
288  this->mesh_,
289  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
290  )
291  );
292  }
293 }
294 
295 
296 template<class modelType>
297 bool Foam::phaseSystem::foundSubModel(const phasePair& key) const
298 {
299  const word name(IOobject::groupName(modelType::typeName, key.name()));
300 
301  if (key.ordered())
302  {
303  if (mesh().foundObject<modelType>(name))
304  {
305  return true;
306  }
307  else
308  {
309  return false;
310  }
311  }
312  else
313  {
314  if
315  (
316  mesh().foundObject<modelType>(name)
317  ||
318  mesh().foundObject<modelType>
319  (
320  IOobject::groupName(modelType::typeName, key.otherName())
321  )
322  )
323  {
324  return true;
325  }
326  else
327  {
328  return false;
329  }
330  }
331 }
332 
333 
334 template<class modelType>
335 const modelType& Foam::phaseSystem::lookupSubModel(const phasePair& key) const
336 {
337  const word name(IOobject::groupName(modelType::typeName, key.name()));
338 
339  if (key.ordered() || mesh().foundObject<modelType>(name))
340  {
341  return mesh().lookupObject<modelType>(name);
342  }
343  else
344  {
345  return
346  mesh().lookupObject<modelType>
347  (
348  IOobject::groupName(modelType::typeName, key.otherName())
349  );
350  }
351 }
352 
353 
354 template<class modelType>
356 (
357  const phaseModel& dispersed,
358  const phaseModel& continuous
359 ) const
360 {
361  return foundSubModel<modelType>(orderedPhasePair(dispersed, continuous));
362 }
363 
364 
365 template<class modelType>
366 const modelType& Foam::phaseSystem::lookupSubModel
367 (
368  const phaseModel& dispersed,
369  const phaseModel& continuous
370 ) const
371 {
372  return lookupSubModel<modelType>(orderedPhasePair(dispersed, continuous));
373 }
374 
375 
376 template<class modelType>
377 bool Foam::phaseSystem::foundBlendedSubModel(const phasePair& key) const
378 {
379  if
380  (
381  mesh().foundObject<BlendedInterfacialModel<modelType>>
382  (
384  (
385  BlendedInterfacialModel<modelType>::typeName,
386  key.name()
387  )
388  )
389  || mesh().foundObject<BlendedInterfacialModel<modelType>>
390  (
392  (
393  BlendedInterfacialModel<modelType>::typeName,
394  key.otherName()
395  )
396  )
397  )
398  {
399  return true;
400  }
401  else
402  {
403  return false;
404  }
405 }
406 
407 
408 template<class modelType>
410 Foam::phaseSystem::lookupBlendedSubModel(const phasePair& key) const
411 {
412  const word name
413  (
415  (
416  BlendedInterfacialModel<modelType>::typeName,
417  key.name()
418  )
419  );
420 
421  if (mesh().foundObject<BlendedInterfacialModel<modelType>>(name))
422  {
423  return mesh().lookupObject<BlendedInterfacialModel<modelType>>(name);
424  }
425  else
426  {
427  return
428  mesh().lookupObject<BlendedInterfacialModel<modelType>>
429  (
431  (
432  BlendedInterfacialModel<modelType>::typeName,
433  key.otherName()
434  )
435  );
436  }
437 }
438 
439 
440 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
441 
442 namespace Foam
443 {
444 
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 
447 template<class GeoField, class Group>
448 inline void addField
449 (
450  const Group& group,
451  const word& name,
452  tmp<GeoField> field,
453  PtrList<GeoField>& fieldList
454 )
455 {
456  if (fieldList.set(group.index()))
457  {
458  fieldList[group.index()] += field;
459  }
460  else
461  {
462  fieldList.set
463  (
464  group.index(),
465  new GeoField
466  (
467  IOobject::groupName(name, group.name()),
468  field
469  )
470  );
471  }
472 }
473 
474 
475 template<class GeoField, class Group>
476 inline void addField
477 (
478  const Group& group,
479  const word& name,
480  const GeoField& field,
481  PtrList<GeoField>& fieldList
482 )
483 {
484  addField(group, name, tmp<GeoField>(field), fieldList);
485 }
486 
487 
488 template<class GeoField, class Group>
489 inline void addField
490 (
491  const Group& group,
492  const word& name,
493  tmp<GeoField> field,
494  HashPtrTable<GeoField>& fieldTable
495 )
496 {
497  if (fieldTable.found(group.name()))
498  {
499  *fieldTable[group.name()] += field;
500  }
501  else
502  {
503  fieldTable.set
504  (
505  group.name(),
506  new GeoField
507  (
508  IOobject::groupName(name, group.name()),
509  field
510  )
511  );
512  }
513 }
514 
515 
516 template<class GeoField, class Group>
517 inline void addField
518 (
519  const Group& group,
520  const word& name,
521  const GeoField& field,
522  HashPtrTable<GeoField>& fieldTable
523 )
524 {
525  addField(group, name, tmp<GeoField>(field), fieldTable);
526 }
527 
528 
529 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
530 
531 } // End namespace Foam
532 
533 // ************************************************************************* //
const BlendedInterfacialModel< modelType > & lookupBlendedSubModel(const phasePair &key) const
Return a blended sub model between a phase pair.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
#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
label phasei
Definition: pEqn.H:27
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
autoPtr< T > set(const label, const word &key, T *)
Set element to pointer provided and return old element.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:163
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:177
bool foundBlendedSubModel(const phasePair &key) const
Check availability of a blended sub model for a given phase pair.
bool foundObject(const word &name) const
Is the named Type found?
word member() const
Return member (name without the extension)
Definition: IOobject.C:346
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:148
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
static word groupName(Name name, const word &group)
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:141
bool foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
const word & name() const
Name function is needed to disambiguate those inherited.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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.
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:175
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models, const bool correctFixedFluxBCs=true)
Generate pairs and sub-model tables.
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
rDeltaTY field()
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition: phaseSystem.H:110
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812