Quemada.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) 2022-2023 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 "Quemada.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace mixtureViscosityModels
35 {
37 
39  (
41  Quemada,
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
53 )
54 :
56  q_(optionalSubDict(typeName + "Coeffs").lookupOrDefault("q", scalar(2))),
57  muMax_
58  (
59  "muMax",
61  optionalSubDict(typeName + "Coeffs").lookup("muMax")
62  )
63 {}
64 
65 
66 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
67 
70 (
71  const volScalarField& muc,
72  const volVectorField& U
73 ) const
74 {
75  return min
76  (
77  muc*pow(max(1 - mixture_.alphad()/mixture_.alphaMax(), small), -q_),
78  muMax_
79  );
80 }
81 
82 
84 {
86  {
87  const dictionary& dict = optionalSubDict(typeName + "Coeffs");
88 
89  dict.lookup("q") >> q_;
90  dict.lookup("muMax") >> muMax_;
91 
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 
101 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Class to represent a mixture of two constant density phases.
An abstract base class for incompressible mixtureViscosityModels.
Quemada viscosity model for colloidal dispersions.
Definition: Quemada.H:75
virtual tmp< volScalarField > mu(const volScalarField &muc, const volVectorField &U) const
Return the mixture viscosity.
Definition: Quemada.C:70
Quemada(const incompressibleDriftFluxMixture &mixture)
Construct from mixture.
Definition: Quemada.C:51
virtual bool read()
Read phaseProperties dictionary.
Definition: Quemada.C:83
A class for managing temporary objects.
Definition: tmp.H:55
virtual bool read()=0
Read physicalProperties dictionary.
U
Definition: pEqn.H:72
defineTypeNameAndDebug(BinghamPlastic, 0)
addToRunTimeSelectionTable(mixtureViscosityModel, BinghamPlastic, dictionary)
Namespace for OpenFOAM.
const dimensionSet dimDynamicViscosity
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict