RanzMarshall.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 "RanzMarshall.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace heatTransferModels
35 {
36  defineTypeNameAndDebug(RanzMarshall, 0);
37  addToRunTimeSelectionTable(heatTransferModel, RanzMarshall, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phasePair& pair
48 )
49 :
50  heatTransferModel(dict, pair)
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
63 Foam::heatTransferModels::RanzMarshall::K(const scalar residualAlpha) const
64 {
65  volScalarField Nu(scalar(2) + 0.6*sqrt(pair_.Re())*cbrt(pair_.Pr()));
66 
67  return
68  6.0
69  *max(pair_.dispersed(), residualAlpha)
70  *pair_.continuous().kappa()
71  *Nu
72  /sqr(pair_.dispersed().d());
73 }
74 
75 
76 // ************************************************************************* //
dictionary dict
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
RanzMarshall(const dictionary &interfaceDict, const volScalarField &alpha1, const phaseModel &phase1, const phaseModel &phase2)
Construct from components.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar cbrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
tmp< volScalarField > K() const
The heat transfer function K used in the enthalpy equation.
defineTypeNameAndDebug(combustionModel, 0)
virtual scalar Nu(const scalar Re, const scalar Pr) const
Nusselt number.
Definition: RanzMarshall.C:59
A class for managing temporary objects.
Definition: PtrList.H:54
Namespace for OpenFOAM.