dispersedSidedPhaseInterface.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) 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 \*---------------------------------------------------------------------------*/
25 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34  (
36  separatorsToTypeName
37  ({
40  }).c_str(),
41  0
42  );
44  (
47  word
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
55 (
56  const phaseInterface& interface,
57  bool strict
58 ) const
59 {
60  return
61  (!strict || isType<dispersedSidedPhaseInterface>(interface))
62  && dispersedPhaseInterface::same(interface, false)
63  && sidedPhaseInterface::same(interface, false);
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 (
71  const phaseModel& dispersed,
72  const phaseModel& continuous,
73  const phaseModel& phase
74 )
75 :
76  phaseInterface(dispersed, continuous),
77  dispersedPhaseInterface(dispersed, continuous),
78  sidedPhaseInterface(phase, phaseInterface::otherPhase(phase))
79 {}
80 
81 
83 (
84  const phaseSystem& fluid,
85  const word& name
86 )
87 :
88  phaseInterface(fluid, name),
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95 
97 {}
98 
99 
100 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101 
103 {
104  return
106  + '_'
108  + '_'
109  + phase().name();
110 }
111 
112 
113 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
static word separator()
Return the separator that delimits this interface's name.
virtual bool same(const phaseInterface &interface, bool strict) const
Return true if the phase interfaces are the same.
Class to represent a certain side of an interface between phases where one phase is considered disper...
dispersedSidedPhaseInterface(const phaseModel &dispersed, const phaseModel &continuous, const phaseModel &phase)
Construct from phases.
virtual bool same(const phaseInterface &interface, bool strict) const
Return true if the phase interfaces are the same.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Class to represent a system of phases.
Definition: phaseSystem.H:74
Class to represent a certain side of an interface between phases.
static word separator()
Return the separator that delimits this interface's name.
virtual bool same(const phaseInterface &interface, bool strict) const
Return true if the phase interfaces are the same.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
defineTypeNameAndDebugWithName(dispersedDisplacedPhaseInterface, separatorsToTypeName({ dispersedPhaseInterface::separator(), displacedPhaseInterface::separator() }).c_str(), 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.