tmdlib is hosted by Hepforge, IPPP Durham
TMDlib  1.0.X
TMDlib_UnintegratedGluonKS.h
Go to the documentation of this file.
1 //------------------------------------------------------------------------------
2 //
3 // Copyright (c) 2012, Sebastian Sapeta
4 //
5 // This is part of KSgluon package.
6 //
7 // KSgluon is free software; you can redistribute it and/or modify
8 // it under the terms of the GNU General Public License as published by
9 // the Free Software Foundation; either version 2 of the License, or
10 // (at your option) any later version.
11 //
12 // KSgluon is distributed in the hope that it will be useful,
13 // but WITHOUT ANY WARRANTY; without even the implied warranty of
14 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 // GNU General Public License for more details.
16 //
17 // You should have received a copy of the GNU General Public License
18 // along with FastJet; if not, write to the Free Software
19 // Foundation, Inc.:
20 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
21 //
22 //------------------------------------------------------------------------------
23 #ifndef __TMDlib_UNINTEGRATEDGLUONKS_H__
24 #define __TMDlib_UNINTEGRATEDGLUONKS_H__
25 
26 #include <string>
27 #include <vector>
28 #include "TMDlib_InterpolationKS.h"
29 using namespace std;
30 
31 //------------------------------------------------------------------------------
34  grid2d, // gluon interpolated from 2d grid
35  grid3d // gluon interpolated from 3d grid
36 };
37 
38 //------------------------------------------------------------------------------
74  public:
80  UnintegratedGluon (string filename, gsl_interp_type int_type);
82  //pair<unsigned int,unsigned int> grid_size() const;
83  vector<int> grid_size() const;
89  vector<pair <double, double> > grid_limits() const {
90  return _updf->grid_limits();}
99  double xg(double logx, double logkt2, double logmu2 = 0.0);
100  virtual ~UnintegratedGluon ();
101  private:
102  // function determining the size of the 2d grid
103  void _get_2dgrid_dimensions(string filename);
104  // function determining the size of the 3d grid
105  void _get_3dgrid_dimensions(string filename);
116  void _set_from_2dgrid(string filename);
117  void _set_from_3dgrid(string filename);
118  // functions setting gbw pdf
119  void _set_gbw();
120  UnintegratedGluonName _determine_grid_type(string filename);
121  unsigned int _n1, _n2, _n3;
122  ofstream *_errorfile;
123  static double _norm_F(const vector<double> &x);
124  BaseInterpolation *_updf;
125  UnintegratedGluonName _name;
126  // pointer to the normalization function for pdf
127  double (*_norm)(const vector<double> &x);
128  gsl_interp_type _int_type;
129  double _some_eps;
130 };
131 
132 
133 #endif // __TMDlib_UNINTEGRATEDGLUONKS_H__
vector< pair< double, double > > grid_limits() const
Definition: TMDlib_UnintegratedGluonKS.h:89
Definition: TMDlib_UnintegratedGluonKS.h:73
Definition: TMDlib_InterpolationKS.h:47
string filename
Definition: ksPDF.cc:18
UnintegratedGluonName
enums for gluon names
Definition: TMDlib_UnintegratedGluonKS.h:33
Definition: TMDlib_UnintegratedGluonKS.h:35
Definition: TMDlib_UnintegratedGluonKS.h:34