AttouFerschneider.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2018 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 "AttouFerschneider.H"
27 #include "phasePair.H"
28 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace dragModels
36 {
37  defineTypeNameAndDebug(AttouFerschneider, 0);
38  addToRunTimeSelectionTable(dragModel, AttouFerschneider, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 Foam::dragModels::AttouFerschneider::KGasLiquid
47 (
48  const phaseModel& gas,
49  const phaseModel& liquid
50 ) const
51 {
52  const phaseModel& solid = gas.fluid().phases()[solidName_];
53 
54  const volScalarField oneMinusGas(max(1 - gas, liquid.residualAlpha()));
55  const volScalarField cbrtR
56  (
57  cbrt(max(solid, solid.residualAlpha())/oneMinusGas)
58  );
59  const volScalarField magURel(mag(gas.U() - liquid.U()));
60 
61  return
62  E2_*gas.mu()*sqr(oneMinusGas/solid.d())*sqr(cbrtR)
63  /max(gas, gas.residualAlpha())
64  + E2_*gas.rho()*magURel*(1 - gas)/solid.d()*cbrtR;
65 }
66 
67 
69 Foam::dragModels::AttouFerschneider::KGasSolid
70 (
71  const phaseModel& gas,
72  const phaseModel& solid
73 ) const
74 {
75  const volScalarField oneMinusGas(max(1 - gas, solid.residualAlpha()));
76  const volScalarField cbrtR
77  (
78  cbrt(max(solid, solid.residualAlpha())/oneMinusGas)
79  );
80 
81  return
82  E1_*gas.mu()*sqr(oneMinusGas/solid.d())*sqr(cbrtR)
83  /max(gas, gas.residualAlpha())
84  + E2_*gas.rho()*mag(gas.U())*(1 - gas)/solid.d()*cbrtR;
85 }
86 
87 
89 Foam::dragModels::AttouFerschneider::KLiquidSolid
90 (
91  const phaseModel& liquid,
92  const phaseModel& solid
93 ) const
94 {
95  const phaseModel& gas = liquid.fluid().phases()[gasName_];
96 
97  return
98  E1_*liquid.mu()*sqr(max(solid, solid.residualAlpha())/solid.d())
99  /max(liquid, liquid.residualAlpha())
100  + E2_*liquid.rho()*mag(gas.U())*solid/solid.d();
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
107 (
108  const dictionary& dict,
109  const phasePair& pair,
110  const bool registerObject
111 )
112 :
113  dragModel(dict, pair, registerObject),
114  gasName_(dict.lookup("gas")),
115  liquidName_(dict.lookup("liquid")),
116  solidName_(dict.lookup("solid")),
117  E1_("E1", dimless, dict),
118  E2_("E2", dimless, dict)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
123 
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
132 {
134  << "Not implemented."
135  << "Drag coefficient is not defined for the AttouFerschneider model."
136  << exit(FatalError);
137 
138  return tmp<volScalarField>(nullptr);
139 }
140 
141 
144 {
145  switch (Pair<word>::compare(pair_, phasePairKey(gasName_, liquidName_)))
146  {
147  case 1:
148  return KGasLiquid(pair_.phase1(), pair_.phase2());
149  case -1:
150  return KGasLiquid(pair_.phase2(), pair_.phase1());
151  }
152 
153  switch (Pair<word>::compare(pair_, phasePairKey(gasName_, solidName_)))
154  {
155  case 1:
156  return KGasSolid(pair_.phase1(), pair_.phase2());
157  case -1:
158  return KGasSolid(pair_.phase2(), pair_.phase1());
159  }
160 
161  switch (Pair<word>::compare(pair_, phasePairKey(liquidName_, solidName_)))
162  {
163  case 1:
164  return KLiquidSolid(pair_.phase1(), pair_.phase2());
165  case -1:
166  return KLiquidSolid(pair_.phase2(), pair_.phase1());
167  }
168 
170  << "The pair does not contain two of out of the gas, liquid and solid "
171  << "phase models."
172  << exit(FatalError);
173 
174  return tmp<volScalarField>(nullptr);
175 }
176 
177 
180 {
181  return fvc::interpolate(K());
182 }
183 
184 
185 // ************************************************************************* //
dictionary dict
const phaseModel & phase2() const
Return phase 2.
Definition: phasePairI.H:34
virtual tmp< volScalarField > K() const
The drag coefficient used in the momentum equation.
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
const phaseModel & phase1() const
Return phase 1.
Definition: phasePairI.H:28
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
virtual tmp< surfaceScalarField > Kf() const
The drag coefficient used in the face-momentum equations.
virtual ~AttouFerschneider()
Destructor.
CGAL::Exact_predicates_exact_constructions_kernel K
const phasePair & pair_
Phase pair.
Macros for easy insertion into run-time selection tables.
AttouFerschneider(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar cbrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
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 dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
Namespace for OpenFOAM.