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-2022 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 template<class ModelType>
33 (
34  const phaseInterface& interface
35 ) const
36 {
37  if (interface.phase1().stationary() || interface.phase2().stationary())
38  {
40  << "A " << ModelType::typeName << " was specified for pair "
41  << interface.name() << ", but one of these phases is stationary. "
42  << "Mass transfer is not supported on stationary phases"
43  << exit(FatalError);
44  }
45 }
46 
47 
48 namespace Foam
49 {
50 
51 template<class Type>
52 class wordListAndType
53 {
54  public:
55 
56  wordList wl;
57  Type t;
58 
59  wordListAndType()
60  {}
61 
62  wordListAndType(Istream& is)
63  :
64  wl(is),
65  t(is)
66  {}
67 };
68 
69 template<class Type>
70 inline Istream& operator>>(Istream& is, wordListAndType<Type>& wlat)
71 {
72  return is >> wlat.wl >> wlat.t;
73 }
74 
75 template<class Type>
76 inline Ostream& operator<<(Ostream& os, const wordListAndType<Type>& wlat)
77 {
78  return os << wlat.wl << wlat.t;
79 }
80 
81 template<class Type>
82 inline bool operator==
83 (
84  const wordListAndType<Type>& a,
85  const wordListAndType<Type>& b
86 )
87 {
88  return a.wl == b.wl && a.t == b.t;
89 }
90 
91 template<class Type>
92 inline bool operator!=
93 (
94  const wordListAndType<Type>& a,
95  const wordListAndType<Type>& b
96 )
97 {
98  return !(a == b);
99 }
100 
101 }
102 
103 
104 template<class Type>
106 {
107  bool found = false;
108  dictionary dict(name);
109 
110  // If it is a dictionary then merge it in
111  if (this->isDict(name))
112  {
113  found = true;
114  dict.merge(this->subDict(name));
115  }
116 
117  // Function to add old-format list/table entries
118  auto add = [&](const word& sidePhaseName)
119  {
120  const word nameSidePhaseName =
121  IOobject::groupName(name, sidePhaseName);
122 
123  if (!this->found(nameSidePhaseName)) return;
124 
125  found = true;
126 
127  List<wordListAndType<Type>> wlats(this->lookup(nameSidePhaseName));
128 
129  forAll(wlats, i)
130  {
131  word keyword =
132  phaseInterface::oldNamePartsToName(*this, wlats[i].wl);
133 
134  if (sidePhaseName != word::null)
135  {
136  keyword.append
137  (
138  "_"
140  + "_"
141  + sidePhaseName
142  );
143  }
144 
145  dict.add(keyword, wlats[i].t);
146  }
147  };
148 
149  // If a dictionary entry was not found then try a list/table entry
150  if (!found)
151  {
152  add(word::null);
153  }
154 
155  // Add sided models
156  forAll(phases(), phasei)
157  {
158  add(phases()[phasei].name());
159  }
160 
161  // Barf if nothing was found
162  if (!found)
163  {
164  return this->subDict(name);
165  }
166 
167  return dict;
168 }
169 
170 
171 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
172 
173 template<class Type, template<class> class PatchField, class GeoMesh>
175 (
176  const word& name,
177  const dimensionSet& dims,
178  PtrList<GeometricField<Type, PatchField, GeoMesh>>& fieldList
179 ) const
180 {
181  forAll(this->phaseModels_, phasei)
182  {
183  if (fieldList.set(phasei))
184  {
185  continue;
186  }
187 
188  const phaseModel& phase = this->phaseModels_[phasei];
189 
190  fieldList.set
191  (
192  phasei,
193  new GeometricField<Type, PatchField, GeoMesh>
194  (
195  IOobject
196  (
197  IOobject::groupName(name, phase.name()),
198  this->mesh_.time().timeName(),
199  this->mesh_
200  ),
201  this->mesh_,
202  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
203  )
204  );
205  }
206 }
207 
208 
209 template<class Type, template<class> class PatchField, class GeoMesh>
211 (
212  const word& name,
213  const dimensionSet& dims,
214  HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>& fieldTable
215 ) const
216 {
217  forAll(this->phaseModels_, phasei)
218  {
219  const phaseModel& phase = this->phaseModels_[phasei];
220 
221  if (fieldTable.set(phase.name()))
222  {
223  continue;
224  }
225 
226  fieldTable.set
227  (
228  phase.name(),
229  new GeometricField<Type, PatchField, GeoMesh>
230  (
231  IOobject
232  (
233  IOobject::groupName(name, phase.name()),
234  this->mesh_.time().timeName(),
235  this->mesh_
236  ),
237  this->mesh_,
238  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
239  )
240  );
241  }
242 }
243 
244 
245 template<class ModelType>
247 {
248  word name = ModelType::typeName;
249 
250  // Extract the innermost part of the template
251  const word::size_type i0 = name.find_last_of('<');
252  if (i0 != word::npos)
253  {
254  const word::size_type i1 = name.find_first_of('>', i0 + 1);
255  if (i1 != word::npos)
256  {
257  name = name(i0 + 1, i1 - i0 - 1);
258  }
259  }
260 
261  // Strip "Model" off the end of the name
262  if (name(name.size() - 5, 5) == "Model")
263  {
264  name = name(name.size() - 5);
265  }
266 
267  return name;
268 }
269 
270 
271 template<class ModelType, class ... InterfaceTypes>
273 (
274  const dictionary& dict,
275  const phaseInterface& interface,
276  PtrList<phaseInterface>& interfaces,
277  PtrList<ModelType>& models
278 ) const
279 {
280  // Construct sub-dictionaries and associated interfaces
281  hashedWordList names;
282  PtrList<dictionary> dicts;
283  forAllConstIter(dictionary, dict, iter)
284  {
285  // Get the model sub dictionary and its associated interface
286  const dictionary& modelDict = iter().dict();
287  autoPtr<phaseInterface> modelInterfacePtr =
288  phaseInterface::New(*this, iter().keyword());
289 
290  // Cast the interface down to the first specified type possible
291  autoPtr<phaseInterface> interfacePtr;
292  List<bool>
293  ({
294  interfacePtr.empty()
295  && isA<InterfaceTypes>(modelInterfacePtr())
296  && (
297  interfacePtr.set
298  (
299  new InterfaceTypes
300  (
301  refCast<InterfaceTypes>(modelInterfacePtr())
302  )
303  ),
304  true
305  )
306  ...
307  });
308  if (!interfacePtr.valid())
309  {
311  << "Interface " << modelInterfacePtr->name()
312  << " is not of suitable type for construction of a "
313  << ModelType::typeName
314  << exit(FatalError);
315  }
316 
317  // If constructing for a specific interface then combine with this
318  // interface. This ensures interface information propagates through
319  // hierarchical model generation.
320  if (notNull(interface))
321  {
322  interfacePtr = phaseInterface::New(interface, interfacePtr());
323  }
324 
325  // Find an existing dictionary to add to or create a new one
326  const word name = interfacePtr->name();
327  if (!names.found(name))
328  {
329  names.append(name);
330  dicts.append(new dictionary(name));
331  interfaces.append(interfacePtr.ptr());
332  models.append(nullptr);
333  }
334 
335  // Add the model dictionary
336  dicts[names[name]].add
337  (
338  modelInterfacePtr->name(),
339  modelDict
340  );
341  }
342 
343  // Construct the models
344  forAll(interfaces, i)
345  {
346  models.set(i, ModelType::New(dicts[i], interfaces[i]));
347  }
348 }
349 
350 
351 template<class ModelType>
353 (
354  const dictionary& dict,
355  HashTable
356  <
357  autoPtr<ModelType>,
358  phaseInterfaceKey,
359  phaseInterfaceKey::hash
360  >& models
361 ) const
362 {
363  // Construct lists of interfaces and models
364  PtrList<phaseInterface> listInterfaces;
365  PtrList<ModelType> listModels;
366  generateInterfacialModels<ModelType, phaseInterface>
367  (
368  dict,
369  NullObjectRef<phaseInterface>(),
370  listInterfaces,
371  listModels
372  );
373 
374  // Transfer to a keyed table
375  forAll(listInterfaces, i)
376  {
377  models.insert(listInterfaces[i], listModels.set(i, nullptr));
378  }
379 }
380 
381 
382 
383 template<class ModelType>
385 (
386  HashTable
387  <
388  autoPtr<ModelType>,
389  phaseInterfaceKey,
390  phaseInterfaceKey::hash
391  >& models
392 ) const
393 {
395  (
396  interfacialDict<dictionary>(modelName<ModelType>()),
397  models
398  );
399 }
400 
401 
402 template<class ValueType>
404 (
405  const dictionary& dict,
406  HashTable<ValueType, phaseInterfaceKey, phaseInterfaceKey::hash>& values
407 ) const
408 {
409  forAllConstIter(dictionary, dict, iter)
410  {
411  autoPtr<phaseInterface> interfacePtr =
412  phaseInterface::New(*this, iter().keyword());
413 
414  const ValueType value(pTraits<ValueType>(iter().stream()));
415 
416  values.insert(interfacePtr(), value);
417  }
418 }
419 
420 
421 template<class ValueType>
423 (
424  const word& valueName,
425  HashTable<ValueType, phaseInterfaceKey, phaseInterfaceKey::hash>& values
426 ) const
427 {
428  generateInterfacialValues(interfacialDict<ValueType>(valueName), values);
429 }
430 
431 
432 template<class ModelType>
434 (
435  const dictionary& dict
436 )
437 {
438  if (dict.size() != 1)
439  {
441  << "Too many matching entries for construction of a "
442  << ModelType::typeName << nl << dict.toc()
443  << exit(FatalError);
444  }
445 
446  if (!dict.first()->isDict())
447  {
449  << "Non-sub-dictionary entries found for specification of a "
450  << ModelType::typeName
451  << exit(FatalError);
452  }
453 
454  return dict.first()->dict();
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 // ************************************************************************* //
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static autoPtr< phaseInterface > New(const phaseSystem &fluid, const word &name)
Select given fluid and name.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
autoPtr< T > set(const label, const word &key, T *)
Set element to pointer provided and return old element.
friend Istream & operator>>(Istream &, dictionary &)
Read dictionary from Istream.
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
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
bool foundObject(const word &name) const
Is the named Type found?
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:956
phaseModelList phaseModels_
Phase models.
Definition: phaseSystem.H:133
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1153
static const dictionary & modelSubDict(const dictionary &dict)
Return the dictionary from which to construct a low-level.
static word separator()
Return the separator that delimits this interface&#39;s name.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:440
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:35
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
void generateInterfacialModels(const dictionary &dict, const phaseInterface &interface, PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models) const
Generate interfacial-model lists.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
static const word null
An empty word.
Definition: word.H:77
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:126
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
void generateInterfacialValues(const dictionary &dict, HashTable< ValueType, phaseInterfaceKey, phaseInterfaceKey::hash > &values) const
Generate interfacial-model tables.
static const char nl
Definition: Ostream.H:260
static word oldNamePartsToName(const phaseSystem &fluid, const wordList &oldNameParts)
Convert old-format interface name parts to an interface name. Used.
dictionary interfacialDict(const word &name) const
Return the dictionary containing interfacial model or value.
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
List< word > wordList
A List of words.
Definition: fileName.H:54
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 & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
void validateMassTransfer(const phaseInterface &interface) const
Check that mass transfer is supported across the given interface.
rDeltaTY field()
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:28
word modelName() const
Return the model name. This is the same as the model&#39;s typename.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864