Changeset e4c3fef in git for modules/FastJetFinder.cc
- Timestamp:
- Apr 17, 2014, 10:54:17 AM (11 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- b96d99b
- Parents:
- 2e6a81b
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/FastJetFinder.cc
r2e6a81b re4c3fef 47 47 #include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh" 48 48 49 #include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh" 49 #include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh" 50 50 #include "fastjet/contribs/Nsubjettiness/Njettiness.hh" 51 51 #include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh" … … 76 76 void FastJetFinder::Init() 77 77 { 78 78 79 79 JetDefinition::Plugin *plugin = NULL; 80 80 JetDefinition::Recombiner *recomb = NULL; 81 81 NjettinessPlugin *njet_plugin = NULL; 82 82 83 83 // read eta ranges 84 84 … … 110 110 111 111 //-- N(sub)jettiness parameters -- 112 112 113 113 fComputeNsubjettiness = GetBool("ComputeNsubjettiness", false); 114 114 fBeta = GetDouble("Beta", 1.0); 115 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 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 119 119 // --- Jet Area Parameters --- 120 120 fAreaAlgorithm = GetInt("AreaAlgorithm", 0); 121 121 fComputeRho = GetBool("ComputeRho", false); 122 122 123 123 // - ghost based areas - 124 124 fGhostEtaMax = GetDouble("GhostEtaMax", 5.0); … … 128 128 fPtScatter = GetDouble("PtScatter", 0.1); 129 129 fMeanGhostPt = GetDouble("MeanGhostPt", 1.0E-100); 130 130 131 131 // - voronoi based areas - 132 132 fEffectiveRfact = GetDouble("EffectiveRfact", 1.0); … … 188 188 break; 189 189 } 190 191 190 192 191 fPlugin = plugin; 193 192 fRecomb = recomb; 194 193 fNjettinessPlugin = njet_plugin; 195 194 196 195 ClusterSequence::print_banner(); 197 196 … … 299 298 inputList.clear(); 300 299 inputList = sequence->constituents(*itOutputList); 301 300 302 301 for(itInputList = inputList.begin(); itInputList != inputList.end(); ++itInputList) 303 302 { … … 308 307 if(deta > detaMax) detaMax = deta; 309 308 if(dphi > dphiMax) dphiMax = dphi; 310 309 311 310 time += TMath::Sqrt(constituent->Momentum.E())*(constituent->Position.T()); 312 311 weightTime += TMath::Sqrt(constituent->Momentum.E()); 313 312 314 313 candidate->AddCandidate(constituent); 315 314 } 316 315 317 316 avTime = time/weightTime; 318 317 … … 325 324 326 325 // --- compute N-subjettiness with N = 1,2,3,4,5 ---- 327 326 328 327 if(fComputeNsubjettiness) 329 328 { 330 329 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 330 331 switch(fAxisMode) 332 { 333 default: 334 case 1: 335 axisMode = Njettiness::wta_kt_axes; 336 break; 337 case 2: 338 axisMode = Njettiness::onepass_wta_kt_axes; 339 break; 340 case 3: 341 axisMode = Njettiness::kt_axes; 342 break; 343 case 4: 344 axisMode = Njettiness::onepass_kt_axes; 345 break; 346 } 347 337 348 Njettiness::MeasureMode measureMode = Njettiness::unnormalized_measure; 338 349 339 350 Nsubjettiness nSub1(1, axisMode, measureMode, fBeta); 340 351 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);352 Nsubjettiness nSub3(3, axisMode, measureMode, fBeta); 353 Nsubjettiness nSub4(4, axisMode, measureMode, fBeta); 354 Nsubjettiness nSub5(5, axisMode, measureMode, fBeta); 355 356 candidate->Tau[0] = nSub1(*itOutputList); 357 candidate->Tau[1] = nSub2(*itOutputList); 358 candidate->Tau[2] = nSub3(*itOutputList); 359 candidate->Tau[3] = nSub4(*itOutputList); 360 candidate->Tau[4] = nSub5(*itOutputList); 350 361 } 351 362
Note:
See TracChangeset
for help on using the changeset viewer.