radialModel.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 | Copyright (C) 2011-2015 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  // Private member functions
54 
55  //- Disallow default bitwise copy construct
56  radialModel(const radialModel&);
57 
58  //- Disallow default bitwise assignment
59  void operator=(const radialModel&);
60 
61 
62 protected:
63 
64  // Protected data
65 
66  const dictionary& dict_;
67 
68 
69 public:
70 
71  //- Runtime type information
72  TypeName("radialModel");
73 
74  // Declare runtime constructor selection table
76  (
77  autoPtr,
78  radialModel,
79  dictionary,
80  (
81  const dictionary& dict
82  ),
83  (dict)
84  );
85 
86 
87  // Constructors
88 
89  //- Construct from components
90  radialModel(const dictionary& dict);
91 
92 
93  // Selectors
94 
95  static autoPtr<radialModel> New
96  (
97  const dictionary& dict
98  );
99 
100 
101  //- Destructor
102  virtual ~radialModel();
103 
104 
105  // Member Functions
106 
107  //- Radial distribution function
108  virtual tmp<volScalarField> g0
109  (
110  const volScalarField& alpha,
111  const dimensionedScalar& alphaMinFriction,
113  ) const = 0;
114 
115  //- Derivative of the radial distribution function
116  virtual tmp<volScalarField> g0prime
117  (
118  const volScalarField& alpha,
119  const dimensionedScalar& alphaMinFriction,
121  ) const = 0;
123  virtual bool read()
124  {
125  return true;
126  }
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
TypeName("radialModel")
Runtime type information.
declareRunTimeSelectionTable(autoPtr, radialModel, dictionary,(const dictionary &dict),(dict))
static autoPtr< radialModel > New(const dictionary &dict)
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Macros to ease declaration of run-time selection tables.
virtual tmp< volScalarField > g0prime(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const =0
Derivative of the radial distribution function.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.