|
|
|||
File indexing completed on 2026-06-02 08:17:12
0001 // 0002 // APFEL++ 2017 0003 // 0004 // Author: Valerio Bertone: valerio.bertone@cern.ch 0005 // 0006 0007 #pragma once 0008 0009 #include "apfel/grid.h" 0010 #include "apfel/dglap.h" 0011 #include "apfel/constants.h" 0012 0013 namespace apfel 0014 { 0015 /** 0016 * @brief Structure that contains all the precomputed quantities 0017 * needed to perform the DGLAP evolution in QCDxQED, 0018 * i.e. perturbative coefficients of splitting functions and 0019 * matching conditions, and the heavy-quark/lepton thresholds. 0020 */ 0021 struct DglapObjectsQCDQED 0022 { 0023 double Threshold; 0024 PartonSpecies Species; 0025 std::vector<int> ActiveFlavours; 0026 Set<Operator> UnitySet; 0027 std::map<std::pair<int, int>, Set<Operator>> SplittingFunctions; 0028 std::map<std::pair<int, int>, Set<Operator>> MatchingConditions; 0029 std::map<std::pair<int, int>, Set<Distribution>> InhomogeneousTerms; 0030 }; 0031 0032 /** 0033 * @name DGLAP object initializers in QCDxQED 0034 * Collection of functions that initialise DglapObjectsQED structure 0035 * for the different kinds of QCDxQED evolution currently available. 0036 */ 0037 ///@{ 0038 /** 0039 * @brief The InitializeDglapObjectsQCDQED function precomputes the 0040 * perturbative coefficients of space-like unpolarised splitting 0041 * functions and matching conditions for QCDxQED evolution and store 0042 * them into a 'DglapObjectsQCDQED' structure. 0043 * @param g: the x-space grid 0044 * @param QuarkThresholds: the quark thresholds 0045 * @param LeptonThresholds: the lepton thresholds 0046 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0047 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0048 * @param n3lo: whether N3LO corrections to splitting and matching functions are computer (default: false) 0049 * @param IMod: the vector of switches to vary the parameterisation of the approximated N3LO splitting functions (only relevant for n3lo = true) (default: all zero's) 0050 * @return A map of DglapObjectsQCDQED objects, one for each possible nf 0051 * @note This function assumes that masses and thresholds coincide. 0052 */ 0053 std::map<int, DglapObjectsQCDQED> InitializeDglapObjectsQCDQED(Grid const& g, 0054 std::vector<double> const& QuarkThresholds, 0055 std::vector<double> const& LeptonThresholds, 0056 bool const& OpEvol = false, 0057 double const& IntEps = 1e-5, 0058 bool const& n3lo = false, 0059 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0060 0061 /** 0062 * @brief The InitializeDglapObjectsQCDQED function precomputes the 0063 * perturbative coefficients of space-like unpolarised splitting 0064 * functions and matching conditions for QCDxQED evolution and store 0065 * them into a 'DglapObjectsQCDQED' structure. This variant uses 0066 * libome for the matching condictions at threshold. 0067 * @param g: the x-space grid 0068 * @param QuarkThresholds: the quark thresholds 0069 * @param LeptonThresholds: the lepton thresholds 0070 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0071 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0072 * @param IMod: the vector of switches to vary the parameterisation of the approximated N3LO splitting functions (only relevant at N3LO) (default: all zero's) 0073 * @return A map of DglapObjectsQCDQED objects, one for each possible nf 0074 * @note This function assumes that masses and thresholds coincide. 0075 */ 0076 std::map<int, DglapObjectsQCDQED> InitializeDglapObjectsQCDQEDOme(Grid const& g, 0077 std::vector<double> const& QuarkThresholds, 0078 std::vector<double> const& LeptonThresholds, 0079 bool const& OpEvol = false, 0080 double const& IntEps = 1e-5, 0081 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0082 0083 /** 0084 * @brief The InitializeDglapObjectsPhoton function precomputes the 0085 * perturbative coefficients of space-like unpolarised splitting 0086 * functions, matching conditions, and inhomogeneous terms for QCD 0087 * evolution of the photon and store them into a 0088 * 'DglapObjectsQCDQED' structure. 0089 * @param g: the x-space grid 0090 * @param Thresholds: the quark thresholds 0091 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0092 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0093 * @param DISgamma: the flag to enable the DISgamma scheme (default: false) 0094 * @return A map of DglapObjectsQCDQED objects, one for each possible nf 0095 * @note This function assumes that masses and thresholds 0096 * coincide. Also, N3LO corrections are not included. 0097 */ 0098 std::map<int, DglapObjectsQCDQED> InitializeDglapObjectsPhoton(Grid const& g, 0099 std::vector<double> const& Thresholds, 0100 bool const& OpEvol = false, 0101 double const& IntEps = 1e-5, 0102 bool const& DISgamma = false); 0103 ///@} 0104 0105 /** 0106 * @name DGLAP builders for QCDxQED evolution 0107 * Collection of functions that build a Dglap object used to perform 0108 * the DGLAP QCDxQED evolution of distributions or operators. 0109 */ 0110 ///@{ 0111 /** 0112 * @brief The SplittingFunctionsQCDQED function constructs a 0113 * function of the number of active flavours and the scale variable 0114 * t = 2log(mu) that returns the full set of splitting functions. 0115 * @param DglapObj: structure with the coefficients of the perturbative objects 0116 * @param PerturbativeOrder: the perturbative order of the evolution 0117 * @param Alphas: the function returning the strong coupling 0118 * @param Alphaem: the function returning the eletromagnetic coupling 0119 * @return A standard function returning a set of operators. 0120 */ 0121 std::function<Set<Operator>(int const&, double const&)> SplittingFunctionsQCDQED(std::map<int, DglapObjectsQCDQED> const& DglapObj, 0122 int const& PerturbativeOrder, 0123 std::function<double(double const&)> const& Alphas, 0124 std::function<double(double const&)> const& Alphaem); 0125 0126 /** 0127 * @brief The MatchingConditionsQCDQED function constructs a 0128 * function of whether the matching is to be done upwards or 0129 * downwards and the number of active flavours that returns the full 0130 * set of matching functions. 0131 * @param DglapObj: structure with the coefficients of the perturbative objects 0132 * @param PerturbativeOrder: the perturbative order of the evolution 0133 * @param AlphasTh: the values of the strong coupling below and above threshold 0134 * @param AlphaemTh: the values of the eletromagnetic coupling below and above threshold 0135 * @return A standard function returning a set of operators. 0136 */ 0137 std::function<Set<Operator>(bool const&, int const&)> MatchingConditionsQCDQED(std::map<int, DglapObjectsQCDQED> const& DglapObj, 0138 int const& PerturbativeOrder, 0139 std::map<int, std::pair<double, double>> const& AlphasTh, 0140 std::map<int, std::pair<double, double>> const& AlphaemTh); 0141 0142 /** 0143 * @brief The InhomogeneousTermsQCDQED function constructs a 0144 * function of the number of active flavours and the scale variable 0145 * t = 2log(mu) that returns the full set of inhomogeneous terms. 0146 * @param DglapObj: structure with the coefficients of the perturbative objects 0147 * @param PerturbativeOrder: the perturbative order of the evolution 0148 * @param Alphas: the function returning the strong coupling 0149 * @param Alphaem: the function returning the eletromagnetic coupling 0150 * @return A standard function returning a set of distributions. 0151 */ 0152 std::function<Set<Distribution>(int const&, double const&)> InhomogeneousTermsQCDQED(std::map<int, DglapObjectsQCDQED> const& DglapObj, 0153 int const& PerturbativeOrder, 0154 std::function<double(double const&)> const& Alphas, 0155 std::function<double(double const&)> const& Alphaem); 0156 0157 /** 0158 * @brief The BuildDglap function builds the actual dglap object 0159 * that performs the DGLAP QCDxQED evolution for distributions. 0160 * @param DglapObj: structure with the coefficients of the perturbative objects 0161 * @param InDistFunc: the distributions at the reference scale 0162 * @param MuRef: the reference scale 0163 * @param PerturbativeOrder: the perturbative order of the evolution 0164 * @param Alphas: the function returning the strong coupling 0165 * @param Alphaem: the function returning the eletromagnetic coupling 0166 * @param nsteps: the number of steps of the ODE solver (default: 10). 0167 * @return A unique pointer to a Dglap object 0168 */ 0169 std::unique_ptr<Dglap<Distribution>> BuildDglap(std::map<int, DglapObjectsQCDQED> const& DglapObj, 0170 std::function<std::map<int, double>(double const&, double const&)> const& InDistFunc, 0171 double const& MuRef, 0172 int const& PerturbativeOrder, 0173 std::function<double(double const&)> const& Alphas, 0174 std::function<double(double const&)> const& Alphaem, 0175 int const& nsteps = 10); 0176 0177 /** 0178 * @brief The BuildDglap function builds the actual dglap object 0179 * that performs the DGLAP QCDxQED evolution for operators. 0180 * @param DglapObj: structure with the coefficients of the perturbative objects 0181 * @param MuRef: the reference scale 0182 * @param PerturbativeOrder: the perturbative order of the evolution 0183 * @param Alphas: the function returning the strong coupling 0184 * @param Alphaem: the function returning the eletromagnetic coupling 0185 * @param nsteps: the number of steps of the ODE solver (default: 10). 0186 * @return A unique pointer to a Dglap object 0187 */ 0188 std::unique_ptr<Dglap<Operator>> BuildDglap(std::map<int, DglapObjectsQCDQED> const& DglapObj, 0189 double const& MuRef, 0190 int const& PerturbativeOrder, 0191 std::function<double(double const&)> const& Alphas, 0192 std::function<double(double const&)> const& Alphaem, 0193 int const& nsteps = 10); 0194 ///@} 0195 }
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|