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-2016 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"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace dragModels
38 {
39  defineTypeNameAndDebug(segregated, 0);
40  addToRunTimeSelectionTable(dragModel, segregated, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const dictionary& dict,
50  const phasePair& pair,
51  const bool registerObject
52 )
53 :
54  dragModel(dict, pair, registerObject),
55  m_("m", dimless, dict),
56  n_("n", dimless, dict)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
71  << "Not implemented."
72  << "Drag coefficient not defined for the segregated model."
73  << exit(FatalError);
74 
75  return pair_.phase1();
76 }
77 
78 
80 {
81  const fvMesh& mesh(pair_.phase1().mesh());
82 
85 
86  const volScalarField& rho1(pair_.phase1().rho());
87  const volScalarField& rho2(pair_.phase2().rho());
88 
89  tmp<volScalarField> tnu1(pair_.phase1().nu());
90  tmp<volScalarField> tnu2(pair_.phase2().nu());
91 
92  const volScalarField& nu1(tnu1());
93  const volScalarField& nu2(tnu2());
94 
96  (
97  IOobject
98  (
99  "L",
100  mesh.time().timeName(),
101  mesh
102  ),
103  mesh,
104  dimensionedScalar("L", dimLength, 0),
106  );
107  L.primitiveFieldRef() = cbrt(mesh.V());
108  L.correctBoundaryConditions();
109 
111  (
112  alpha1
113  /max
114  (
115  alpha1 + alpha2,
117  )
118  );
119  volScalarField magGradI
120  (
121  max
122  (
123  mag(fvc::grad(I)),
125  )
126  );
127 
128  volScalarField muI
129  (
130  rho1*nu1*rho2*nu2
131  /(rho1*nu1 + rho2*nu2)
132  );
133  volScalarField muAlphaI
134  (
135  alpha1*rho1*nu1*alpha2*rho2*nu2
136  /(alpha1*rho1*nu1 + alpha2*rho2*nu2)
137  );
138 
139  volScalarField ReI
140  (
141  pair_.rho()
142  *pair_.magUr()
143  /(magGradI*muI)
144  );
145 
146  volScalarField lambda(m_*ReI + n_*muAlphaI/muI);
147 
148  return lambda*sqr(magGradI)*muI;
149 }
150 
151 
153 {
154  return fvc::interpolate(K());
155 }
156 
157 
158 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
static const char *const typeName
Definition: Field.H:94
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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:114
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
virtual ~segregated()
Destructor.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
virtual tmp< surfaceScalarField > Kf() const
The drag function Kf used in the face-momentum equations.
const phasePair & pair_
Phase pair.
Definition: dragModel.H:62
Macros for easy insertion into run-time selection tables.
rho1
Definition: pEqn.H:114
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
static const Identity< scalar > I
Definition: Identity.H:93
Calculate the gradient of the given field.
rho2
Definition: pEqn.H:115
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
alpha2
Definition: alphaEqn.H:112
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > K() const
The drag function used in the momentum equation.
volScalarField & alpha1
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.
const phaseModel & phase1() const
Definition: phasePairI.H:28
segregated(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from components.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
tmp< volScalarField > rho() const
Average density.
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< volScalarField > magUr() const
Relative velocity magnitude.
const phaseModel & phase2() const
Definition: phasePairI.H:34
Namespace for OpenFOAM.