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-2023 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 "phaseSystem.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace dragModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 Foam::dragModels::AttouFerschneider::KGasLiquid
46 (
47  const phaseModel& gas,
48  const phaseModel& liquid
49 ) const
50 {
51  const phaseModel& solid = gas.fluid().phases()[solidName_];
52 
53  const volScalarField oneMinusGas(max(1 - gas, liquid.residualAlpha()));
54  const volScalarField cbrtR
55  (
56  cbrt(max(solid, solid.residualAlpha())/oneMinusGas)
57  );
58  const volScalarField magURel(mag(gas.U() - liquid.U()));
59 
60  return
61  E2_*gas.thermo().mu()*sqr(oneMinusGas/solid.d())*sqr(cbrtR)
62  /max(gas, gas.residualAlpha())
63  + E2_*gas.rho()*magURel*(1 - gas)/solid.d()*cbrtR;
64 }
65 
66 
68 Foam::dragModels::AttouFerschneider::KGasSolid
69 (
70  const phaseModel& gas,
71  const phaseModel& solid
72 ) const
73 {
74  const volScalarField oneMinusGas(max(1 - gas, solid.residualAlpha()));
75  const volScalarField cbrtR
76  (
77  cbrt(max(solid, solid.residualAlpha())/oneMinusGas)
78  );
79 
80  return
81  E1_*gas.thermo().mu()*sqr(oneMinusGas/solid.d())*sqr(cbrtR)
82  /max(gas, gas.residualAlpha())
83  + E2_*gas.rho()*mag(gas.U())*(1 - gas)/solid.d()*cbrtR;
84 }
85 
86 
88 Foam::dragModels::AttouFerschneider::KLiquidSolid
89 (
90  const phaseModel& liquid,
91  const phaseModel& solid
92 ) const
93 {
94  const phaseModel& gas = liquid.fluid().phases()[gasName_];
95 
96  return
97  E1_*liquid.thermo().mu()
98  *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 phaseInterface& interface,
110  const bool registerObject
111 )
112 :
113  dragModel(dict, interface, registerObject),
114  interface_(interface),
115  gasName_(dict.lookup("gas")),
116  liquidName_(dict.lookup("liquid")),
117  solidName_(dict.lookup("solid")),
118  E1_("E1", dimless, dict),
119  E2_("E2", dimless, dict)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
124 
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
133 {
134  const phaseModel& gas = interface_.fluid().phases()[gasName_];
135  const phaseModel& liquid = interface_.fluid().phases()[liquidName_];
136  const phaseModel& solid = interface_.fluid().phases()[solidName_];
137 
138  if (interface_.contains(gas) && interface_.contains(liquid))
139  {
140  return KGasLiquid(gas, liquid);
141  }
142  if (interface_.contains(gas) && interface_.contains(solid))
143  {
144  return KGasSolid(gas, solid);
145  }
146  if (interface_.contains(liquid) && interface_.contains(solid))
147  {
148  return KLiquidSolid(liquid, solid);
149  }
150 
152  << "The interface " << interface_.name() << " does not contain two "
153  << "out of the gas, liquid and solid phase models."
154  << exit(FatalError);
155 
156  return tmp<volScalarField>(nullptr);
157 }
158 
159 
162 {
163  return fvc::interpolate(K());
164 }
165 
166 
167 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Model for drag between phases.
Definition: dragModel.H:55
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.
AttouFerschneider(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from a dictionary and an interface.
virtual tmp< volScalarField > K() const
The drag coefficient used in the momentum equation.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:127
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
K
Definition: pEqn.H:75
defineTypeNameAndDebug(aerosolDrag, 0)
addToRunTimeSelectionTable(dragModel, aerosolDrag, dictionary)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict