diameterModel.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::diameterModel
26 
27 Description
28  Abstract base-class for dispersed-phase particle diameter models.
29 
30 SourceFiles
31  diameterModel.C
32  diameterModelNew.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef diameterModel_H
37 #define diameterModel_H
38 
39 #include "dictionary.H"
40 #include "phaseModel.H"
41 #include "runTimeSelectionTables.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class diameterModel Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 class diameterModel
53 {
54  // Private Data
55 
56  //- The phase diameter properties dictionary
57  dictionary diameterProperties_;
58 
59  //- The phase that this model applies
60  const phaseModel& phase_;
61 
62  //- Optionally stored diameter field
64 
65  //- Optionally stored surface area per unit volume field
67 
68 
69 protected:
70 
71  // Access
72 
73  //- Get a reference to the stored diameter field
75 
76  //- Get a reference to the stored surface area per unit volume field
78 
79  //- Get the diameter field
80  virtual tmp<volScalarField> calcD() const = 0;
81 
82  //- Get the surface area per unit volume field
83  virtual tmp<volScalarField> calcA() const = 0;
84 
85 
86 public:
87 
88  //- Runtime type information
89  TypeName("diameterModel");
90 
91 
92  // Declare runtime construction
93 
95  (
96  autoPtr,
98  dictionary,
99  (
101  const phaseModel& phase
102  ),
103  (diameterProperties, phase)
104  );
105 
106 
107  // Constructors
108 
110  (
111  const dictionary& diameterProperties,
112  const phaseModel& phase
113  );
114 
115 
116  //- Destructor
117  virtual ~diameterModel();
118 
119 
120  // Selectors
121 
123  (
124  const dictionary& diameterProperties,
125  const phaseModel& phase
126  );
127 
128 
129  // Member Functions
130 
131  //- Return the phase diameter properties dictionary
132  const dictionary& diameterProperties() const
133  {
134  return diameterProperties_;
135  }
136 
137  //- Return the phase
138  const phaseModel& phase() const
139  {
140  return phase_;
141  }
142 
143  //- Return the diameter
144  tmp<volScalarField> d() const;
145 
146  //- Return the surface area per unit volume
147  tmp<volScalarField> a() const;
148 
149  //- Correct the diameter field
150  virtual void correct();
151 
152  //- Read phaseProperties dictionary
153  virtual bool read(const dictionary& phaseProperties);
154 };
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace Foam
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 #endif
164 
165 // ************************************************************************* //
diameterModel(const dictionary &diameterProperties, const phaseModel &phase)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual ~diameterModel()
Destructor.
Helper class to manage multi-specie phase properties.
static autoPtr< diameterModel > New(const dictionary &diameterProperties, const phaseModel &phase)
tmp< volScalarField > a() const
Return the surface area per unit volume.
volScalarField & aRef()
Get a reference to the stored surface area per unit volume field.
TypeName("diameterModel")
Runtime type information.
const phaseModel & phase() const
Return the phase.
virtual tmp< volScalarField > calcA() const =0
Get the surface area per unit volume field.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:51
virtual tmp< volScalarField > calcD() const =0
Get the diameter field.
volScalarField & dRef()
Get a reference to the stored diameter field.
tmp< volScalarField > d() const
Return the diameter.
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
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
declareRunTimeSelectionTable(autoPtr, diameterModel, dictionary,(const dictionary &diameterProperties, const phaseModel &phase),(diameterProperties, phase))
virtual void correct()
Correct the diameter field.
Namespace for OpenFOAM.