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-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"
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 )
70 {
71  dictTable modelDicts(lookup(modelName));
72 
73  generatePairs(modelDicts);
74 
75  createSubModels(modelDicts, models);
76 }
77 
78 
79 template<class modelType>
81 (
82  const word& modelName,
83  HashTable
84  <
85  autoPtr<BlendedInterfacialModel<modelType>>,
86  phasePairKey,
87  phasePairKey::hash
88  >& models,
89  const bool correctFixedFluxBCs
90 )
91 {
92  typedef
93  HashTable<autoPtr<modelType>, phasePairKey, phasePairKey::hash>
94  modelTypeTable;
95 
96  modelTypeTable tempModels;
97  generatePairsAndSubModels(modelName, tempModels);
98 
99  const blendingMethod& blending
100  (
101  blendingMethods_.found(modelName)
102  ? blendingMethods_[modelName]
103  : blendingMethods_.found(member(modelName))
104  ? blendingMethods_[member(modelName)]
105  : blendingMethods_["default"]
106  );
107 
108  autoPtr<modelType> noModel(nullptr);
109 
110  forAllConstIter(typename modelTypeTable, tempModels, iter)
111  {
112  if (!iter().valid())
113  {
114  continue;
115  }
116 
117  const phasePairKey key(iter.key().first(), iter.key().second());
118  const phasePairKey key1In2(key.first(), key.second(), true);
119  const phasePairKey key2In1(key.second(), key.first(), true);
120 
121  models.insert
122  (
123  key,
124  autoPtr<BlendedInterfacialModel<modelType>>
125  (
126  new BlendedInterfacialModel<modelType>
127  (
128  phaseModels_[key.first()],
129  phaseModels_[key.second()],
130  blending,
131  tempModels.found(key ) ? tempModels[key ] : noModel,
132  tempModels.found(key1In2) ? tempModels[key1In2] : noModel,
133  tempModels.found(key2In1) ? tempModels[key2In1] : noModel,
134  correctFixedFluxBCs
135  )
136  )
137  );
138 
139  if (!phasePairs_.found(key))
140  {
142  (
143  key,
144  autoPtr<phasePair>
145  (
146  new phasePair
147  (
148  phaseModels_[key.first()],
149  phaseModels_[key.second()]
150  )
151  )
152  );
153  }
154  }
155 }
156 
157 
158 template<class modelType>
160 (
161  const word& modelName,
162  HashTable
163  <
164  Pair<autoPtr<modelType>>,
165  phasePairKey,
166  phasePairKey::hash
167  >& models,
168  const bool correctFixedFluxBCs
169 )
170 {
171  typedef
172  HashTable<autoPtr<modelType>, phasePairKey, phasePairKey::hash>
173  modelTypeTable;
174 
176  {
177  const phaseModel& phase = phaseModels_[phasei];
178 
179  modelTypeTable tempModels;
181  (
182  IOobject::groupName(modelName, phase.name()),
183  tempModels,
184  correctFixedFluxBCs
185  );
186 
187  forAllIter(typename modelTypeTable, tempModels, tempModelIter)
188  {
189  const phasePairKey& key(tempModelIter.key());
190 
191  if (!models.found(key))
192  {
193  models.insert
194  (
195  key,
196  Pair<autoPtr<modelType>>()
197  );
198  }
199 
200  const phasePair& pair = phasePairs_[key];
201 
202  if (!pair.contains(phase))
203  {
205  << "A two-sided " << modelType::typeName << " was "
206  << "specified for the " << phase.name() << " side of the "
207  << pair << " pair, but that phase is not part of that pair."
208  << exit(FatalError);
209  }
210 
211  models[key][pair.index(phase)] = tempModelIter().ptr();
212  }
213  }
214 }
215 
216 
217 template<class GeoField>
219 (
220  const phaseModel& phase,
221  const word& fieldName,
222  tmp<GeoField> field,
223  PtrList<GeoField>& fieldList
224 ) const
225 {
226  if (fieldList.set(phase.index()))
227  {
228  fieldList[phase.index()] += field;
229  }
230  else
231  {
232  fieldList.set
233  (
234  phase.index(),
235  new GeoField
236  (
237  IOobject::groupName(fieldName, phase.name()),
238  field
239  )
240  );
241  }
242 }
243 
244 
245 template<class GeoField>
247 (
248  const phaseModel& phase,
249  const word& fieldName,
250  const GeoField& field,
251  PtrList<GeoField>& fieldList
252 ) const
253 {
254  addField(phase, fieldName, tmp<GeoField>(field), fieldList);
255 }
256 
257 
258 template<class GeoField>
260 (
261  const phaseModel& phase,
262  const word& fieldName,
263  tmp<GeoField> field,
264  HashPtrTable<GeoField>& fieldTable
265 ) const
266 {
267  if (fieldTable.found(phase.name()))
268  {
269  *fieldTable[phase.name()] += field;
270  }
271  else
272  {
273  fieldTable.set
274  (
275  phase.name(),
276  new GeoField
277  (
278  IOobject::groupName(fieldName, phase.name()),
279  field
280  )
281  );
282  }
283 }
284 
285 
286 template<class GeoField>
288 (
289  const phaseModel& phase,
290  const word& fieldName,
291  const GeoField& field,
292  HashPtrTable<GeoField>& fieldTable
293 ) const
294 {
295  addField(phase, fieldName, tmp<GeoField>(field), fieldTable);
296 }
297 
298 
299 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
300 
301 template<class Type, template<class> class PatchField, class GeoMesh>
303 (
304  const word& name,
305  const dimensionSet& dims,
306  PtrList<GeometricField<Type, PatchField, GeoMesh>>& fieldList
307 ) const
308 {
309  forAll(this->phaseModels_, phasei)
310  {
311  if (fieldList.set(phasei))
312  {
313  continue;
314  }
315 
316  const phaseModel& phase = this->phaseModels_[phasei];
317 
318  fieldList.set
319  (
320  phasei,
321  new GeometricField<Type, PatchField, GeoMesh>
322  (
323  IOobject
324  (
325  IOobject::groupName(name, phase.name()),
326  this->mesh_.time().timeName(),
327  this->mesh_
328  ),
329  this->mesh_,
330  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
331  )
332  );
333  }
334 }
335 
336 
337 template<class Type, template<class> class PatchField, class GeoMesh>
339 (
340  const word& name,
341  const dimensionSet& dims,
342  HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>& fieldTable
343 ) const
344 {
345  forAll(this->phaseModels_, phasei)
346  {
347  const phaseModel& phase = this->phaseModels_[phasei];
348 
349  if (fieldTable.set(phase.name()))
350  {
351  continue;
352  }
353 
354  fieldTable.set
355  (
356  phase.name(),
357  new GeometricField<Type, PatchField, GeoMesh>
358  (
359  IOobject
360  (
361  IOobject::groupName(name, phase.name()),
362  this->mesh_.time().timeName(),
363  this->mesh_
364  ),
365  this->mesh_,
366  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
367  )
368  );
369  }
370 }
371 
372 
373 template<class modelType>
374 const modelType& Foam::phaseSystem::lookupSubModel(const phasePair& key) const
375 {
376  const word name(IOobject::groupName(modelType::typeName, key.name()));
377 
378  if (key.ordered() || mesh().foundObject<modelType>(name))
379  {
380  return mesh().lookupObject<modelType>(name);
381  }
382  else
383  {
384  return
385  mesh().lookupObject<modelType>
386  (
387  IOobject::groupName(modelType::typeName, key.otherName())
388  );
389  }
390 }
391 
392 
393 template<class modelType>
394 const modelType& Foam::phaseSystem::lookupSubModel
395 (
396  const phaseModel& dispersed,
397  const phaseModel& continuous
398 ) const
399 {
400  return lookupSubModel<modelType>(orderedPhasePair(dispersed, continuous));
401 }
402 
403 
404 template<class modelType>
406 Foam::phaseSystem::lookupBlendedSubModel(const phasePair& key) const
407 {
408  const word name
409  (
411  (
412  BlendedInterfacialModel<modelType>::typeName,
413  key.name()
414  )
415  );
416 
417  if (mesh().foundObject<BlendedInterfacialModel<modelType>>(name))
418  {
419  return mesh().lookupObject<BlendedInterfacialModel<modelType>>(name);
420  }
421  else
422  {
423  return
424  mesh().lookupObject<BlendedInterfacialModel<modelType>>
425  (
427  (
428  BlendedInterfacialModel<modelType>::typeName,
429  key.otherName()
430  )
431  );
432  }
433 }
434 
435 
436 // ************************************************************************* //
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:428
bool found(const word &) const
Search DictionaryBase for given keyword.
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:453
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:142
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:209
bool foundObject(const word &name) const
Is the named Type found?
word member() const
Return member (name without the extension)
Definition: IOobject.C:385
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
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:243
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:127
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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
static word groupName(Name name, const word &group)
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:124
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.
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
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,.
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:154
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition: phaseSystem.H:96
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576