1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "HerschelBulkley.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "surfaceFields.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 namespace viscosityModels
37 defineTypeNameAndDebug(HerschelBulkley, 0);
39 addToRunTimeSelectionTable
49 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
51 Foam::tmp<Foam::volScalarField>
52 Foam::viscosityModels::HerschelBulkley::calcNu() const
54 dimensionedScalar tone("tone", dimTime, 1.0);
55 dimensionedScalar rtone("rtone", dimless/dimTime, 1.0);
56 tmp<volScalarField> sr(strainRate());
57 return (min(nu0_,(tau0_ + k_* rtone *( pow(tone * sr(), n_)
58 + pow(tone*tau0_/nu0_,n_))) / (max(sr(), dimensionedScalar
59 ("VSMALL", dimless/dimTime, VSMALL)))));
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 Foam::viscosityModels::HerschelBulkley::HerschelBulkley
68 const dictionary& viscosityProperties,
69 const volVectorField& U,
70 const surfaceScalarField& phi
73 viscosityModel(name, viscosityProperties, U, phi),
74 HerschelBulkleyCoeffs_(viscosityProperties.subDict(typeName + "Coeffs")),
75 k_(HerschelBulkleyCoeffs_.lookup("k")),
76 n_(HerschelBulkleyCoeffs_.lookup("n")),
77 tau0_(HerschelBulkleyCoeffs_.lookup("tau0")),
78 nu0_(HerschelBulkleyCoeffs_.lookup("nu0")),
94 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 bool Foam::viscosityModels::HerschelBulkley::read
98 const dictionary& viscosityProperties
101 viscosityModel::read(viscosityProperties);
103 HerschelBulkleyCoeffs_ = viscosityProperties.subDict(typeName + "Coeffs");
105 HerschelBulkleyCoeffs_.lookup("k") >> k_;
106 HerschelBulkleyCoeffs_.lookup("n") >> n_;
107 HerschelBulkleyCoeffs_.lookup("tau0") >> tau0_;
108 HerschelBulkleyCoeffs_.lookup("nu0") >> nu0_;
114 // ************************************************************************* //