|
|
|||
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/dglapnonlinear.h" 0012 0013 namespace apfel 0014 { 0015 /** 0016 * @brief Structure that contains all the precomputed quantities 0017 * needed to perform the DGLAP evolution, i.e. perturbative 0018 * coefficients of splitting functions and matching conditions, and 0019 * the quark thresholds. 0020 */ 0021 struct DglapObjects 0022 { 0023 double Threshold; 0024 Set<Operator> UnitySet; 0025 std::map<int, Set<Operator>> SplittingFunctions; 0026 std::map<int, Set<Operator>> MatchingConditions; 0027 }; 0028 0029 /** 0030 * @name DGLAP object initializers in QCD 0031 * Collection of functions that initialise DglapObjects structure 0032 * for the different kinds of QCD evolution currently available. 0033 */ 0034 ///@{ 0035 /** 0036 * @brief The InitializeDglapObjectsQCD function precomputes the 0037 * perturbative coefficients of space-like unpolarised splitting 0038 * functions and matching conditions and store them into a 0039 * 'DglapObjects' structure. 0040 * @param g: the x-space grid 0041 * @param Masses: the quark masses 0042 * @param Thresholds: the quark thresholds 0043 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0044 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0045 * @param n3lo: whether N3LO corrections to splitting and matching functions are computer (default: false) 0046 * @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) 0047 * @return A map of DglapObject objects, one for each possible nf 0048 */ 0049 std::map<int, DglapObjects> InitializeDglapObjectsQCD(Grid const& g, 0050 std::vector<double> const& Masses, 0051 std::vector<double> const& Thresholds, 0052 bool const& OpEvol = false, 0053 double const& IntEps = 1e-5, 0054 bool const& n3lo = false, 0055 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0056 0057 /** 0058 * @brief The InitializeDglapObjectsQCD function precomputes the 0059 * perturbative coefficients of space-like unpolarised splitting 0060 * functions and matching conditions and store them into a 0061 * 'DglapObjects' structure. 0062 * @param g: the x-space grid 0063 * @param Thresholds: the quark thresholds 0064 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0065 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0066 * @param n3lo: whether N3LO corrections to splitting and matching functions are computer (default: false) 0067 * @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) 0068 * @return A map of DglapObject objects, one for each possible nf 0069 * @note This function assumes that masses and thresholds coincide. 0070 */ 0071 std::map<int, DglapObjects> InitializeDglapObjectsQCD(Grid const& g, 0072 std::vector<double> const& Thresholds, 0073 bool const& OpEvol = false, 0074 double const& IntEps = 1e-5, 0075 bool const& n3lo = false, 0076 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0077 0078 /** 0079 * @brief The InitializeDglapObjectsQCDome function precomputes the 0080 * perturbative coefficients of space-like unpolarised splitting 0081 * functions and matching conditions and store them into a 0082 * 'DglapObjects' structure. This variant uses libome for the 0083 * matching condictions at threshold. 0084 * @param g: the x-space grid 0085 * @param Masses: the quark masses 0086 * @param Thresholds: the quark thresholds 0087 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0088 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0089 * @param IMod: the vector of switches to vary the parameterisation of the approximated N3LO splitting functions (only relevant at N3LO = true) (default: all zero's) 0090 * @return A map of DglapObject objects, one for each possible nf 0091 */ 0092 std::map<int, DglapObjects> InitializeDglapObjectsQCDOme(Grid const& g, 0093 std::vector<double> const& Masses, 0094 std::vector<double> const& Thresholds, 0095 bool const& OpEvol = false, 0096 double const& IntEps = 1e-5, 0097 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0098 0099 /** 0100 * @brief The InitializeDglapObjectsQCDome function precomputes the 0101 * perturbative coefficients of space-like unpolarised splitting 0102 * functions and matching conditions and store them into a 0103 * 'DglapObjects' structure. This variant uses libome for the 0104 * matching condictions at threshold. 0105 * @param g: the x-space grid 0106 * @param Masses: the quark masses 0107 * @param Thresholds: the quark thresholds 0108 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0109 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0110 * @param IMod: the vector of switches to vary the parameterisation of the approximated N3LO splitting functions (only relevant at N3LO = true) (default: all zero's) 0111 * @return A map of DglapObject objects, one for each possible nf 0112 * @note This function assumes that masses and thresholds coincide. 0113 */ 0114 std::map<int, DglapObjects> InitializeDglapObjectsQCDOme(Grid const& g, 0115 std::vector<double> const& Thresholds, 0116 bool const& OpEvol = false, 0117 double const& IntEps = 1e-5, 0118 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0119 0120 /** 0121 * @brief The InitializeDglapObjectsQCD function precomputes the 0122 * perturbative coefficients of space-like unpolarised splitting 0123 * functions and matching conditions and store them into a 0124 * 'DglapObjects' structure. This function prepares the evolution in 0125 * the physical basis. 0126 * @param g: the x-space grid 0127 * @param Masses: the quark masses 0128 * @param Thresholds: the quark thresholds 0129 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0130 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0131 * @param n3lo: whether N3LO corrections to splitting and matching functions are computer (default: false) 0132 * @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) 0133 * @return A map of DglapObject objects, one for each possible nf 0134 */ 0135 std::map<int, DglapObjects> InitializeDglapObjectsQCDPhys(Grid const& g, 0136 std::vector<double> const& Masses, 0137 std::vector<double> const& Thresholds, 0138 bool const& OpEvol = false, 0139 double const& IntEps = 1e-5, 0140 bool const& n3lo = false, 0141 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0142 0143 /** 0144 * @brief The InitializeDglapObjectsQCD function precomputes the 0145 * perturbative coefficients of space-like unpolarised splitting 0146 * functions and matching conditions and store them into a 0147 * 'DglapObjects' structure. This function prepares the evolution in 0148 * the physical basis. 0149 * @param g: the x-space grid 0150 * @param Thresholds: the quark thresholds 0151 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0152 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0153 * @param n3lo: whether N3LO corrections to splitting and matching functions are computer (default: false) 0154 * @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) 0155 * @return A map of DglapObject objects, one for each possible nf 0156 * @note This function assumes that masses and thresholds coincide. 0157 */ 0158 std::map<int, DglapObjects> InitializeDglapObjectsQCDPhys(Grid const& g, 0159 std::vector<double> const& Thresholds, 0160 bool const& OpEvol = false, 0161 double const& IntEps = 1e-5, 0162 bool const& n3lo = false, 0163 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0164 0165 /** 0166 * @brief The InitializeDglapObjectsQCD function precomputes the 0167 * perturbative coefficients of space-like unpolarised splitting 0168 * functions and matching conditions and store them into a 0169 * 'DglapObjects' structure. This function prepares the evolution in 0170 * the physical basis. This variant uses libome for the matching 0171 * condictions at threshold. 0172 * @param g: the x-space grid 0173 * @param Masses: the quark masses 0174 * @param Thresholds: the quark thresholds 0175 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0176 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0177 * @param IMod: the vector of switches to vary the parameterisation of the approximated N3LO splitting functions (only relevant at N3LO) (default: all zero's) 0178 * @return A map of DglapObject objects, one for each possible nf 0179 */ 0180 std::map<int, DglapObjects> InitializeDglapObjectsQCDPhysOme(Grid const& g, 0181 std::vector<double> const& Masses, 0182 std::vector<double> const& Thresholds, 0183 bool const& OpEvol = false, 0184 double const& IntEps = 1e-5, 0185 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0186 0187 /** 0188 * @brief The InitializeDglapObjectsQCD function precomputes the 0189 * perturbative coefficients of space-like unpolarised splitting 0190 * functions and matching conditions and store them into a 0191 * 'DglapObjects' structure. This function prepares the evolution in 0192 * the physical basis. This variant uses libome for the matching 0193 * condictions at threshold. 0194 * @param g: the x-space grid 0195 * @param Thresholds: the quark thresholds 0196 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0197 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0198 * @param IMod: the vector of switches to vary the parameterisation of the approximated N3LO splitting functions (only relevant at N3LO) (default: all zero's) 0199 * @return A map of DglapObject objects, one for each possible nf 0200 * @note This function assumes that masses and thresholds coincide. 0201 */ 0202 std::map<int, DglapObjects> InitializeDglapObjectsQCDPhysOme(Grid const& g, 0203 std::vector<double> const& Thresholds, 0204 bool const& OpEvol = false, 0205 double const& IntEps = 1e-5, 0206 std::vector<int> const& IMod = {0, 0, 0, 0, 0, 0, 0}); 0207 0208 /** 0209 * @brief The InitializeDglapObjectsQCD function precomputes the 0210 * perturbative coefficients of space-like unpolarised splitting 0211 * functions and matching conditions with MSbar masses and store 0212 * them into a 'DglapObjects' structure. 0213 * @param g: the x-space grid 0214 * @param Masses: the quark masses 0215 * @param Thresholds: the quark thresholds 0216 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0217 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0218 * @return A map of DglapObject objects, one for each possible nf 0219 */ 0220 std::map<int, DglapObjects> InitializeDglapObjectsQCDMSbarMass(Grid const& g, 0221 std::vector<double> const& Masses, 0222 std::vector<double> const& Thresholds, 0223 bool const& OpEvol = false, 0224 double const& IntEps = 1e-5); 0225 0226 /** 0227 * @brief The InitializeDglapObjectsQCD function precomputes the 0228 * perturbative coefficients of space-like unpolarised splitting 0229 * functions and matching conditions with MSbar masses and store 0230 * them into a 'DglapObjects' structure. 0231 * @param g: the x-space grid 0232 * @param Thresholds: the quark masses 0233 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0234 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0235 * @return A map of DglapObject objects, one for each possible nf 0236 * @note This function assumes that masses and thresholds coincide. 0237 */ 0238 std::map<int, DglapObjects> InitializeDglapObjectsQCDMSbarMass(Grid const& g, 0239 std::vector<double> const& Thresholds, 0240 bool const& OpEvol = false, 0241 double const& IntEps = 1e-5); 0242 0243 /** 0244 * @brief The InitializeDglapObjectsQCDpol function precomputes the 0245 * perturbative coefficients of space-like longitudinally polarised 0246 * splitting functions and matching conditions (assumed to be equal 0247 * to the unpolarised ones) and store them into a 'DglapObjects' 0248 * structure. 0249 * @param g: the x-space grid 0250 * @param Masses: the quark masses 0251 * @param Thresholds: the quark thresholds 0252 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0253 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0254 * @return A map of DglapObject objects, one for each possible nf 0255 */ 0256 std::map<int, DglapObjects> InitializeDglapObjectsQCDpol(Grid const& g, 0257 std::vector<double> const& Masses, 0258 std::vector<double> const& Thresholds, 0259 bool const& OpEvol = false, 0260 double const& IntEps = 1e-5); 0261 0262 /** 0263 * @brief The InitializeDglapObjectsQCDpol function precomputes the 0264 * perturbative coefficients of space-like longitudinally polarised 0265 * splitting functions and matching conditions (assumed to be equal 0266 * to the unpolarised ones) and store them into a 'DglapObjects' 0267 * structure. 0268 * @param g: the x-space grid 0269 * @param Thresholds: the quark thresholds 0270 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0271 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0272 * @return A map of DglapObject objects, one for each possible nf 0273 * @note This function assumes that masses and thresholds coincide. 0274 */ 0275 std::map<int, DglapObjects> InitializeDglapObjectsQCDpol(Grid const& g, 0276 std::vector<double> const& Thresholds, 0277 bool const& OpEvol = false, 0278 double const& IntEps = 1e-5); 0279 0280 /** 0281 * @brief The InitializeDglapObjectsQCDpolPhys function precomputes 0282 * the perturbative coefficients of space-like longitudinally 0283 * polarised splitting functions and matching conditions (assumed to 0284 * be equal to the unpolarised ones) and store them into a 0285 * 'DglapObjects' structure. This function prepares the evolution in 0286 * the physical basis. 0287 * @param g: the x-space grid 0288 * @param Masses: the quark masses 0289 * @param Thresholds: the quark thresholds 0290 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0291 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0292 * @return A map of DglapObject objects, one for each possible nf 0293 */ 0294 std::map<int, DglapObjects> InitializeDglapObjectsQCDpolPhys(Grid const& g, 0295 std::vector<double> const& Masses, 0296 std::vector<double> const& Thresholds, 0297 bool const& OpEvol = false, 0298 double const& IntEps = 1e-5); 0299 0300 /** 0301 * @brief The InitializeDglapObjectsQCDpolPhys function precomputes 0302 * the perturbative coefficients of space-like longitudinally 0303 * polarised splitting functions and matching conditions (assumed to 0304 * be equal to the unpolarised ones) and store them into a 0305 * 'DglapObjects' structure. This function prepares the evolution in 0306 * the physical basis. 0307 * @param g: the x-space grid 0308 * @param Thresholds: the quark thresholds 0309 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0310 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0311 * @return A map of DglapObject objects, one for each possible nf 0312 * @note This function assumes that masses and thresholds coincide. 0313 */ 0314 std::map<int, DglapObjects> InitializeDglapObjectsQCDpolPhys(Grid const& g, 0315 std::vector<double> const& Thresholds, 0316 bool const& OpEvol = false, 0317 double const& IntEps = 1e-5); 0318 0319 /** 0320 * @brief The InitializeDglapObjectsQCDT function precomputes the 0321 * perturbative coefficients of time-like unpolarised splitting 0322 * functions and matching conditions and store them into a 0323 * 'DglapObjects' structure. 0324 * @param g: the x-space grid 0325 * @param Masses: the quark masses 0326 * @param Thresholds: the quark thresholds 0327 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0328 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0329 * @return A map of DglapObject objects, one for each possible nf 0330 */ 0331 std::map<int, DglapObjects> InitializeDglapObjectsQCDT(Grid const& g, 0332 std::vector<double> const& Masses, 0333 std::vector<double> const& Thresholds, 0334 bool const& OpEvol = false, 0335 double const& IntEps = 1e-5); 0336 0337 /** 0338 * @brief The InitializeDglapObjectsQCDT function precomputes the 0339 * perturbative coefficients of time-like unpolarised splitting 0340 * functions and matching conditions and store them into a 0341 * 'DglapObjects' structure. 0342 * @param g: the x-space grid 0343 * @param Thresholds: the quark thresholds 0344 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0345 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0346 * @return A map of DglapObject objects, one for each possible nf 0347 * @note This function assumes that masses and thresholds coincide. 0348 */ 0349 std::map<int, DglapObjects> InitializeDglapObjectsQCDT(Grid const& g, 0350 std::vector<double> const& Thresholds, 0351 bool const& OpEvol = false, 0352 double const& IntEps = 1e-5); 0353 0354 /** 0355 * @brief The InitializeDglapObjectsQCDTPhys function precomputes 0356 * the perturbative coefficients of time-like unpolarised splitting 0357 * functions and matching conditions and store them into a 0358 * 'DglapObjects' structure.This function prepares the evolution in 0359 * the physical basis. 0360 * @param g: the x-space grid 0361 * @param Masses: the quark masses 0362 * @param Thresholds: the quark thresholds 0363 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0364 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0365 * @return A map of DglapObject objects, one for each possible nf 0366 */ 0367 std::map<int, DglapObjects> InitializeDglapObjectsQCDTPhys(Grid const& g, 0368 std::vector<double> const& Masses, 0369 std::vector<double> const& Thresholds, 0370 bool const& OpEvol = false, 0371 double const& IntEps = 1e-5); 0372 0373 /** 0374 * @brief The InitializeDglapObjectsQCDTPhys function precomputes 0375 * the perturbative coefficients of time-like unpolarised splitting 0376 * functions and matching conditions and store them into a 0377 * 'DglapObjects' structure.This function prepares the evolution in 0378 * the physical basis. 0379 * @param g: the x-space grid 0380 * @param Thresholds: the quark thresholds 0381 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0382 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0383 * @return A map of DglapObject objects, one for each possible nf 0384 * @note This function assumes that masses and thresholds coincide. 0385 */ 0386 std::map<int, DglapObjects> InitializeDglapObjectsQCDTPhys(Grid const& g, 0387 std::vector<double> const& Thresholds, 0388 bool const& OpEvol = false, 0389 double const& IntEps = 1e-5); 0390 0391 /** 0392 * @brief The InitializeDglapObjectsQCDTpol function precomputes the 0393 * perturbative coefficients of time-like longitudinally polarised 0394 * splitting functions and matching conditions and store them into a 0395 * 'DglapObjects' structure. Limited to LO for now. 0396 * @param g: the x-space grid 0397 * @param Masses: the quark masses 0398 * @param Thresholds: the quark thresholds 0399 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0400 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0401 * @return A map of DglapObject objects, one for each possible nf 0402 */ 0403 std::map<int, DglapObjects> InitializeDglapObjectsQCDTpol(Grid const& g, 0404 std::vector<double> const& Masses, 0405 std::vector<double> const& Thresholds, 0406 bool const& OpEvol = false, 0407 double const& IntEps = 1e-5); 0408 0409 /** 0410 * @brief The InitializeDglapObjectsQCDTpol function precomputes the 0411 * perturbative coefficients of time-like longitudinally polarised 0412 * splitting functions and matching conditions and store them into a 0413 * 'DglapObjects' structure. Limited to LO for now. 0414 * @param g: the x-space grid 0415 * @param Thresholds: the quark thresholds 0416 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0417 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0418 * @return A map of DglapObject objects, one for each possible nf 0419 * @note This function assumes that masses and thresholds coincide. 0420 */ 0421 std::map<int, DglapObjects> InitializeDglapObjectsQCDTpol(Grid const& g, 0422 std::vector<double> const& Thresholds, 0423 bool const& OpEvol = false, 0424 double const& IntEps = 1e-5); 0425 0426 /** 0427 * @brief The InitializeDglapObjectsQCDtrans function precomputes 0428 * the perturbative coefficients of space-like transverity splitting 0429 * functions and matching conditions and store them into a 0430 * 'DglapObjects' structure. 0431 * @param g: the x-space grid 0432 * @param Masses: the quark masses 0433 * @param Thresholds: the quark thresholds 0434 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0435 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0436 * @return A map of DglapObject objects, one for each possible nf 0437 */ 0438 std::map<int, DglapObjects> InitializeDglapObjectsQCDtrans(Grid const& g, 0439 std::vector<double> const& Masses, 0440 std::vector<double> const& Thresholds, 0441 bool const& OpEvol = false, 0442 double const& IntEps = 1e-5); 0443 0444 /** 0445 * @brief The InitializeDglapObjectsQCDtrans function precomputes 0446 * the perturbative coefficients of space-like transverity splitting 0447 * functions and matching conditions and store them into a 0448 * 'DglapObjects' structure. 0449 * @param g: the x-space grid 0450 * @param Thresholds: the quark thresholds 0451 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0452 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0453 * @return A map of DglapObject objects, one for each possible nf 0454 * @note This function assumes that masses and thresholds coincide. 0455 */ 0456 std::map<int, DglapObjects> InitializeDglapObjectsQCDtrans(Grid const& g, 0457 std::vector<double> const& Thresholds, 0458 bool const& OpEvol = false, 0459 double const& IntEps = 1e-5); 0460 0461 /** 0462 * @brief The InitializeDglapObjectsQCDtrans function precomputes 0463 * the perturbative coefficients of timelike-like transverity splitting 0464 * functions and matching conditions and store them into a 0465 * 'DglapObjects' structure. 0466 * @param g: the x-space grid 0467 * @param Masses: the quark masses 0468 * @param Thresholds: the quark thresholds 0469 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0470 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0471 * @return A map of DglapObject objects, one for each possible nf 0472 */ 0473 std::map<int, DglapObjects> InitializeDglapObjectsQCDTtrans(Grid const& g, 0474 std::vector<double> const& Masses, 0475 std::vector<double> const& Thresholds, 0476 bool const& OpEvol = false, 0477 double const& IntEps = 1e-5); 0478 0479 /** 0480 * @brief The InitializeDglapObjectsQCDtrans function precomputes 0481 * the perturbative coefficients of time-like transverity splitting 0482 * functions and matching conditions and store them into a 0483 * 'DglapObjects' structure. 0484 * @param g: the x-space grid 0485 * @param Thresholds: the quark thresholds 0486 * @param OpEvol: the switch for the computation of the evolution operator (default: false) 0487 * @param IntEps: the integration accuracy (default: 10<SUP>-5</SUP>) 0488 * @return A map of DglapObject objects, one for each possible nf 0489 * @note This function assumes that masses and thresholds coincide. 0490 */ 0491 std::map<int, DglapObjects> InitializeDglapObjectsQCDTtrans(Grid const& g, 0492 std::vector<double> const& Thresholds, 0493 bool const& OpEvol = false, 0494 double const& IntEps = 1e-5); 0495 ///@} 0496 0497 /** 0498 * @name DGLAP builders in QCD 0499 * Collection of functions that build a Dglap object used to perform 0500 * the DGLAP evolution of distributions or operators in QCD. 0501 */ 0502 ///@{ 0503 /** 0504 * @brief The BuildDglap function builds the actual dglap object 0505 * that performs the DGLAP evolution for distributions. 0506 * @param DglapObj: structure with the coefficients of the perturbative objects 0507 * @param InDistFunc: the distributions at the reference scale 0508 * @param MuRef: the reference scale 0509 * @param PerturbativeOrder: the perturbative order of the evolution 0510 * @param Alphas: the function returning the strong coupling 0511 * @param xi: the scale-variation parameter (default: 1) 0512 * @param nsteps: the number of steps of the ODE solver (default: 10). 0513 * @return A unique pointer to a Dglap object 0514 */ 0515 std::unique_ptr<Dglap<Distribution>> BuildDglap(std::map<int, DglapObjects> const& DglapObj, 0516 std::function<std::map<int, double>(double const&, double const&)> const& InDistFunc, 0517 double const& MuRef, 0518 int const& PerturbativeOrder, 0519 std::function<double(double const&)> const& Alphas, 0520 double const& xi = 1, 0521 int const& nsteps = 10); 0522 0523 /** 0524 * @brief The BuildDglap function builds the actual dglap object 0525 * that performs the DGLAP evolution for operators. 0526 * @param DglapObj: structure with the coefficients of the perturbative objects 0527 * @param MuRef: the reference scale 0528 * @param PerturbativeOrder: the perturbative order of the evolution 0529 * @param Alphas: the function returning the strong coupling 0530 * @param xi: the scale-variation parameter (default: 1) 0531 * @param nsteps: the number of steps of the ODE solver (default: 10). 0532 * @return A unique pointer to a Dglap object 0533 */ 0534 std::unique_ptr<Dglap<Operator>> BuildDglap(std::map<int, DglapObjects> const& DglapObj, 0535 double const& MuRef, 0536 int const& PerturbativeOrder, 0537 std::function<double(double const&)> const& Alphas, 0538 double const& xi = 1, 0539 int const& nsteps = 10); 0540 0541 /** 0542 * @brief The BuildDglap function builds the actual dglap object 0543 * that performs the DGLAP evolution for distributions. 0544 * @param DglapObj: DglapObjects-valued function that returns the 0545 * structure with the coefficients of the perturbative objects as 0546 * functions of the scale. In this function, the DglapObjects are 0547 * still assumed to have a perturbative expansion where each 0548 * coefficient is mu dependent. Therefore, one can still decide the 0549 * perturbative order at this stage. 0550 * @param Thresholds: the quark thresholds 0551 * @param InDistFunc: the distributions at the reference scale 0552 * @param MuRef: the reference scale 0553 * @param PerturbativeOrder: the perturbative order of the evolution 0554 * @param Alphas: the function returning the strong coupling 0555 * @param nsteps: the number of steps of the ODE solver (default: 10). 0556 * @return A unique pointer to a Dglap object 0557 */ 0558 std::unique_ptr<Dglap<Distribution>> BuildDglap(std::function<DglapObjects(double const&)> const& DglapObj, 0559 std::vector<double> const& Thresholds, 0560 std::function<std::map<int, double>(double const&, double const&)> const& InDistFunc, 0561 double const& MuRef, 0562 int const& PerturbativeOrder, 0563 std::function<double(double const&)> const& Alphas, 0564 int const& nsteps = 10); 0565 0566 /** 0567 * @brief The BuildDglap function builds the actual dglap object 0568 * that performs the DGLAP evolution for distributions. 0569 * @param DglapObj: DglapObjects-valued function that returns the 0570 * structure with the perturbative objects as functions of the 0571 * scale. In this function, the DglapObjects do not longer admit a 0572 * truncated perturbative expasion (as for example in 0573 * resummation). The relevant objects are to be stored in the key 0574 * zero of the map of DglapObjects. 0575 * @param Thresholds: the quark thresholds 0576 * @param InDistFunc: the distributions at the reference scale 0577 * @param MuRef: the reference scale 0578 * @param nsteps: the number of steps of the ODE solver (default: 10). 0579 * @return A unique pointer to a Dglap object 0580 */ 0581 std::unique_ptr<Dglap<Distribution>> BuildDglap(std::function<DglapObjects(double const&)> const& DglapObj, 0582 std::vector<double> const& Thresholds, 0583 std::function<std::map<int, double>(double const&, double const&)> const& InDistFunc, 0584 double const& MuRef, 0585 int const& nsteps = 10); 0586 0587 /** 0588 * @brief The BuildDglap function builds the actual dglap object 0589 * that performs a non-liner DGLAP evolution. 0590 * @param DglapObj: structure with the coefficients of the perturbative objects 0591 * @param TranformationFuncs: set of functions that trasform distributions in the l.h.s. of DGLAP 0592 * @param InDistFunc: the distributions at the reference scale 0593 * @param MuRef: the reference scale 0594 * @param PerturbativeOrder: the perturbative order of the evolution 0595 * @param Alphas: the function returning the strong coupling 0596 * @param xi: the scale-variation parameter (default: 1) 0597 * @param nsteps: the number of steps of the ODE solver (default: 10). 0598 * @return A unique pointer to a Dglap object 0599 */ 0600 std::unique_ptr<DglapNonLinear> BuildDglapNonLinear(std::map<int, DglapObjects> const& DglapObj, 0601 std::function<std::map<int, std::function<double(std::map<int, double> const&)>>(double const&)> const& TranformationFuncs, 0602 std::function<std::map<int, double>(double const&, double const&)> const& InDistFunc, 0603 double const& MuRef, 0604 int const& PerturbativeOrder, 0605 std::function<double(double const&)> const& Alphas, 0606 double const& xi = 1, 0607 int const& nsteps = 10); 0608 ///@} 0609 }
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
|