Changeset 8abab33 in git
- Timestamp:
- Oct 16, 2015, 2:03:01 PM (9 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- 5e4c40d
- Parents:
- 0bbf248
- Location:
- modules
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/FastJetFinder.cc
r0bbf248 r8abab33 78 78 79 79 FastJetFinder::FastJetFinder() : 80 fPlugin(0), fRecomb(0), fNjettinessPlugin(0), fDefinition(0), fAreaDefinition(0), fItInputArray(0) 80 fPlugin(0), fRecomb(0), fAxesDef(0), fMeasureDef(0), fNjettinessPlugin(0), 81 fDefinition(0), fAreaDefinition(0), fItInputArray(0) 81 82 { 82 83 … … 124 125 fRcutOff = GetDouble("RcutOff", 0.8); // used only if Njettiness is used as jet clustering algo (case 8) 125 126 fN = GetInt("N", 2); // used only if Njettiness is used as jet clustering algo (case 8) 127 128 fMeasureDef = new UnnormalizedMeasure(fBeta); 129 130 switch(fAxisMode) 131 { 132 default: 133 case 1: 134 fAxesDef = new WTA_KT_Axes(); 135 break; 136 case 2: 137 fAxesDef = new OnePass_WTA_KT_Axes(); 138 break; 139 case 3: 140 fAxesDef = new KT_Axes(); 141 break; 142 case 4: 143 fAxesDef = new OnePass_KT_Axes(); 144 } 126 145 127 146 //-- Trimming parameters -- … … 219 238 } 220 239 240 241 221 242 fPlugin = plugin; 222 243 fRecomb = recomb; … … 271 292 if(fRecomb) delete static_cast<JetDefinition::Recombiner*>(fRecomb); 272 293 if(fNjettinessPlugin) delete static_cast<JetDefinition::Plugin*>(fNjettinessPlugin); 294 if(fAxesDef) delete fAxesDef; 295 if(fMeasureDef) delete fMeasureDef; 273 296 } 274 297 … … 399 422 //------------------------------------ 400 423 401 424 if(fComputeTrimming) 402 425 { 403 426 … … 416 439 candidate->NSubJetsTrimmed = subjets.size(); 417 440 418 for (size_t i = 0; i < subjets.size() and i < 4; i++){ 419 if(subjets.at(i).pt() < 0) continue ; 420 candidate->TrimmedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m()); 441 for (size_t i = 0; i < subjets.size() and i < 4; i++) 442 { 443 if(subjets.at(i).pt() < 0) continue ; 444 candidate->TrimmedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m()); 421 445 } 422 446 } … … 443 467 candidate->NSubJetsPruned = subjets.size(); 444 468 445 for (size_t i = 0; i < subjets.size() and i < 4; i++){ 446 if(subjets.at(i).pt() < 0) continue ; 447 candidate->PrunedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m()); 469 for (size_t i = 0; i < subjets.size() and i < 4; i++) 470 { 471 if(subjets.at(i).pt() < 0) continue ; 472 candidate->PrunedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m()); 448 473 } 449 474 … … 457 482 { 458 483 459 contrib::SoftDrop 484 contrib::SoftDrop softDrop(fBetaSoftDrop,fSymmetryCutSoftDrop,fR0SoftDrop); 460 485 fastjet::PseudoJet softdrop_jet = softDrop(*itOutputList); 461 486 … … 469 494 candidate->NSubJetsSoftDropped = softdrop_jet.pieces().size(); 470 495 471 for (size_t i = 0; i < subjets.size() and i < 4; i++){ 472 if(subjets.at(i).pt() < 0) continue ; 473 candidate->SoftDroppedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m()); 496 for (size_t i = 0; i < subjets.size() and i < 4; i++) 497 { 498 if(subjets.at(i).pt() < 0) continue ; 499 candidate->SoftDroppedP4[i+1].SetPtEtaPhiM(subjets.at(i).pt(), subjets.at(i).eta(), subjets.at(i).phi(), subjets.at(i).m()); 474 500 } 475 501 } … … 479 505 if(fComputeNsubjettiness) 480 506 { 481 Njettiness::AxesMode axisMode; 482 483 switch(fAxisMode) 484 { 485 default: 486 case 1: 487 axisMode = Njettiness::wta_kt_axes; 488 break; 489 case 2: 490 axisMode = Njettiness::onepass_wta_kt_axes; 491 break; 492 case 3: 493 axisMode = Njettiness::kt_axes; 494 break; 495 case 4: 496 axisMode = Njettiness::onepass_kt_axes; 497 break; 498 } 499 500 Njettiness::MeasureMode measureMode = Njettiness::unnormalized_measure; 501 502 Nsubjettiness nSub1(1, axisMode, measureMode, fBeta); 503 Nsubjettiness nSub2(2, axisMode, measureMode, fBeta); 504 Nsubjettiness nSub3(3, axisMode, measureMode, fBeta); 505 Nsubjettiness nSub4(4, axisMode, measureMode, fBeta); 506 Nsubjettiness nSub5(5, axisMode, measureMode, fBeta); 507 507 508 Nsubjettiness nSub1(1, *fAxesDef, *fMeasureDef); 509 Nsubjettiness nSub2(2, *fAxesDef, *fMeasureDef); 510 Nsubjettiness nSub3(3, *fAxesDef, *fMeasureDef); 511 Nsubjettiness nSub4(4, *fAxesDef, *fMeasureDef); 512 Nsubjettiness nSub5(5, *fAxesDef, *fMeasureDef); 513 508 514 candidate->Tau[0] = nSub1(*itOutputList); 509 515 candidate->Tau[1] = nSub2(*itOutputList); … … 511 517 candidate->Tau[3] = nSub4(*itOutputList); 512 518 candidate->Tau[4] = nSub5(*itOutputList); 519 513 520 } 514 521 -
modules/FastJetFinder.h
r0bbf248 r8abab33 41 41 namespace contrib { 42 42 class NjettinessPlugin; 43 class AxesDefinition; 44 class MeasureDefinition; 43 45 } 44 46 } … … 78 80 79 81 Bool_t fComputeNsubjettiness; 82 fastjet::contrib::AxesDefinition *fAxesDef; 83 fastjet::contrib::MeasureDefinition *fMeasureDef; 80 84 Double_t fBeta; 81 85 Int_t fAxisMode;
Note:
See TracChangeset
for help on using the changeset viewer.