12#include <yaml-cpp/yaml.h> 
   46    ConvolutionTable(YAML::Node 
const& table, std::vector<double> 
const& cutParam = {}, std::vector<std::shared_ptr<Cut>> 
const& cuts = {}, 
double const& acc = 1e-7);
 
   59    ConvolutionTable(std::string 
const& infile, std::vector<double> 
const& cutParam = {}, std::vector<std::shared_ptr<Cut>> 
const& cuts = {}, 
double const& acc = 1e-7);
 
   67    std::map<double, double> 
ConvoluteDY(std::function<
double(
double const&, 
double const&, 
double const&)> 
const& fNP) 
const;
 
   76    std::map<double, double> 
ConvoluteSIDIS(std::function<
double(
double const&, 
double const&, 
double const&)> 
const& fNP,
 
   77                                            std::function<
double(
double const&, 
double const&, 
double const&)> 
const& DNP) 
const;
 
   86    virtual std::vector<double> 
GetPredictions(std::function<
double(
double const&, 
double const&, 
double const&)> 
const& fNP1,
 
   87                                               std::function<
double(
double const&, 
double const&, 
double const&)> 
const& fNP2) 
const;
 
   98    virtual std::vector<double> 
GetPredictions(std::function<
double(
double const&, 
double const&, 
double const&, 
int const&)> 
const& fNP) 
const;
 
  115    virtual  std::vector<double> 
GetPredictions(std::function<
double(
double const&, 
double const&, 
double const&, 
int const&)> 
const& fNP,
 
  116                                                std::function<
double(
double const&, 
double const&, 
double const&, 
int const&)> 
const& dNP) 
const;
 
  142    std::vector<double>                                                   
const _qTv;      
 
  143    std::vector<std::vector<double>>                                      
const _qTmap;    
 
  148    std::vector<double>                                                   
const _Qg;       
 
  152    std::map<double,std::vector<std::vector<double>>>                           
_PSRed;    
 
  153    std::map<double,std::vector<std::vector<double>>>                           
_dPSRed;   
 
  154    std::map<double,std::vector<std::vector<std::vector<double>>>>              
_WDY;      
 
  155    std::map<double,std::vector<std::vector<std::vector<std::vector<double>>>>> 
_WSIDIS;   
 
  159    std::vector<std::shared_ptr<Cut>>                                           
_cuts;     
 
  168    virtual void SetInputFFs(std::function<std::map<int, double>(
double const &, 
double const &)> 
const &InDistFunc) {};
 
  169    virtual void SetInputFFs(std::function<apfel::Set<apfel::Distribution>(
double const&)> 
const& InDistFunc) {};
 
  170    virtual std::vector<double> 
GetPredictions(std::function<
double(
double const&, 
double const&, 
double const&)> 
const&)
 const { 
return {}; };
 
Class that implements the methods fot the numerical convolution of the interpolation tables with user...
Definition: convolutiontable.h:25
 
std::vector< double > _xig
Grid in xi;.
Definition: convolutiontable.h:149
 
std::string GetName() const
Definition: convolutiontable.h:123
 
std::map< double, std::vector< std::vector< std::vector< std::vector< double > > > > > _WSIDIS
The weights for SIDIS.
Definition: convolutiontable.h:155
 
std::vector< double > GetcutParam() const
Definition: convolutiontable.h:128
 
virtual std::vector< double > GetPredictions(std::function< double(double const &, double const &, double const &)> const &) const
Definition: convolutiontable.h:170
 
std::vector< double > const _cutParam
The parameters needed to compute the ratio qT / Q.
Definition: convolutiontable.h:157
 
std::vector< double > _zg
Grid in xi;.
Definition: convolutiontable.h:151
 
std::vector< std::shared_ptr< Cut > > _cuts
Cut objects.
Definition: convolutiontable.h:159
 
int const _proc
Index of the process (0: DY, 1: SIDIS)
Definition: convolutiontable.h:139
 
