- Timestamp:
- Apr 16, 2014, 4:32:53 PM (11 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 82575a3
- Parents:
- a0431dc
- Location:
- modules
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/FastJetFinder.cc
ra0431dc r9687203 47 47 #include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh" 48 48 49 #include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh" 50 #include "fastjet/contribs/Nsubjettiness/Njettiness.hh" 51 #include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh" 52 #include "fastjet/contribs/Nsubjettiness/WinnerTakeAllRecombiner.hh" 53 49 54 using namespace std; 50 55 using namespace fastjet; 56 using namespace fastjet::contrib; 57 51 58 52 59 //------------------------------------------------------------------------------ 53 60 54 61 FastJetFinder::FastJetFinder() : 55 fPlugin(0), fDefinition(0), fAreaDefinition(0), fItInputArray(0) 62 fPlugin(0), fDefinition(0), fAreaDefinition(0), fItInputArray(0), fRecomb(0), fNjettinessPlugin(0) 56 63 { 57 64 … … 69 76 void FastJetFinder::Init() 70 77 { 78 71 79 JetDefinition::Plugin *plugin = NULL; 72 80 JetDefinition::Recombiner *recomb = NULL; 81 NjettinessPlugin *njet_plugin = NULL; 82 73 83 // read eta ranges 74 84 … … 99 109 fJetPTMin = GetDouble("JetPTMin", 10.0); 100 110 111 //-- N(sub)jettiness parameters -- 112 113 fComputeNsubjettiness = GetBool("ComputeNsubjettiness", false); 114 fBeta = GetDouble("Beta", 1.0); 115 fAxisMode = GetInt("AxisMode", 1); 116 fRcutOff = GetDouble("RcutOff", 0.8); //used only if Njettiness is used as jet clustering algo (case 8) 117 fN = GetInt("N", 2); //used only if Njettiness is used as jet clustering algo (case 8) 118 101 119 // --- Jet Area Parameters --- 102 120 fAreaAlgorithm = GetInt("AreaAlgorithm", 0); 103 121 fComputeRho = GetBool("ComputeRho", false); 122 104 123 // - ghost based areas - 105 124 fGhostEtaMax = GetDouble("GhostEtaMax", 5.0); … … 109 128 fPtScatter = GetDouble("PtScatter", 0.1); 110 129 fMeanGhostPt = GetDouble("MeanGhostPt", 1.0E-100); 130 111 131 // - voronoi based areas - 112 132 fEffectiveRfact = GetDouble("EffectiveRfact", 1.0); … … 159 179 fDefinition = new fastjet::JetDefinition(fastjet::antikt_algorithm, fParameterR); 160 180 break; 161 } 181 case 7: 182 recomb = new fastjet::contrib::WinnerTakeAllRecombiner(); 183 fDefinition = new fastjet::JetDefinition(fastjet::antikt_algorithm, fParameterR, recomb, Best); 184 break; 185 case 8: 186 njet_plugin = new fastjet::contrib::NjettinessPlugin(fN, Njettiness::wta_kt_axes, Njettiness::unnormalized_cutoff_measure, fBeta, fRcutOff); 187 fDefinition = new fastjet::JetDefinition(njet_plugin); 188 break; 189 } 190 162 191 163 192 fPlugin = plugin; 164 193 fRecomb = recomb; 194 fNjettinessPlugin = njet_plugin; 195 165 196 ClusterSequence::print_banner(); 166 197 … … 184 215 if(fAreaDefinition) delete fAreaDefinition; 185 216 if(fPlugin) delete static_cast<JetDefinition::Plugin*>(fPlugin); 217 if(fRecomb) delete static_cast<JetDefinition::Recombiner*>(fRecomb); 218 if(fNjettinessPlugin) delete static_cast<JetDefinition::Plugin*>(fNjettinessPlugin); 186 219 } 187 220 … … 248 281 outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin)); 249 282 283 250 284 // loop over all jets and export them 251 285 detaMax = 0.0; … … 265 299 inputList.clear(); 266 300 inputList = sequence->constituents(*itOutputList); 301 267 302 for(itInputList = inputList.begin(); itInputList != inputList.end(); ++itInputList) 268 303 { … … 289 324 candidate->DeltaPhi = dphiMax; 290 325 326 // --- compute N-subjettiness with N = 1,2,3,4,5 ---- 327 328 if(fComputeNsubjettiness) 329 { 330 Njettiness::AxesMode axisMode; 331 332 if (fAxisMode == 1) axisMode = Njettiness::wta_kt_axes; 333 if (fAxisMode == 2) axisMode = Njettiness::onepass_wta_kt_axes; 334 if (fAxisMode == 3) axisMode = Njettiness::kt_axes; 335 if (fAxisMode == 4) axisMode = Njettiness::onepass_kt_axes; 336 337 Njettiness::MeasureMode measureMode = Njettiness::unnormalized_measure; 338 339 Nsubjettiness nSub1(1, axisMode, measureMode, fBeta); 340 Nsubjettiness nSub2(2, axisMode, measureMode, fBeta); 341 Nsubjettiness nSub3(3, axisMode, measureMode, fBeta); 342 Nsubjettiness nSub4(4, axisMode, measureMode, fBeta); 343 Nsubjettiness nSub5(5, axisMode, measureMode, fBeta); 344 345 candidate -> Tau1 = nSub1(*itOutputList); 346 candidate -> Tau2 = nSub2(*itOutputList); 347 candidate -> Tau3 = nSub3(*itOutputList); 348 candidate -> Tau4 = nSub4(*itOutputList); 349 candidate -> Tau5 = nSub5(*itOutputList); 350 } 351 352 291 353 fOutputArray->Add(candidate); 292 354 } -
modules/FastJetFinder.h
ra0431dc r9687203 41 41 42 42 void *fPlugin; //! 43 void *fRecomb; //! 44 void *fNjettinessPlugin; //! 45 43 46 fastjet::JetDefinition *fDefinition; //! 44 47 … … 55 58 Double_t fOverlapThreshold; 56 59 60 //-- N (sub)jettiness parameters -- 61 62 Bool_t fComputeNsubjettiness; 63 Double_t fBeta; 64 Int_t fAxisMode; 65 Double_t fRcutOff; 66 Int_t fN ; 67 57 68 // --- FastJet Area method -------- 58 69 -
modules/TreeWriter.cc
ra0431dc r9687203 554 554 entry->FracPt[3] = candidate->FracPt[3]; 555 555 entry->FracPt[4] = candidate->FracPt[4]; 556 556 557 //--- N-subjettiness variables ---- 558 559 entry->Tau1 = candidate->Tau1; 560 entry->Tau2 = candidate->Tau2; 561 entry->Tau3 = candidate->Tau3; 562 entry->Tau4 = candidate->Tau4; 563 entry->Tau5 = candidate->Tau5; 564 557 565 FillParticles(candidate, &entry->Particles); 558 566 }
Note:
See TracChangeset
for help on using the changeset viewer.