uniformInterpolationTable.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-2019 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::uniformInterpolationTable
26 
27 Description
28  Table with uniform interval in independent variable, with linear
29  interpolation
30 
31  Example usage (scalar): values specified within a dictionary:
32 
33  \verbatim
34  {
35  x0 0; // lower limit
36  dx 0.2; // fixed interval
37  log10 true; // take log(10) when interpolating?
38  data // list of dependent data values
39  (
40  7870 // value at x0
41  7870 // value at x0 + dx
42  ...
43  7870 // value at x0 + n*dx
44  );
45  }
46  \endverbatim
47 
48 SourceFiles
49  uniformInterpolationTable.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef uniformInterpolationTable_H
54 #define uniformInterpolationTable_H
55 
56 #include "List.H"
57 #include "Switch.H"
58 #include "IOobject.H"
59 #include "objectRegistry.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class uniformInterpolationTable Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 template<class Type>
72 :
73  public IOobject,
74  public List<Type>
75 {
76  // Private Data
77 
78  // Control parameters
79 
80  //- Lower limit
81  scalar x0_;
82 
83  //- Fixed interval
84  scalar dx_;
85 
86  //- Flag to indicate that x data are given in log10(x) form
87  Switch log10_;
88 
89  //- Bound x values
90  Switch bound_;
91 
92 
93  // Private Member Functions
94 
95  //- Check that the table is valid
96  void checkTable() const;
97 
98 
99 public:
100 
101  // Constructors
102 
103  //- Construct from IOobject and readFields flag.
104  // Creates a null object if readFields = false
105  uniformInterpolationTable(const IOobject&, const bool readFields);
106 
107  //- Construct from name, objectRegistry and dictionary.
108  // If initialiseOnly flag is set, control parameters are read from
109  // the dictionary, but not the data table
111  (
112  const word& tableName,
113  const objectRegistry&,
114  const dictionary&,
115  const bool initialiseOnly = false
116  );
117 
118  //- Copy constructor
120 
121 
122  //- Destructor
124 
125 
126  // Member Functions
127 
128  // Access
129 
130  //- Return the lower limit
131  inline scalar x0() const;
132 
133  //- Return the fixed interval
134  inline scalar dx() const;
135 
136  //- Return the log10(x) flag
137  inline const Switch& log10() const;
138 
139  //- Return the bound flag
140  inline const Switch& bound() const;
141 
142 
143  // Edit
144 
145  //- Return the lower limit
146  inline scalar& x0();
147 
148  //- Return the fixed interval
149  inline scalar& dx();
150 
151  //- Return the log10(x) flag
152  inline Switch& log10();
153 
154  //- Return the bound flag
155  inline Switch& bound();
156 
157 
158  // Evaluation
159 
160  //- Return the minimum x value
161  inline scalar xMin() const;
162 
163  //- Return the maximum x value
164  inline scalar xMax() const;
165 
166  //- Interpolate
167  Type interpolate(scalar x) const;
168 
169  //- Interpolate - takes log10 flag into account
170  Type interpolateLog10(scalar x) const;
171 
172 
173  // Override ancestor size() function and [] operator
174 
175  //- Return the size of the table
176  using List<Type>::size;
177 
178  //- Use List[] operator for read/write access
179  using List<Type>::operator[];
180 
181 
182  // I-O
183 
184  //- Write
185  void write() const;
186 
187 
188  // Member Operators
189 
190  //- Disallow default bitwise assignment
191  void operator=(const uniformInterpolationTable&) = delete;
192 };
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace Foam
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 #ifdef NoRepository
206  #include "uniformInterpolationTable.C"
207 #endif
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 #endif
212 
213 // ************************************************************************* //
scalar xMin() const
Return the minimum x value.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Table with uniform interval in independent variable, with linear interpolation.
scalar x0() const
Return the lower limit.
uniformInterpolationTable(const IOobject &, const bool readFields)
Construct from IOobject and readFields flag.
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
A class for handling words, derived from string.
Definition: word.H:59
const Switch & bound() const
Return the bound flag.
scalar xMax() const
Return the maximum x value.
scalar dx() const
Return the fixed interval.
const Switch & log10() const
Return the log10(x) flag.
Registry of regIOobjects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void operator=(const uniformInterpolationTable &)=delete
Disallow default bitwise assignment.
Namespace for OpenFOAM.
Type interpolate(scalar x) const
Interpolate.