phaseInterfaceI.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-2023 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 "phaseInterface.H"
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
31 {
32  return phase1();
33 }
34 
35 
37 {
38  return phase2();
39 }
40 
41 
43 {
44  return phase1().thermo();
45 }
46 
47 
49 {
50  return phase2().thermo();
51 }
52 
53 
55 {
56  return phase1().rho();
57 }
58 
59 
61 {
62  return phase2().rho();
63 }
64 
65 
67 {
68  return phase1_;
69 }
70 
71 
73 {
74  return phase2_;
75 }
76 
77 
78 inline bool Foam::phaseInterface::contains(const phaseModel& phase) const
79 {
80  return &phase1_ == &phase || &phase2_ == &phase;
81 }
82 
83 
85 (
86  const phaseModel& phase
87 ) const
88 {
89  if (&phase1_ == &phase)
90  {
91  return phase2_;
92  }
93  else if (&phase2_ == &phase)
94  {
95  return phase1_;
96  }
97  else
98  {
100  << "this phaseInterface does not contain phase " << phase.name()
101  << exit(FatalError);
102 
103  return phase;
104  }
105 }
106 
107 
109 {
110  if (&phase1_ == &phase)
111  {
112  return 0;
113  }
114  else if (&phase2_ == &phase)
115  {
116  return 1;
117  }
118  else
119  {
121  << "this phaseInterface does not contain phase " << phase.name()
122  << exit(FatalError);
123 
124  return -1;
125  }
126 }
127 
128 
130 {
131  return phase1().fluid();
132 }
133 
134 
136 {
137  return phase1().mesh();
138 }
139 
140 
143 {
144  return g_;
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Iterators * * * * * * * * * * * * * * * * //
149 
150 inline Foam::phaseInterface::const_iterator::const_iterator
151 (
152  const phaseInterface& pair,
153  const label index
154 )
155 :
156  pair_(pair),
157  index_(index)
158 {}
159 
160 
161 inline Foam::phaseInterface::const_iterator::const_iterator
162 (
163  const phaseInterface& pair
164 )
165 :
166  const_iterator(pair, 0)
167 {}
168 
169 
171 {
172  return index_;
173 }
174 
175 
176 inline bool Foam::phaseInterface::const_iterator::operator==
177 (
178  const const_iterator& iter
179 ) const
180 {
181  return (this->index_ == iter.index_);
182 }
183 
184 
185 inline bool Foam::phaseInterface::const_iterator::operator!=
186 (
187  const const_iterator& iter
188 ) const
189 {
190  return !(this->operator==(iter));
191 }
192 
193 
194 inline const Foam::phaseModel&
196 {
197  if (index_ == 0)
198  {
199  return pair_.phase1_;
200  }
201  else
202  {
203  return pair_.phase2_;
204  }
205 }
206 
207 
208 inline const Foam::phaseModel&
210 {
211  return operator*();
212 }
213 
214 
215 inline const Foam::phaseModel&
217 {
218  if (index_ == 0)
219  {
220  return pair_.phase2_;
221  }
222  else
223  {
224  return pair_.phase1_;
225  }
226 }
227 
228 
231 {
232  index_++;
233  return *this;
234 }
235 
236 
239 {
240  const_iterator old = *this;
241  this->operator++();
242  return old;
243 }
244 
245 
247 {
248  return const_iterator(*this);
249 }
250 
251 
253 {
254  return const_iterator(*this, 2);
255 }
256 
257 
259 {
260  return const_iterator(*this);
261 }
262 
263 
265 {
266  return const_iterator(*this, 2);
267 }
268 
269 
270 // ************************************************************************* //
Generic GeometricField class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const phaseModel & operator*() const
const phaseModel & otherPhase() const
const phaseModel & operator()() const
label index() const
Return the current index.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const uniformDimensionedVectorField & g() const
Return gravitational acceleration.
const_iterator end() const
const_iterator set to beyond the end of the pair
const phaseSystem & fluid() const
Return the phase system.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the other phase relative to the given phase.
virtual const rhoThermo & thermo2() const
Return the thermo for phase 2.
const_iterator begin() const
const_iterator set to the beginning of the pair
const_iterator cend() const
const_iterator set to beyond the end of the pair
virtual const volScalarField & rho2() const
Return the density of phase 2.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
const_iterator cbegin() const
const_iterator set to the beginning of the pair
const fvMesh & mesh() const
Return the mesh.
bool contains(const phaseModel &phase) const
Return true if this phaseInterface contains the given phase.
virtual const volScalarField & rho1() const
Return the density of phase 1.
virtual const volScalarField & alpha2() const
Return the volume fraction of phase 2.
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.
virtual const rhoThermo & thermo1() const
Return the thermo for phase 1.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
tmp< fvMatrix< Type > > operator*(const volScalarField::Internal &, const fvMatrix< Type > &)
error FatalError