Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2012 OpenFOAM Foundation
9 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "RanzMarshall.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34namespace Foam
36namespace multiphaseEuler
37{
38namespace heatTransferModels
39{
41
43 (
47 );
49}
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const dictionary& interfaceDict,
59 const phaseModel& phase1,
60 const phaseModel& phase2
61)
63 heatTransferModel(interfaceDict, alpha1, phase1, phase2)
64{}
65
66
67// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68
70{}
71
72
73// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74
77(
78 const volScalarField& Ur
79) const
80{
81 volScalarField Re(max(Ur*phase1_.d()/phase2_.nu(), scalar(1.0e-3)));
83 phase2_.rho()*phase2_.nu()*phase2_.Cp()/phase2_.kappa();
84 volScalarField Nu(scalar(2) + 0.6*sqrt(Re)*cbrt(Prb));
85
86 return 6.0*phase2_.kappa()*Nu/sqr(phase1_.d());
87}
88
89
90// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
phaseModel & phase1
phaseModel & phase2
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
heatTransferModel(const dictionary &dict, const volScalarField &alpha1, const phaseModel &phase1, const phaseModel &phase2)
tmp< volScalarField > K(const volScalarField &Ur) const
The heat-transfer function K used in the enthalpy eq.
RanzMarshall(const dictionary &interfaceDict, const volScalarField &alpha1, const phaseModel &phase1, const phaseModel &phase2)
Construct from components.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sqrt(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.