All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
radialModel.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::radialModel
26 
27 SourceFiles
28  radialModel.C
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef radialModel_H
33 #define radialModel_H
34 
35 #include "dictionary.H"
36 #include "volFields.H"
37 #include "dimensionedTypes.H"
38 #include "runTimeSelectionTables.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace kineticTheoryModels
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class radialModel Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 class radialModel
52 {
53 protected:
54 
55  // Protected data
56 
57  const dictionary& dict_;
58 
59 
60 public:
61 
62  //- Runtime type information
63  TypeName("radialModel");
64 
65  // Declare runtime constructor selection table
67  (
68  autoPtr,
70  dictionary,
71  (
72  const dictionary& dict
73  ),
74  (dict)
75  );
76 
77 
78  // Constructors
79 
80  //- Construct from components
81  radialModel(const dictionary& dict);
82 
83  //- Disallow default bitwise copy construction
84  radialModel(const radialModel&) = delete;
85 
86 
87  // Selectors
88 
90  (
91  const dictionary& dict
92  );
93 
94 
95  //- Destructor
96  virtual ~radialModel();
97 
98 
99  // Member Functions
100 
101  //- Radial distribution function
102  virtual tmp<volScalarField> g0
103  (
104  const volScalarField& alpha,
105  const dimensionedScalar& alphaMinFriction,
107  ) const = 0;
108 
109  //- Derivative of the radial distribution function
111  (
112  const volScalarField& alpha,
113  const dimensionedScalar& alphaMinFriction,
115  ) const = 0;
117  virtual bool read()
118  {
119  return true;
120  }
121 
122 
123  // Member Operators
124 
125  //- Disallow default bitwise assignment
126  void operator=(const radialModel&) = delete;
127 };
128 
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 } // End namespace kineticTheoryModels
133 } // End namespace Foam
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 #endif
138 
139 // ************************************************************************* //
virtual tmp< volScalarField > g0(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const =0
Radial distribution function.
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void operator=(const radialModel &)=delete
Disallow default bitwise assignment.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
TypeName("radialModel")
Runtime type information.
declareRunTimeSelectionTable(autoPtr, radialModel, dictionary,(const dictionary &dict),(dict))
static autoPtr< radialModel > New(const dictionary &dict)
radialModel(const dictionary &dict)
Construct from components.
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
virtual tmp< volScalarField > g0prime(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const =0
Derivative of the radial distribution function.
Namespace for OpenFOAM.