CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
FuncGrid.h
Go to the documentation of this file.
1 /*******************************************************************
2  * Copyright (C) 2010 by Federico Marulli and Alfonso Veropalumbo *
3  * federico.marulli3@unibo.it *
4  * *
5  * This program is free software; you can redistribute it and/or *
6  * modify it under the terms of the GNU General Public License as *
7  * published by the Free Software Foundation; either version 2 of *
8  * the License, or (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful,*
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public *
16  * License along with this program; if not, write to the Free *
17  * Software Foundation, Inc., *
18  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
19  *******************************************************************/
20 
31 #ifndef __FUNCGRID__
32 #define __FUNCGRID__
33 
34 #include "GSLwrapper.h"
35 #include "CUBAwrapper.h"
36 
37 
38 // =====================================================================================
39 
40 
41 namespace cbl {
42 
43  namespace glob {
44 
54  class FuncGrid
55  {
56 
57  private:
58 
60  std::vector<double> m_x;
61 
63  std::vector<double> m_y;
64 
66  size_t m_size;
67 
69  std::string m_interpType;
70 
72  std::shared_ptr<gsl_spline> m_spline;
73 
75  const gsl_interp_type *m_type;
76 
78  std::shared_ptr<gsl_interp_accel> m_acc;
79 
81  double m_xmin;
82 
84  double m_xmax;
85 
88 
89 
90  public:
91 
96 
100  FuncGrid () = default;
101 
111  FuncGrid (const std::vector<double> x, const std::vector<double> y, const std::string interpType, const BinType bin_type=BinType::_linear_);
112 
117  ~FuncGrid () = default;
118 
120 
121 
126 
131  std::vector<double> x () const { return m_x; }
132 
139  double x (const int i) const { return m_x[i]; }
140 
145  std::vector<double> y () const { return m_y; }
146 
153  double y (const int i) const { return m_y[i]; }
154 
159  double size () const { return m_size; }
160 
165  double xmin () const { return m_xmin; }
166 
171  double xmax () const { return m_xmax; }
172 
174 
175 
180 
185  void free ();
186 
193  double operator () (const double xx) const;
194 
201  std::vector<double> eval_func (const std::vector<double> xx) const;
202 
208  double D1v (const double xx) const;
209 
215  double D2v (const double xx) const;
216 
227  double integrate_qag (const double a, const double b, const double rel_err=1.e-2, const double abs_err=1.e-6, const int limit_size=1000, const int rule=6);
228 
239  double integrate_cquad (const double a, const double b, const double rel_err=1.e-2, const double abs_err=1.e-6, const int nevals=100);
240 
254  double integrate_qaws (const double a, const double b, const double alpha=0, const double beta=0, const int mu=0, const int nu=0, const double rel_err=1.e-2, const double abs_err=1.e-6, const int limit_size=1000);
255 
265  double root (const double x_low, const double x_up, const double fx0=0, const double rel_err=1.e-2, const double abs_err=1.e-6);
266 
277  double root_D1v (const double x_low, const double x_up, const double fx0=0, const double rel_err=1.e-2, const double abs_err=1.e-6);
278 
289  double root_D2v (const double x_low, const double x_up, const double fx0=0, const double rel_err=1.e-2, const double abs_err=1.e-6);
290  };
291 
302  {
303 
304  private:
305 
307  std::vector<double> m_x;
308 
310  std::vector<double> m_y;
311 
313  std::shared_ptr<double> m_fxy;
314 
316  size_t m_size_x;
317 
319  size_t m_size_y;
320 
322  std::string m_interpType;
323 
325  const gsl_interp2d_type *m_type;
326 
328  std::shared_ptr<gsl_interp_accel> m_acc_x;
329 
331  std::shared_ptr<gsl_interp_accel> m_acc_y;
332 
334  std::shared_ptr<gsl_spline2d> m_spline;
335 
337  double m_xmin;
338 
340  double m_xmax;
341 
343  double m_ymin;
344 
346  double m_ymax;
347 
348  public:
349 
354 
358  FuncGrid2D () = default;
359 
369  FuncGrid2D (const std::vector<double> x, const std::vector<double> y, const std::vector<std::vector<double>> fxy, const std::string interpType);
370 
375  ~FuncGrid2D () = default;
376 
378 
379 
384 
389  std::vector<double> x () const { return m_x; }
390 
397  double x (const int i) const { return m_x[i]; }
398 
403  std::vector<double> y () const { return m_y; }
404 
411  double y (const int i) const { return m_y[i]; }
412 
417  double size_x () const { return m_size_x; }
418 
423  double xmin () const { return m_xmin; }
424 
429  double xmax () const { return m_xmax; }
430 
435  double size_y () const { return m_size_y; }
436 
441  double ymin () const { return m_ymin; }
442 
447  double ymax () const { return m_ymax; }
448 
450 
451 
456 
461  void free ();
462 
471  double operator () (const double xx, const double yy) const;
472 
479  std::vector<double> eval_func (const std::vector<std::vector<double>> xx) const;
480 
496  double IntegrateVegas (const double xmin, const double xmax, const double ymin, const double ymax) const;
497 
513  double IntegrateSuave (const double xmin, const double xmax, const double ymin, const double ymax) const;
514 
530  double IntegrateDivonne (const double xmin, const double xmax, const double ymin, const double ymax) const;
531 
547  double IntegrateCuhre (const double xmin, const double xmax, const double ymin, const double ymax) const;
548 
549  };
550 
551 
552  }
553 }
554 
555 #endif
class CUBAwrapper that wrap CUBA routines for multidimensional integration
functions that wrap GSL routines for integration, root finding and minimization
The class FuncGrid2D.
Definition: FuncGrid.h:302
double xmax() const
get the private member FuncGrid::m_xmax
Definition: FuncGrid.h:429
double ymax() const
get the private member FuncGrid::m_ymax
Definition: FuncGrid.h:447
double xmin() const
get the private member FuncGrid::m_xmin
Definition: FuncGrid.h:423
std::vector< double > x() const
get the private member FuncGrid::m_x
Definition: FuncGrid.h:389
std::shared_ptr< gsl_interp_accel > m_acc_x
GSL accelerator object.
Definition: FuncGrid.h:328
const gsl_interp2d_type * m_type
GSL object used to set the interpolation type.
Definition: FuncGrid.h:325
double IntegrateSuave(const double xmin, const double xmax, const double ymin, const double ymax) const
evaluate the 2D integral of the interpolated function using the Suave routine from CUBA libraries
Definition: FuncGrid.cpp:371
FuncGrid2D()=default
default constructor
std::vector< double > m_x
x values
Definition: FuncGrid.h:307
double operator()(const double xx, const double yy) const
overloading of the () operator
Definition: FuncGrid.cpp:319
double m_xmin
minimum x value
Definition: FuncGrid.h:337
double IntegrateCuhre(const double xmin, const double xmax, const double ymin, const double ymax) const
evaluate the 2D integral of the interpolated function using the Cuhre routine from CUBA libraries
Definition: FuncGrid.cpp:403
double m_xmax
maximum x value
Definition: FuncGrid.h:340
std::vector< double > m_y
y values
Definition: FuncGrid.h:310
std::shared_ptr< gsl_spline2d > m_spline
GSL object used to interpolate.
Definition: FuncGrid.h:334
double y(const int i) const
get the i-th element of the private member FuncGrid::m_y
Definition: FuncGrid.h:411
std::string m_interpType
method used to interpolate
Definition: FuncGrid.h:322
void free()
free the GSL objects
Definition: FuncGrid.cpp:305
double size_y() const
get the private member FuncGrid::m_size
Definition: FuncGrid.h:435
double IntegrateVegas(const double xmin, const double xmax, const double ymin, const double ymax) const
evaluate the 2D integral of the interpolated function using the Vegas routine from CUBA libraries
Definition: FuncGrid.cpp:355
size_t m_size_x
size of the x vector
Definition: FuncGrid.h:316
size_t m_size_y
size of the x vector
Definition: FuncGrid.h:319
std::shared_ptr< gsl_interp_accel > m_acc_y
GSL accelerator object.
Definition: FuncGrid.h:331
double m_ymin
minimum x value
Definition: FuncGrid.h:343
double IntegrateDivonne(const double xmin, const double xmax, const double ymin, const double ymax) const
evaluate the 2D integral of the interpolated function using the Divonne routine from CUBA libraries
Definition: FuncGrid.cpp:387
std::vector< double > y() const
get the private member FuncGrid::m_y
Definition: FuncGrid.h:403
double size_x() const
get the private member FuncGrid::m_size
Definition: FuncGrid.h:417
~FuncGrid2D()=default
default destructor
std::shared_ptr< double > m_fxy
y values
Definition: FuncGrid.h:313
double x(const int i) const
get the i-th element of the private member FuncGrid::m_x
Definition: FuncGrid.h:397
double m_ymax
maximum x value
Definition: FuncGrid.h:346
double ymin() const
get the private member FuncGrid::m_ymin
Definition: FuncGrid.h:441
std::vector< double > eval_func(const std::vector< std::vector< double >> xx) const
evaluate the function at the xx points
Definition: FuncGrid.cpp:341
The class FuncGrid.
Definition: FuncGrid.h:55
~FuncGrid()=default
default destructor
BinType m_binType
bin type
Definition: FuncGrid.h:87
double size() const
get the private member FuncGrid::m_size
Definition: FuncGrid.h:159
std::vector< double > m_y
y values
Definition: FuncGrid.h:63
const gsl_interp_type * m_type
GSL object used to set the interpolation type.
Definition: FuncGrid.h:75
void free()
free the GSL objects
Definition: FuncGrid.cpp:110
std::shared_ptr< gsl_spline > m_spline
GSL object used to interpolate.
Definition: FuncGrid.h:72
std::vector< double > eval_func(const std::vector< double > xx) const
evaluate the function at the xx points
Definition: FuncGrid.cpp:153
double integrate_cquad(const double a, const double b, const double rel_err=1.e-2, const double abs_err=1.e-6, const int nevals=100)
compute the definite integral with GSL cquad method
Definition: FuncGrid.cpp:198
double integrate_qaws(const double a, const double b, const double alpha=0, const double beta=0, const int mu=0, const int nu=0, const double rel_err=1.e-2, const double abs_err=1.e-6, const int limit_size=1000)
compute the definite integral with GSL qaws method
Definition: FuncGrid.cpp:209
double root(const double x_low, const double x_up, const double fx0=0, const double rel_err=1.e-2, const double abs_err=1.e-6)
find roots with GSL brent method
Definition: FuncGrid.cpp:220
double m_xmax
maximum x value
Definition: FuncGrid.h:84
double D2v(const double xx) const
compute the second derivative at xx
Definition: FuncGrid.cpp:178
double operator()(const double xx) const
overloading of the () operator
Definition: FuncGrid.cpp:120
std::shared_ptr< gsl_interp_accel > m_acc
GSL accelerator object.
Definition: FuncGrid.h:78
double D1v(const double xx) const
compute the first derivative at xx
Definition: FuncGrid.cpp:167
std::vector< double > y() const
get the private member FuncGrid::m_y
Definition: FuncGrid.h:145
size_t m_size
size of the x,y vectors, i.e. the grid size
Definition: FuncGrid.h:66
double root_D1v(const double x_low, const double x_up, const double fx0=0, const double rel_err=1.e-2, const double abs_err=1.e-6)
find roots with GSL brent method for the first derivative
Definition: FuncGrid.cpp:231
double integrate_qag(const double a, const double b, const double rel_err=1.e-2, const double abs_err=1.e-6, const int limit_size=1000, const int rule=6)
compute the definite integral with GSL qag method
Definition: FuncGrid.cpp:187
double root_D2v(const double x_low, const double x_up, const double fx0=0, const double rel_err=1.e-2, const double abs_err=1.e-6)
find roots with GSL brent method for the second derivative
Definition: FuncGrid.cpp:242
double xmax() const
get the private member FuncGrid::m_xmax
Definition: FuncGrid.h:171
double xmin() const
get the private member FuncGrid::m_xmin
Definition: FuncGrid.h:165
std::vector< double > x() const
get the private member FuncGrid::m_x
Definition: FuncGrid.h:131
FuncGrid()=default
default constructor
double m_xmin
minimum x value
Definition: FuncGrid.h:81
std::string m_interpType
method used to interpolate
Definition: FuncGrid.h:69
std::vector< double > m_x
x values
Definition: FuncGrid.h:60
double y(const int i) const
get the i-th element of the private member FuncGrid::m_y
Definition: FuncGrid.h:153
double x(const int i) const
get the i-th element of the private member FuncGrid::m_x
Definition: FuncGrid.h:139
static const double alpha
: the fine-structure constant
Definition: Constants.h:233
The global namespace of the CosmoBolognaLib
Definition: CAMB.h:38
BinType
the binning type
Definition: Kernel.h:505
@ _linear_
linear binning