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 using namespace Foam::constant;
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace dragModels
39 {
40  defineTypeNameAndDebug(aerosolDrag, 0);
41  addToRunTimeSelectionTable(dragModel, aerosolDrag, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const phasePair& pair,
52  const bool registerObject
53 )
54 :
55  dragModel(dict, pair, registerObject),
56  A1_(dict.lookupOrDefault<scalar>("A1", 2.514)),
57  A2_(dict.lookupOrDefault<scalar>("A2", 0.8)),
58  A3_(dict.lookupOrDefault<scalar>("A3", 0.55)),
59  dm_(dimensionedScalar::lookupOrDefault("dm", dict, dimLength, 364e-12))
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
72 {
74  (
76  *pair_.continuous().thermo().T()
77  /sqrt(2.0)
79  /pair_.continuous().thermo().p()
80  /sqr(dm_)
81  );
82 
83  const volScalarField dp(pair_.dispersed().d());
84 
85  const volScalarField Cc
86  (
87  1 + lambda/dp*(A1_ + A2_*exp(-A3_*dp/lambda))
88  );
89 
90  return 24/Cc;
91 }
92 
93 
94 // ************************************************************************* //
Collection of constants.
dictionary dict
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual volScalarField & p()
Pressure [Pa].
Definition: fluidThermo.C:79
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:482
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.
Macros for easy insertion into run-time selection tables.
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:120
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
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.