virtual std::vector< double > GetPredictions(std::function< double(double const &, double const &, double const &, int const &)> const &fNP, std::function< double(double const &, double const &, double const &, int const &)> const &dNP) const
This function returns a vector of predictions given two user-defined non-perturbative function.
 
int GetProcess() const
Definition: convolutiontable.h:124
 
ConvolutionTable(std::string const &infile, std::vector< double > const &cutParam={}, std::vector< std::shared_ptr< Cut > > const &cuts={}, double const &acc=1e-7)
The "ConvolutionTable" constructor.
 
std::map< double, double > ConvoluteSIDIS(std::function< double(double const &, double const &, double const &)> const &fNP, std::function< double(double const &, double const &, double const &)> const &DNP) const
This function convolutes a SIDIS input convolution table with two user-defined non-perturbative funct...
 
std::map< double, std::vector< std::vector< double > > > _dPSRed
The derivative of the phase-space reduction factors.
Definition: convolutiontable.h:153
 
double const _Vs
Center of mass energy.
Definition: convolutiontable.h:140
 
std::vector< std::vector< double > > GetqTBins() const
Definition: convolutiontable.h:126
 
virtual std::vector< double > GetPredictions(std::function< double(double const &, double const &, double const &, int const &)> const &fNP) const
This function returns a vector of predictions with a single user-defined non-perturbative function.
 
ConvolutionTable(YAML::Node const &table, std::vector< double > const &cutParam={}, std::vector< std::shared_ptr< Cut > > const &cuts={}, double const &acc=1e-7)
The "ConvolutionTable" constructor.
 
virtual void SetInputFFs(std::function< apfel::Set< apfel::Distribution >(double const &)> const &InDistFunc)
Definition: convolutiontable.h:169
 
double const _prefact
Overall prefactor.
Definition: convolutiontable.h:145
 
std::valarray< bool > _cutmask
Mask of points that pass the cuts.
Definition: convolutiontable.h:160
 
std::vector< double > _xbg
Grid in xi;.
Definition: convolutiontable.h:150
 
std::valarray< bool > GetCutMask() const
This function returns the mask of points that pass all the cuts.
Definition: convolutiontable.h:135
 
std::string const _name
Name of the table.
Definition: convolutiontable.h:135
 
double GetCME() const
Definition: convolutiontable.h:125
 
double const _prefact2
Second overall prefactor (do we really need it?)
Definition: convolutiontable.h:146
 
virtual void SetInputFFs(std::function< std::map< int, double >(double const &, double const &)> const &InDistFunc)
Definition: convolutiontable.h:168
 
std::vector< double > const _zOgata
Unscaled Ogata coordinate.
Definition: convolutiontable.h:147
 
ConvolutionTable()
The "ConvolutionTable" constructor.
 
double _acc
The Ogata-quadrature accuracy.
Definition: convolutiontable.h:158
 
std::map< double, double > ConvoluteDY(std::function< double(double const &, double const &, double const &)> const &fNP) const
This function convolutes a Drell-Yan input convolution table with a user-defined non-perturbative fun...
 
std::vector< double > const _Qg
Grid in Q.
Definition: convolutiontable.h:148
 
std::vector< double > const _qTv
Vector of qT bin-bounds.
Definition: convolutiontable.h:142
 
std::vector< double > const _qTfact
Bin-by-bin factors.
Definition: convolutiontable.h:144
 
bool const _IntqT
Whether the bin are integrated in qT or not.
Definition: convolutiontable.h:141
 
std::vector< std::vector< double > > const _qTmap
Vector of bounds for each qT bin.
Definition: convolutiontable.h:143
 
std::map< double, std::vector< std::vector< double > > > _PSRed
The phase-space reduction factors.
Definition: convolutiontable.h:152
 
virtual std::vector< double > GetPredictions(std::function< double(double const &, double const &, double const &)> const &fNP1, std::function< double(double const &, double const &, double const &)> const &fNP2) const
This function returns a vector of predictions given two user-defined non-perturbative functions.
 
std::map< double, std::vector< std::vector< std::vector< double > > > > _WDY
The weights for Drell-Yan.
Definition: convolutiontable.h:154