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-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
37namespace Foam
38{
39namespace dragModels
40{
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 (
100 (
101 "L",
102 mesh.time().timeName(),
103 mesh
104 ),
105 mesh,
108 );
109 L.primitiveFieldRef() = cbrt(mesh.V());
110 L.correctBoundaryConditions();
111
112 const volScalarField I
113 (
114 alpha1
115 /max
116 (
117 alpha1 + alpha2,
118 pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha()
119 )
120 );
121 const volScalarField magGradI
122 (
123 max
124 (
125 mag(fvc::grad(I)),
126 (pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha())/L
127 )
128 );
129
130 const volScalarField muI
131 (
132 rho1*nu1*rho2*nu2
133 /(rho1*nu1 + rho2*nu2)
134 );
135
136 const volScalarField limitedAlpha1
137 (
138 max(alpha1, pair_.phase1().residualAlpha())
139 );
140
141 const volScalarField limitedAlpha2
142 (
143 max(alpha2, pair_.phase2().residualAlpha())
144 );
145
146 const volScalarField muAlphaI
147 (
148 alpha1*rho1*nu1*alpha2*rho2*nu2
149 /(limitedAlpha1*rho1*nu1 + limitedAlpha2*rho2*nu2)
150 );
151
152 const volScalarField ReI
153 (
154 pair_.rho()
155 *pair_.magUr()
156 /(magGradI*limitedAlpha1*limitedAlpha2*muI)
157 );
158
159 const volScalarField lambda(m_*ReI + n_*muAlphaI/muI);
160
161 return lambda*sqr(magGradI)*muI;
162}
163
164
166{
167 return fvc::interpolate(K());
168}
169
170
171// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
volScalarField & rho2
const volScalarField & alpha2
volScalarField & rho1
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Segregated drag model for use in regions with no obvious dispersed phase.
Definition: segregated.H:65
virtual tmp< surfaceScalarField > Kf() const
The drag function Kf used in the face-momentum equations.
Definition: segregated.C:165
virtual tmp< volScalarField > K() const
The drag function used in the momentum equation.
Definition: segregated.C:81
virtual ~segregated()
Destructor.
Definition: segregated.C:64
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
Definition: segregated.C:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
A class for managing temporary objects.
Definition: tmp.H:65
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Calculate the gradient of the given field.
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.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
static const Identity< scalar > I
Definition: Identity.H:94
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
const vector L(dict.get< vector >("L"))