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 | Copyright (C) 2014-2015 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 "phasePair.H"
28 #include "fvcGrad.H"
29 #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 phasePair& pair,
50  const bool registerObject
51 )
52 :
53  dragModel(dict, pair, registerObject),
54  m_("m", dimless, dict),
55  n_("n", dimless, dict)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
68 {
70  << "Not implemented."
71  << "Drag coefficient not defined for the segregated model."
72  << exit(FatalError);
73 
74  return pair_.phase1();
75 }
76 
77 
79 {
80  const fvMesh& mesh(pair_.phase1().mesh());
81 
82  const volScalarField& alpha1(pair_.phase1());
83  const volScalarField& alpha2(pair_.phase2());
84 
85  const volScalarField& rho1(pair_.phase1().rho());
86  const volScalarField& rho2(pair_.phase2().rho());
87 
88  tmp<volScalarField> tnu1(pair_.phase1().nu());
89  tmp<volScalarField> tnu2(pair_.phase2().nu());
90 
91  const volScalarField& nu1(tnu1());
92  const volScalarField& nu2(tnu2());
93 
95  (
96  IOobject
97  (
98  "L",
99  mesh.time().timeName(),
100  mesh
101  ),
102  mesh,
103  dimensionedScalar("L", dimLength, 0),
104  zeroGradientFvPatchField<scalar>::typeName
105  );
106  L.internalField() = cbrt(mesh.V());
107  L.correctBoundaryConditions();
108 
110  (
111  alpha1
112  /max
113  (
114  alpha1 + alpha2,
115  pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha()
116  )
117  );
118  volScalarField magGradI
119  (
120  max
121  (
122  mag(fvc::grad(I)),
123  (pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha())/L
124  )
125  );
126 
127  volScalarField muI
128  (
129  rho1*nu1*rho2*nu2
130  /(rho1*nu1 + rho2*nu2)
131  );
132  volScalarField muAlphaI
133  (
134  alpha1*rho1*nu1*alpha2*rho2*nu2
135  /(
136  max(alpha1, pair_.phase1().residualAlpha())*rho1*nu1
137  + max(alpha2, pair_.phase2().residualAlpha())*rho2*nu2
138  )
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 // ************************************************************************* //
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
Foam::dragModels::segregated::K
virtual tmp< volScalarField > K() const
The drag function used in the momentum equation.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
rho2
rho2
Definition: pEqn.H:125
Foam::I
static const sphericalTensor I(1)
alpha2
alpha2
Definition: alphaEqn.H:112
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dragModels::segregated::CdRe
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
surfaceInterpolate.H
Surface Interpolation.
rho1
rho1
Definition: pEqn.H:124
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
alpha1
volScalarField & alpha1
Definition: createFields.H:15
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::dragModels::segregated::Kf
virtual tmp< surfaceScalarField > Kf() const
The drag function Kf used in the face-momentum equations.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
fvcGrad.H
Calculate the gradient of the given field.
Foam::dragModels::segregated::~segregated
virtual ~segregated()
Destructor.
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
lambda
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
Foam::dragModels::segregated::segregated
segregated(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from components.