AttouFerschneider.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) 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 "AttouFerschneider.H"
29#include "phasePair.H"
30#include "phaseSystem.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace dragModels
38{
41}
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
48Foam::dragModels::AttouFerschneider::KGasLiquid
49(
50 const phaseModel& gas,
51 const phaseModel& liquid
52) const
53{
54 const phaseModel& solid = gas.fluid().phases()[solidName_];
55
56 const volScalarField oneMinusGas(max(1 - gas, liquid.residualAlpha()));
57 const volScalarField cbrtR
58 (
59 cbrt(max(solid, solid.residualAlpha())/oneMinusGas)
60 );
61 const volScalarField magURel(mag(gas.U() - liquid.U()));
62
63 return
64 E2_*gas.mu()*sqr(oneMinusGas/solid.d())*sqr(cbrtR)
65 /max(gas, gas.residualAlpha())
66 + E2_*gas.rho()*magURel*(1 - gas)/solid.d()*cbrtR;
67}
68
69
71Foam::dragModels::AttouFerschneider::KGasSolid
72(
73 const phaseModel& gas,
74 const phaseModel& solid
75) const
76{
77 const volScalarField oneMinusGas(max(1 - gas, solid.residualAlpha()));
78 const volScalarField cbrtR
79 (
80 cbrt(max(solid, solid.residualAlpha())/oneMinusGas)
81 );
82
83 return
84 E1_*gas.mu()*sqr(oneMinusGas/solid.d())*sqr(cbrtR)
85 /max(gas, gas.residualAlpha())
86 + E2_*gas.rho()*mag(gas.U())*(1 - gas)/solid.d()*cbrtR;
87}
88
89
91Foam::dragModels::AttouFerschneider::KLiquidSolid
92(
93 const phaseModel& liquid,
94 const phaseModel& solid
95) const
96{
97 const phaseModel& gas = liquid.fluid().phases()[gasName_];
98
99 return
100 E1_*liquid.mu()*sqr(max(solid, solid.residualAlpha())/solid.d())
101 /max(liquid, liquid.residualAlpha())
102 + E2_*liquid.rho()*mag(gas.U())*solid/solid.d();
103}
104
105
106// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107
109(
110 const dictionary& dict,
111 const phasePair& pair,
112 const bool registerObject
113)
114:
115 dragModel(dict, pair, registerObject),
116 gasName_(dict.lookup("gas")),
117 liquidName_(dict.lookup("liquid")),
118 solidName_(dict.lookup("solid")),
119 E1_("E1", dimless, dict),
120 E2_("E2", dimless, dict)
121{}
122
123
124// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125
128{
130 << "Not implemented."
131 << "Drag coefficient is not defined for the AttouFerschneider model."
132 << exit(FatalError);
133
134 return tmp<volScalarField>(nullptr);
135}
136
137
140{
141 switch (Pair<word>::compare(pair_, phasePairKey(gasName_, liquidName_)))
142 {
143 case 1:
144 return KGasLiquid(pair_.phase1(), pair_.phase2());
145 case -1:
146 return KGasLiquid(pair_.phase2(), pair_.phase1());
147 }
148
149 switch (Pair<word>::compare(pair_, phasePairKey(gasName_, solidName_)))
150 {
151 case 1:
152 return KGasSolid(pair_.phase1(), pair_.phase2());
153 case -1:
154 return KGasSolid(pair_.phase2(), pair_.phase1());
155 }
156
157 switch (Pair<word>::compare(pair_, phasePairKey(liquidName_, solidName_)))
158 {
159 case 1:
160 return KLiquidSolid(pair_.phase1(), pair_.phase2());
161 case -1:
162 return KLiquidSolid(pair_.phase2(), pair_.phase1());
163 }
164
166 << "The pair does not contain two of out of the gas, liquid and solid "
167 << "phase models."
168 << exit(FatalError);
169
170 return tmp<volScalarField>(nullptr);
171}
172
173
176{
177 return fvc::interpolate(K());
178}
179
180
181// ************************************************************************* //
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.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Attou and Ferschneider's Drag model for film flow through packed beds. The implementation follows the...
virtual tmp< surfaceScalarField > Kf() const
The drag coefficient used in the face-momentum equations.
virtual tmp< volScalarField > K() const
The drag coefficient used in the momentum equation.
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
Lookup type of boundary radiation properties.
Definition: lookup.H:66
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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.
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)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
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