CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/HLTrigger/HLTanalyzers/src/HLTJets.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <sstream>
00003 #include <istream>
00004 #include <fstream>
00005 #include <iomanip>
00006 #include <string>
00007 #include <cmath>
00008 #include <functional>
00009 #include <stdlib.h>
00010 #include <string.h>
00011 
00012 #include "HLTrigger/HLTanalyzers/interface/HLTJets.h"
00013 
00014 HLTJets::HLTJets() {
00015     evtCounter=0;
00016     
00017     //set parameter defaults 
00018     _Monte=false;
00019     _Debug=false;
00020     _CalJetMin=0.;
00021     _GenJetMin=0.;
00022 }
00023 
00024 /*  Setup the analysis to put the branch-variables into the tree. */
00025 void HLTJets::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
00026     
00027     edm::ParameterSet myJetParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
00028     std::vector<std::string> parameterNames = myJetParams.getParameterNames() ;
00029     
00030     for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
00031          iParam != parameterNames.end(); iParam++ ){
00032         if  ( (*iParam) == "Monte" ) _Monte =  myJetParams.getParameter<bool>( *iParam );
00033         else if ( (*iParam) == "Debug" ) _Debug =  myJetParams.getParameter<bool>( *iParam );
00034         else if ( (*iParam) == "CalJetMin" ) _CalJetMin =  myJetParams.getParameter<double>( *iParam );
00035         else if ( (*iParam) == "GenJetMin" ) _GenJetMin =  myJetParams.getParameter<double>( *iParam );
00036     }
00037 
00038     const int kMaxRecoPFJet = 10000;
00039     jpfrecopt=new float[kMaxRecoPFJet];
00040     jpfrecophi=new float[kMaxRecoPFJet];
00041     jpfrecoeta=new float[kMaxRecoPFJet];
00042     jpfreconeutralHadronFraction=new float[kMaxRecoPFJet];
00043     jpfreconeutralEMFraction=new float[kMaxRecoPFJet];
00044     jpfrecochargedHadronFraction=new float[kMaxRecoPFJet];
00045     jpfrecochargedEMFraction=new float[kMaxRecoPFJet];
00046     jpfreconeutralMultiplicity=new int[kMaxRecoPFJet];
00047     jpfrecochargedMultiplicity=new int[kMaxRecoPFJet];
00048     
00049     const int kMaxJetCal = 10000;
00050     jcalpt = new float[kMaxJetCal];
00051     jcalphi = new float[kMaxJetCal];
00052     jcaleta = new float[kMaxJetCal];
00053     jcale = new float[kMaxJetCal];
00054     jcalemf = new float[kMaxJetCal]; 
00055     jcaln90 = new float[kMaxJetCal]; 
00056     
00057     jcorcalpt = new float[kMaxJetCal]; 
00058     jcorcalphi = new float[kMaxJetCal]; 
00059     jcorcaleta = new float[kMaxJetCal]; 
00060     jcorcale = new float[kMaxJetCal]; 
00061     jcorcalemf = new float[kMaxJetCal]; 
00062     jcorcaln90 = new float[kMaxJetCal]; 
00063     
00064     const int kMaxJetgen = 10000;
00065     jgenpt = new float[kMaxJetgen];
00066     jgenphi = new float[kMaxJetgen];
00067     jgeneta = new float[kMaxJetgen];
00068     jgene = new float[kMaxJetgen];
00069     const int kMaxTower = 10000;
00070     towet = new float[kMaxTower];
00071     toweta = new float[kMaxTower];
00072     towphi = new float[kMaxTower];
00073     towen = new float[kMaxTower];
00074     towem = new float[kMaxTower];
00075     towhd = new float[kMaxTower];
00076     towoe = new float[kMaxTower];
00077     const int kMaxTau = 500;
00078     l2tauemiso = new float[kMaxTau];
00079     l25tauPt = new float[kMaxTau];
00080     l3tautckiso = new int[kMaxTau];
00081     tauEta = new float[kMaxTau];
00082     tauPt = new float[kMaxTau];
00083     tauPhi = new float[kMaxTau];
00084     
00085     const int kMaxPFTau = 500;
00086     ohpfTauEta         =  new float[kMaxPFTau];
00087     ohpfTauPhi         =  new float[kMaxPFTau];
00088     ohpfTauPt          =  new float[kMaxPFTau];
00089     ohpfTauJetPt       =  new float[kMaxPFTau];
00090     ohpfTauLeadTrackPt =  new float[kMaxPFTau];
00091     ohpfTauLeadPionPt  =  new float[kMaxPFTau];
00092     ohpfTauTrkIso      =  new float[kMaxPFTau];
00093     ohpfTauGammaIso    =  new float[kMaxPFTau];
00094     
00095     recopfTauEta         =  new float[kMaxPFTau];
00096     recopfTauPhi         =  new float[kMaxPFTau];
00097     recopfTauPt          =  new float[kMaxPFTau];
00098     recopfTauJetPt       =  new float[kMaxPFTau];
00099     recopfTauLeadTrackPt =  new float[kMaxPFTau];
00100     recopfTauLeadPionPt  =  new float[kMaxPFTau];
00101     recopfTauTrkIso      =  new int[kMaxPFTau];
00102     recopfTauGammaIso    =  new int[kMaxPFTau];
00103     recopfTauDiscrByTancOnePercent     =  new float[kMaxPFTau];
00104     recopfTauDiscrByTancHalfPercent    =  new float[kMaxPFTau];
00105     recopfTauDiscrByTancQuarterPercent =  new float[kMaxPFTau]; 
00106     recopfTauDiscrByTancTenthPercent   =  new float[kMaxPFTau];
00107     recopfTauDiscrByIso        =  new float[kMaxPFTau]; 
00108     recopfTauDiscrAgainstMuon  =  new float[kMaxPFTau];
00109     recopfTauDiscrAgainstElec  =  new float[kMaxPFTau];
00110     
00111     pfMHT   = -100;    
00112     const int kMaxPFJet = 500;
00113     pfJetEta         = new float[kMaxPFJet];
00114     pfJetPhi         = new float[kMaxPFJet];
00115     pfJetPt         = new float[kMaxPFJet];
00116     
00117     // Jet- MEt-specific branches of the tree 
00118     HltTree->Branch("NrecoJetCal",&njetcal,"NrecoJetCal/I");
00119     HltTree->Branch("NrecoJetGen",&njetgen,"NrecoJetGen/I");
00120     HltTree->Branch("NrecoTowCal",&ntowcal,"NrecoTowCal/I");
00121     HltTree->Branch("recoJetCalPt",jcalpt,"recoJetCalPt[NrecoJetCal]/F");
00122     HltTree->Branch("recoJetCalPhi",jcalphi,"recoJetCalPhi[NrecoJetCal]/F");
00123     HltTree->Branch("recoJetCalEta",jcaleta,"recoJetCalEta[NrecoJetCal]/F");
00124     HltTree->Branch("recoJetCalE",jcale,"recoJetCalE[NrecoJetCal]/F");
00125     HltTree->Branch("recoJetCalEMF",jcalemf,"recoJetCalEMF[NrecoJetCal]/F");
00126     HltTree->Branch("recoJetCalN90",jcaln90,"recoJetCalN90[NrecoJetCal]/F");
00127     
00128     HltTree->Branch("recoJetGenPt",jgenpt,"recoJetGenPt[NrecoJetGen]/F");
00129     HltTree->Branch("recoJetGenPhi",jgenphi,"recoJetGenPhi[NrecoJetGen]/F");
00130     HltTree->Branch("recoJetGenEta",jgeneta,"recoJetGenEta[NrecoJetGen]/F");
00131     HltTree->Branch("recoJetGenE",jgene,"recoJetGenE[NrecoJetGen]/F");
00132     HltTree->Branch("recoTowEt",towet,"recoTowEt[NrecoTowCal]/F");
00133     HltTree->Branch("recoTowEta",toweta,"recoTowEta[NrecoTowCal]/F");
00134     HltTree->Branch("recoTowPhi",towphi,"recoTowPhi[NrecoTowCal]/F");
00135     HltTree->Branch("recoTowE",towen,"recoTowE[NrecoTowCal]/F");
00136     HltTree->Branch("recoTowEm",towem,"recoTowEm[NrecoTowCal]/F");
00137     HltTree->Branch("recoTowHad",towhd,"recoTowHad[NrecoTowCal]/F");
00138     HltTree->Branch("recoTowOE",towoe,"recoTowOE[NrecoTowCal]/F");
00139     HltTree->Branch("recoMetCal",&mcalmet,"recoMetCal/F");
00140     HltTree->Branch("recoMetCalPhi",&mcalphi,"recoMetCalPhi/F");
00141     HltTree->Branch("recoMetCalSum",&mcalsum,"recoMetCalSum/F");
00142     HltTree->Branch("recoMetGen",&mgenmet,"recoMetGen/F");
00143     HltTree->Branch("recoMetGenPhi",&mgenphi,"recoMetGenPhi/F");
00144     HltTree->Branch("recoMetGenSum",&mgensum,"recoMetGenSum/F");
00145     HltTree->Branch("recoHTCal",&htcalet,"recoHTCal/F");
00146     HltTree->Branch("recoHTCalPhi",&htcalphi,"recoHTCalPhi/F");
00147     HltTree->Branch("recoHTCalSum",&htcalsum,"recoHTCalSum/F");
00148     //for(int ieta=0;ieta<NETA;ieta++){std::cout << " ieta " << ieta << " eta min " << CaloTowerEtaBoundries[ieta] <<std::endl;}
00149     
00150     HltTree->Branch("NrecoJetCorCal",&ncorjetcal,"NrecoJetCorCal/I"); 
00151     HltTree->Branch("recoJetCorCalPt",jcorcalpt,"recoJetCorCalPt[NrecoJetCorCal]/F"); 
00152     HltTree->Branch("recoJetCorCalPhi",jcorcalphi,"recoJetCorCalPhi[NrecoJetCorCal]/F"); 
00153     HltTree->Branch("recoJetCorCalEta",jcorcaleta,"recoJetCorCalEta[NrecoJetCorCal]/F"); 
00154     HltTree->Branch("recoJetCorCalE",jcorcale,"recoJetCorCalE[NrecoJetCorCal]/F"); 
00155     HltTree->Branch("recoJetCorCalEMF",jcorcalemf,"recoJetCorCalEMF[NrecoJetCorCal]/F");
00156     HltTree->Branch("recoJetCorCalN90",jcorcaln90,"recoJetCorCalN90[NrecoJetCorCal]/F");
00157     
00158     // Taus
00159     HltTree->Branch("NohTau",&nohtau,"NohTau/I");
00160     HltTree->Branch("ohTauEta",tauEta,"ohTauEta[NohTau]/F");
00161     HltTree->Branch("ohTauPhi",tauPhi,"ohTauPhi[NohTau]/F");
00162     HltTree->Branch("ohTauPt",tauPt,"ohTauPt[NohTau]/F");
00163     HltTree->Branch("ohTauEiso",l2tauemiso,"ohTauEiso[NohTau]/F");
00164     HltTree->Branch("ohTauL25Tpt",l25tauPt,"ohTauL25Tpt[NohTau]/F");
00165     HltTree->Branch("ohTauL3Tiso",l3tautckiso,"ohTauL3Tiso[NohTau]/I");
00166 
00167     //ohpfTaus
00168     HltTree->Branch("NohpfTau",&nohPFTau,"NohpfTau/I");
00169     HltTree->Branch("ohpfTauPt",ohpfTauPt,"ohpfTauPt[NohpfTau]/F");
00170     HltTree->Branch("ohpfTauEta",ohpfTauEta,"ohpfTauEta[NohpfTau]/F");
00171     HltTree->Branch("ohpfTauPhi",ohpfTauPhi,"ohpfTauPhi[NohpfTau]/F");
00172     HltTree->Branch("ohpfTauLeadTrackPt",ohpfTauLeadTrackPt,"ohpfTauLeadTrackPt[NohpfTau]/F");
00173     HltTree->Branch("ohpfTauLeadPionPt",ohpfTauLeadPionPt,"ohpfTauLeadPionPt[NohpfTau]/F");
00174     HltTree->Branch("ohpfTauTrkIso",ohpfTauTrkIso,"ohpfTauTrkIso[NohpfTau]/F");
00175     HltTree->Branch("ohpfTauGammaIso",ohpfTauGammaIso,"ohpfTauGammaIso[NohpfTau]/F");
00176     HltTree->Branch("ohpfTauJetPt",ohpfTauJetPt,"ohpfTauJetPt[NohpfTau]/F");    
00177    
00178    //Reco PFTaus
00179     nRecoPFTau = 0;
00180     HltTree->Branch("NRecoPFTau",&nRecoPFTau,"NRecoPFTau/I");
00181     HltTree->Branch("recopfTauPt",recopfTauPt,"recopfTauPt[NRecoPFTau]/F");
00182     HltTree->Branch("recopfTauEta",recopfTauEta,"recopfTauEta[NRecoPFTau]/F");
00183     HltTree->Branch("recopfTauPhi",recopfTauPhi,"recopfTauPhi[NRecoPFTau]/F");
00184     HltTree->Branch("recopfTauLeadTrackPt",recopfTauLeadTrackPt,"recopfTauLeadTrackPt[NRecoPFTau]/F");
00185     HltTree->Branch("recopfTauLeadPionPt",recopfTauLeadPionPt,"recopfTauLeadPionPt[NRecoPFTau]/F");
00186     HltTree->Branch("recopfTauTrkIso",recopfTauTrkIso,"recopfTauTrkIso[NRecoPFTau]/I");
00187     HltTree->Branch("recopfTauGammaIso",recopfTauGammaIso,"recopfTauGammaIso[NRecoPFTau]/I");
00188     HltTree->Branch("recopfTauJetPt",recopfTauJetPt,"recopfTauJetPt[NRecoPFTau]/F");   
00189     HltTree->Branch("recopfTauDiscrByTancOnePercent",recopfTauDiscrByTancOnePercent,"recopfTauDiscrByTancOnePercent[NRecoPFTau]/F");   
00190     HltTree->Branch("recopfTauDiscrByTancHalfPercent",recopfTauDiscrByTancHalfPercent,"recopfTauDiscrByTancHalfPercent[NRecoPFTau]/F");   
00191     HltTree->Branch("recopfTauDiscrByTancQuarterPercent",recopfTauDiscrByTancQuarterPercent,"recopfTauDiscrByTancQuarterPercent[NRecoPFTau]/F");   
00192     HltTree->Branch("recopfTauDiscrByTancTenthPercent",recopfTauDiscrByTancTenthPercent,"recopfTauDiscrByTancTenthPercent[NRecoPFTau]/F");       
00193     HltTree->Branch("recopfTauDiscrByIso",recopfTauDiscrByIso,"recopfTauDiscrByIso[NRecoPFTau]/F");   
00194     HltTree->Branch("recopfTauDiscrAgainstMuon",recopfTauDiscrAgainstMuon,"recopfTauDiscrAgainstMuon[NRecoPFTau]/F");
00195     HltTree->Branch("recopfTauDiscrAgainstElec",recopfTauDiscrAgainstElec,"recopfTauDiscrAgainstElec[NRecoPFTau]/F");
00196   
00197     //PFJets
00198     HltTree->Branch("pfMHT",&pfMHT,"pfMHT/F");
00199     HltTree->Branch("NohPFJet",&nohPFJet,"NohPFJet/I");
00200     HltTree->Branch("pfJetPt",pfJetPt,"pfJetPt[NohPFJet]/F");
00201     HltTree->Branch("pfJetEta",pfJetEta,"pfJetEta[NohPFJet]/F");
00202     HltTree->Branch("pfJetPhi",pfJetPhi,"pfJetPhi[NohPFJet]/F");
00203 
00204     //RECO PFJets
00205     HltTree->Branch("nrpj",&nrpj,"nrpj/I");
00206     HltTree->Branch("recopfJetpt",                    jpfrecopt,                     "recopfJetpt[nrpj]/F");
00207     HltTree->Branch("recopfJetphi",                   jpfrecophi,                    "recopfJetphi[nrpj]/F");
00208     HltTree->Branch("recopfJeteta",                   jpfrecoeta,                    "recopfJeteta[nrpj]/F");
00209     HltTree->Branch("recopfJetneutralHadronFraction", jpfreconeutralHadronFraction,  "recopfJetneutralHadronFraction[nrpj]/F");
00210     HltTree->Branch("recopfJetneutralEMFraction",     jpfreconeutralEMFraction,      "recopfJetneutralEMFraction[nrpj]/F");
00211     HltTree->Branch("recopfJetchargedHadronFraction", jpfrecochargedHadronFraction,  "recopfJetchargedHadronFraction[nrpj]/F");
00212     HltTree->Branch("recopfJetchargedEMFraction",     jpfrecochargedEMFraction,      "recopfJetchargedEMFraction[nrpj]/F");
00213     HltTree->Branch("recopfJetneutralMultiplicity",   jpfreconeutralMultiplicity,    "recopfJetneutralMultiplicity[nrpj]/I");
00214     HltTree->Branch("recopfJetchargedMultiplicity",   jpfrecochargedMultiplicity,    "recopfJetchargedMultiplicity[nrpj]/I"); 
00215     
00216 }
00217 
00218 /* **Analyze the event** */
00219 void HLTJets::analyze(const edm::Handle<reco::CaloJetCollection>      & calojets,
00220                       const edm::Handle<reco::CaloJetCollection>      & calocorjets,
00221                       const edm::Handle<reco::GenJetCollection>       & genjets,
00222                       const edm::Handle<reco::CaloMETCollection>      & recmets,
00223                       const edm::Handle<reco::GenMETCollection>       & genmets,
00224                       const edm::Handle<reco::METCollection>          & ht,
00225                       const edm::Handle<reco::HLTTauCollection>       & taujets,
00226                       const edm::Handle<reco::PFTauCollection>        & pfTaus,
00227                       const edm::Handle<reco::PFJetCollection>        & pfJets,
00228                       const edm::Handle<reco::PFTauCollection>        & recoPfTaus,  
00229                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCOnePercent,
00230                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCHalfPercent,
00231                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCQuarterPercent,
00232                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCTenthPercent,                      
00233                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByIsolation,
00234                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrAgainstElec,
00235                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrAgainstMuon,
00236                       const edm::Handle<reco::PFJetCollection>        & recoPFJets,
00237                       const edm::Handle<CaloTowerCollection>          & caloTowers,
00238                       double thresholdForSavingTowers, 
00239                       double                minPtCH,
00240                       double               minPtGamma,
00241                       TTree * HltTree) {
00242     
00243     if (_Debug) std::cout << " Beginning HLTJets " << std::endl;
00244     
00245     //initialize branch variables
00246     njetcal=0; ncorjetcal=0; njetgen=0;ntowcal=0;
00247     mcalmet=0.; mcalphi=0.;
00248     mgenmet=0.; mgenphi=0.;
00249     htcalet=0.,htcalphi=0.,htcalsum=0.;
00250     
00251     if (calojets.isValid()) {
00252         reco::CaloJetCollection mycalojets;
00253         mycalojets=*calojets;
00254         std::sort(mycalojets.begin(),mycalojets.end(),PtGreater());
00255         typedef reco::CaloJetCollection::const_iterator cjiter;
00256         int jcal=0;
00257         for ( cjiter i=mycalojets.begin(); i!=mycalojets.end(); i++) {
00258             
00259             if (i->pt()>_CalJetMin){
00260                 jcalpt[jcal] = i->pt();
00261                 jcalphi[jcal] = i->phi();
00262                 jcaleta[jcal] = i->eta();
00263                 jcale[jcal] = i->energy();
00264                 jcalemf[jcal] = i->emEnergyFraction();
00265                 jcaln90[jcal] = i->n90();
00266                 jcal++;
00267             }
00268             
00269         }
00270         njetcal = jcal;
00271     }
00272     else {njetcal = 0;}
00273     
00274     if (calocorjets.isValid()) {
00275         reco::CaloJetCollection mycalocorjets;
00276         mycalocorjets=*calocorjets;
00277         std::sort(mycalocorjets.begin(),mycalocorjets.end(),PtGreater());
00278         typedef reco::CaloJetCollection::const_iterator ccorjiter;
00279         int jcorcal=0;
00280         for ( ccorjiter i=mycalocorjets.begin(); i!=mycalocorjets.end(); i++) {
00281             
00282             if (i->pt()>_CalJetMin){
00283                 jcorcalpt[jcorcal] = i->pt();
00284                 jcorcalphi[jcorcal] = i->phi();
00285                 jcorcaleta[jcorcal] = i->eta();
00286                 jcorcale[jcorcal] = i->energy();
00287                 jcorcalemf[jcorcal] = i->emEnergyFraction();
00288                 jcorcaln90[jcorcal] = i->n90();
00289                 jcorcal++;
00290             }
00291             
00292         }
00293         ncorjetcal = jcorcal;
00294     }
00295     else {ncorjetcal = 0;}
00296     
00297     if (caloTowers.isValid()) {
00298         //    ntowcal = caloTowers->size();
00299         int jtow = 0;
00300         for ( CaloTowerCollection::const_iterator tower=caloTowers->begin(); tower!=caloTowers->end(); tower++) {
00301             if(tower->energy() > thresholdForSavingTowers)
00302             {
00303                 towet[jtow] = tower->et();
00304                 toweta[jtow] = tower->eta();
00305                 towphi[jtow] = tower->phi();
00306                 towen[jtow] = tower->energy();
00307                 towem[jtow] = tower->emEnergy();
00308                 towhd[jtow] = tower->hadEnergy();
00309                 towoe[jtow] = tower->outerEnergy();
00310                 jtow++;
00311             }
00312         }
00313         ntowcal = jtow;
00314     }
00315     else {ntowcal = 0;}
00316     
00317     if (recmets.isValid()) {
00318         typedef reco::CaloMETCollection::const_iterator cmiter;
00319         for ( cmiter i=recmets->begin(); i!=recmets->end(); i++) {
00320             mcalmet = i->pt();
00321             mcalphi = i->phi();
00322             mcalsum = i->sumEt();
00323         }
00324     }
00325     
00326     if (ht.isValid()) {
00327         typedef reco::METCollection::const_iterator iter;
00328         for ( iter i=ht->begin(); i!=ht->end(); i++) {
00329             htcalet = i->pt();
00330             htcalphi = i->phi();
00331             htcalsum = i->sumEt();
00332         }
00333     }
00334     
00335     if (_Monte){
00336         
00337         if (genjets.isValid()) {
00338             reco::GenJetCollection mygenjets;
00339             mygenjets=*genjets;
00340             std::sort(mygenjets.begin(),mygenjets.end(),PtGreater());
00341             typedef reco::GenJetCollection::const_iterator gjiter;
00342             int jgen=0;
00343             for ( gjiter i=mygenjets.begin(); i!=mygenjets.end(); i++) {
00344                 
00345                 if (i->pt()>_GenJetMin){
00346                     jgenpt[jgen] = i->pt();
00347                     jgenphi[jgen] = i->phi();
00348                     jgeneta[jgen] = i->eta();
00349                     jgene[jgen] = i->energy();
00350                     jgen++;
00351                 }
00352                 
00353             }
00354             njetgen = jgen;
00355         }
00356         else {njetgen = 0;}
00357         
00358         if (genmets.isValid()) {
00359             typedef reco::GenMETCollection::const_iterator gmiter;
00360             for ( gmiter i=genmets->begin(); i!=genmets->end(); i++) {
00361                 mgenmet = i->pt();
00362                 mgenphi = i->phi();
00363                 mgensum = i->sumEt();
00364             }
00365         }
00366         
00367     }
00368     
00369     
00371     
00372     if (taujets.isValid()) {      
00373         nohtau = taujets->size();
00374         reco::HLTTauCollection mytaujets;
00375         mytaujets=*taujets;
00376         std::sort(mytaujets.begin(),mytaujets.end(),GetPtGreater());
00377         typedef reco::HLTTauCollection::const_iterator tauit;
00378         int itau=0;
00379         for(tauit i=mytaujets.begin(); i!=mytaujets.end(); i++){
00380             //Ask for Eta,Phi and Et of the tau:
00381             tauEta[itau] = i->getEta();
00382             tauPhi[itau] = i->getPhi();
00383             tauPt[itau] = i->getPt();
00384             //Ask for L2 EMIsolation cut: Nominal cut : < 5
00385             l2tauemiso[itau] = i->getEMIsolationValue();
00386             //Get L25 LeadTrackPt : Nominal cut : > 20 GeV
00387             l25tauPt[itau] = i->getL25LeadTrackPtValue();
00388             //Get TrackIsolation response (returns 0 = failed or 1= passed)
00389             l3tautckiso[itau] = i->getL3TrackIsolationResponse();
00390             //MET : > 65
00391             itau++;
00392         }      
00393     }
00394     else {nohtau = 0;}
00395 
00396     
00398     if(pfTaus.isValid()) {
00399         //float minTrkPt = minPtCH;
00400         //float minGammaPt = minPtGamma;
00401         nohPFTau  = pfTaus->size();
00402         reco::PFTauCollection taus = *pfTaus;
00403         std::sort(taus.begin(),taus.end(),GetPFPtGreater());
00404         typedef reco::PFTauCollection::const_iterator pftauit;
00405         int ipftau=0;
00406         float pfMHTx = 0;
00407         float pfMHTy = 0;
00408         for(pftauit i=taus.begin(); i!=taus.end(); i++){
00409             //Ask for Eta,Phi and Et of the tau:
00410             ohpfTauEta[ipftau] = i->eta();
00411             ohpfTauPhi[ipftau] = i->phi();
00412             ohpfTauPt[ipftau] = i->pt();
00413             ohpfTauJetPt[ipftau] = i->pfTauTagInfoRef()->pfjetRef()->pt();
00414 
00415             pfMHTx = pfMHTx + i->pfTauTagInfoRef()->pfjetRef()->px();
00416             pfMHTy = pfMHTy + i->pfTauTagInfoRef()->pfjetRef()->py();
00417             
00418   /*
00419             if( (i->leadPFCand()).isNonnull())
00420                 pfTauLeadPionPt[ipftau] = i->leadPFCand()->pt();            
00421 */
00422             if( (i->leadPFNeutralCand()).isNonnull())
00423                 ohpfTauLeadPionPt[ipftau] = i->leadPFNeutralCand()->pt();        
00424             if((i->leadPFChargedHadrCand()).isNonnull())
00425                 ohpfTauLeadTrackPt[ipftau] = i->leadPFChargedHadrCand()->pt();
00426             float maxPtTrkIso = 0;
00427             for (unsigned int iTrk = 0; iTrk < i->isolationPFChargedHadrCands().size(); iTrk++)
00428             {
00429                 if(i->isolationPFChargedHadrCands()[iTrk]->pt() > maxPtTrkIso) maxPtTrkIso = i->isolationPFChargedHadrCands()[iTrk]->pt();
00430             }
00431                 
00432             ohpfTauTrkIso[ipftau] = maxPtTrkIso;
00433             float maxPtGammaIso = 0;
00434             for (unsigned int iGamma = 0; iGamma < i->isolationPFGammaCands().size(); iGamma++)
00435             {
00436                 if(i->isolationPFGammaCands()[iGamma]->pt() > maxPtGammaIso) maxPtGammaIso = i->isolationPFGammaCands()[iGamma]->pt();
00437             }                        
00438             ohpfTauGammaIso[ipftau] = maxPtGammaIso;
00439             ipftau++;
00440         } 
00441         pfMHT = sqrt(pfMHTx*pfMHTx + pfMHTy*pfMHTy);
00442         
00443     }
00444     
00446       
00447     if(recoPfTaus.isValid()) {
00448         float minTrkPt = minPtCH;
00449         float minGammaPt = minPtGamma;
00450         nRecoPFTau  = pfTaus->size();
00451         reco::PFTauCollection taus = *recoPfTaus;
00452         std::sort(taus.begin(),taus.end(),GetPFPtGreater());
00453         typedef reco::PFTauCollection::const_iterator pftauit;
00454         int ipftau=0;
00455         
00456         for(pftauit i=taus.begin(); i!=taus.end(); i++){
00457             //Ask for Eta,Phi and Et of the tau:
00458             recopfTauEta[ipftau] = i->eta();
00459             recopfTauPhi[ipftau] = i->phi();
00460             recopfTauPt[ipftau]  = i->pt();
00461             recopfTauJetPt[ipftau] = i->jetRef()->pt();
00462 
00463             if( (i->leadPFNeutralCand()).isNonnull())
00464                 recopfTauLeadPionPt[ipftau] = i->leadPFNeutralCand()->pt();        
00465             if((i->leadPFChargedHadrCand()).isNonnull())
00466                 recopfTauLeadTrackPt[ipftau] = i->leadPFChargedHadrCand()->pt();
00467             int myTrks=0;
00468             for (unsigned int iTrk = 0; iTrk < i->isolationPFChargedHadrCands().size(); iTrk++)
00469             {
00470                 if(i->isolationPFChargedHadrCands()[iTrk]->pt() > minTrkPt) myTrks++;
00471             }
00472                
00473             recopfTauTrkIso[ipftau] = myTrks;
00474             int myGammas=0;
00475             for (unsigned int iGamma = 0; iGamma < i->isolationPFGammaCands().size(); iGamma++)
00476             {
00477                 if(i->isolationPFGammaCands()[iGamma]->pt() > minGammaPt) myGammas++;
00478             }                        
00479             recopfTauGammaIso[ipftau] = myGammas;
00480             
00481             const reco::PFTauRef thisTauRef(recoPfTaus,ipftau);
00482             
00483             if(theRecoPFTauDiscrByTanCOnePercent.isValid()){
00484             recopfTauDiscrByTancOnePercent[ipftau] = (*theRecoPFTauDiscrByTanCOnePercent)[thisTauRef];}
00485             if(theRecoPFTauDiscrByIsolation.isValid()){ 
00486             recopfTauDiscrByIso[ipftau] = (*theRecoPFTauDiscrByIsolation)[thisTauRef];} 
00487             if(theRecoPFTauDiscrAgainstMuon.isValid()){
00488             recopfTauDiscrAgainstMuon[ipftau] = (*theRecoPFTauDiscrAgainstMuon)[thisTauRef];}
00489             if(theRecoPFTauDiscrAgainstElec.isValid()){
00490             recopfTauDiscrAgainstElec[ipftau] = (*theRecoPFTauDiscrAgainstElec)[thisTauRef];}
00491             if(theRecoPFTauDiscrByTanCHalfPercent.isValid()){
00492             recopfTauDiscrByTancHalfPercent[ipftau] = (*theRecoPFTauDiscrByTanCHalfPercent)[thisTauRef];}
00493             if(theRecoPFTauDiscrByTanCQuarterPercent.isValid()){
00494             recopfTauDiscrByTancQuarterPercent[ipftau] = (*theRecoPFTauDiscrByTanCQuarterPercent)[thisTauRef];}
00495             if(theRecoPFTauDiscrByTanCTenthPercent.isValid()){
00496             recopfTauDiscrByTancTenthPercent[ipftau] = (*theRecoPFTauDiscrByTanCTenthPercent)[thisTauRef];}
00497 
00498             ipftau++;
00499         }        
00500     }
00501    
00503     if(pfJets.isValid()) {
00504         nohPFJet  = pfJets->size();
00505         reco::PFJetCollection Jets = *pfJets;
00506         std::sort(Jets.begin(),Jets.end(),GetPFPtGreater());
00507         typedef reco::PFJetCollection::const_iterator pfJetit;
00508         int ipfJet=0;
00509         float pfMHTx = 0;
00510         float pfMHTy = 0;
00511         for(pfJetit i=Jets.begin(); i!=Jets.end(); i++){
00512             //Ask for Eta,Phi and Et of the Jet:
00513             pfJetEta[ipfJet] = i->eta();
00514             pfJetPhi[ipfJet] = i->phi();
00515             pfJetPt[ipfJet] = i->pt();           
00516             
00517             pfMHTx = pfMHTx + i->px();
00518             pfMHTy = pfMHTy + i->py();
00519             ipfJet++;   
00520         } 
00521         pfMHT = sqrt(pfMHTx*pfMHTx + pfMHTy*pfMHTy);
00522         
00523     }
00524     
00526     nrpj = 0;
00527     if(recoPFJets.isValid()){
00528             nrpj = recoPFJets->size();
00529             reco::PFJetCollection Jets = *recoPFJets;
00530             std::sort(Jets.begin(),Jets.end(),GetPFPtGreater());
00531             typedef reco::PFJetCollection::const_iterator pfJetit;
00532             int ipfJet=0;
00533             for(pfJetit i=Jets.begin(); i!=Jets.end(); i++){
00534                     //Ask for Eta,Phi and Et of the Jet:
00535                     jpfrecoeta[ipfJet] = i->eta();
00536                     jpfrecophi[ipfJet] = i->phi();
00537                     jpfrecopt[ipfJet] = i->pt();           
00538                     jpfreconeutralHadronFraction[ipfJet] = i->neutralHadronEnergyFraction ();
00539                     jpfrecochargedHadronFraction[ipfJet] = i->chargedHadronEnergyFraction ();
00540                     jpfreconeutralMultiplicity[ipfJet] = i->neutralMultiplicity ();
00541                     jpfrecochargedMultiplicity[ipfJet] = i->chargedMultiplicity ();
00542                     jpfreconeutralEMFraction[ipfJet] = i->neutralEmEnergyFraction ();
00543                     jpfrecochargedEMFraction[ipfJet] = i->chargedEmEnergyFraction ();
00544 
00545                     ipfJet++;   
00546             } 
00547 
00548     }
00549     
00550 }