aerosolDrag.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) 2019-2020 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 "aerosolDrag.H"
27 #include "phasePair.H"
28 #include "swarmCorrection.H"
30 #include "constants.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace dragModels
37 {
38  defineTypeNameAndDebug(aerosolDrag, 0);
39  addToRunTimeSelectionTable(dragModel, aerosolDrag, dictionary);
40 }
41 }
42 
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const phasePair& pair,
53  const bool registerObject
54 )
55 :
56  dragModel(dict, pair, registerObject),
57  A1_(dict.lookupOrDefault<scalar>("A1", 2.514)),
58  A2_(dict.lookupOrDefault<scalar>("A2", 0.8)),
59  A3_(dict.lookupOrDefault<scalar>("A3", 0.55)),
60  sigma_("sigma", dimLength, dict)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
74  const volScalarField& T = pair_.continuous().thermo().T();
75  const volScalarField& p = pair_.continuous().thermo().p();
76  tmp<volScalarField> td(pair_.dispersed().d());
77  const volScalarField& d = td();
78 
79  const volScalarField lambda(k*T/(sqrt(2.0)*pi*p*sqr(sigma_)));
80 
81  return 24/(1 + lambda/d*(A1_ + A2_*exp(-A3_*d/lambda)));
82 }
83 
84 
85 // ************************************************************************* //
dictionary dict
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
aerosolDrag(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
virtual volScalarField & p()=0
Pressure [Pa].
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual const phaseModel & continuous() const
Continuous phase.
dimensionedScalar exp(const dimensionedScalar &ds)
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual const volScalarField & T() const =0
Temperature [K].
const dimensionedScalar k
Boltzmann constant.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
virtual ~aerosolDrag()
Destructor.
const phasePair & pair_
Phase pair.
Namespace for OpenFOAM.