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-2021 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 template<class modelType>
223 void Foam::phaseSystem::validateMassTransfer(const phasePair& key) const
224 {
225  if (key.phase1().stationary() || key.phase2().stationary())
226  {
228  << "A " << modelType::typeName << " was specified for pair "
229  << key.name() << ", but one of these phases is stationary. "
230  << "Mass transfer is not supported on stationary phases"
231  << exit(FatalError);
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
237 
238 template<class Type, template<class> class PatchField, class GeoMesh>
240 (
241  const word& name,
242  const dimensionSet& dims,
243  PtrList<GeometricField<Type, PatchField, GeoMesh>>& fieldList
244 ) const
245 {
246  forAll(this->phaseModels_, phasei)
247  {
248  if (fieldList.set(phasei))
249  {
250  continue;
251  }
252 
253  const phaseModel& phase = this->phaseModels_[phasei];
254 
255  fieldList.set
256  (
257  phasei,
258  new GeometricField<Type, PatchField, GeoMesh>
259  (
260  IOobject
261  (
262  IOobject::groupName(name, phase.name()),
263  this->mesh_.time().timeName(),
264  this->mesh_
265  ),
266  this->mesh_,
267  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
268  )
269  );
270  }
271 }
272 
273 
274 template<class Type, template<class> class PatchField, class GeoMesh>
276 (
277  const word& name,
278  const dimensionSet& dims,
279  HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>& fieldTable
280 ) const
281 {
282  forAll(this->phaseModels_, phasei)
283  {
284  const phaseModel& phase = this->phaseModels_[phasei];
285 
286  if (fieldTable.set(phase.name()))
287  {
288  continue;
289  }
290 
291  fieldTable.set
292  (
293  phase.name(),
294  new GeometricField<Type, PatchField, GeoMesh>
295  (
296  IOobject
297  (
298  IOobject::groupName(name, phase.name()),
299  this->mesh_.time().timeName(),
300  this->mesh_
301  ),
302  this->mesh_,
303  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
304  )
305  );
306  }
307 }
308 
309 
310 template<class modelType>
311 bool Foam::phaseSystem::foundSubModel(const phasePair& key) const
312 {
313  const word name(IOobject::groupName(modelType::typeName, key.name()));
314 
315  if (key.ordered())
316  {
317  if (mesh().foundObject<modelType>(name))
318  {
319  return true;
320  }
321  else
322  {
323  return false;
324  }
325  }
326  else
327  {
328  if
329  (
330  mesh().foundObject<modelType>(name)
331  ||
332  mesh().foundObject<modelType>
333  (
334  IOobject::groupName(modelType::typeName, key.otherName())
335  )
336  )
337  {
338  return true;
339  }
340  else
341  {
342  return false;
343  }
344  }
345 }
346 
347 
348 template<class modelType>
349 const modelType& Foam::phaseSystem::lookupSubModel(const phasePair& key) const
350 {
351  const word name(IOobject::groupName(modelType::typeName, key.name()));
352 
353  if (key.ordered() || mesh().foundObject<modelType>(name))
354  {
355  return mesh().lookupObject<modelType>(name);
356  }
357  else
358  {
359  return
360  mesh().lookupObject<modelType>
361  (
362  IOobject::groupName(modelType::typeName, key.otherName())
363  );
364  }
365 }
366 
367 
368 template<class modelType>
370 (
371  const phaseModel& dispersed,
372  const phaseModel& continuous
373 ) const
374 {
375  return foundSubModel<modelType>(orderedPhasePair(dispersed, continuous));
376 }
377 
378 
379 template<class modelType>
380 const modelType& Foam::phaseSystem::lookupSubModel
381 (
382  const phaseModel& dispersed,
383  const phaseModel& continuous
384 ) const
385 {
386  return lookupSubModel<modelType>(orderedPhasePair(dispersed, continuous));
387 }
388 
389 
390 template<class modelType>
391 bool Foam::phaseSystem::foundBlendedSubModel(const phasePair& key) const
392 {
393  if
394  (
395  mesh().foundObject<BlendedInterfacialModel<modelType>>
396  (
398  (
399  BlendedInterfacialModel<modelType>::typeName,
400  key.name()
401  )
402  )
403  || mesh().foundObject<BlendedInterfacialModel<modelType>>
404  (
406  (
407  BlendedInterfacialModel<modelType>::typeName,
408  key.otherName()
409  )
410  )
411  )
412  {
413  return true;
414  }
415  else
416  {
417  return false;
418  }
419 }
420 
421 
422 template<class modelType>
424 Foam::phaseSystem::lookupBlendedSubModel(const phasePair& key) const
425 {
426  const word name
427  (
429  (
430  BlendedInterfacialModel<modelType>::typeName,
431  key.name()
432  )
433  );
434 
435  if (mesh().foundObject<BlendedInterfacialModel<modelType>>(name))
436  {
437  return mesh().lookupObject<BlendedInterfacialModel<modelType>>(name);
438  }
439  else
440  {
441  return
442  mesh().lookupObject<BlendedInterfacialModel<modelType>>
443  (
445  (
446  BlendedInterfacialModel<modelType>::typeName,
447  key.otherName()
448  )
449  );
450  }
451 }
452 
453 
454 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
455 
456 namespace Foam
457 {
458 
459 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
460 
461 template<class GeoField, class Group>
462 inline void addField
463 (
464  const Group& group,
465  const word& name,
466  tmp<GeoField> field,
467  PtrList<GeoField>& fieldList
468 )
469 {
470  if (fieldList.set(group.index()))
471  {
472  fieldList[group.index()] += field;
473  }
474  else
475  {
476  fieldList.set
477  (
478  group.index(),
479  new GeoField
480  (
481  IOobject::groupName(name, group.name()),
482  field
483  )
484  );
485  }
486 }
487 
488 
489 template<class GeoField, class Group>
490 inline void addField
491 (
492  const Group& group,
493  const word& name,
494  const GeoField& field,
495  PtrList<GeoField>& fieldList
496 )
497 {
498  addField(group, name, tmp<GeoField>(field), fieldList);
499 }
500 
501 
502 template<class GeoField, class Group>
503 inline void addField
504 (
505  const Group& group,
506  const word& name,
507  tmp<GeoField> field,
508  HashPtrTable<GeoField>& fieldTable
509 )
510 {
511  if (fieldTable.found(group.name()))
512  {
513  *fieldTable[group.name()] += field;
514  }
515  else
516  {
517  fieldTable.set
518  (
519  group.name(),
520  new GeoField
521  (
522  IOobject::groupName(name, group.name()),
523  field
524  )
525  );
526  }
527 }
528 
529 
530 template<class GeoField, class Group>
531 inline void addField
532 (
533  const Group& group,
534  const word& name,
535  const GeoField& field,
536  HashPtrTable<GeoField>& fieldTable
537 )
538 {
539  addField(group, name, tmp<GeoField>(field), fieldTable);
540 }
541 
542 
543 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
544 
545 } // End namespace Foam
546 
547 // ************************************************************************* //
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:323
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:166
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:167
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:336
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
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:151
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
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:144
bool foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
void validateMassTransfer(const phasePair &key) const
Check that mass transfer is supported across the given interface.
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:178
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:113
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844