|
Nanga Parbat 1.1.0
A TMD fitting framework
|
Classes | |
| class | ChiSquare |
| The "ChiSquare" class computes the χ2's given a set of "DataHandler" objects and the corresponding "ConvolutionTable" objects. The computation depends of the non-perturbative functions given as input. More... | |
| class | ConvolutionTable |
| Class that implements the methods fot the numerical convolution of the interpolation tables with user-defined non-perturbative functions. More... | |
| class | Cut |
| Purely virtual mother class that implements the main feautures of a cut function on a given dataset. Given a DataHandler object as an input, it returns a mask of bools that selects the points that pass the specific cut. More... | |
| class | CutFactory |
| class | DataHandler |
| The "DataHandler" class provides a common interface to all datasets. It provides methods to get kinematics, central values, uncertainties, etc. More... | |
| class | DWS |
| Davies-Webber-Stirling parameterisation derived from the "Parameterisation" mother class. More... | |
| class | FastInterface |
| Class that implements the methods of the computation of the interpolation tables. More... | |
| class | FcnCeres |
| The "FcnCeres" class to be used with ceres-solver to minimise the χ2. More... | |
| class | FcnMinuit |
| The "FcnMinuit" class is a derived class of of the "FCNBase" mother class of Minuit that returns the χ2 to be minimised. More... | |
| class | FcnMinuitGrad |
| The "FcnMinuit" class is a derived class of of the "FCNGradientBase" mother class of Minuit that returns the χ2 to be minimised along with it derivatives w.r.t. the free parameters. More... | |
| struct | FourDGrid |
| Structure that contains the 4D grid in Q, x, z and qT. More... | |
| class | MAP21test |
| MAP 2021 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | MeanReplica |
| Parameterisation derived from the "Parameterisation" mother class to compute the mean replica of a Monte Carlo set. More... | |
| class | Parameterisation |
| Mother class that implements the main feautures of a functional parameterisation of non-perturbative functions. More... | |
| class | PV17 |
| Pavia 2017 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | PV19 |
| Pavia 2019 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | PV19b |
| Pavia 2019 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | PV19x |
| Pavia 2019 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | PV20Sivers |
| Pavia 2017 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | QCut |
| Derivation of the class Cut to impose cut on the variable "Q". More... | |
| class | QGG13 |
| Pavia 2019 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | QGG6 |
| Pavia 2019 parameterisation derived from the "Parameterisation" mother class. More... | |
| class | StructGrid |
| Class for the interpolation of a single TMD grid. More... | |
| struct | ThreeDGrid |
| Structure that contains the 3D grid in Q, x, and qT. More... | |
| class | TMDGrid |
| Class for the interpolation of a single TMD grid. More... | |
| class | TrainingCut |
| Derivation of the class Cut to impose Cross Validation cut. More... | |
| class | ZCut |
| Derivation of the class Cut to impose cut on the variable "z". More... | |
Functions | |
| double | bstarCSS (double const &b, double const &Q) |
| b* prescription a la Collins-Soper-Sterman More... | |
| double | bstarmin (double const &b, double const &Q) |
| b* prescription with bmin More... | |
| YAML::Emitter & | operator<< (YAML::Emitter &os, ChiSquare const &chi2) |
| Method which prints ChiSquare feautures with cout <<. More... | |
| FourDGrid | Inter4DGrid (std::string const &pf) |
| Function that returns FourDGrid object better tuned for F_UUT, according to the input string. More... | |
| void | ProduceStructGrid (std::string const &GridsDirectory, std::string const &GridTMDPDFfolder, std::string const &GridTMDFFfolder, std::string const &Output, std::string const &repID="none", std::string const &structype="FUUT") |
| This function encapsulates and streamlines the production of an interpolation grid for a structure function starting from TMD grids. More... | |
| std::unique_ptr< YAML::Emitter > | EmitStructGrid (std::string const &GridsDirectory, std::string const &GridTMDPDFfolder, std::string const &GridTMDFFfolder, int const &repnumber, std::string const &pf, FourDGrid const &fdg, int const &qToQcut=5) |
| Function that produces the structure function interpolation grid in momentum space. This is supposed to resamble an LHAPDF grid. We use plain YAML format. More... | |
| std::unique_ptr< YAML::Emitter > | EmitStructGridDirect (std::string const &FitDirectory, int const &repnumber, std::string const &pf, FourDGrid const &fdg, int const &qToQcut) |
| Function that produces the structure function interpolation grid in momentum space. This is supposed to resamble an LHAPDF grid. We use plain YAML format. More... | |
| std::unique_ptr< YAML::Emitter > | EmitStructInfo (std::string const &GridsDirectory, std::string const &GridTMDPDFfolder, std::string const &GridTMDFFfolder, YAML::Node const &config, int const &NumMembers, std::string const &pf, FourDGrid const &fdg, std::vector< int > const &Flavors={-5, -4, -3, -2, -1, 1, 2, 3, 4, 5}, std::string const &SetDesc="Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)", std::string const &Target="proton", std::string const &Hadron="Pip", std::string const &Authors="A. Bacchetta, F. Delcarro, C. Pisano, M. Radici, A. Signori", std::string const &Reference="arXiv:1703.10157", std::string const &SetIndex="000000", std::string const &Format="TMDlib2", std::string const &DataVersion="1", std::string const &ErrorType="Monte Carlo") |
| Function that produces the info file of the TMD set. This is suppose to resamble an LHAPDF info file for the TMDs. We use plain YAML format. More... | |
| ThreeDGrid | Inter3DGrid (std::string const &pf) |
| Function that returns ThreeDGrid object better tuned for PDFs or FFs, according to the input string. More... | |
| void | ProduceTMDGrid (std::string const &ReportFolder, std::string const &Output, std::string const &distype="pdf") |
| This function encapsulates and streamlines the production of an interpolation grid starting from the report produced by a NangaParbat fit. More... | |
| std::unique_ptr< YAML::Emitter > | EmitTMDGrid (YAML::Node const &config, std::string const ¶meterisation, std::vector< double > const ¶ms, std::string const &pf, ThreeDGrid const &tdg) |
| Function that produces the TMD interpolation grid in momentum space. This is supposed to resamble an LHAPDF grid for the TMDs. We use plain YAML format. More... | |
| std::unique_ptr< YAML::Emitter > | EmitTMDInfo (YAML::Node const &config, int const &NumMembers, std::string const &pf, ThreeDGrid const &tdg, std::vector< int > const &Flavors={-5, -4, -3, -2, -1, 1, 2, 3, 4, 5}, std::string const &SetDesc="Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)", std::string const &Authors="M.A.P. Collaboration", std::string const &Reference="arXiv:xxxx.xxxx", std::string const &SetIndex="000000", std::string const &SetName="PV17nll", std::string const &Format="TMDlib2", std::string const &DataVersion="1", std::string const &ErrorType="Monte Carlo", std::string const &FlavorScheme="LHAPDF style") |
| Function that produces the info file of the TMD set. This is suppose to resamble an LHAPDF info file for the TMDs. We use plain YAML format. More... | |
| std::ostream & | operator<< (std::ostream &os, DataHandler const &DH) |
| bool | dir_exists (std::string const &dir) |
| Function that checks whether the input folder exists in the current folder. More... | |
| TMDGrid * | mkTMD (std::string const &name, int const &mem=0) |
| Factory that returns a "TMDGrid" object. More... | |
| TMDGrid * | mkTMD (std::string const &name, std::string const &folder, int const &mem=0) |
| Factory that returns a "TMDGrid" object. More... | |
| std::vector< TMDGrid * > | mkTMDs (std::string const &name) |
| Factory that returns a vector of "TMDGrid" objects. More... | |
| StructGrid * | mkSF (std::string const &name, int const &mem=0) |
| Factory that returns a "StructGrid" object. More... | |
| StructGrid * | mkSF (std::string const &name, std::string const &folder, int const &mem=0) |
| Factory that returns a "StructGrid" object. More... | |
| std::vector< StructGrid * > | mkSFs (std::string const &name) |
| Factory that returns a vector of "StructGrid" objects. More... | |
| std::function< double(double const &, double const &, double const &, double const &)> | Convolution (TMDGrid const *TMD1, TMDGrid const *TMD2, std::function< std::vector< double >(double const &)> const &Charges, double const &kTCutOff=1, double const &IntEps=1e-5) |
| Function that performs the convolution of two TMD distributions in kT space. More... | |
| std::function< double(double const &, double const &, double const &, double const &)> | Convolution (TMDGrid const *TMD, std::function< std::vector< double >(double const &)> const &Charges, double const &kTCutOff=1, double const &IntEps=1e-5) |
| Function that performs the convolution of two TMD distributions in kT space assuming that first and second distributions are equal. More... | |
| std::vector< double > | GenerateGrid (int const &n, double const &min, double const &max, int const &ext=0, bool const &lgt=false) |
| Utility function to generate an interpolation grid. More... | |
| apfel::matrix< double > | CholeskyDecomposition (apfel::matrix< double > const V) |
| Cholesky decomposition of the covariance matrix. More... | |
| std::vector< double > | SolveLowerSystem (apfel::matrix< double > L, std::vector< double > y) |
| Solve lower-diagonal system of equations by forward substitution. More... | |
| std::vector< double > | SolveUpperSystem (apfel::matrix< double > U, std::vector< double > y) |
| Solve upper-diagonal system of equations by backward substitution. More... | |
| std::vector< double > | SolveSymmetricSystem (apfel::matrix< double > A, std::vector< double > rho) |
| Solve symmetric system of equations. More... | |
| std::vector< std::string > | list_dir (std::string const &path) |
| Function that lists elements in a directory. More... | |
| bool | MinuitMinimiser (ChiSquare const &chi2, YAML::Node const ¶meters, gsl_rng *rng=NULL) |
| The "MinuitMinimiser" function using Minuit2 as implemented in ROOT. More... | |
| bool | CeresMinimiser (ChiSquare const &chi2, YAML::Node const ¶meters, gsl_rng *rng=NULL) |
| The "CeresMinimiser" function using ceres-solver. More... | |
| bool | NoMinimiser (ChiSquare const &chi2, YAML::Node const ¶meters) |
| The "NoMinimiser" function simply returns predictions. More... | |
| bool | MinuitScan (ChiSquare const &chi2, YAML::Node const ¶meters, std::string const &outfolder) |
| The "MinuitScan" function performs a scan around the parameters using Minuit2. More... | |
| Parameterisation * | GetParametersation (std::string const &name) |
| Utility function that returns a pointer to a NangaParbat::Parameterisation object pointing to a specific parameterisation. More... | |
| std::string | num_to_string (int const &i, int const &len=4) |
| Function that, given an integer 'i', outputs a string version of 'i' of length 'len' preceeded by as many zeros as needed to reach lenght 'len'. More... | |
| std::string | PreprocessE288 (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the E288 datasets. More... | |
| std::string | PreprocessE605 (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the E605 datasets. More... | |
| std::string | PreprocessPHENIX200 (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the STAR dataset at 510 GeV. More... | |
| std::string | PreprocessSTAR510 (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the STAR dataset at 510 GeV. More... | |
| std::string | PreprocessCDFRunI (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the CDF RunI dataset. More... | |
| std::string | PreprocessCDFRunII (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the CDF RunII dataset. More... | |
| std::string | PreprocessD0RunI (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the D0 RunI dataset. More... | |
| std::string | PreprocessD0RunIImu (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the D0 RunII dataset for muons. More... | |
| std::string | PreprocessD0RunII (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the D0 RunII dataset. More... | |
| std::string | PreprocessLHCb7TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the LHCb dataset at 7 TeV. More... | |
| std::string | PreprocessLHCb8TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the LHCb dataset at 8 TeV. More... | |
| std::string | PreprocessLHCb13TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the LHCb dataset at 13 TeV. More... | |
| std::string | PreprocessCMS7TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the CMS dataset at 7 TeV. More... | |
| std::string | PreprocessCMS8TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the CMS dataset at 8 TeV. More... | |
| std::string | PreprocessATLAS7TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the ATLAS dataset at 7 TeV. More... | |
| std::string | PreprocessATLAS8TeV (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the ATLAS dataset at 8 TeV. More... | |
| std::string | PreprocessEICPseudodata (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the EIC pseudodata. More... | |
| std::string | PreprocessHERMES (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the HERMES multiplicities. More... | |
| std::string | PreprocessCOMPASS (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the COMPASS multiplicities. More... | |
| std::string | PreprocessE537 (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the E537 data. More... | |
| std::string | PreprocessE537_xF (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the E537 data. More... | |
| std::string | PreprocessE615 (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the E615 data. More... | |
| std::string | PreprocessE615_xF (std::string const &RawDataPath, std::string const &ProcessedDataPath, bool const &PDFError=true) |
| Preprocessing of the E615 data. More... | |
| std::string | to_string_with_precision (const double a_value, const int n=3) |
| Function that, given a number 'a_value', outputs a string version of 'a_value' rounding the double to show only 'n' significant figures. More... | |
Variables | |
| const std::map< std::string, std::function< double(double const &, double const &)> > | bstarMap |
| Map of available b* star functions. More... | |
| const std::map< int, std::string > | PtOrderMap {{0, "LL"}, {1, "NLL"}, {-1, "NLL'"}, {2, "NNLL"}, {-2, "NNLL'"}, {3, "NNNLL"}} |
| Map of perturbative orders. More... | |
| const std::map< std::string, Parameterisation * > | AvPars |
| Map of currently available parameterisations. Each of them must correspond to a header file containing a class deriving from the NangaParbat::Parameterisation mother class. More... | |
| double NangaParbat::bstarCSS | ( | double const & | b, |
| double const & | Q | ||
| ) |
b* prescription a la Collins-Soper-Sterman
| b | the impact parameter |
| Q | the hard scale |
| double NangaParbat::bstarmin | ( | double const & | b, |
| double const & | Q | ||
| ) |
b* prescription with bmin
| b | the impact parameter |
| Q | the hard scale |
| bool NangaParbat::CeresMinimiser | ( | ChiSquare const & | chi2, |
| YAML::Node const & | parameters, | ||
| gsl_rng * | rng = NULL |
||
| ) |
The "CeresMinimiser" function using ceres-solver.
| chi2 | the "ChiSquare" object that returns the values of all chi2's |
| parameters | the "YAML::Node" object that contains the parameters to be minimised along with the relevant information |
| rng | GSL random number object |
| apfel::matrix< double > NangaParbat::CholeskyDecomposition | ( | apfel::matrix< double > const | V | ) |
Cholesky decomposition of the covariance matrix.
| V | the covariance matrix |
| std::function< double(double const &, double const &, double const &, double const &)> NangaParbat::Convolution | ( | TMDGrid const * | TMD, |
| std::function< std::vector< double >(double const &)> const & | Charges, | ||
| double const & | kTCutOff = 1, |
||
| double const & | IntEps = 1e-5 |
||
| ) |
Function that performs the convolution of two TMD distributions in kT space assuming that first and second distributions are equal.
| TMD | distribution |
| Charges | to be used as weights of the partonic combinations |
| kTCutOff | cutoff on the integration in kT relative to Q (default: 1) |
| IntEps | integration relative accuracy (default: 1e-5) |
| std::function< double(double const &, double const &, double const &, double const &)> NangaParbat::Convolution | ( | TMDGrid const * | TMD1, |
| TMDGrid const * | TMD2, | ||
| std::function< std::vector< double >(double const &)> const & | Charges, | ||
| double const & | kTCutOff = 1, |
||
| double const & | IntEps = 1e-5 |
||
| ) |
Function that performs the convolution of two TMD distributions in kT space.
| TMD1 | first distribution |
| TMD2 | second distribution |
| Charges | to be used as weights of the partonic combinations |
| kTCutOff | cutoff on the integration in kT relative to Q (default: 1) |
| IntEps | integration relative accuracy (default: 1e-5) |
| bool NangaParbat::dir_exists | ( | std::string const & | dir | ) |
Function that checks whether the input folder exists in the current folder.
| dir | directory name |
| std::unique_ptr< YAML::Emitter > NangaParbat::EmitStructGrid | ( | std::string const & | GridsDirectory, |
| std::string const & | GridTMDPDFfolder, | ||
| std::string const & | GridTMDFFfolder, | ||
| int const & | repnumber, | ||
| std::string const & | pf, | ||
| FourDGrid const & | fdg, | ||
| int const & | qToQcut = 5 |
||
| ) |
Function that produces the structure function interpolation grid in momentum space. This is supposed to resamble an LHAPDF grid. We use plain YAML format.
| GridsDirectory | path to main folder |
| GridTMDPDFfolder | name of TMDPDF grids (subfolder of main folder) |
| GridTMDFFfolder | name of TMDFF grids (subfolder of main folder) |
| repnumber | replica number |
| pf | whether F_UUT or others (not implemented yet) |
| fdg | 4D grid used |
| qToQcut | cut for the convolution integral (default: 5) |
| std::unique_ptr< YAML::Emitter > NangaParbat::EmitStructGridDirect | ( | std::string const & | FitDirectory, |
| int const & | repnumber, | ||
| std::string const & | pf, | ||
| FourDGrid const & | fdg, | ||
| int const & | qToQcut | ||
| ) |
Function that produces the structure function interpolation grid in momentum space. This is supposed to resamble an LHAPDF grid. We use plain YAML format.
| FitDirectory | path to main folder, output of NangaParbat fit |
| repnumber | replica number |
| pf | whether F_UUT or others (not implemented yet) |
| fdg | 4D grid used |
| qToQcut | cut for the convolution integral |
| std::unique_ptr< YAML::Emitter > NangaParbat::EmitStructInfo | ( | std::string const & | GridsDirectory, |
| std::string const & | GridTMDPDFfolder, | ||
| std::string const & | GridTMDFFfolder, | ||
| YAML::Node const & | config, | ||
| int const & | NumMembers, | ||
| std::string const & | pf, | ||
| FourDGrid const & | fdg, | ||
| std::vector< int > const & | Flavors = {-5, -4, -3, -2, -1, 1, 2, 3, 4, 5}, |
||
| std::string const & | SetDesc = "Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)", |
||
| std::string const & | Target = "proton", |
||
| std::string const & | Hadron = "Pip", |
||
| std::string const & | Authors = "A. Bacchetta, F. Delcarro, C. Pisano, M. Radici, A. Signori", |
||
| std::string const & | Reference = "arXiv:1703.10157", |
||
| std::string const & | SetIndex = "000000", |
||
| std::string const & | Format = "TMDlib2", |
||
| std::string const & | DataVersion = "1", |
||
| std::string const & | ErrorType = "Monte Carlo" |
||
| ) |
Function that produces the info file of the TMD set. This is suppose to resamble an LHAPDF info file for the TMDs. We use plain YAML format.
| GridsDirectory | path to main folder |
| GridTMDPDFfolder | name of TMDPDF grids (subfolder of main folder) |
| GridTMDFFfolder | name of TMDFF grids (subfolder of main folder) |
| config | the YAML node with the theory settings |
| NumMembers | number of members |
| pf | whether F_UUT or others (not implemented yet) |
| fdg | 4D grid used |
| Flavors | vector of flavours |
| SetDesc | grid description |
| Target | target type |
| Hadron | hadron type |
| Authors | list of authors |
| Reference | reference |
| SetIndex | ID index of the set |
| Format | TMDlib format |
| DataVersion | version of the grid |
| ErrorType | error-type descriptor |
| std::unique_ptr< YAML::Emitter > NangaParbat::EmitTMDGrid | ( | YAML::Node const & | config, |
| std::string const & | parameterisation, | ||
| std::vector< double > const & | params, | ||
| std::string const & | pf, | ||
| ThreeDGrid const & | tdg | ||
| ) |
Function that produces the TMD interpolation grid in momentum space. This is supposed to resamble an LHAPDF grid for the TMDs. We use plain YAML format.
| config | the YAML node with the theory settings |
| parameterisation | the parameterisation type |
| params | the vector of parameters to be used for the tabulation |
| pf | whether PDFs ("pdf") of FFs ("ff") |
| tdg | 3D grid used |
| std::unique_ptr< YAML::Emitter > NangaParbat::EmitTMDInfo | ( | YAML::Node const & | config, |
| int const & | NumMembers, | ||
| std::string const & | pf, | ||
| ThreeDGrid const & | tdg, | ||
| std::vector< int > const & | Flavors = {-5, -4, -3, -2, -1, 1, 2, 3, 4, 5}, |
||
| std::string const & | SetDesc = "Set produced with NangaParbat + APFEL++ (please cite arXiv:1912.07550 and arXiv:1708.00911)", |
||
| std::string const & | Authors = "M.A.P. Collaboration", |
||
| std::string const & | Reference = "arXiv:xxxx.xxxx", |
||
| std::string const & | SetIndex = "000000", |
||
| std::string const & | SetName = "PV17nll", |
||
| std::string const & | Format = "TMDlib2", |
||
| std::string const & | DataVersion = "1", |
||
| std::string const & | ErrorType = "Monte Carlo", |
||
| std::string const & | FlavorScheme = "LHAPDF style" |
||
| ) |
Function that produces the info file of the TMD set. This is suppose to resamble an LHAPDF info file for the TMDs. We use plain YAML format.
| config | the YAML not with the theory settings |
| NumMembers | number of members |
| pf | whether PDFs ("pdf") of FFs ("ff") |
| tdg | 3D grid used |
| Flavors | vector of flavours |
| SetDesc | grid description |
| Authors | list of authors |
| Reference | reference |
| SetIndex | ID index of the set |
| SetName | name of the set |
| Format | TMDlib format |
| DataVersion | version of the grid |
| ErrorType | error-type descriptor |
| FlavorScheme | flavour scheme |
| std::vector< double > NangaParbat::GenerateGrid | ( | int const & | n, |
| double const & | min, | ||
| double const & | max, | ||
| int const & | ext = 0, |
||
| bool const & | lgt = false |
||
| ) |
Utility function to generate an interpolation grid.
| n | the number of nodes of the grid |
| min | the lower bound |
| max | the upper bound |
| ext | the number of extra nodes (default: 0) |
| lgt | whether the grid is logarithmically spaced (default: false) |
| Parameterisation * NangaParbat::GetParametersation | ( | std::string const & | name | ) |
Utility function that returns a pointer to a NangaParbat::Parameterisation object pointing to a specific parameterisation.
| name | name of the parameterisation |
| ThreeDGrid NangaParbat::Inter3DGrid | ( | std::string const & | pf | ) |
Function that returns ThreeDGrid object better tuned for PDFs or FFs, according to the input string.
| pf | whether PDFs ("pdf") of FFs ("ff") |
| FourDGrid NangaParbat::Inter4DGrid | ( | std::string const & | pf | ) |
Function that returns FourDGrid object better tuned for F_UUT, according to the input string.
| pf | whether F_UUT or others (not implemented yet) |
| std::vector< std::string > NangaParbat::list_dir | ( | std::string const & | path | ) |
Function that lists elements in a directory.
| path | path to the directory |
| bool NangaParbat::MinuitMinimiser | ( | ChiSquare const & | chi2, |
| YAML::Node const & | parameters, | ||
| gsl_rng * | rng = NULL |
||
| ) |
The "MinuitMinimiser" function using Minuit2 as implemented in ROOT.
| chi2 | the "ChiSquare" object that returns the values of all chi2's |
| parameters | the "YAML::Node" object that contains the parameters to be minimised along with the relevant information |
| rng | GSL random number object |
| bool NangaParbat::MinuitScan | ( | ChiSquare const & | chi2, |
| YAML::Node const & | parameters, | ||
| std::string const & | outfolder | ||
| ) |
The "MinuitScan" function performs a scan around the parameters using Minuit2.
| chi2 | the "ChiSquare" object that returns the values of all chi2's |
| parameters | the "YAML::Node" object that contains the parameters |
| outfolder | folder for the output |
| StructGrid * NangaParbat::mkSF | ( | std::string const & | name, |
| int const & | mem = 0 |
||
| ) |
Factory that returns a "StructGrid" object.
| name | name of the structure function set (assumed to be in the current folder for now) |
| mem | member to be read (defaul: 0, i.e. central member) |
| StructGrid * NangaParbat::mkSF | ( | std::string const & | name, |
| std::string const & | folder, | ||
| int const & | mem = 0 |
||
| ) |
Factory that returns a "StructGrid" object.
| name | name of the structure function set |
| folder | name of the folder where the structure function set is |
| mem | member to be read (defaul: 0, i.e. central member) |
| std::vector< StructGrid * > NangaParbat::mkSFs | ( | std::string const & | name | ) |
Factory that returns a vector of "StructGrid" objects.
| name | name of the structure function set (assumed to be in the current folder for now) |
| TMDGrid * NangaParbat::mkTMD | ( | std::string const & | name, |
| int const & | mem = 0 |
||
| ) |
Factory that returns a "TMDGrid" object.
| name | name of the TMD set (assumed to be in the current folder for now) |
| mem | member to be read (defaul: 0, i.e. central member) |
| TMDGrid * NangaParbat::mkTMD | ( | std::string const & | name, |
| std::string const & | folder, | ||
| int const & | mem = 0 |
||
| ) |
Factory that returns a "TMDGrid" object.
| name | name of the TMD set |
| folder | name of the folder where the TMD set is |
| mem | member to be read (defaul: 0, i.e. central member) |
| std::vector< TMDGrid * > NangaParbat::mkTMDs | ( | std::string const & | name | ) |
Factory that returns a vector of "TMDGrid" objects.
| name | name of the TMD set (assumed to be in the current folder for now) |
| bool NangaParbat::NoMinimiser | ( | ChiSquare const & | chi2, |
| YAML::Node const & | parameters | ||
| ) |
The "NoMinimiser" function simply returns predictions.
| chi2 | the "ChiSquare" object that returns the values of all chi2's |
| parameters | the "YAML::Node" object that contains the parameters to be minimised along with the relevant information |
| std::string NangaParbat::num_to_string | ( | int const & | i, |
| int const & | len = 4 |
||
| ) |
Function that, given an integer 'i', outputs a string version of 'i' of length 'len' preceeded by as many zeros as needed to reach lenght 'len'.
| i | number to transform in string |
| len | length of the string (default: 4) |
| std::ostream & NangaParbat::operator<< | ( | std::ostream & | os, |
| DataHandler const & | DH | ||
| ) |
| YAML::Emitter & NangaParbat::operator<< | ( | YAML::Emitter & | os, |
| ChiSquare const & | chi2 | ||
| ) |
Method which prints ChiSquare feautures with cout <<.
| std::string NangaParbat::PreprocessATLAS7TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the ATLAS dataset at 7 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessATLAS8TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the ATLAS dataset at 8 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessCDFRunI | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the CDF RunI dataset.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessCDFRunII | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the CDF RunII dataset.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessCMS7TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the CMS dataset at 7 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessCMS8TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the CMS dataset at 8 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessCOMPASS | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the COMPASS multiplicities.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessD0RunI | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the D0 RunI dataset.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessD0RunII | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the D0 RunII dataset.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessD0RunIImu | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the D0 RunII dataset for muons.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessE288 | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the E288 datasets.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessE537 | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the E537 data.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessE537_xF | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the E537 data.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessE605 | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the E605 datasets.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessE615 | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the E615 data.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessE615_xF | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the E615 data.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessEICPseudodata | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the EIC pseudodata.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessHERMES | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the HERMES multiplicities.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessLHCb13TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the LHCb dataset at 13 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessLHCb7TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the LHCb dataset at 7 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessLHCb8TeV | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the LHCb dataset at 8 TeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessPHENIX200 | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the STAR dataset at 510 GeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| std::string NangaParbat::PreprocessSTAR510 | ( | std::string const & | RawDataPath, |
| std::string const & | ProcessedDataPath, | ||
| bool const & | PDFError = true |
||
| ) |
Preprocessing of the STAR dataset at 510 GeV.
| RawDataPath | the path to the raw-data folder |
| ProcessedDataPath | the path to the folder where the processed data will be stored |
| PDFError | whether PDF uncertainties have to be included in the processed datafiles (default: true) |
| void NangaParbat::ProduceStructGrid | ( | std::string const & | GridsDirectory, |
| std::string const & | GridTMDPDFfolder, | ||
| std::string const & | GridTMDFFfolder, | ||
| std::string const & | Output, | ||
| std::string const & | repID = "none", |
||
| std::string const & | structype = "FUUT" |
||
| ) |
This function encapsulates and streamlines the production of an interpolation grid for a structure function starting from TMD grids.
| GridsDirectory | path to main folder |
| GridTMDPDFfolder | name of TMDPDF grids (subfolder of main folder) |
| GridTMDFFfolder | name of TMDFF grids (subfolder of main folder) |
| Output | name of the output grid |
| repID | number of the replica |
| structype | whether F_UUT or others (not implemented yet) |
| void NangaParbat::ProduceTMDGrid | ( | std::string const & | ReportFolder, |
| std::string const & | Output, | ||
| std::string const & | distype = "pdf" |
||
| ) |
This function encapsulates and streamlines the production of an interpolation grid starting from the report produced by a NangaParbat fit.
| ReportFolder | path to the report folder |
| Output | name of the output grid |
| distype | whether PDFs ("pdf") of FFs ("ff") |
| std::vector< double > NangaParbat::SolveLowerSystem | ( | apfel::matrix< double > | L, |
| std::vector< double > | y | ||
| ) |
Solve lower-diagonal system of equations by forward substitution.
| L | lower-diagonal matrix |
| y | vector of constants |
| std::vector< double > NangaParbat::SolveSymmetricSystem | ( | apfel::matrix< double > | A, |
| std::vector< double > | rho | ||
| ) |
Solve symmetric system of equations.
| A | symmetric matrix |
| rho | vector of constants |
| std::vector< double > NangaParbat::SolveUpperSystem | ( | apfel::matrix< double > | U, |
| std::vector< double > | y | ||
| ) |
Solve upper-diagonal system of equations by backward substitution.
| U | upper-diagonal matrix |
| y | vector of constants |
| std::string NangaParbat::to_string_with_precision | ( | const double | a_value, |
| const int | n = 3 |
||
| ) |
Function that, given a number 'a_value', outputs a string version of 'a_value' rounding the double to show only 'n' significant figures.
| a_value | number to transform in string |
| n | length of the string (default: 3) |
| const std::map<std::string, Parameterisation*> NangaParbat::AvPars |
Map of currently available parameterisations. Each of them must correspond to a header file containing a class deriving from the NangaParbat::Parameterisation mother class.
| const std::map<std::string, std::function<double(double const&, double const&)> > NangaParbat::bstarMap |
Map of available b* star functions.
| const std::map<int, std::string> NangaParbat::PtOrderMap {{0, "LL"}, {1, "NLL"}, {-1, "NLL'"}, {2, "NNLL"}, {-2, "NNLL'"}, {3, "NNNLL"}} |
Map of perturbative orders.