phaseInterface.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) 2021-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 Class
25  Foam::phaseInterface
26 
27 Description
28  Class to represent an interface between phases. Derivations can further
29  specify the configuration of that interface; e.g., representing dispersal,
30  displacement or sidedness.
31 
32 SourceFiles
33  phaseInterface.C
34  phaseInterfaceI.H
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef phaseInterface_H
39 #define phaseInterface_H
40 
41 #include "phaseModel.H"
42 #include "compressibleTwoPhases.H"
43 #include "runTimeSelectionTables.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 class phaseInterfaceKey;
51 
52 /*---------------------------------------------------------------------------*\
53  Class phaseInterface Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class phaseInterface
57 :
59 {
60 private:
61 
62  // Private Static Data
63 
64  //- List of head separators. This is the first separator in an
65  // interface name, and it delimits the names of phases on either side
66  // of the interface. Other separators are modifiers (e.g.,
67  // displacedBy, inThe, ...).
68  static wordList headSeparators_;
69 
70  //- Map from old-format separators to separators
71  static HashTable<word> oldSeparatorToSeparator_;
72 
73 
74 public:
75 
76  // Public Static Functions
77 
78  //- Get a reference to phase1 after sorting the phases by index
79  static const phaseModel& getPhase1
80  (
81  const phaseModel& phase1,
82  const phaseModel& phase2
83  );
84 
85  //- Get a reference to phase2 after sorting the phases by index
86  static const phaseModel& getPhase2
87  (
88  const phaseModel& phase1,
89  const phaseModel& phase2
90  );
91 
92  //- Add a head separator to the list
93  static bool addHeadSeparator(const word& separator);
94 
95  //- Add a old separator to separator to the table
96  static bool addOldSeparatorToSeparator
97  (
98  const word& oldSeparator,
99  const word& separator
100  );
101 
102  //- Split an interface name and return all its parts
104  (
105  const phaseSystem& fluid,
106  const word& name
107  );
108 
109  //- Split an interface name and return its separators
111  (
112  const phaseSystem& fluid,
113  const word& name
114  );
115 
116  //- Convert a list of separators into a type name
118  (
119  const wordList& separators
120  );
121 
122  //- Convert an interface name into a type name. Essentially just
123  // replaces valid phase names with a "<phase>" placeholder.
124  static word nameToTypeName
125  (
126  const phaseSystem& fluid,
127  const word& name
128  );
129 
130  //- Convert interface name parts to an interface name
131  static word namePartsToName
132  (
133  const phaseSystem& fluid,
134  const wordList& nameParts
135  );
136 
137  //- Convert old-format interface name parts to an interface name. Used
138  // in phaseSystem to provide backwards compatible input.
139  static word oldNamePartsToName
140  (
141  const phaseSystem& fluid,
142  const wordList& oldNameParts
143  );
144 
145  //- Return references to the phases associated with a given name, and a
146  // list of valid separators
148  (
149  const phaseSystem& fluid,
150  const word& name,
151  const wordList& separators
152  );
153 
154 
155 private:
156 
157  // Private Data
158 
159  //- Phase 1
160  const phaseModel& phase1_;
161 
162  //- Phase 2
163  const phaseModel& phase2_;
164 
165 
166 protected:
167 
168  // Protected Member Functions
169 
170  //- Return true if the phase interfaces are the same
171  virtual bool same(const phaseInterface& interface, bool strict) const;
172 
173 
174 public:
175 
176  //- Runtime type information
177  TypeName("phaseInterface");
178 
179 
180  // Declare runtime construction
181 
183  (
184  autoPtr,
186  word,
187  (
188  const phaseSystem& fluid,
189  const word& name
190  ),
191  (fluid, name)
192  );
193 
194 
195  // Constructors
196 
197  //- Construct from phases
199  (
200  const phaseModel& phase1,
201  const phaseModel& phase2
202  );
203 
204  //- Construct from phases
206  (
208  );
209 
210  //- Construct from fluid and name
212  (
213  const phaseSystem& fluid,
214  const word& name
215  );
216 
217  //- Construct from fluid and key
219  (
220  const phaseSystem& fluid,
221  const phaseInterfaceKey& name
222  );
223 
224  //- Clone function
225  virtual autoPtr<phaseInterface> clone() const;
226 
227 
228  //- Destructor
229  virtual ~phaseInterface();
230 
231 
232  // Selectors
233 
234  //- Select given fluid and name
236  (
237  const phaseSystem& fluid,
238  const word& name
239  );
240 
241  //- Select by combining two interfaces
243  (
244  const phaseInterface& interface1,
245  const phaseInterface& interface2
246  );
247 
248  //- Class used for construction of PtrLists of phaseInterfaces
249  class iNew
250  {
251  const phaseSystem& fluid_;
252 
253  public:
254 
255  iNew(const phaseSystem& fluid)
256  :
257  fluid_(fluid)
258  {}
259 
261  {
262  return phaseInterface::New(fluid_, word(is));
263  }
264  };
265 
266 
267  // Static Member Functions
268 
269  //- Return the separator that delimits this interface's name
270  static word separator()
271  {
272  return word::null;
273  }
274 
275 
276  // Member Functions
277 
278  //- Name
279  virtual word name() const;
280 
281  //- Cast to derived type for use in a model
282  template<class ModelType, class Derived>
283  const Derived& modelCast() const
284  {
285  if (!isA<Derived>(*this))
286  {
288  << "Constructing " << ModelType::typeName
289  << " for interface " << name()
290  << " which is not of the required type "
291  << Derived::typeName << exit(FatalError);
292  }
293 
294  return refCast<const Derived>(*this);
295  }
296 
297 
298  // Access
299 
300  //- Return the volume fraction of phase 1
301  virtual inline const volScalarField& alpha1() const;
302 
303  //- Return the volume fraction of phase 2
304  virtual inline const volScalarField& alpha2() const;
305 
306  //- Return the thermo for phase 1
307  virtual inline const rhoFluidThermo& thermo1() const;
308 
309  //- Return the thermo for phase 2
310  virtual inline const rhoFluidThermo& thermo2() const;
311 
312  //- Return the density of phase 1
313  virtual inline const volScalarField& rho1() const;
314 
315  //- Return the density of phase 2
316  virtual inline const volScalarField& rho2() const;
317 
318  //- Return phase 1
319  inline const phaseModel& phase1() const;
320 
321  //- Return phase 2
322  inline const phaseModel& phase2() const;
323 
324  //- Return true if this phaseInterface contains the given phase
325  inline bool contains(const phaseModel& phase) const;
326 
327  //- Return true if the phase interfaces are the same
328  inline bool operator==(const phaseInterface& interface) const;
329 
330  //- Return true if the phase interfaces are different
331  inline bool operator!=(const phaseInterface& interface) const;
332 
333  //- Return the other phase relative to the given phase
334  // Generates a FatalError if this phaseInterface does not contain
335  // the given phase
336  inline const phaseModel& otherPhase(const phaseModel& phase) const;
337 
338  //- Return the index of the given phase. Generates a FatalError if
339  // this phaseInterface does not contain the given phase
340  inline label index(const phaseModel& phase) const;
341 
342  //- Access the phase for the given index. Generates a FatalError if
343  // the index is not 0 or 1.
344  inline const phaseModel& operator[](const label index) const;
345 
346  //- Return the phase system
347  inline const phaseSystem& fluid() const;
348 
349  //- Return the mesh
350  inline const fvMesh& mesh() const;
351 
352 
353  // Properties
354 
355  //- Average density
356  tmp<volScalarField> rho() const;
357 
358  //- Relative velocity magnitude
359  tmp<volScalarField> magUr() const;
360 
361  //- Surface tension coefficient
362  tmp<volScalarField> sigma() const;
363 
364 
365  // Iterators
366 
367  //- STL const_iterator
368  class const_iterator
369  {
370  // Private Data
371 
372  //- Reference to the interface for which this is an iterator
373  const phaseInterface& interface_;
374 
375  //- Current index
376  label index_;
377 
378  //- Construct an iterator with the given index
379  inline const_iterator(const phaseInterface&, const label index);
380 
381  public:
382 
383  friend class phaseInterface;
384 
385  // Constructors
386 
387  //- Construct from interface, moving to its 'begin' position
388  inline explicit const_iterator(const phaseInterface&);
389 
390 
391  // Access
392 
393  //- Return the current index
394  inline label index() const;
395 
396 
397  // Member Operators
398 
399  inline bool operator==(const const_iterator&) const;
400 
401  inline bool operator!=(const const_iterator&) const;
402 
403  inline const phaseModel& operator*() const;
404  inline const phaseModel& operator()() const;
405 
406  inline const phaseModel& otherPhase() const;
407 
408  inline const_iterator& operator++();
409  inline const_iterator operator++(int);
410  };
411 
412 
413  //- const_iterator set to the beginning of the interface
414  inline const_iterator cbegin() const;
415 
416  //- const_iterator set to beyond the end of the interface
417  inline const_iterator cend() const;
418 
419  //- const_iterator set to the beginning of the interface
420  inline const_iterator begin() const;
421 
422  //- const_iterator set to beyond the end of the interface
423  inline const_iterator end() const;
424 };
425 
426 
427 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
428 
429 } // End namespace Foam
430 
431 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
432 
433 #include "phaseInterfaceI.H"
434 
435 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
436 
437 #endif
438 
439 // ************************************************************************* //
Generic GeometricField class.
An STL-conforming hash table.
Definition: HashTable.H:127
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Word-pair based class used for keying interface models in hash tables.
const phaseModel & operator*() const
bool operator!=(const const_iterator &) const
const phaseModel & otherPhase() const
const phaseModel & operator()() const
bool operator==(const const_iterator &) const
label index() const
Return the current index.
Class used for construction of PtrLists of phaseInterfaces.
iNew(const phaseSystem &fluid)
autoPtr< phaseInterface > operator()(Istream &is) const
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.
static const phaseModel & getPhase1(const phaseModel &phase1, const phaseModel &phase2)
Get a reference to phase1 after sorting the phases by index.
static Tuple2< const phaseModel &, const phaseModel & > identifyPhases(const phaseSystem &fluid, const word &name, const wordList &separators)
Return references to the phases associated with a given name, and a.
const_iterator end() const
const_iterator set to beyond the end of the interface
const phaseSystem & fluid() const
Return the phase system.
TypeName("phaseInterface")
Runtime type information.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the other phase relative to the given phase.
static word separator()
Return the separator that delimits this interface's name.
virtual autoPtr< phaseInterface > clone() const
Clone function.
const phaseModel & operator[](const label index) const
Access the phase for the given index. Generates a FatalError if.
virtual word name() const
Name.
static wordList nameToNameParts(const phaseSystem &fluid, const word &name)
Split an interface name and return all its parts.
static word nameToTypeName(const phaseSystem &fluid, const word &name)
Convert an interface name into a type name. Essentially just.
virtual const rhoFluidThermo & thermo2() const
Return the thermo for phase 2.
const_iterator begin() const
const_iterator set to the beginning of the interface
virtual const rhoFluidThermo & thermo1() const
Return the thermo for phase 1.
const_iterator cend() const
const_iterator set to beyond the end of the interface
static bool addHeadSeparator(const word &separator)
Add a head separator to the list.
virtual ~phaseInterface()
Destructor.
tmp< volScalarField > magUr() const
Relative velocity magnitude.
bool operator==(const phaseInterface &interface) const
Return true if the phase interfaces are the same.
virtual bool same(const phaseInterface &interface, bool strict) const
Return true if the phase interfaces are the same.
static word separatorsToTypeName(const wordList &separators)
Convert a list of separators into a type name.
static const phaseModel & getPhase2(const phaseModel &phase1, const phaseModel &phase2)
Get a reference to phase2 after sorting the phases by index.
virtual const volScalarField & rho2() const
Return the density of phase 2.
declareRunTimeSelectionTable(autoPtr, phaseInterface, word,(const phaseSystem &fluid, const word &name),(fluid, name))
tmp< volScalarField > sigma() const
Surface tension coefficient.
static wordList nameToSeparators(const phaseSystem &fluid, const word &name)
Split an interface name and return its separators.
phaseInterface(const phaseModel &phase1, const phaseModel &phase2)
Construct from phases.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
static bool addOldSeparatorToSeparator(const word &oldSeparator, const word &separator)
Add a old separator to separator to the table.
const_iterator cbegin() const
const_iterator set to the beginning of the interface
const fvMesh & mesh() const
Return the mesh.
bool contains(const phaseModel &phase) const
Return true if this phaseInterface contains the given phase.
tmp< volScalarField > rho() const
Average density.
virtual const volScalarField & rho1() const
Return the density of phase 1.
const Derived & modelCast() const
Cast to derived type for use in a model.
static word namePartsToName(const phaseSystem &fluid, const wordList &nameParts)
Convert interface name parts to an interface name.
virtual const volScalarField & alpha2() const
Return the volume fraction of phase 2.
bool operator!=(const phaseInterface &interface) const
Return true if the phase interfaces are different.
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 const volScalarField & alpha1() const
Return the volume fraction of phase 1.
Class to represent a system of phases.
Definition: phaseSystem.H:74
Base-class for fluid thermodynamic properties based on density.
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:334
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
error FatalError
Macros to ease declaration of run-time selection tables.