dispersedPhaseInterface.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 
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
36 
39  (
40  "in",
42  );
43 }
44 
45 namespace Foam
46 {
48  (
50  separatorsToTypeName({separator()}).c_str(),
51  0
52  );
54 }
55 
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
60 (
61  const phaseInterface& interface,
62  bool strict
63 ) const
64 {
65  return
66  (!strict || isType<dispersedPhaseInterface>(interface))
67  && (strict || isA<dispersedPhaseInterface>(interface))
68  && &dispersed_
69  == &refCast<const dispersedPhaseInterface>(interface).dispersed_
70  && phaseInterface::same(interface, false);
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
77 (
78  const phaseModel& dispersed,
79  const phaseModel& continuous
80 )
81 :
82  phaseInterface(dispersed, continuous),
83  dispersed_(dispersed)
84 {}
85 
86 
88 (
89  const phaseSystem& fluid,
90  const word& name
91 )
92 :
93  phaseInterface(fluid, name),
94  dispersed_(identifyPhases(fluid, name, {separator()}).first())
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
107 {
108  return dispersed().name() + '_' + separator() + '_' + continuous().name();
109 }
110 
111 
113 {
114  return dispersed_;
115 }
116 
117 
119 {
120  return otherPhase(dispersed_);
121 }
122 
123 
125 {
126  return &dispersed_ == &phase1() ? +1 : -1;
127 }
128 
129 
131 {
132  return dispersed().U() - continuous().U();
133 }
134 
135 
137 {
138  return magUr()*dispersed().d()/continuous().fluidThermo().nu();
139 }
140 
141 
143 {
144  return
145  continuous().fluidThermo().nu()
146  *continuous().thermo().Cp()
147  *continuous().rho()
148  /continuous().thermo().kappa();
149 }
150 
151 
153 {
154  return Eo(dispersed().d());
155 }
156 
157 
159 (
160  const volScalarField& d
161 ) const
162 {
165 
166  return
167  mag(dispersed().rho() - continuous().rho())
168  *mag(g)
169  *sqr(d)
170  /sigma();
171 }
172 
173 
175 {
178 
179  return
180  mag(g)
181  *continuous().fluidThermo().nu()
182  *pow3
183  (
184  continuous().fluidThermo().nu()
185  *continuous().rho()
186  /sigma()
187  );
188 }
189 
190 
192 {
193  return Re()*pow(Mo(), 0.23);
194 }
195 
196 
197 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
tmp< volScalarField > Eo() const
Eotvos number.
static word separator()
Return the separator that delimits this interface's name.
tmp< volScalarField > Mo() const
Morton Number.
virtual ~dispersedPhaseInterface()
Destructor.
dispersedPhaseInterface(const phaseModel &dispersed, const phaseModel &continuous)
Construct from phases.
scalar sign() const
Return the sign. +1 if the dispersed phase is first. -1 if the.
tmp< volScalarField > Pr() const
Prandtl number.
virtual bool same(const phaseInterface &interface, bool strict) const
Return true if the phase interfaces are the same.
const phaseModel & continuous() const
Continuous phase.
tmp< volScalarField > Re() const
Reynolds number.
tmp< volVectorField > Ur() const
Relative velocity.
const phaseModel & dispersed() const
Dispersed phase.
tmp< volScalarField > Ta() const
Takahashi Number.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:56
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
static bool addHeadSeparator(const word &separator)
Add a head separator to the list.
virtual bool same(const phaseInterface &interface, bool strict) const
Return true if the phase interfaces are the same.
static bool addOldSeparatorToSeparator(const word &oldSeparator, const word &separator)
Add a old separator to separator to the table.
Class to represent a system of phases.
Definition: phaseSystem.H:74
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
defineTypeNameAndDebugWithName(dispersedDisplacedPhaseInterface, separatorsToTypeName({ dispersedPhaseInterface::separator(), displacedPhaseInterface::separator() }).c_str(), 0)
bool dispersedPhaseInterfaceAddedOldSeparatorToSeparator
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
bool dispersedPhaseInterfaceAddedHeadSeparator
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97