EOL BSpline Library  v2.2
All Classes Files Functions Variables Enumerations Enumerator Pages
BSplineBase.h
Go to the documentation of this file.
1 /************************************************************************
2  * Copyright 2022 University Corporation for Atmospheric Research.
3  * All rights reserved.
4  *
5  * Use of this code is subject to the standard BSD license:
6  *
7  * http://www.opensource.org/licenses/bsd-license.html
8  *
9  * See the COPYRIGHT file in the source distribution for the license text,
10  * or see this web page:
11  *
12  * http://www.eol.ucar.edu/homes/granger/bspline/doc/
13  *
14  *************************************************************************/
15 #ifndef BSPLINEBASE_H_
16 #define BSPLINEBASE_H_
17 
18 #include <memory>
19 
26 template <class T> class BSpline;
27 
28 /*
29  * Opaque member structure to hide the matrix implementation.
30  */
31 template <class T> struct BSplineBaseP;
32 
166 template <class T> class BSplineBase
167 {
168 public:
169  // Datum type
170  typedef T datum_type;
171 
173  static const char* Version();
174 
180  static bool Debug(int on = -1);
181 
186  {
192  BC_ZERO_SECOND = 2
193  };
194 
195 public:
202  BSplineBase(const T* x, int nx, double wl, int bc_type = BC_ZERO_SECOND,
203  int num_nodes = 0);
204 
206  BSplineBase(const BSplineBase& bb);
207 
209  BSplineBase& operator=(const BSplineBase& right);
210 
211  // Declare the destructor so it is not implicitly instantiated. It can
212  // only be instantiated with the full implementation in BSplineBase.cpp.
213  ~BSplineBase();
214 
240  bool setDomain(const T* x, int nx, double wl, int bc_type = BC_ZERO_SECOND,
241  int num_nodes = 0);
242 
250  BSpline<T>* apply(const T* y);
251 
257  const T* nodes(int* nnodes);
258 
262  int
264  {
265  return M + 1;
266  }
267 
271  int
272  nX()
273  {
274  return NX;
275  }
276 
278  T
280  {
281  return xmin;
282  }
283 
285  T
287  {
288  return xmin + (M * DX);
289  }
290 
295  double Alpha(double wavelength);
296 
300  double
302  {
303  return alpha;
304  }
305 
314  bool
315  ok()
316  {
317  return OK;
318  }
319 
320 protected:
321  typedef BSplineBaseP<T> Base;
322 
323  // Provided
324  double waveLength; // Cutoff wavelength (l sub c)
325  int NX;
326  int K; // Degree of derivative constraint (currently fixed at 2)
327  int BC; // Boundary conditions type (0,1,2)
328 
329  // Derived
330  T xmax;
331  T xmin;
332  int M; // Number of intervals (M+1 nodes)
333  double DX; // Interval length in same units as X
334  double alpha;
335  bool OK;
336  // Hide more complicated state members from the public interface.
337  std::unique_ptr<Base> base;
338 
339  bool Setup(int num_nodes = 0);
340  void calculateQ();
341  double qDelta(int m1, int m2);
342  double Beta(int m);
343  void addP();
344  bool factor();
345  double Basis(int m, T x);
346  double DBasis(int m, T x);
347 
348  static const double BoundaryConditions[3][4];
349  static const double PI;
350 
351  double Ratiod(int&, double&, double&);
352 };
353 
354 #endif /*BSPLINEBASE_H_*/
Definition: BSplineBase.h:167
static const char * Version()
Return a string describing the bspline library version.
Definition: BSplineBase.cpp:127
double Alpha()
Definition: BSplineBase.h:301
double Basis(int m, T x)
Definition: BSplineBase.cpp:319
double Beta(int m)
Definition: BSplineBase.cpp:290
BoundaryConditionTypes
Definition: BSplineBase.h:186
@ BC_ZERO_FIRST
Set the first derivative of the spline to zero at the endpoints.
Definition: BSplineBase.h:190
@ BC_ZERO_ENDPOINTS
Set the endpoints of the spline to zero.
Definition: BSplineBase.h:188
@ BC_ZERO_SECOND
Set the second derivative to zero.
Definition: BSplineBase.h:192
static bool Debug(int on=-1)
Definition: BSplineBase.cpp:115
T Xmin()
Minimum x value found.
Definition: BSplineBase.h:279
static const double BoundaryConditions[3][4]
Definition: BSplineBase.h:348
int nX()
Definition: BSplineBase.h:272
double qDelta(int m1, int m2)
Definition: BSplineBase.cpp:382
bool Setup(int num_nodes=0)
Setup number of nodes and deltax for the given domain and cutoff wavelength.
Definition: BSplineBase.cpp:561
BSpline< T > * apply(const T *y)
Definition: BSplineBase.cpp:308
bool ok()
Definition: BSplineBase.h:315
const T * nodes(int *nnodes)
Definition: BSplineBase.cpp:694
BSplineBase & operator=(const BSplineBase &right)
Assignment operator.
Definition: BSplineBase.cpp:142
double DBasis(int m, T x)
Definition: BSplineBase.cpp:348
int nNodes()
Definition: BSplineBase.h:263
BSplineBase(const T *x, int nx, double wl, int bc_type=BC_ZERO_SECOND, int num_nodes=0)
Definition: BSplineBase.cpp:170
bool setDomain(const T *x, int nx, double wl, int bc_type=BC_ZERO_SECOND, int num_nodes=0)
Definition: BSplineBase.cpp:189
T Xmax()
Maximum x value found.
Definition: BSplineBase.h:286
Definition: BSpline.h:32
Private state structure, hiding use of matrix template classes.
Definition: BSplineBase.cpp:88