Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
32template<class modelType>
34(
35 const dictTable& modelDicts,
37 <
41 >& models
42)
43{
44 forAllConstIter(dictTable, modelDicts, iter)
45 {
46 const phasePairKey& key = iter.key();
47
48 models.insert
49 (
50 key,
51 modelType::New
52 (
53 *iter,
54 phasePairs_[key]
55 )
56 );
57 }
58}
59
60
61template<class modelType>
63(
64 const word& modelName,
65 HashTable
66 <
67 autoPtr<modelType>,
68 phasePairKey,
70 >& models
71)
72{
73 dictTable modelDicts(lookup(modelName));
74
75 generatePairs(modelDicts);
76
77 createSubModels(modelDicts, models);
78}
79
80
81template<class modelType>
83(
84 const word& modelName,
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 const autoPtr<modelType>& modelPtr = iter.val();
115
116 if (!modelPtr)
117 {
118 continue;
119 }
120
121 const phasePairKey key(iter.key().first(), iter.key().second());
122 const phasePairKey key1In2(key.first(), key.second(), true);
123 const phasePairKey key2In1(key.second(), key.first(), true);
124
125 models.insert
126 (
127 key,
128 autoPtr<BlendedInterfacialModel<modelType>>
129 (
130 new BlendedInterfacialModel<modelType>
131 (
132 phaseModels_[key.first()],
133 phaseModels_[key.second()],
134 blending,
135 tempModels.found(key ) ? tempModels[key ] : noModel,
136 tempModels.found(key1In2) ? tempModels[key1In2] : noModel,
137 tempModels.found(key2In1) ? tempModels[key2In1] : noModel,
138 correctFixedFluxBCs
139 )
140 )
141 );
142
143 if (!phasePairs_.found(key))
144 {
145 phasePairs_.insert
146 (
147 key,
148 autoPtr<phasePair>
149 (
150 new phasePair
151 (
152 phaseModels_[key.first()],
153 phaseModels_[key.second()]
154 )
155 )
156 );
157 }
158 }
159}
160
161
162template<class modelType>
164(
165 const word& modelName,
167 <
171 >& models,
172 const bool correctFixedFluxBCs
173)
174{
175 typedef
177 modelTypeTable;
178
179 forAll(phaseModels_, phasei)
180 {
181 const phaseModel& phase = phaseModels_[phasei];
182
183 modelTypeTable tempModels;
184 generatePairsAndSubModels
185 (
186 IOobject::groupName(modelName, phase.name()),
187 tempModels,
188 correctFixedFluxBCs
189 );
190
191 forAllIter(typename modelTypeTable, tempModels, tempModelIter)
192 {
193 const phasePairKey& key(tempModelIter.key());
194
195 if (!models.found(key))
196 {
197 models.insert
198 (
199 key,
201 );
202 }
203
204 const phasePair& pair = phasePairs_[key];
205
206 if (!pair.contains(phase))
207 {
209 << "A two-sided " << modelType::typeName << " was "
210 << "specified for the " << phase.name() << " side of the "
211 << pair << " pair, but that phase is not part of that pair."
212 << exit(FatalError);
213 }
214
215 models[key][pair.index(phase)].reset(tempModelIter().ptr());
216 }
217 }
218}
219
220
221template<class GeoField>
223(
224 const phaseModel& phase,
225 const word& fieldName,
227 PtrList<GeoField>& fieldList
228) const
229{
230 if (fieldList.set(phase.index()))
231 {
232 fieldList[phase.index()] += field;
233 }
234 else
235 {
236 fieldList.set
237 (
238 phase.index(),
239 new GeoField
240 (
241 IOobject::groupName(fieldName, phase.name()),
242 field
244 );
245 }
246}
247
248
249template<class GeoField>
251(
252 const phaseModel& phase,
253 const word& fieldName,
254 const GeoField& field,
255 PtrList<GeoField>& fieldList
256) const
257{
258 addField(phase, fieldName, tmp<GeoField>(field), fieldList);
259}
260
261
262template<class GeoField>
264(
265 const phaseModel& phase,
266 const word& fieldName,
268 HashPtrTable<GeoField>& fieldTable
269) const
270{
271 if (fieldTable.found(phase.name()))
272 {
273 *fieldTable[phase.name()] += field;
274 }
275 else
276 {
277 fieldTable.set
278 (
279 phase.name(),
280 new GeoField
281 (
282 IOobject::groupName(fieldName, phase.name()),
283 field
285 );
286 }
287}
288
289
290template<class GeoField>
292(
293 const phaseModel& phase,
294 const word& fieldName,
295 const GeoField& field,
296 HashPtrTable<GeoField>& fieldTable
297) const
298{
299 addField(phase, fieldName, tmp<GeoField>(field), fieldTable);
300}
301
302
303// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
304
305template<class Type, template<class> class PatchField, class GeoMesh>
307(
308 const word& name,
309 const dimensionSet& dims,
311) const
312{
313 forAll(this->phaseModels_, phasei)
314 {
315 if (fieldList.set(phasei))
316 {
317 continue;
318 }
319
320 const phaseModel& phase = this->phaseModels_[phasei];
321
322 fieldList.set
323 (
324 phasei,
325 new GeometricField<Type, PatchField, GeoMesh>
326 (
327 IOobject
328 (
329 IOobject::groupName(name, phase.name()),
330 this->mesh_.time().timeName(),
331 this->mesh_
332 ),
333 this->mesh_,
334 dimensioned<Type>(dims, Zero)
336 );
337 }
338}
339
340
341template<class Type, template<class> class PatchField, class GeoMesh>
343(
344 const word& name,
345 const dimensionSet& dims,
346 HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>& fieldTable
347) const
348{
349 forAll(this->phaseModels_, phasei)
350 {
351 const phaseModel& phase = this->phaseModels_[phasei];
352
353 if (fieldTable.set(phase.name()))
354 {
355 continue;
356 }
357
358 fieldTable.set
359 (
360 phase.name(),
362 (
364 (
366 this->mesh_.time().timeName(),
367 this->mesh_
368 ),
369 this->mesh_,
372 );
373 }
374}
375
376
377template<class modelType>
378bool Foam::phaseSystem::foundSubModel(const phasePair& key) const
379{
380 const word name(IOobject::groupName(modelType::typeName, key.name()));
381
382 if (key.ordered())
383 {
384 if (mesh().foundObject<modelType>(name))
385 {
386 return true;
387 }
388 else
389 {
390 return false;
391 }
392 }
393 else
394 {
395 if
396 (
397 mesh().foundObject<modelType>(name)
398 ||
399 mesh().foundObject<modelType>
400 (
401 IOobject::groupName(modelType::typeName, key.otherName())
402 )
403 )
404 {
405 return true;
406 }
407 else
408 {
409 return false;
410 }
411 }
412}
413
414
415template<class modelType>
416const modelType& Foam::phaseSystem::lookupSubModel(const phasePair& key) const
417{
418 const word name(IOobject::groupName(modelType::typeName, key.name()));
419
420 if (key.ordered() || mesh().foundObject<modelType>(name))
421 {
422 return mesh().lookupObject<modelType>(name);
423 }
424 else
425 {
426 return
427 mesh().lookupObject<modelType>
428 (
429 IOobject::groupName(modelType::typeName, key.otherName())
430 );
431 }
432}
433
434
435template<class modelType>
437(
438 const phaseModel& dispersed,
439 const phaseModel& continuous
440) const
441{
442 return foundSubModel<modelType>(orderedPhasePair(dispersed, continuous));
443}
444
445
446template<class modelType>
448(
449 const phaseModel& dispersed,
450 const phaseModel& continuous
451) const
452{
453 return lookupSubModel<modelType>(orderedPhasePair(dispersed, continuous));
454}
455
456
457template<class modelType>
459{
460 if
461 (
463 (
465 (
467 key.name()
468 )
469 )
471 (
473 (
475 key.otherName()
476 )
477 )
478 )
479 {
480 return true;
481 }
482 else
483 {
484 return false;
485 }
486}
487
488
489template<class modelType>
491Foam::phaseSystem::lookupBlendedSubModel(const phasePair& key) const
492{
493 const word name
494 (
496 (
497 BlendedInterfacialModel<modelType>::typeName,
498 key.name()
499 )
500 );
501
502 if (mesh().foundObject<BlendedInterfacialModel<modelType>>(name))
503 {
504 return mesh().lookupObject<BlendedInterfacialModel<modelType>>(name);
505 }
506 else
507 {
508 return
510 (
512 (
514 key.otherName()
515 )
516 );
517 }
518}
519
520
521// ************************************************************************* //
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
Generic GeometricField class.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word member(const word &name)
Return member (name without the extension).
Definition IOobject.C:313
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Generic dimensioned Type class.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const word & name() const
Definition phaseModel.H:166
An ordered or unorder pair of phase names. Typically specified as follows.
phasePairKey::hasher hash
Alternative name for functor.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
Definition phasePairI.H:65
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Definition phasePairI.H:35
void createSubModels(const dictTable &modelDicts, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
void generatePairs(const dictTable &modelDicts)
Generate pairs.
Definition phaseSystem.C:67
phaseModelList phaseModels_
Phase models.
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.
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
const BlendedInterfacialModel< modelType > & lookupBlendedSubModel(const phasePair &key) const
Return a blended sub model between a phase pair.
bool foundBlendedSubModel(const phasePair &key) const
Check availability of a blended sub model for a given phase pair.
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
bool foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
phasePairTable phasePairs_
Phase pairs.
const fvMesh & mesh() const
Return the mesh.
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.
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Definition phaseSystem.H:98
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
const word & name() const
Definition phase.H:114
Lookup type of boundary radiation properties.
Definition lookup.H:60
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
label phasei
Definition pEqn.H:27
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:334
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:356