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 | Copyright (C) 2015-2017 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 )
90 {
91  typedef
92  HashTable<autoPtr<modelType>, phasePairKey, phasePairKey::hash>
93  modelTypeTable;
94 
95  modelTypeTable tempModels;
96  generatePairsAndSubModels(modelName, tempModels);
97 
98  const blendingMethod& blending
99  (
100  blendingMethods_.found(modelName)
101  ? blendingMethods_[modelName]
102  : blendingMethods_["default"]
103  );
104 
105  autoPtr<modelType> noModel(nullptr);
106 
107  forAllConstIter(typename modelTypeTable, tempModels, iter)
108  {
109  if (!iter().valid())
110  {
111  continue;
112  }
113 
114  const phasePairKey key(iter.key().first(), iter.key().second());
115  const phasePairKey key1In2(key.first(), key.second(), true);
116  const phasePairKey key2In1(key.second(), key.first(), true);
117 
118  models.insert
119  (
120  key,
121  autoPtr<BlendedInterfacialModel<modelType>>
122  (
123  new BlendedInterfacialModel<modelType>
124  (
125  phaseModels_[key.first()],
126  phaseModels_[key.second()],
127  blending,
128  tempModels.found(key ) ? tempModels[key ] : noModel,
129  tempModels.found(key1In2) ? tempModels[key1In2] : noModel,
130  tempModels.found(key2In1) ? tempModels[key2In1] : noModel
131  )
132  )
133  );
134 
135  if (!phasePairs_.found(key))
136  {
138  (
139  key,
140  autoPtr<phasePair>
141  (
142  new phasePair
143  (
144  phaseModels_[key.first()],
145  phaseModels_[key.second()]
146  )
147  )
148  );
149  }
150  }
151 }
152 
153 
154 template<class modelType>
156 (
157  const word& modelName,
158  HashTable
159  <
160  HashTable<autoPtr<modelType>>,
161  phasePairKey,
162  phasePairKey::hash
163  >& models
164 )
165 {
166  typedef
167  HashTable<autoPtr<modelType>, phasePairKey, phasePairKey::hash>
168  modelTypeTable;
169 
171  {
172  modelTypeTable tempModels;
174  (
176  tempModels
177  );
178 
179  forAllConstIter(typename modelTypeTable, tempModels, tempModelIter)
180  {
181  const phasePairKey key(tempModelIter.key());
182 
183  if (!models.found(key))
184  {
185  models.insert
186  (
187  key,
188  HashTable<autoPtr<modelType>>()
189  );
190  }
191 
192  models[tempModelIter.key()].insert
193  (
195  *tempModelIter
196  );
197  }
198  }
199 }
200 
201 template<class modelType>
202 const modelType& Foam::phaseSystem::lookupSubModel(const phasePair& key) const
203 {
204  return
205  mesh().lookupObject<modelType>
206  (
207  IOobject::groupName(modelType::typeName, key.name())
208  );
209 }
210 
211 
212 template<class modelType>
213 const modelType& Foam::phaseSystem::lookupSubModel
214 (
215  const phaseModel& dispersed,
216  const phaseModel& continuous
217 ) const
218 {
219  return lookupSubModel<modelType>(orderedPhasePair(dispersed, continuous));
220 }
221 
222 
223 // ************************************************************************* //
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.
label phasei
Definition: pEqn.H:27
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:167
const modelType & lookupSubModel(const phasePair &key) const
Access a sub model between a phase pair.
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:164
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 word & name() const
Name function is needed to disambiguate those inherited.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
blendingMethodTable blendingMethods_
Blending methods.
Definition: phaseSystem.H:179
const fvMesh & mesh() const
Constant access the mesh.
Definition: phaseSystemI.H:28
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition: phaseSystem.H:129
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576