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-2024 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 Type>
33 {
34  bool found = false;
36 
37  // If it is a dictionary then use it
38  if (this->isDict(name))
39  {
40  found = true;
41  dict = this->subDict(name);
42  }
43 
44  // Function to add old-format list/table entries
45  auto add = [&](const word& sidePhaseName)
46  {
47  const word nameSidePhaseName =
48  IOobject::groupName(name, sidePhaseName);
49 
50  if (!this->found(nameSidePhaseName)) return;
51 
52  found = true;
53 
54  List<wordListAndType<Type>> wlats(this->lookup(nameSidePhaseName));
55 
56  forAll(wlats, i)
57  {
58  word keyword =
59  phaseInterface::oldNamePartsToName(*this, wlats[i].wl);
60 
61  if (sidePhaseName != word::null)
62  {
63  keyword.append
64  (
65  "_"
67  + "_"
68  + sidePhaseName
69  );
70  }
71 
72  dict.add(keyword, wlats[i].t);
73  }
74  };
75 
76  // If a dictionary entry was not found then try a list/table entry
77  if (!found)
78  {
79  add(word::null);
80  }
81 
82  // Add sided models
83  forAll(phases(), phasei)
84  {
85  add(phases()[phasei].name());
86  }
87 
88  // Barf if nothing was found
89  if (!found)
90  {
91  return this->subDict(name);
92  }
93 
94  return dict;
95 }
96 
97 
98 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
99 
100 template<class ModelType>
102 {
103  word name = ModelType::typeName;
104 
105  // Extract the innermost part of the template
106  const word::size_type i0 = name.find_last_of('<');
107  if (i0 != word::npos)
108  {
109  const word::size_type i1 = name.find_first_of('>', i0 + 1);
110  if (i1 != word::npos)
111  {
112  name = name(i0 + 1, i1 - i0 - 1);
113  }
114  }
115 
116  // Strip "Model" off the end of the name
117  if (name(name.size() - 5, 5) == "Model")
118  {
119  name = name(name.size() - 5);
120  }
121 
122  return name;
123 }
124 
125 
126 template<class ModelType, class ... InterfaceTypes>
128 (
129  const dictionary& dict,
130  const phaseInterface& interface,
131  PtrList<phaseInterface>& interfaces,
132  PtrList<ModelType>& models
133 ) const
134 {
135  // Construct sub-dictionaries and associated interfaces
136  hashedWordList names;
137  PtrList<dictionary> dicts;
139  {
140  // Get the model sub dictionary and its associated interface
141  const dictionary& modelDict = iter().dict();
142  autoPtr<phaseInterface> modelInterfacePtr =
143  phaseInterface::New(*this, iter().keyword());
144 
145  // Cast the interface down to the first specified type possible
146  autoPtr<phaseInterface> interfacePtr;
147  List<bool>
148  ({
149  interfacePtr.empty()
150  && isA<InterfaceTypes>(modelInterfacePtr())
151  && (
152  interfacePtr.set
153  (
154  new InterfaceTypes
155  (
156  refCast<InterfaceTypes>(modelInterfacePtr())
157  )
158  ),
159  true
160  )
161  ...
162  });
163  if (!interfacePtr.valid())
164  {
166  << "Interface " << modelInterfacePtr->name()
167  << " is not of suitable type for construction of a "
168  << ModelType::typeName
169  << exit(FatalError);
170  }
171 
172  // If constructing for a specific interface then combine with this
173  // interface. This ensures interface information propagates through
174  // hierarchical model generation.
175  if (notNull(interface))
176  {
177  interfacePtr = phaseInterface::New(interface, interfacePtr());
178  }
179 
180  // Find an existing dictionary to add to or create a new one
181  const word name = interfacePtr->name();
182  if (!names.found(name))
183  {
184  names.append(name);
185  dicts.append(new dictionary(dict.name()));
186  interfaces.append(interfacePtr.ptr());
187  models.append(nullptr);
188  }
189 
190  // Add the model dictionary
191  dicts[names[name]].add
192  (
193  modelInterfacePtr->name(),
194  modelDict
195  );
196  }
197 
198  // Construct the models
199  forAll(interfaces, i)
200  {
201  models.set(i, ModelType::New(dicts[i], interfaces[i]));
202  }
203 }
204 
205 
206 template<class ModelType>
208 (
209  const dictionary& dict,
210  HashTable
211  <
215  >& models
216 ) const
217 {
218  // Construct lists of interfaces and models
219  PtrList<phaseInterface> listInterfaces;
220  PtrList<ModelType> listModels;
221  generateInterfacialModels<ModelType, phaseInterface>
222  (
223  dict,
224  NullObjectRef<phaseInterface>(),
225  listInterfaces,
226  listModels
227  );
228 
229  // Transfer to a keyed table
230  forAll(listInterfaces, i)
231  {
232  models.insert(listInterfaces[i], listModels.set(i, nullptr));
233  }
234 }
235 
236 
237 template<class ModelType>
239 (
240  HashTable
241  <
245  >& models
246 ) const
247 {
248  generateInterfacialModels
249  (
250  interfacialDict<dictionary>(modelName<ModelType>()),
251  models
252  );
253 }
254 
255 
256 template<class ValueType>
258 (
259  const dictionary& dict,
261 ) const
262 {
264  {
265  autoPtr<phaseInterface> interfacePtr =
266  phaseInterface::New(*this, iter().keyword());
267 
268  const ValueType value(pTraits<ValueType>(iter().stream()));
269 
270  values.insert(interfacePtr(), value);
271  }
272 }
273 
274 
275 template<class ValueType>
277 (
278  const word& valueName,
280 ) const
281 {
282  generateInterfacialValues(interfacialDict<ValueType>(valueName), values);
283 }
284 
285 
286 template<class ModelType>
288 (
289  const dictionary& dict
290 )
291 {
292  if (dict.size() != 1)
293  {
295  << "Too many matching entries for construction of a "
296  << ModelType::typeName << nl << dict.toc()
297  << exit(FatalError);
298  }
299 
300  if (!dict.first()->isDict())
301  {
303  << "Non-sub-dictionary entries found for specification of a "
304  << ModelType::typeName
305  << exit(FatalError);
306  }
307 
308  return dict.first()->dict();
309 }
310 
311 
312 template<class ModelType>
314 (
315  const phaseInterface& interface
316 ) const
317 {
318  if (interface.phase1().stationary() || interface.phase2().stationary())
319  {
321  << "A " << ModelType::typeName << " was specified for pair "
322  << interface.name() << ", but one of these phases is stationary. "
323  << "Mass transfer is not supported on stationary phases"
324  << exit(FatalError);
325  }
326 }
327 
328 
329 template<class ModelType>
331 (
332  const phaseInterface& interface
333 ) const
334 {
335  return
336  mesh().foundObject<ModelType>
337  (
338  IOobject::groupName(ModelType::typeName, interface.name())
339  );
340 }
341 
342 
343 template<class ModelType>
345 (
346  const phaseInterface& interface
347 ) const
348 {
349  return
350  mesh().lookupObject<ModelType>
351  (
352  IOobject::groupName(ModelType::typeName, interface.name())
353  );
354 }
355 
356 
357 // ************************************************************************* //
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
An STL-conforming hash table.
Definition: HashTable.H:127
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
static word groupName(Name name, const word &group)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
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
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:111
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:797
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1014
wordList toc() const
Return the table of contents.
Definition: dictionary.C:976
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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.
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:102
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.
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 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:334
word valueName(const direction argument)
Return the name of the value entry for the given argument.
Definition: Product2I.H:81
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
error FatalError
static const char nl
Definition: Ostream.H:266
dictionary dict