segregated.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "segregated.H"
29 #include "phasePair.H"
30 #include "fvcGrad.H"
31 #include "surfaceInterpolate.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace dragModels
40 {
41  defineTypeNameAndDebug(segregated, 0);
42  addToRunTimeSelectionTable(dragModel, segregated, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const phasePair& pair,
53  const bool registerObject
54 )
55 :
56  dragModel(dict, pair, registerObject),
57  m_("m", dimless, dict),
58  n_("n", dimless, dict)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 {
73  << "Not implemented."
74  << "Drag coefficient not defined for the segregated model."
75  << exit(FatalError);
76 
77  return pair_.phase1();
78 }
79 
80 
82 {
83  const fvMesh& mesh(pair_.phase1().mesh());
84 
85  const volScalarField& alpha1(pair_.phase1());
86  const volScalarField& alpha2(pair_.phase2());
87 
88  const volScalarField& rho1(pair_.phase1().rho());
89  const volScalarField& rho2(pair_.phase2().rho());
90 
91  tmp<volScalarField> tnu1(pair_.phase1().nu());
92  tmp<volScalarField> tnu2(pair_.phase2().nu());
93 
94  const volScalarField& nu1(tnu1());
95  const volScalarField& nu2(tnu2());
96 
98  (
99  IOobject
100  (
101  "L",
102  mesh.time().timeName(),
103  mesh
104  ),
105  mesh,
107  zeroGradientFvPatchField<scalar>::typeName
108  );
109  L.primitiveFieldRef() = cbrt(mesh.V());
110  L.correctBoundaryConditions();
111 
113  (
114  alpha1
115  /max
116  (
117  alpha1 + alpha2,
118  pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha()
119  )
120  );
121  volScalarField magGradI
122  (
123  max
124  (
125  mag(fvc::grad(I)),
126  (pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha())/L
127  )
128  );
129 
130  volScalarField muI
131  (
132  rho1*nu1*rho2*nu2
133  /(rho1*nu1 + rho2*nu2)
134  );
135  volScalarField muAlphaI
136  (
137  alpha1*rho1*nu1*alpha2*rho2*nu2
138  /(alpha1*rho1*nu1 + alpha2*rho2*nu2)
139  );
140 
141  volScalarField ReI
142  (
143  pair_.rho()
144  *pair_.magUr()
145  /(magGradI*muI)
146  );
147 
148  volScalarField lambda(m_*ReI + n_*muAlphaI/muI);
149 
150  return lambda*sqr(magGradI)*muI;
151 }
152 
153 
155 {
156  return fvc::interpolate(K());
157 }
158 
159 
160 // ************************************************************************* //
L
const vector L(dict.get< vector >("L"))
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::dragModels::defineTypeNameAndDebug
defineTypeNameAndDebug(blended, 0)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::dragModels::segregated::~segregated
virtual ~segregated()
Definition: segregated.C:57
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:52
Foam::dragModels::segregated::K
virtual tmp< volScalarField > K() const
Definition: segregated.C:74
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
rho2
volScalarField & rho2
Definition: setRegionFluidFields.H:30
Foam::dragModels::segregated::Kf
virtual tmp< surfaceScalarField > Kf() const
Definition: segregated.C:158
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
rho1
volScalarField & rho1
Definition: setRegionFluidFields.H:27
Foam::dragModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(dragModel, blended, dictionary)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::dragModels::segregated::CdRe
virtual tmp< volScalarField > CdRe() const
Definition: segregated.C:63
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
fvcGrad.H
Calculate the gradient of the given field.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:148
zeroGradientFvPatchFields.H
Foam::dragModels::segregated::segregated
segregated(const dictionary &dict, const phasePair &pair, const bool registerObject)
Definition: segregated.C:43
Foam::I
static const Identity< scalar > I
Definition: Identity.H:89
Foam::dimless
const dimensionSet dimless
Definition: dimensionSets.C:182
segregated.H