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-2022 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 {
38  defineTypeNameAndDebug(segregated, 0);
39  addToRunTimeSelectionTable(dragModel, segregated, dictionary);
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 
74 
77 
78  tmp<volScalarField> tnu1(interface_.phase1().thermo().nu());
79  tmp<volScalarField> tnu2(interface_.phase2().thermo().nu());
80 
81  const volScalarField& nu1(tnu1());
82  const volScalarField& nu2(tnu2());
83 
85  (
86  IOobject
87  (
88  "L",
89  mesh.time().timeName(),
90  mesh
91  ),
92  mesh,
95  );
96  L.primitiveFieldRef() = cbrt(mesh.V());
97  L.correctBoundaryConditions();
98 
99  const dimensionedScalar residualAlpha
100  (
101  (
104  )/2
105  );
106 
107  const volScalarField I1(alpha1/max(alpha1 + alpha2, residualAlpha));
108  const volScalarField I2(alpha2/max(alpha1 + alpha2, residualAlpha));
109  const volScalarField magGradI
110  (
111  max
112  (
113  (rho2*mag(fvc::grad(I1)) + rho1*mag(fvc::grad(I2)))/(rho1 + rho2),
114  residualAlpha/2/L
115  )
116  );
117 
118  const volScalarField muI(rho1*nu1*rho2*nu2/(rho1*nu1 + rho2*nu2));
119 
120  const volScalarField limitedAlpha1
121  (
123  );
124 
125  const volScalarField limitedAlpha2
126  (
128  );
129 
130  const volScalarField muAlphaI
131  (
132  limitedAlpha1*rho1*nu1*limitedAlpha2*rho2*nu2
133  /(limitedAlpha1*rho1*nu1 + limitedAlpha2*rho2*nu2)
134  );
135 
136  const volScalarField ReI
137  (
139  /(magGradI*limitedAlpha1*limitedAlpha2*muI)
140  );
141 
142  const volScalarField lambda(m_*ReI + n_*muAlphaI/muI);
143 
144  return lambda*sqr(magGradI)*muI;
145 }
146 
147 
149 {
150  return fvc::interpolate(K());
151 }
152 
153 
154 // ************************************************************************* //
const volScalarField & rho1
virtual tmp< volScalarField > K() const
The drag function used in the momentum equation.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static const char *const typeName
Definition: Field.H:105
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< volScalarField > rho() const
Average density.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:167
virtual tmp< volScalarField > rho() const =0
Return the density field.
dimensionedScalar lambda(viscosity->lookup("lambda"))
volScalarField & alpha1(mixture.alpha1())
const phaseModel & phase1() const
Return phase 1.
const dimensionSet dimless
virtual ~segregated()
Destructor.
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
alpha2
Definition: alphaEqn.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the gradient of the given field.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
const phaseModel & phase2() const
Return phase 2.
tmp< volScalarField > magUr() const
Relative velocity magnitude.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const volScalarField & rho2
const dispersedPhaseInterface interface_
Interface.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:73
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
segregated(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from a dictionary and an interface.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< surfaceScalarField > Kf() const
The drag function Kf used in the face-momentum equations.
Namespace for OpenFOAM.