JohnsonJacksonFrictionalStress.H
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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson
28
29Description
30
31SourceFiles
32 JohnsonJacksonFrictionalStress.C
33
34\*---------------------------------------------------------------------------*/
35
36#ifndef JohnsonJackson_H
37#define JohnsonJackson_H
38
39#include "frictionalStressModel.H"
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace Foam
44{
45namespace kineticTheoryModels
46{
47namespace frictionalStressModels
48{
49
50/*---------------------------------------------------------------------------*\
51 Class JohnsonJackson Declaration
52\*---------------------------------------------------------------------------*/
53
54class JohnsonJackson
55:
56 public frictionalStressModel
57{
58 // Private data
59
60 dictionary coeffDict_;
61
62 //- Material constant for frictional normal stress
64
65 //- Material constant for frictional normal stress
67
68 //- Material constant for frictional normal stress
70
71 //- Angle of internal friction
73
74 //- Lower limit for (alphaMax - alpha1)
75 dimensionedScalar alphaDeltaMin_;
76
77
78public:
79
80 //- Runtime type information
81 TypeName("JohnsonJackson");
82
83
84 // Constructors
85
86 //- Construct from components
88
89
90 //- Destructor
91 virtual ~JohnsonJackson();
92
93
94 // Member functions
97 (
98 const phaseModel& phase,
99 const dimensionedScalar& alphaMinFriction,
101 ) const;
104 (
105 const phaseModel& phase,
106 const dimensionedScalar& alphaMinFriction,
108 ) const;
110 virtual tmp<volScalarField> nu
111 (
112 const phaseModel& phase,
113 const dimensionedScalar& alphaMinFriction,
115 const volScalarField& pf,
116 const volSymmTensorField& D
117 ) const;
119 virtual bool read();
120};
121
122
123// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124
125} // End namespace frictionalStressModels
126} // End namespace kineticTheoryModels
127} // End namespace Foam
128
129// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130
131#endif
132
133// ************************************************************************* //
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
JohnsonJackson(const dictionary &dict)
Construct from components.
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
TypeName("JohnsonJackson")
Runtime type information.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & nu
dictionary dict
const dimensionedScalar & D
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73