Changeset 6965fe1 in git
- Timestamp:
- Dec 5, 2017, 8:33:29 PM (7 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 8d4d51b
- Parents:
- ba75867 (diff), dda357d (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the(diff)
links above to see all the changes relative to each parent. - Files:
-
- 9 added
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
Makefile
rba75867 r6965fe1 1257 1257 external/fastjet/Voronoi.$(SrcSuf) \ 1258 1258 external/fastjet/internal/Voronoi.hh 1259 tmp/external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.$(ObjSuf): \ 1260 external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.$(ObjSuf) 1259 1261 tmp/external/fastjet/contribs/Nsubjettiness/AxesDefinition.$(ObjSuf): \ 1260 1262 external/fastjet/contribs/Nsubjettiness/AxesDefinition.$(SrcSuf) … … 1417 1419 external/fastjet/plugins/CDFCones/fastjet/CDFMidPointPlugin.hh \ 1418 1420 external/fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh \ 1421 external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.hh \ 1419 1422 external/fastjet/contribs/Nsubjettiness/Nsubjettiness.hh \ 1420 1423 external/fastjet/contribs/Nsubjettiness/Njettiness.hh \ … … 1444 1447 external/fastjet/plugins/CDFCones/fastjet/CDFMidPointPlugin.hh \ 1445 1448 external/fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh \ 1449 external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.hh \ 1446 1450 external/fastjet/contribs/Nsubjettiness/Nsubjettiness.hh \ 1447 1451 external/fastjet/contribs/Nsubjettiness/Njettiness.hh \ … … 1501 1505 tmp/external/fastjet/TilingExtent.$(ObjSuf) \ 1502 1506 tmp/external/fastjet/Voronoi.$(ObjSuf) \ 1507 tmp/external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.$(ObjSuf) \ 1503 1508 tmp/external/fastjet/contribs/Nsubjettiness/AxesDefinition.$(ObjSuf) \ 1504 1509 tmp/external/fastjet/contribs/Nsubjettiness/ExtraRecombiners.$(ObjSuf) \ … … 1894 1899 @touch $@ 1895 1900 1901 1902 external/fastjet/contribs/ValenciaPlugin/ValenciaPlugin.hh: \ 1903 external/fastjet/JetDefinition.hh \ 1904 external/fastjet/ClusterSequence.hh 1905 @touch $@ 1906 1896 1907 external/fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh: \ 1897 1908 external/fastjet/ClusterSequence.hh \ -
cards/delphes_card_CMS.tcl
rba75867 r6965fe1 169 169 170 170 # resolution formula for electrons 171 # based on arXiv:1 405.6569171 # based on arXiv:1502.02701 172 172 set ResolutionFormula { (abs(eta) <= 0.5) * (pt > 0.1) * sqrt(0.03^2 + pt^2*1.3e-3^2) + 173 173 (abs(eta) > 0.5 && abs(eta) <= 1.5) * (pt > 0.1) * sqrt(0.05^2 + pt^2*1.7e-3^2) + -
modules/FastJetFinder.cc
rba75867 r6965fe1 66 66 #include "fastjet/contribs/Nsubjettiness/ExtraRecombiners.hh" 67 67 68 #include "fastjet/contribs/ValenciaPlugin/ValenciaPlugin.hh" 69 68 70 #include "fastjet/tools/Filter.hh" 69 71 #include "fastjet/tools/Pruner.hh" … … 79 81 FastJetFinder::FastJetFinder() : 80 82 fPlugin(0), fRecomb(0), fAxesDef(0), fMeasureDef(0), fNjettinessPlugin(0), 81 fDefinition(0), fAreaDefinition(0), fItInputArray(0) 83 fDefinition(0), fAreaDefinition(0), fItInputArray(0), fValenciaPlugin(0) 82 84 { 83 85 … … 118 120 fJetPTMin = GetDouble("JetPTMin", 10.0); 119 121 122 120 123 //-- N(sub)jettiness parameters -- 121 124 … … 125 128 fRcutOff = GetDouble("RcutOff", 0.8); // used only if Njettiness is used as jet clustering algo (case 8) 126 129 fN = GetInt("N", 2); // used only if Njettiness is used as jet clustering algo (case 8) 127 130 131 //-- Exclusive clustering for e+e- collisions -- 132 133 fNJets = GetInt("NJets",2); 134 fExclusiveClustering = GetBool("ExclusiveClustering", false); 135 136 //-- Valencia Linear Collider algorithm 137 fGamma = GetDouble("Gamma", 1.0); 138 //fBeta parameter see above 139 128 140 fMeasureDef = new NormalizedMeasure(fBeta, fParameterR); 129 141 … … 236 248 fDefinition = new JetDefinition(fNjettinessPlugin); 237 249 break; 250 case 9: 251 fValenciaPlugin = new ValenciaPlugin(fParameterR, fBeta, fGamma); 252 fDefinition = new JetDefinition(fValenciaPlugin); 253 break; 254 238 255 } 239 256 … … 294 311 if(fAxesDef) delete fAxesDef; 295 312 if(fMeasureDef) delete fMeasureDef; 313 if(fValenciaPlugin) delete static_cast<JetDefinition::Plugin*>(fValenciaPlugin); 314 296 315 } 297 316 … … 357 376 358 377 outputList.clear(); 359 outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin)); 360 378 379 if(fExclusiveClustering) 380 { 381 outputList = sorted_by_pt(sequence->exclusive_jets( fNJets )); 382 } 383 else 384 { 385 outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin)); 386 } 361 387 362 388 // loop over all jets and export them -
modules/FastJetFinder.h
rba75867 r6965fe1 41 41 namespace contrib { 42 42 class NjettinessPlugin; 43 class ValenciaPlugin; 43 44 class AxesDefinition; 44 45 class MeasureDefinition; … … 61 62 void *fPlugin; //! 62 63 void *fRecomb; //! 64 65 fastjet::contrib::AxesDefinition *fAxesDef; 66 fastjet::contrib::MeasureDefinition *fMeasureDef; 67 63 68 fastjet::contrib::NjettinessPlugin *fNjettinessPlugin; //! 64 69 fastjet::contrib::ValenciaPlugin *fValenciaPlugin; //! 65 70 fastjet::JetDefinition *fDefinition; //! 66 71 67 72 Int_t fJetAlgorithm; 68 73 Double_t fParameterR; 74 69 75 Double_t fJetPTMin; 70 76 Double_t fConeRadius; … … 77 83 Double_t fOverlapThreshold; 78 84 85 //-- Exclusive clustering for e+e- collisions -- 86 87 Int_t fNJets; 88 Bool_t fExclusiveClustering; 89 90 //-- Valencia Linear Collider algorithm 91 Double_t fGamma; 92 79 93 //-- N (sub)jettiness parameters -- 80 94 81 95 Bool_t fComputeNsubjettiness; 82 fastjet::contrib::AxesDefinition *fAxesDef;83 fastjet::contrib::MeasureDefinition *fMeasureDef;84 96 Double_t fBeta; 85 97 Int_t fAxisMode;
Note:
See TracChangeset
for help on using the changeset viewer.