Moraga.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) 2014-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 "Moraga.H"
27 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace liftModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phaseInterface& interface
48 )
49 :
50  dispersedLiftModel(dict, interface)
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
63 {
64  volScalarField Re(interface_.Re());
65 
66  volScalarField sqrSr
67  (
68  sqr(interface_.dispersed().d())
69  /interface_.continuous().thermo().nu()
70  *mag(fvc::grad(interface_.continuous().U()))
71  );
72 
73  if
74  (
75  min(Re).value() < 1200.0
76  || max(Re).value() > 18800.0
77  || min(sqrSr).value() < 0.0016
78  || max(sqrSr).value() > 0.04
79  )
80  {
82  << "Re and/or Sr are out of the range of applicability of the "
83  << "Moraga model. Clamping to range bounds"
84  << endl;
85  }
86 
87  Re.min(1200.0);
88  Re.max(18800.0);
89 
90  sqrSr.min(0.0016);
91  sqrSr.max(0.04);
92 
93  return 0.2*exp(- Re*sqrSr/3.6e5 - 0.12)*exp(Re*sqrSr/3.0e7);
94 }
95 
96 
97 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
void max(const dimensioned< Type > &)
void min(const dimensioned< Type > &)
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Model for the lift force between two phases.
Definition: liftModel.H:53
Lift model of Moraga et al.
Definition: Moraga.H:67
virtual tmp< volScalarField > Cl() const
Lift coefficient.
Definition: Moraga.C:62
Moraga(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
Definition: Moraga.C:45
virtual ~Moraga()
Destructor.
Definition: Moraga.C:56
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the gradient of the given field.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
addToRunTimeSelectionTable(liftModel, constantLiftCoefficient, dictionary)
defineTypeNameAndDebug(constantLiftCoefficient, 0)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict