slurry.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) 2014-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 "slurry.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace mixtureViscosityModels
34 {
35  defineTypeNameAndDebug(slurry, 0);
36 
38  (
39  mixtureViscosityModel,
40  slurry,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const incompressibleTwoPhaseInteractingMixture& mixture
52 )
53 :
54  mixtureViscosityModel(mixture)
55 {}
56 
57 
58 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
59 
62 (
63  const volScalarField& muc,
64  const volVectorField& U
65 ) const
66 {
67  const volScalarField& alphad = mixture_.alphad();
68 
69  return
70  (
71  muc*(1.0 + 2.5*alphad + 10.05*sqr(alphad) + 0.00273*exp(16.6*alphad))
72  );
73 }
74 
75 
77 {
79 }
80 
81 
82 // ************************************************************************* //
const volScalarField & alphad() const
Return const-access to the dispersed phase-fraction.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
Macros for easy insertion into run-time selection tables.
virtual bool read()=0
Read physicalProperties dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
dimensionedScalar exp(const dimensionedScalar &ds)
slurry(const incompressibleTwoPhaseInteractingMixture &mixture)
Construct from mixture.
const incompressibleTwoPhaseInteractingMixture & mixture_
Mixture properties.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > mu(const volScalarField &muc, const volVectorField &U) const
Return the mixture viscosity.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual bool read()
Read phaseProperties dictionary.
Namespace for OpenFOAM.