Loading...
Searching...
No Matches
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"
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 * * * * * * * * * * * * * * * //
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 tmp<volScalarField> trho1(pair_.phase1().rho());
89 const tmp<volScalarField> trho2(pair_.phase2().rho());
90
91 const volScalarField& rho1 = trho1();
92 const volScalarField& rho2 = trho2();
93
94 tmp<volScalarField> tnu1(pair_.phase1().nu());
95 tmp<volScalarField> tnu2(pair_.phase2().nu());
96
97 const volScalarField& nu1(tnu1());
98 const volScalarField& nu2(tnu2());
99
101 (
103 (
104 "L",
105 mesh.time().timeName(),
106 mesh
107 ),
108 mesh,
111 );
112 L.primitiveFieldRef() = cbrt(mesh.V());
113 L.correctBoundaryConditions();
114
115 const volScalarField I
116 (
117 alpha1
118 /max
119 (
120 alpha1 + alpha2,
121 pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha()
122 )
123 );
124 const volScalarField magGradI
125 (
126 max
127 (
128 mag(fvc::grad(I)),
129 (pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha())/L
130 )
131 );
132
133 const volScalarField muI
134 (
135 rho1*nu1*rho2*nu2
136 /(rho1*nu1 + rho2*nu2)
137 );
138
139 const volScalarField limitedAlpha1
140 (
141 max(alpha1, pair_.phase1().residualAlpha())
142 );
143
144 const volScalarField limitedAlpha2
145 (
146 max(alpha2, pair_.phase2().residualAlpha())
147 );
148
149 const volScalarField muAlphaI
150 (
151 alpha1*rho1*nu1*alpha2*rho2*nu2
152 /(limitedAlpha1*rho1*nu1 + limitedAlpha2*rho2*nu2)
153 );
154
155 const volScalarField ReI
156 (
157 pair_.rho()
158 *pair_.magUr()
159 /(magGradI*limitedAlpha1*limitedAlpha2*muI)
160 );
162 const volScalarField lambda(m_*ReI + n_*muAlphaI/muI);
163
164 return lambda*sqr(magGradI)*muI;
165}
166
167
169{
170 return fvc::interpolate(K());
171}
172
173
174// ************************************************************************* //
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
bool registerObject() const noexcept
Should objects created with this IOobject be registered?
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dragModel(const phasePair &pair, const bool registerObject)
Definition dragModel.C:43
const phasePair & pair_
Phase pair.
Definition dragModel.H:63
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:161
segregated(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from components.
Definition segregated.C:43
virtual tmp< volScalarField > K() const
The drag function used in the momentum equation.
Definition segregated.C:74
virtual ~segregated()
Destructor.
Definition segregated.C:57
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
Definition segregated.C:63
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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:47
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:40
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static const Identity< scalar > I
Definition Identity.H:100
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
const vector L(dict.get< vector >("L"))