1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2011 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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "HerschelBulkley.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "surfaceFields.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace viscosityModels
36 defineTypeNameAndDebug(HerschelBulkley, 0);
38 addToRunTimeSelectionTable
48 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50 Foam::tmp<Foam::volScalarField>
51 Foam::viscosityModels::HerschelBulkley::calcNu() const
53 dimensionedScalar tone("tone", dimTime, 1.0);
54 dimensionedScalar rtone("rtone", dimless/dimTime, 1.0);
56 tmp<volScalarField> sr(strainRate());
63 // (tau0_ + k_*rtone*(pow(tone*sr(), n_) - pow(tone*tau0_/nu0_, n_)))
64 // /max(sr(), dimensionedScalar("VSMALL", dimless/dimTime, VSMALL))
73 (tau0_ + k_*rtone*pow(tone*sr(), n_))
74 /(max(sr(), dimensionedScalar ("VSMALL", dimless/dimTime, VSMALL)))
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 Foam::viscosityModels::HerschelBulkley::HerschelBulkley
85 const dictionary& viscosityProperties,
86 const volVectorField& U,
87 const surfaceScalarField& phi
90 viscosityModel(name, viscosityProperties, U, phi),
91 HerschelBulkleyCoeffs_(viscosityProperties.subDict(typeName + "Coeffs")),
92 k_(HerschelBulkleyCoeffs_.lookup("k")),
93 n_(HerschelBulkleyCoeffs_.lookup("n")),
94 tau0_(HerschelBulkleyCoeffs_.lookup("tau0")),
95 nu0_(HerschelBulkleyCoeffs_.lookup("nu0")),
101 U_.time().timeName(),
111 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
113 bool Foam::viscosityModels::HerschelBulkley::read
115 const dictionary& viscosityProperties
118 viscosityModel::read(viscosityProperties);
120 HerschelBulkleyCoeffs_ = viscosityProperties.subDict(typeName + "Coeffs");
122 HerschelBulkleyCoeffs_.lookup("k") >> k_;
123 HerschelBulkleyCoeffs_.lookup("n") >> n_;
124 HerschelBulkleyCoeffs_.lookup("tau0") >> tau0_;
125 HerschelBulkleyCoeffs_.lookup("nu0") >> nu0_;
131 // ************************************************************************* //