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-2023 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 "phaseSystem.H"
27 #include "sidedPhaseInterface.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 template<class Type>
36 {
37  public:
38 
40  Type t;
41 
43  {}
44 
46  :
47  wl(is),
48  t(is)
49  {}
50 };
51 
52 template<class Type>
54 {
55  return is >> wlat.wl >> wlat.t;
56 }
57 
58 template<class Type>
60 {
61  return os << wlat.wl << wlat.t;
62 }
63 
64 template<class Type>
65 inline bool operator==
66 (
67  const wordListAndType<Type>& a,
69 )
70 {
71  return a.wl == b.wl && a.t == b.t;
72 }
73 
74 template<class Type>
75 inline bool operator!=
76 (
77  const wordListAndType<Type>& a,
79 )
80 {
81  return !(a == b);
82 }
83 
84 }
85 
86 
87 template<class Type>
89 {
90  bool found = false;
92 
93  // If it is a dictionary then merge it in
94  if (this->isDict(name))
95  {
96  found = true;
97  dict.merge(this->subDict(name));
98  }
99 
100  // Function to add old-format list/table entries
101  auto add = [&](const word& sidePhaseName)
102  {
103  const word nameSidePhaseName =
104  IOobject::groupName(name, sidePhaseName);
105 
106  if (!this->found(nameSidePhaseName)) return;
107 
108  found = true;
109 
110  List<wordListAndType<Type>> wlats(this->lookup(nameSidePhaseName));
111 
112  forAll(wlats, i)
113  {
114  word keyword =
115  phaseInterface::oldNamePartsToName(*this, wlats[i].wl);
116 
117  if (sidePhaseName != word::null)
118  {
119  keyword.append
120  (
121  "_"
123  + "_"
124  + sidePhaseName
125  );
126  }
127 
128  dict.add(keyword, wlats[i].t);
129  }
130  };
131 
132  // If a dictionary entry was not found then try a list/table entry
133  if (!found)
134  {
135  add(word::null);
136  }
137 
138  // Add sided models
139  forAll(phases(), phasei)
140  {
141  add(phases()[phasei].name());
142  }
143 
144  // Barf if nothing was found
145  if (!found)
146  {
147  return this->subDict(name);
148  }
149 
150  return dict;
151 }
152 
153 
154 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
155 
156 template<class Type, template<class> class PatchField, class GeoMesh>
158 (
159  const word& name,
160  const dimensionSet& dims,
162 ) const
163 {
164  forAll(this->phaseModels_, phasei)
165  {
166  if (fieldList.set(phasei))
167  {
168  continue;
169  }
170 
171  const phaseModel& phase = this->phaseModels_[phasei];
172 
173  fieldList.set
174  (
175  phasei,
177  (
178  IOobject
179  (
180  IOobject::groupName(name, phase.name()),
181  this->mesh_.time().name(),
182  this->mesh_
183  ),
184  this->mesh_,
186  )
187  );
188  }
189 }
190 
191 
192 template<class Type, template<class> class PatchField, class GeoMesh>
194 (
195  const word& name,
196  const dimensionSet& dims,
198 ) const
199 {
200  forAll(this->phaseModels_, phasei)
201  {
202  const phaseModel& phase = this->phaseModels_[phasei];
203 
204  if (fieldTable.set(phase.name()))
205  {
206  continue;
207  }
208 
209  fieldTable.set
210  (
211  phase.name(),
213  (
214  IOobject
215  (
216  IOobject::groupName(name, phase.name()),
217  this->mesh_.time().name(),
218  this->mesh_
219  ),
220  this->mesh_,
222  )
223  );
224  }
225 }
226 
227 
228 template<class ModelType>
230 {
231  word name = ModelType::typeName;
232 
233  // Extract the innermost part of the template
234  const word::size_type i0 = name.find_last_of('<');
235  if (i0 != word::npos)
236  {
237  const word::size_type i1 = name.find_first_of('>', i0 + 1);
238  if (i1 != word::npos)
239  {
240  name = name(i0 + 1, i1 - i0 - 1);
241  }
242  }
243 
244  // Strip "Model" off the end of the name
245  if (name(name.size() - 5, 5) == "Model")
246  {
247  name = name(name.size() - 5);
248  }
249 
250  return name;
251 }
252 
253 
254 template<class ModelType, class ... InterfaceTypes>
256 (
257  const dictionary& dict,
258  const phaseInterface& interface,
259  PtrList<phaseInterface>& interfaces,
260  PtrList<ModelType>& models
261 ) const
262 {
263  // Construct sub-dictionaries and associated interfaces
264  hashedWordList names;
265  PtrList<dictionary> dicts;
267  {
268  // Get the model sub dictionary and its associated interface
269  const dictionary& modelDict = iter().dict();
270  autoPtr<phaseInterface> modelInterfacePtr =
271  phaseInterface::New(*this, iter().keyword());
272 
273  // Cast the interface down to the first specified type possible
274  autoPtr<phaseInterface> interfacePtr;
275  List<bool>
276  ({
277  interfacePtr.empty()
278  && isA<InterfaceTypes>(modelInterfacePtr())
279  && (
280  interfacePtr.set
281  (
282  new InterfaceTypes
283  (
284  refCast<InterfaceTypes>(modelInterfacePtr())
285  )
286  ),
287  true
288  )
289  ...
290  });
291  if (!interfacePtr.valid())
292  {
294  << "Interface " << modelInterfacePtr->name()
295  << " is not of suitable type for construction of a "
296  << ModelType::typeName
297  << exit(FatalError);
298  }
299 
300  // If constructing for a specific interface then combine with this
301  // interface. This ensures interface information propagates through
302  // hierarchical model generation.
303  if (notNull(interface))
304  {
305  interfacePtr = phaseInterface::New(interface, interfacePtr());
306  }
307 
308  // Find an existing dictionary to add to or create a new one
309  const word name = interfacePtr->name();
310  if (!names.found(name))
311  {
312  names.append(name);
313  dicts.append(new dictionary(name));
314  interfaces.append(interfacePtr.ptr());
315  models.append(nullptr);
316  }
317 
318  // Add the model dictionary
319  dicts[names[name]].add
320  (
321  modelInterfacePtr->name(),
322  modelDict
323  );
324  }
325 
326  // Construct the models
327  forAll(interfaces, i)
328  {
329  models.set(i, ModelType::New(dicts[i], interfaces[i]));
330  }
331 }
332 
333 
334 template<class ModelType>
336 (
337  const dictionary& dict,
338  HashTable
339  <
343  >& models
344 ) const
345 {
346  // Construct lists of interfaces and models
347  PtrList<phaseInterface> listInterfaces;
348  PtrList<ModelType> listModels;
349  generateInterfacialModels<ModelType, phaseInterface>
350  (
351  dict,
352  NullObjectRef<phaseInterface>(),
353  listInterfaces,
354  listModels
355  );
356 
357  // Transfer to a keyed table
358  forAll(listInterfaces, i)
359  {
360  models.insert(listInterfaces[i], listModels.set(i, nullptr));
361  }
362 }
363 
364 
365 
366 template<class ModelType>
368 (
369  HashTable
370  <
374  >& models
375 ) const
376 {
378  (
379  interfacialDict<dictionary>(modelName<ModelType>()),
380  models
381  );
382 }
383 
384 
385 template<class ValueType>
387 (
388  const dictionary& dict,
390 ) const
391 {
393  {
394  autoPtr<phaseInterface> interfacePtr =
395  phaseInterface::New(*this, iter().keyword());
396 
397  const ValueType value(pTraits<ValueType>(iter().stream()));
398 
399  values.insert(interfacePtr(), value);
400  }
401 }
402 
403 
404 template<class ValueType>
406 (
407  const word& valueName,
409 ) const
410 {
411  generateInterfacialValues(interfacialDict<ValueType>(valueName), values);
412 }
413 
414 
415 template<class ModelType>
417 (
418  const dictionary& dict
419 )
420 {
421  if (dict.size() != 1)
422  {
424  << "Too many matching entries for construction of a "
425  << ModelType::typeName << nl << dict.toc()
426  << exit(FatalError);
427  }
428 
429  if (!dict.first()->isDict())
430  {
432  << "Non-sub-dictionary entries found for specification of a "
433  << ModelType::typeName
434  << exit(FatalError);
435  }
436 
437  return dict.first()->dict();
438 }
439 
440 
441 template<class ModelType>
443 (
444  const phaseInterface& interface
445 ) const
446 {
447  if (interface.phase1().stationary() || interface.phase2().stationary())
448  {
450  << "A " << ModelType::typeName << " was specified for pair "
451  << interface.name() << ", but one of these phases is stationary. "
452  << "Mass transfer is not supported on stationary phases"
453  << exit(FatalError);
454  }
455 }
456 
457 
458 template<class ModelType>
460 (
461  const phaseInterface& interface
462 ) const
463 {
464  return
465  mesh().foundObject<ModelType>
466  (
467  IOobject::groupName(ModelType::typeName, interface.name())
468  );
469 }
470 
471 
472 template<class ModelType>
474 (
475  const phaseInterface& interface
476 ) const
477 {
478  return
479  mesh().lookupObject<ModelType>
480  (
481  IOobject::groupName(ModelType::typeName, interface.name())
482  );
483 }
484 
485 
486 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
487 
488 namespace Foam
489 {
490 
491 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
492 
493 template<class GeoField, class Group>
494 inline void addField
495 (
496  const Group& group,
497  const word& name,
498  tmp<GeoField> field,
499  PtrList<GeoField>& fieldList
500 )
501 {
502  if (fieldList.set(group.index()))
503  {
504  fieldList[group.index()] += field;
505  }
506  else
507  {
508  fieldList.set
509  (
510  group.index(),
511  new GeoField
512  (
513  IOobject::groupName(name, group.name()),
514  field
515  )
516  );
517  }
518 }
519 
520 
521 template<class GeoField, class Group>
522 inline void addField
523 (
524  const Group& group,
525  const word& name,
526  const GeoField& field,
527  PtrList<GeoField>& fieldList
528 )
529 {
530  addField(group, name, tmp<GeoField>(field), fieldList);
531 }
532 
533 
534 template<class GeoField, class Group>
535 inline void addField
536 (
537  const Group& group,
538  const word& name,
539  tmp<GeoField> field,
540  HashPtrTable<GeoField>& fieldTable
541 )
542 {
543  if (fieldTable.found(group.name()))
544  {
545  *fieldTable[group.name()] += field;
546  }
547  else
548  {
549  fieldTable.set
550  (
551  group.name(),
552  new GeoField
553  (
554  IOobject::groupName(name, group.name()),
555  field
556  )
557  );
558  }
559 }
560 
561 
562 template<class GeoField, class Group>
563 inline void addField
564 (
565  const Group& group,
566  const word& name,
567  const GeoField& field,
568  HashPtrTable<GeoField>& fieldTable
569 )
570 {
571  addField(group, name, tmp<GeoField>(field), fieldTable);
572 }
573 
574 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
575 
576 } // End namespace Foam
577 
578 // ************************************************************************* //
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
An STL-conforming hash table.
Definition: HashTable.H:127
bool set(const Key &, const T &newElmt)
Assign a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
T * first()
Return the first entry.
Definition: UILList.H:109
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
bool empty() const
Return true if the autoPtr is empty (ie, no pointer set)
Definition: autoPtrI.H:76
T * ptr()
Return object pointer for reuse.
Definition: autoPtrI.H:90
void set(T *)
Set pointer to that given.
Definition: autoPtrI.H:99
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:952
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1169
wordList toc() const
Return the table of contents.
Definition: dictionary.C:1131
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
Dimension set for the base types.
Definition: dimensionSet.H:122
Generic dimensioned Type class.
A wordList with hashed indices for faster lookup by name.
bool found(const word &) const
Does the list contain the specified name.
void append(const word &)
Append an element at the end of the list.
Traits class for primitives.
Definition: pTraits.H:53
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
static word oldNamePartsToName(const phaseSystem &fluid, const wordList &oldNameParts)
Convert old-format interface name parts to an interface name. Used.
virtual word name() const
Name.
static autoPtr< phaseInterface > New(const phaseSystem &fluid, const word &name)
Select given fluid and name.
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual bool stationary() const =0
Return whether the phase is stationary.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
dictionary interfacialDict(const word &name) const
Return the dictionary containing interfacial model or value.
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
void validateMassTransfer(const phaseInterface &interface) const
Check that mass transfer is supported across the given interface.
void generateInterfacialModels(const dictionary &dict, const phaseInterface &interface, PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models) const
Generate interfacial-model lists.
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
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.
word modelName() const
Return the model name. This is the same as the model's typename.
void generateInterfacialValues(const dictionary &dict, HashTable< ValueType, phaseInterfaceKey, phaseInterfaceKey::hash > &values) const
Generate interfacial-model tables.
void generateInterfacialModels(HashTable< autoPtr< ModelType >, phaseInterfaceKey, phaseInterfaceKey::hash > &models) const
Generate interfacial-model tables.
static const dictionary & modelSubDict(const dictionary &dict)
Return the dictionary from which to construct a low-level.
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
static word separator()
Return the separator that delimits this interface's name.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
volScalarField & b
Definition: createFields.H:27
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const char *const group
Group name for atomic constants.
Namespace for OpenFOAM.
void addField(const Group &group, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:58
Istream & operator>>(Istream &, directionInfo &)
error FatalError
Ostream & operator<<(Ostream &, const ensightPart &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
dictionary dict