Fork me on GitHub

Ignore:
Timestamp:
Jun 30, 2015, 2:38:29 PM (9 years ago)
Author:
Michele Selvaggi <michele.selvaggi@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
c627b07
Parents:
6153fb0
Message:

fixed BTagging without LHEF

File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/JetFlavorAssociation.cc

    r6153fb0 r1180bc1  
    151151  fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
    152152  fItPartonInputArray = fPartonInputArray->MakeIterator();
     153  fPartonFilter = new ExRootFilter(fPartonInputArray);
    153154
    154155  fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/allParticles"));
    155156  fItParticleInputArray = fParticleInputArray->MakeIterator();
    156157
    157   fParticleLHEFInputArray = ImportArray(GetString("ParticleLHEFInputArray", "Delphes/allParticlesLHEF"));
    158   fItParticleLHEFInputArray = fParticleLHEFInputArray->MakeIterator();
    159 
    160   fPartonFilter = new ExRootFilter(fPartonInputArray);
    161   fParticleLHEFFilter = new ExRootFilter(fParticleLHEFInputArray);
     158  try
     159  {
     160    fParticleLHEFInputArray = ImportArray(GetString("ParticleLHEFInputArray", "Delphes/allParticlesLHEF"));
     161  }
     162  catch(runtime_error &e)
     163  {
     164    fParticleLHEFInputArray = 0;
     165  }
     166
     167  if(fParticleLHEFInputArray)
     168  {
     169    fItParticleLHEFInputArray = fParticleLHEFInputArray->MakeIterator();
     170    fParticleLHEFFilter = new ExRootFilter(fParticleLHEFInputArray);
     171  }
    162172
    163173  fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
     
    183193
    184194  Candidate *jet;
    185   TObjArray *partonArray;
    186   TObjArray *partonLHEFArray;
    187 
    188   // select quark && gluons
     195  TObjArray *partonArray = 0;
     196  TObjArray *partonLHEFArray = 0;
     197
     198  // select quark and gluons
    189199  fPartonFilter->Reset();
    190200  partonArray = fPartonFilter->GetSubArray(fPartonClassifier, 0); // get the filtered parton array
    191 
    192201  if(partonArray == 0) return;
    193   TIter itPartonArray(partonArray);
    194 
    195   fParticleLHEFFilter->Reset();
    196   partonLHEFArray = fParticleLHEFFilter->GetSubArray(fParticleLHEFClassifier, 0); // get the filtered parton array
    197 
    198   if(partonLHEFArray == 0) return;
    199   TIter itPartonLHEFArray(partonLHEFArray);
    200 
     202
     203  if(fParticleLHEFInputArray)
     204  {
     205    fParticleLHEFFilter->Reset();
     206    partonLHEFArray = fParticleLHEFFilter->GetSubArray(fParticleLHEFClassifier, 0); // get the filtered parton array
     207  }
    201208  // loop over all input jets
    202209  fItJetInputArray->Reset();
     
    204211  {
    205212    // get standard flavor
    206     GetAlgoFlavor(jet, itPartonArray, itPartonLHEFArray);
    207     GetPhysicsFlavor(jet, itPartonArray, itPartonLHEFArray);
     213    GetAlgoFlavor(jet, partonArray, partonLHEFArray);
     214    if(fParticleLHEFInputArray) GetPhysicsFlavor(jet, partonArray, partonLHEFArray);
    208215  }
    209216}
     
    213220// https://cmssdt.cern.ch/SDT/lxr/source/PhysicsTools/JetMCAlgos/plugins/JetPartonMatcher.cc?v=CMSSW_7_3_0_pre1
    214221
    215 void JetFlavorAssociation::GetAlgoFlavor(Candidate *jet, TIter &itPartonArray, TIter &itPartonLHEFArray)
     222void JetFlavorAssociation::GetAlgoFlavor(Candidate *jet, TObjArray *partonArray, TObjArray *partonLHEFArray)
    216223{
    217224  float maxPt = 0;
    218   bool isGoodParton = true;
    219225  int daughterCounter = 0;
    220226  Candidate *parton, *partonLHEF;
    221227  Candidate *tempParton = 0, *tempPartonHighestPt = 0;
    222228  int pdgCode, pdgCodeMax = -1;
     229 
     230  TIter itPartonArray(partonArray);
     231  TIter itPartonLHEFArray(partonLHEFArray);
    223232
    224233  itPartonArray.Reset();
     
    232241      if(pdgCodeMax < pdgCode) pdgCodeMax = pdgCode;
    233242    }
    234 
    235     isGoodParton = true;
     243 
     244    if(!fParticleLHEFInputArray) continue;
     245 
    236246    itPartonLHEFArray.Reset();
    237247    while((partonLHEF = static_cast<Candidate *>(itPartonLHEFArray.Next())))
     
    240250         parton->PID == partonLHEF->PID &&
    241251         partonLHEF->Charge == parton->Charge)
    242       {
    243          isGoodParton = false;
     252      {     
    244253         break;
    245254      }
    246 
    247       if(!isGoodParton) continue;
    248255
    249256      // check the daugheter
     
    285292//------------------------------------------------------------------------------
    286293
    287 void JetFlavorAssociation::GetPhysicsFlavor(Candidate *jet, TIter &itPartonArray, TIter &itPartonLHEFArray)
     294void JetFlavorAssociation::GetPhysicsFlavor(Candidate *jet, TObjArray *partonArray, TObjArray *partonLHEFArray)
    288295{
    289296  int partonCounter = 0;
     
    298305  vector<Candidate *>::iterator itContaminations;
    299306
     307  TIter itPartonArray(partonArray);
     308  TIter itPartonLHEFArray(partonLHEFArray);
     309
    300310  contaminations.clear();
    301311
Note: See TracChangeset for help on using the changeset viewer.