Frossling.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) 2015-2022 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 "Frossling.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diffusiveMassTransferModels
34 {
35  defineTypeNameAndDebug(Frossling, 0);
37  (
38  diffusiveMassTransferModel,
39  Frossling,
40  dictionary
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const phaseInterface& interface
52 )
53 :
54  diffusiveMassTransferModel(dict, interface),
55  interface_
56  (
57  interface.modelCast
58  <
59  diffusiveMassTransferModel,
60  dispersedPhaseInterface
61  >()
62  ),
63  Le_("Le", dimless, dict)
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
77 {
78  const volScalarField Sh
79  (
80  2 + 0.552*sqrt(interface_.Re())*cbrt(Le_*interface_.Pr())
81  );
82 
83  return 6*interface_.dispersed()*Sh/sqr(interface_.dispersed().d());
84 }
85 
86 
87 // ************************************************************************* //
dictionary dict
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
dimensionedScalar cbrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
Frossling(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
virtual tmp< volScalarField > K() const
The implicit mass transfer coefficient.