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-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 "Moraga.H"
27 #include "phasePair.H"
28 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace liftModels
36 {
37  defineTypeNameAndDebug(Moraga, 0);
38  addToRunTimeSelectionTable(liftModel, Moraga, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
48  const phasePair& pair
49 )
50 :
51  liftModel(dict, pair)
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
64 {
66 
67  volScalarField sqrSr
68  (
69  sqr(pair_.dispersed().d())
70  /pair_.continuous().nu()
72  );
73 
74  if
75  (
76  min(Re).value() < 1200.0
77  || max(Re).value() > 18800.0
78  || min(sqrSr).value() < 0.0016
79  || max(sqrSr).value() > 0.04
80  )
81  {
83  << "Re and/or Sr are out of the range of applicability of the "
84  << "Moraga model. Clamping to range bounds"
85  << endl;
86  }
87 
88  Re.min(1200.0);
89  Re.max(18800.0);
90 
91  sqrSr.min(0.0016);
92  sqrSr.max(0.04);
93 
94  return 0.2*exp(- Re*sqrSr/3.6e5 - 0.12)*exp(Re*sqrSr/3.0e7);
95 }
96 
97 
98 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const phasePair & pair_
Phase pair.
Definition: liftModel.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< volScalarField > d() const
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const volVectorField & U() const
Definition: phaseModel.H:170
virtual const phaseModel & continuous() const
Continuous phase.
dimensionedScalar exp(const dimensionedScalar &ds)
Calculate the gradient of the given field.
Moraga(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
tmp< volScalarField > Re() const
Reynolds number.
virtual const phaseModel & dispersed() const
Dispersed phase.
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define WarningInFunction
Report a warning using Foam::Warning.
virtual ~Moraga()
Destructor.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< volScalarField > Cl() const
Lift coefficient.
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.