generateInterfacialModels.H
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) 2025 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 Namespace
25  Foam
26 
27 Description
28  Functions for generating tables of interfacial models
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef generateInterfacialModels_H
33 #define generateInterfacialModels_H
34 
35 #include "phaseSystem.H"
36 #include "TypeSet.H"
40 #include "sidedPhaseInterface.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 template<class ModelType>
51 {
52  word name = ModelType::typeName;
53 
54  // Extract the innermost part of the template
55  const word::size_type i0 = name.find_last_of('<');
56  if (i0 != word::npos)
57  {
58  const word::size_type i1 = name.find_first_of('>', i0 + 1);
59  if (i1 != word::npos)
60  {
61  name = name(i0 + 1, i1 - i0 - 1);
62  }
63  }
64 
65  // Strip "Model" off the end of the name
66  if (name(name.size() - 5, 5) == "Model")
67  {
68  name = name(name.size() - 5);
69  }
70 
71  return name;
72 }
73 
74 
75 template<class ModelType>
77 {
78  if (dict.size() != 1)
79  {
81  << "Too many matching entries for construction of a "
82  << ModelType::typeName << nl << dict.toc()
83  << exit(FatalError);
84  }
85 
86  if (!dict.first()->isDict())
87  {
89  << "Non-sub-dictionary entries found for specification of a "
90  << ModelType::typeName
91  << exit(FatalError);
92  }
93 
94  return dict.first()->dict();
95 }
96 
97 
98 template<class ModelType, typename = void>
100 :
101  public std::false_type
102 {
104 
106 
107  typedef TypeSet
108  <
114 };
115 
116 template<class ModelType>
117 struct ModelPhaseInterfaceTypes<ModelType, VoidT<typename ModelType::modelType>>
118 :
119  public std::true_type
120 {
121  typedef typename
123  <
124  typename ModelPhaseInterfaceTypes
125  <
126  typename ModelType::modelType
127  >::required,
128  typename ModelType::requiredPhaseInterfaces
130 
131  typedef typename
133  <
134  typename ModelPhaseInterfaceTypes
135  <
136  typename ModelType::modelType
137  >::allowed,
138  typename ModelType::allowedPhaseInterfaces
140 
141  typedef typename
143  <
145  allowed
147 };
148 
149 
150 template<class ModelType>
152 {
153  return
156 }
157 
158 
159 template<class InterfaceTypeSet>
161 
162 template<class InterfaceType, class ... InterfaceTypes>
163 struct PhaseInterfaceInfo<TypeSet<InterfaceType, InterfaceTypes ...>>
164 {
165  string operator()(const char* andOr) const
166  {
167  const string subS =
168  PhaseInterfaceInfo<TypeSet<InterfaceTypes ...>>()(andOr);
169 
170  string type(InterfaceType::typeName_());
171  type.removeTrailing("PhaseInterface");
172 
173  const string s =
174  type + " (i.e., contains '_" + InterfaceType::separator() + "_')";
175 
176  return
177  sizeof ... (InterfaceTypes) == 0
178  ? s.c_str()
179  : sizeof ... (InterfaceTypes) == 1
180  ? s + ' ' + andOr + ' ' + subS
181  : s + ", " + subS;
182  }
183 };
184 
185 template<>
187 {
188  inline string operator()(const char* andOr) const
189  {
190  return string();
191  }
192 };
193 
194 
195 template<class ModelType>
197 (
198  const phaseSystem& fluid,
199  const dictionary& dict,
200  const wordHashSet& ignoreKeys = wordHashSet()
201 )
202 {
204  {
205  // Get the keyword name and skip if it is ignored
206  const word& modelName = iter().keyword();
207  if (ignoreKeys.found(modelName)) continue;
208 
209  // Construct the associated interface
210  autoPtr<phaseInterface> modelInterfacePtr =
212 
213  // Check the interface against the model's permitted interface types
214  if (!isModelPhaseInterfaceType<ModelType>(modelInterfacePtr()))
215  {
217  << "The interface " << modelName
218  << " is not of suitable type for construction of a "
219  << ModelType::typeName << '.';
220 
221  const string requiredInfo =
223  <
225  >()("and");
226 
227  if (requiredInfo.size())
228  {
230  << " The interface must be of type "
231  << requiredInfo.c_str() << '.';
232  }
233 
234  const string prohibitedInfo =
236  <
238  >()("or");
239 
240  if (prohibitedInfo.size())
241  {
243  << " The interface must NOT be of type "
244  << prohibitedInfo.c_str() << '.';
245  }
246 
248  << exit(FatalIOError);
249  }
250  }
251 }
252 
253 
254 template<class ModelType, class ... InterfaceTypes, class ... Args>
256 (
257  PtrList<phaseInterface>& interfaces,
258  PtrList<ModelType>& models,
259  const phaseSystem& fluid,
260  const dictionary& dict,
261  const wordHashSet& ignoreKeys,
262  const phaseInterface& interface,
263  const Args& ... args
264 )
265 {
266  // Construct sub-dictionaries and associated interfaces
267  hashedWordList names;
268  PtrList<dictionary> dicts;
270  {
271  // Get the keyword name and skip if it is ignored
272  const word& modelName = iter().keyword();
273  if (ignoreKeys.found(modelName)) continue;
274 
275  // Get the model sub dictionary and construct the associated interface
276  const dictionary& modelDict = iter().dict();
277  autoPtr<phaseInterface> modelInterfacePtr =
279 
280  // Check the interface against the model's permitted interface types
281  if
282  (
283  isNull(interface)
284  && !isModelPhaseInterfaceType<ModelType>(modelInterfacePtr())
285  ) continue;
286 
287  // Convert the interface to the first specified type possible
288  autoPtr<phaseInterface> interfacePtr =
289  TypeSet<InterfaceTypes ...>::clone(modelInterfacePtr());
290  if (!interfacePtr.valid())
291  {
293  << "The interface " << modelName
294  << " is not of suitable type for construction of a "
295  << ModelType::typeName
296  << exit(FatalIOError);
297  }
298 
299  // If constructing for a specific interface then combine with this
300  // interface. This ensures interface information propagates through
301  // hierarchical model generation.
302  if (notNull(interface))
303  {
304  interfacePtr = phaseInterface::New(interface, interfacePtr());
305  }
306 
307  // Find an existing dictionary to add to or create a new one
308  const word name = interfacePtr->name();
309  if (!names.found(name))
310  {
311  names.append(name);
312  dicts.append(new dictionary(dict.name()));
313  interfaces.append(interfacePtr.ptr());
314  models.append(nullptr);
315  }
316 
317  // Add the model dictionary
318  dicts[names[name]].add(modelName, modelDict);
319  }
320 
321  // Construct the models
322  forAll(interfaces, i)
323  {
324  models.set(i, ModelType::New(dicts[i], interfaces[i], args ...));
325  }
326 }
327 
328 
329 template<class ModelType, class ... Args>
331 <
332  ModelType,
336 (
337  const phaseSystem& fluid,
338  const dictionary& dict,
339  const wordHashSet& ignoreKeys,
340  const bool ignoreNonModelPhaseInterfaceTypes,
341  const Args& ... args
342 )
343 {
344  // Check that the dictionary has the correct interface types in it
345  if (!ignoreNonModelPhaseInterfaceTypes)
346  {
347  checkInterfacialModelsDict<ModelType>(fluid, dict, ignoreKeys);
348  }
349 
350  // Construct lists of interfaces and models
351  PtrList<phaseInterface> listInterfaces;
352  PtrList<ModelType> listModels;
353  generateInterfacialModels<ModelType, phaseInterface>
354  (
355  listInterfaces,
356  listModels,
357  fluid,
358  dict,
359  ignoreKeys,
360  NullObjectRef<phaseInterface>(),
361  args ...
362  );
363 
364  // Transfer to a keyed table
366  forAll(listInterfaces, i)
367  {
368  models.insert(listInterfaces[i], listModels.set(i, nullptr).ptr());
369  }
370 
371  return models;
372 }
373 
374 
375 template<class ModelType>
377 <
378  ModelType,
382 (
383  const phaseSystem& fluid,
384  const dictionary& dict
385 )
386 {
387  return
388  generateInterfacialModels<ModelType>
389  (
390  fluid,
391  dict,
392  wordHashSet(),
393  false
394  );
395 }
396 
397 
398 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 
400 } // End namespace Foam
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 #endif
405 
406 // ************************************************************************* //
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:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
A HashTable with keys but without contents.
Definition: HashSet.H:62
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:138
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:62
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
T * ptr()
Return object pointer for reuse.
Definition: autoPtrI.H:90
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
wordList toc() const
Return the table of contents.
Definition: dictionary.C:982
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Class to represent an interface between phases which has been displaced to some extent by a third pha...
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.
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 autoPtr< phaseInterface > New(const phaseSystem &fluid, const word &name)
Select given fluid and name.
Class to represent a system of phases.
Definition: phaseSystem.H:74
Class to represent a interface between phases where the two phases are considered to be segregated; t...
Class to represent a certain side of an interface between phases.
A class for handling character strings derived from std::string.
Definition: string.H:79
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
void generateInterfacialModels(PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models, const phaseSystem &fluid, const dictionary &dict, const wordHashSet &ignoreKeys, const phaseInterface &interface, const Args &... args)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool isModelPhaseInterfaceType(const phaseInterface &interface)
void VoidT
Definition: VoidT.H:42
void checkInterfacialModelsDict(const phaseSystem &fluid, const dictionary &dict, const wordHashSet &ignoreKeys=wordHashSet())
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
T clone(const T &t)
Definition: List.H:55
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
error FatalError
static const char nl
Definition: Ostream.H:267
const dictionary & modelSubDict(const dictionary &dict)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
Foam::argList args(argc, argv)
TypeSetRemove< typename ModelPhaseInterfaceTypes< void >::prohibited, allowed >::type prohibited
TypeSetConcatenate< typename ModelPhaseInterfaceTypes< typename ModelType::modelType >::allowed, typename ModelType::allowedPhaseInterfaces >::type allowed
TypeSetConcatenate< typename ModelPhaseInterfaceTypes< typename ModelType::modelType >::required, typename ModelType::requiredPhaseInterfaces >::type required
TypeSet< dispersedPhaseInterface, segregatedPhaseInterface, displacedPhaseInterface, sidedPhaseInterface > prohibited
string operator()(const char *andOr) const
Template meta-programming for operations involving sets of types.
Definition: TypeSet.H:47