segregated.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) 2014-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 "segregated.H"
27 #include "fvcGrad.H"
28 #include "surfaceInterpolate.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace dragModels
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const dictionary& dict,
49  const phaseInterface& interface,
50  const bool registerObject
51 )
52 :
53  dragModel(dict, interface, registerObject),
54  interface_(interface.modelCast<dragModel, segregatedPhaseInterface>()),
55  m_("m", dimless, dict),
56  n_("n", dimless, dict)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
70  const fvMesh& mesh(interface_.phase1().mesh());
71 
72  const volScalarField& alpha1(interface_.phase1());
73  const volScalarField& alpha2(interface_.phase2());
74 
75  const volScalarField::Internal& rho1(interface_.phase1().rho());
76  const volScalarField::Internal& rho2(interface_.phase2().rho());
77 
78  tmp<volScalarField> tnu1(interface_.phase1().thermo().nu());
79  tmp<volScalarField> tnu2(interface_.phase2().thermo().nu());
80 
81  const volScalarField::Internal& nu1(tnu1());
82  const volScalarField::Internal& nu2(tnu2());
83 
84  const volScalarField::Internal L(cbrt(mesh.V()));
85 
86  const dimensionedScalar residualAlpha
87  (
88  (
89  interface_.phase1().residualAlpha()
90  + interface_.phase2().residualAlpha()
91  )/2
92  );
93 
94  const volScalarField I1
95  (
96  alpha1/max(alpha1 + alpha2, residualAlpha)
97  );
98  const volScalarField I2
99  (
100  alpha2/max(alpha1 + alpha2, residualAlpha)
101  );
102  const volScalarField::Internal magGradI
103  (
104  max
105  (
106  (
107  rho2*mag(fvc::grad(I1)()())
108  + rho1*mag(fvc::grad(I2)()())
109  )/(rho1 + rho2),
110  residualAlpha/2/L
111  )
112  );
113 
114  const volScalarField::Internal muI(rho1*nu1*rho2*nu2/(rho1*nu1 + rho2*nu2));
115 
116  const volScalarField::Internal limitedAlpha1
117  (
118  max(alpha1, interface_.phase1().residualAlpha())
119  );
120 
121  const volScalarField::Internal limitedAlpha2
122  (
123  max(alpha2, interface_.phase2().residualAlpha())
124  );
125 
126  const volScalarField::Internal muAlphaI
127  (
128  limitedAlpha1*rho1*nu1*limitedAlpha2*rho2*nu2
129  /(limitedAlpha1*rho1*nu1 + limitedAlpha2*rho2*nu2)
130  );
131 
132  const volScalarField::Internal ReI
133  (
134  (interface_.rho()()*interface_.magUr()()())
135  /(magGradI*limitedAlpha1*limitedAlpha2*muI)
136  );
137 
138  const volScalarField::Internal lambda(m_*ReI + n_*muAlphaI/muI);
139 
141  (
143  (
144  "K",
145  mesh,
146  dimensionedScalar(dimK, 0),
148  )
149  );
150 
151  tK.ref().ref() = lambda*sqr(magGradI)*muI;
152  tK.ref().correctBoundaryConditions();
153 
154  return tK;
155 }
156 
157 
159 {
160  return fvc::interpolate(K());
161 }
162 
163 
164 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Model for drag between phases.
Definition: dragModel.H:55
Segregated drag model for use in regions with no obvious dispersed phase.
Definition: segregated.H:61
virtual tmp< surfaceScalarField > Kf() const
The drag function Kf used in the face-momentum equations.
Definition: segregated.C:158
segregated(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from a dictionary and an interface.
Definition: segregated.C:47
virtual tmp< volScalarField > K() const
The drag function used in the momentum equation.
Definition: segregated.C:68
virtual ~segregated()
Destructor.
Definition: segregated.C:62
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Class to represent a interface between phases where the two phases are considered to be segregated; t...
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Calculate the gradient of the given field.
K
Definition: pEqn.H:75
dimensionedScalar lambda(viscosity->lookup("lambda"))
defineTypeNameAndDebug(aerosolDrag, 0)
addToRunTimeSelectionTable(dragModel, aerosolDrag, dictionary)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict