Fork me on GitHub

Changeset e9c0d73 in git for modules/FastJetFinder.cc


Ignore:
Timestamp:
Dec 5, 2017, 1:11:19 PM (7 years ago)
Author:
Ulrike Schnoor <schnooru@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
0879ed1
Parents:
197fed7
Message:

added exclusive_ymerge for n=2,3,4,5

File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/FastJetFinder.cc

    r197fed7 re9c0d73  
    8080
    8181FastJetFinder::FastJetFinder() :
    82   fPlugin(0), fRecomb(0), fAxesDef(0), fMeasureDef(0), fNjettinessPlugin(0),
    83   fDefinition(0), fAreaDefinition(0), fItInputArray(0), fValenciaPlugin(0)
     82  fPlugin(0), fRecomb(0), fAxesDef(0), fMeasureDef(0), fNjettinessPlugin(0), fValenciaPlugin(0),
     83  fDefinition(0), fAreaDefinition(0), fItInputArray(0)
    8484{
    8585
     
    332332  vector< PseudoJet >::iterator itInputList, itOutputList;
    333333  vector< TEstimatorStruct >::iterator itEstimators;
    334 
     334  Double_t excl_ymerge23 = 0.0;
     335  Double_t excl_ymerge34 = 0.0;
     336  Double_t excl_ymerge45 = 0.0;
     337  Double_t excl_ymerge56 = 0.0;
     338 
    335339  DelphesFactory *factory = GetFactory();
    336340
     
    377381  outputList.clear();
    378382
     383 
    379384  if(fExclusiveClustering)
    380385    {
    381       //not neccessary to sort
    382       outputList = sequence->exclusive_jets( fNJets );
     386      outputList = sorted_by_pt(sequence->exclusive_jets( fNJets ));
     387
     388      excl_ymerge23 = sequence->exclusive_ymerge( 2 );
     389      excl_ymerge34 = sequence->exclusive_ymerge( 3 );
     390      excl_ymerge45 = sequence->exclusive_ymerge( 4 );
     391      excl_ymerge56 = sequence->exclusive_ymerge( 5 );
    383392    }
    384393  else
     
    401410    if(fAreaDefinition) area = itOutputList->area_4vector();
    402411
     412
     413   
    403414    candidate = factory->NewCandidate();
    404415
     
    444455    candidate->NNeutrals = nneutrals;
    445456    candidate->NCharged = ncharged;
     457
     458
     459    //for exclusive clustering, access y_n,n+1 as exclusive_ymerge (fNJets);
     460    candidate->ExclYmerge23 = excl_ymerge23;
     461    candidate->ExclYmerge34 = excl_ymerge34;
     462    candidate->ExclYmerge45 = excl_ymerge45;
     463    candidate->ExclYmerge56 = excl_ymerge56;
    446464   
    447465    //------------------------------------
Note: See TracChangeset for help on using the changeset viewer.