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-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 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
35 
38  (
39  "in",
41  );
42 }
43 
44 namespace Foam
45 {
47  (
49  separatorsToTypeName({separator()}).c_str(),
50  0
51  );
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const phaseModel& dispersed,
61  const phaseModel& continuous
62 )
63 :
64  phaseInterface(dispersed, continuous),
65  dispersed_(dispersed)
66 {}
67 
68 
70 (
71  const phaseSystem& fluid,
72  const word& name
73 )
74 :
75  phaseInterface(fluid, name),
76  dispersed_(identifyPhases(fluid, name, {separator()}).first())
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
87 
89 {
90  return dispersed().name() + '_' + separator() + '_' + continuous().name();
91 }
92 
93 
95 {
96  return dispersed_;
97 }
98 
99 
101 {
102  return otherPhase(dispersed_);
103 }
104 
105 
107 {
108  return dispersed().U() - continuous().U();
109 }
110 
111 
113 {
114  return magUr()*dispersed().d()/continuous().fluidThermo().nu();
115 }
116 
117 
119 {
120  return
121  continuous().fluidThermo().nu()
122  *continuous().thermo().Cp()
123  *continuous().rho()
124  /continuous().thermo().kappa();
125 }
126 
127 
129 {
130  return Eo(dispersed().d());
131 }
132 
133 
135 (
136  const volScalarField& d
137 ) const
138 {
139  return
140  mag(dispersed().rho() - continuous().rho())
141  *mag(g())
142  *sqr(d)
143  /sigma();
144 }
145 
146 
148 {
149  return
150  mag(g())
151  *continuous().fluidThermo().nu()
152  *pow3
153  (
154  continuous().fluidThermo().nu()
155  *continuous().rho()
156  /sigma()
157  );
158 }
159 
160 
162 {
163  return Re()*pow(Mo(), 0.23);
164 }
165 
166 
167 // ************************************************************************* //
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.
virtual word name() const
Name.
dispersedPhaseInterface(const phaseModel &dispersed, const phaseModel &continuous)
Construct from phases.
tmp< volScalarField > Pr() const
Prandtl number.
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:57
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.
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 and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Namespace for OpenFOAM.
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
defineTypeNameAndDebugWithName(dispersedDisplacedPhaseInterface, separatorsToTypeName({ dispersedPhaseInterface::separator(), displacedPhaseInterface::separator() }).c_str(), 0)
bool dispersedPhaseInterfaceAddedOldSeparatorToSeparator
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
bool dispersedPhaseInterfaceAddedHeadSeparator
dimensioned< scalar > mag(const dimensioned< Type > &)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97