frictionalStressModel.H
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) 2011-2020 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 Class
25  Foam::kineticTheoryModels::frictionalStressModel
26 
27 SourceFiles
28  frictionalStressModel.C
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef frictionalStressModel_H
33 #define frictionalStressModel_H
34 
35 #include "dictionary.H"
36 #include "volFields.H"
37 #include "dimensionedTypes.H"
38 #include "runTimeSelectionTables.H"
39 #include "phaseModel.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace kineticTheoryModels
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class frictionalStressModel Declaration
50 \*---------------------------------------------------------------------------*/
51 
53 {
54 protected:
55 
56  // Protected data
57 
58  //- Reference to higher-level dictionary for re-read
59  const dictionary& dict_;
60 
61 
62 public:
63 
64  //- Runtime type information
65  TypeName("frictionalStressModel");
66 
67  // Declare runtime constructor selection table
69  (
70  autoPtr,
72  dictionary,
73  (
74  const dictionary& dict
75  ),
76  (dict)
77  );
78 
79 
80  // Constructors
81 
82  //- Construct from components
83  frictionalStressModel(const dictionary& dict);
84 
85  //- Disallow default bitwise copy construction
87 
88 
89  // Selectors
90 
92  (
93  const dictionary& dict
94  );
95 
96 
97  //- Destructor
98  virtual ~frictionalStressModel();
99 
100 
101  // Member Functions
102 
104  (
105  const phaseModel& phase,
106  const dimensionedScalar& alphaMinFriction,
108  ) const = 0;
109 
111  (
112  const phaseModel& phase,
113  const dimensionedScalar& alphaMinFriction,
115  ) const = 0;
116 
117  virtual tmp<volScalarField> nu
118  (
119  const phaseModel& phase,
120  const dimensionedScalar& alphaMinFriction,
122  const volScalarField& pf,
123  const volSymmTensorField& D
124  ) const = 0;
125 
126  virtual bool read() = 0;
127 
128 
129  // Member Operators
130 
131  //- Disallow default bitwise assignment
132  void operator=(const frictionalStressModel&) = delete;
133 };
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace kineticTheoryModels
139 } // End namespace Foam
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 #endif
144 
145 // ************************************************************************* //
dictionary dict
declareRunTimeSelectionTable(autoPtr, frictionalStressModel, dictionary,(const dictionary &dict),(dict))
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
TypeName("frictionalStressModel")
Runtime type information.
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const =0
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const =0
const dictionary & dict_
Reference to higher-level dictionary for re-read.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
static autoPtr< frictionalStressModel > New(const dictionary &dict)
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
frictionalStressModel(const dictionary &dict)
Construct from components.
void operator=(const frictionalStressModel &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const =0
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Namespace for OpenFOAM.