CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/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 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00014 #include "DataFormats/TrackReco/interface/Track.h"
00015 
00016 HLTJets::HLTJets() {
00017     evtCounter=0;
00018     
00019     //set parameter defaults 
00020     _Monte=false;
00021     _Debug=false;
00022     _CalJetMin=0.;
00023     _GenJetMin=0.;
00024 }
00025 
00026 /*  Setup the analysis to put the branch-variables into the tree. */
00027 void HLTJets::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
00028     
00029     edm::ParameterSet myJetParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
00030     std::vector<std::string> parameterNames = myJetParams.getParameterNames() ;
00031     
00032     for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
00033          iParam != parameterNames.end(); iParam++ ){
00034         if  ( (*iParam) == "Monte" ) _Monte =  myJetParams.getParameter<bool>( *iParam );
00035         else if ( (*iParam) == "Debug" ) _Debug =  myJetParams.getParameter<bool>( *iParam );
00036         else if ( (*iParam) == "CalJetMin" ) _CalJetMin =  myJetParams.getParameter<double>( *iParam );
00037         else if ( (*iParam) == "GenJetMin" ) _GenJetMin =  myJetParams.getParameter<double>( *iParam );
00038     }
00039 
00040     jetID = new reco::helper::JetIDHelper(pSet.getParameter<edm::ParameterSet>("JetIDParams"));
00041 
00042     const int kMaxRecoPFJet = 10000;
00043     jpfrecopt=new float[kMaxRecoPFJet];
00044     jpfrecoe=new float[kMaxRecoPFJet];
00045     jpfrecophi=new float[kMaxRecoPFJet];
00046     jpfrecoeta=new float[kMaxRecoPFJet];
00047     jpfreconeutralHadronFraction=new float[kMaxRecoPFJet];
00048     jpfreconeutralEMFraction=new float[kMaxRecoPFJet];
00049     jpfrecochargedHadronFraction=new float[kMaxRecoPFJet];
00050     jpfrecochargedEMFraction=new float[kMaxRecoPFJet];
00051     jpfreconeutralMultiplicity=new int[kMaxRecoPFJet];
00052     jpfrecochargedMultiplicity=new int[kMaxRecoPFJet];
00053     
00054     const int kMaxJetCal = 10000;
00055     jhcalpt = new float[kMaxJetCal];
00056     jhcalphi = new float[kMaxJetCal];
00057     jhcaleta = new float[kMaxJetCal];
00058     jhcale = new float[kMaxJetCal];
00059     jhcalemf = new float[kMaxJetCal]; 
00060     jhcaln90 = new float[kMaxJetCal]; 
00061     jhcaln90hits = new float[kMaxJetCal]; 
00062     
00063     jhcorcalpt = new float[kMaxJetCal]; 
00064     jhcorcalphi = new float[kMaxJetCal]; 
00065     jhcorcaleta = new float[kMaxJetCal]; 
00066     jhcorcale = new float[kMaxJetCal]; 
00067     jhcorcalemf = new float[kMaxJetCal]; 
00068     jhcorcaln90 = new float[kMaxJetCal]; 
00069     jhcorcaln90hits = new float[kMaxJetCal]; 
00070 
00071     jhcorL1L2L3calpt = new float[kMaxJetCal]; 
00072     jhcorL1L2L3calphi = new float[kMaxJetCal]; 
00073     jhcorL1L2L3caleta = new float[kMaxJetCal]; 
00074     jhcorL1L2L3cale = new float[kMaxJetCal]; 
00075     jhcorL1L2L3calemf = new float[kMaxJetCal]; 
00076     jhcorL1L2L3caln90 = new float[kMaxJetCal]; 
00077     jhcorL1L2L3caln90hits = new float[kMaxJetCal]; 
00078     
00079     jrcalpt = new float[kMaxJetCal];
00080     jrcalphi = new float[kMaxJetCal];
00081     jrcaleta = new float[kMaxJetCal];
00082     jrcale = new float[kMaxJetCal];
00083     jrcalemf = new float[kMaxJetCal]; 
00084     jrcaln90 = new float[kMaxJetCal]; 
00085     jrcaln90hits = new float[kMaxJetCal]; 
00086 
00087     jrcorcalpt = new float[kMaxJetCal];
00088     jrcorcalphi = new float[kMaxJetCal];
00089     jrcorcaleta = new float[kMaxJetCal];
00090     jrcorcale = new float[kMaxJetCal];
00091     jrcorcalemf = new float[kMaxJetCal]; 
00092     jrcorcaln90 = new float[kMaxJetCal]; 
00093     jrcorcaln90hits = new float[kMaxJetCal]; 
00094 
00095     const int kMaxJetgen = 10000;
00096     jgenpt = new float[kMaxJetgen];
00097     jgenphi = new float[kMaxJetgen];
00098     jgeneta = new float[kMaxJetgen];
00099     jgene = new float[kMaxJetgen];
00100     const int kMaxTower = 10000;
00101     towet = new float[kMaxTower];
00102     toweta = new float[kMaxTower];
00103     towphi = new float[kMaxTower];
00104     towen = new float[kMaxTower];
00105     towem = new float[kMaxTower];
00106     towhd = new float[kMaxTower];
00107     towoe = new float[kMaxTower];
00108     const int kMaxTau = 500;
00109     l2tauemiso = new float[kMaxTau];
00110     l25tauPt = new float[kMaxTau];
00111     l3tautckiso = new int[kMaxTau];
00112     tauEta = new float[kMaxTau];
00113     tauPt = new float[kMaxTau];
00114     tauPhi = new float[kMaxTau];
00115     
00116     const int kMaxPFTau = 500;
00117     ohpfTauEta         =  new float[kMaxPFTau];
00118     ohpfTauProngs      =  new int[kMaxPFTau];
00119     ohpfTauPhi         =  new float[kMaxPFTau];
00120     ohpfTauPt          =  new float[kMaxPFTau];
00121     ohpfTauJetPt       =  new float[kMaxPFTau];
00122     ohpfTauLeadTrackPt =  new float[kMaxPFTau];
00123     ohpfTauLeadPionPt  =  new float[kMaxPFTau];
00124     ohpfTauTrkIso      =  new float[kMaxPFTau];
00125     ohpfTauGammaIso    =  new float[kMaxPFTau];
00126 
00127     ohpfTauTightConeProngs      =  new int[kMaxPFTau];
00128     ohpfTauTightConeEta         =  new float[kMaxPFTau];
00129     ohpfTauTightConePhi         =  new float[kMaxPFTau];
00130     ohpfTauTightConePt          =  new float[kMaxPFTau];
00131     ohpfTauTightConeJetPt       =  new float[kMaxPFTau];
00132     ohpfTauTightConeLeadTrackPt =  new float[kMaxPFTau];
00133     ohpfTauTightConeLeadPionPt  =  new float[kMaxPFTau];
00134     ohpfTauTightConeTrkIso      =  new float[kMaxPFTau];
00135     ohpfTauTightConeGammaIso    =  new float[kMaxPFTau];
00136     
00137     recopfTauEta         =  new float[kMaxPFTau];
00138     recopfTauPhi         =  new float[kMaxPFTau];
00139     recopfTauPt          =  new float[kMaxPFTau];
00140     recopfTauJetPt       =  new float[kMaxPFTau];
00141     recopfTauLeadTrackPt =  new float[kMaxPFTau];
00142     recopfTauLeadPionPt  =  new float[kMaxPFTau];
00143     recopfTauTrkIso      =  new int[kMaxPFTau];
00144     recopfTauGammaIso    =  new int[kMaxPFTau];
00145     recopfTauDiscrByTancOnePercent     =  new float[kMaxPFTau];
00146     recopfTauDiscrByTancHalfPercent    =  new float[kMaxPFTau];
00147     recopfTauDiscrByTancQuarterPercent =  new float[kMaxPFTau]; 
00148     recopfTauDiscrByTancTenthPercent   =  new float[kMaxPFTau];
00149     recopfTauDiscrByIso        =  new float[kMaxPFTau]; 
00150     recopfTauDiscrAgainstMuon  =  new float[kMaxPFTau];
00151     recopfTauDiscrAgainstElec  =  new float[kMaxPFTau];
00152     
00153     pfHT    = -100.;
00154     pfMHT   = -100.;    
00155     const int kMaxPFJet = 500;
00156     pfJetEta         = new float[kMaxPFJet];
00157     pfJetPhi         = new float[kMaxPFJet];
00158     pfJetPt         = new float[kMaxPFJet];
00159     pfJetE         = new float[kMaxPFJet];
00160 
00161 
00162     const int kMaxTauIso = 5000;
00163 
00164     // for offlineHPStau 
00165     signalTrToPFTauMatch = new int[kMaxTauIso];// index of reconstructed tau in tau collection
00166     recoPFTauSignalTrDz = new float[kMaxTauIso];
00167     recoPFTauSignalTrPt = new float[kMaxTauIso];
00168 
00169     isoTrToPFTauMatch = new int[kMaxTauIso]; // index of reconstructed tau in tau collection
00170     recoPFTauIsoTrDz = new float[kMaxTauIso];
00171     recoPFTauIsoTrPt = new float[kMaxTauIso];
00172 
00173     // HLT pf taus
00174     hltpftauSignalTrToPFTauMatch = new int[kMaxTauIso]; // index of HLTPF tau in tau collection
00175     HLTPFTauSignalTrDz = new float[kMaxTauIso];
00176     HLTPFTauSignalTrPt = new float[kMaxTauIso];
00177 
00178     hltpftauIsoTrToPFTauMatch = new int[kMaxTauIso]; // index of HLTPF tau in tau collection
00179     HLTPFTauIsoTrDz = new float[kMaxTauIso];
00180     HLTPFTauIsoTrPt = new float[kMaxTauIso];
00181 
00182     // offline pftau isolation and signal cands
00183     HltTree->Branch("NoRecoPFTausSignal",&noRecoPFTausSignal,"NoRecoPFTausSignal/I");
00184     HltTree->Branch("signalTrToPFTauMatch", signalTrToPFTauMatch,"signalTrToPFTauMatch[NoRecoPFTausSignal]/I");
00185     HltTree->Branch("recoPFTauSignalTrDz", recoPFTauSignalTrDz,"recoPFTauSignalTrDz[NoRecoPFTausSignal]/F");
00186     HltTree->Branch("recoPFTauSignalTrPt", recoPFTauSignalTrPt,"recoPFTauSignalTrPt[NoRecoPFTausSignal]/F");
00187 
00188     HltTree->Branch("NoRecoPFTausIso",&noRecoPFTausIso,"NoRecoPFTausIso/I");
00189     HltTree->Branch("isoTrToPFTauMatch", isoTrToPFTauMatch,"isoTrToPFTauMatch[NoRecoPFTausIso]/I");
00190     HltTree->Branch("recoPFTauIsoTrDz", recoPFTauIsoTrDz,"recoPFTauIsoTrDz[NoRecoPFTausIso]/F");
00191     HltTree->Branch("recoPFTauIsoTrPt", recoPFTauIsoTrPt,"recoPFTauIsoTrPt[NoRecoPFTausIso]/F");
00192 
00193     // HLT pftau isolation and signal cands
00194     HltTree->Branch("NoHLTPFTausSignal",&noHLTPFTausSignal,"NoHLTPFTausSignal/I");
00195     HltTree->Branch("hltpftauSignalTrToPFTauMatch",
00196                     hltpftauSignalTrToPFTauMatch,"hltpftauSignalTrToPFTauMatch[NoHLTPFTausSignal]/I");
00197     HltTree->Branch("HLTPFTauSignalTrDz", HLTPFTauSignalTrDz,"HLTPFTauSignalTrDz[NoHLTPFTausSignal]/F");
00198     HltTree->Branch("HLTPFTauSignalTrPt", HLTPFTauSignalTrPt,"HLTPFTauSignalTrPt[NoHLTPFTausSignal]/F");
00199 
00200     HltTree->Branch("NoHLTPFTausIso",&noHLTPFTausIso,"NoHLTPFTausIso/I");
00201     HltTree->Branch("hltpftauIsoTrToPFTauMatch",
00202                     hltpftauIsoTrToPFTauMatch,"hltpftauIsoTrToPFTauMatch[NoHLTPFTausIso]/I");
00203     HltTree->Branch("HLTPFTauIsoTrDz", HLTPFTauIsoTrDz,"HLTPFTauIsoTrDz[NoHLTPFTausIso]/F");
00204     HltTree->Branch("HLTPFTauIsoTrPt", HLTPFTauIsoTrPt,"HLTPFTauIsoTrPt[NoHLTPFTausIso]/F");
00205     
00206     // Jet- MEt-specific branches of the tree 
00207 
00208     HltTree->Branch("NrecoJetGen",&njetgen,"NrecoJetGen/I");
00209     HltTree->Branch("NrecoTowCal",&ntowcal,"NrecoTowCal/I");
00210 
00211     //ccla RECO JETs
00212     HltTree->Branch("NrecoJetCal",&nrjetcal,"NrecoJetCal/I");
00213     HltTree->Branch("recoJetCalPt",jrcalpt,"recoJetCalPt[NrecoJetCal]/F");
00214     HltTree->Branch("recoJetCalPhi",jrcalphi,"recoJetCalPhi[NrecoJetCal]/F");
00215     HltTree->Branch("recoJetCalEta",jrcaleta,"recoJetCalEta[NrecoJetCal]/F");
00216     HltTree->Branch("recoJetCalE",jrcale,"recoJetCalE[NrecoJetCal]/F");
00217     HltTree->Branch("recoJetCalEMF",jrcalemf,"recoJetCalEMF[NrecoJetCal]/F");
00218     HltTree->Branch("recoJetCalN90",jrcaln90,"recoJetCalN90[NrecoJetCal]/F");
00219     HltTree->Branch("recoJetCalN90hits",jrcaln90hits,"recoJetCalN90hits[NrecoJetCal]/F");
00220 
00221     HltTree->Branch("NrecoJetCorCal",&nrcorjetcal,"NrecoJetCorCal/I"); 
00222     HltTree->Branch("recoJetCorCalPt",jrcorcalpt,"recoJetCorCalPt[NrecoJetCorCal]/F"); 
00223     HltTree->Branch("recoJetCorCalPhi",jrcorcalphi,"recoJetCorCalPhi[NrecoJetCorCal]/F"); 
00224     HltTree->Branch("recoJetCorCalEta",jrcorcaleta,"recoJetCorCalEta[NrecoJetCorCal]/F"); 
00225     HltTree->Branch("recoJetCorCalE",jrcorcale,"recoJetCorCalE[NrecoJetCorCal]/F"); 
00226     HltTree->Branch("recoJetCorCalEMF",jrcorcalemf,"recoJetCorCalEMF[NrecoJetCorCal]/F");
00227     HltTree->Branch("recoJetCorCalN90",jrcorcaln90,"recoJetCorCalN90[NrecoJetCorCal]/F");
00228     HltTree->Branch("recoJetCorCalN90hits",jrcorcaln90hits,"recoJetCorCalN90hits[NrecoJetCorCal]/F");
00229  
00230     //ccla HLTJETS
00231     HltTree->Branch("NohJetCal",&nhjetcal,"NohJetCal/I");
00232     HltTree->Branch("ohJetCalPt",jhcalpt,"ohJetCalPt[NohJetCal]/F");
00233     HltTree->Branch("ohJetCalPhi",jhcalphi,"ohJetCalPhi[NohJetCal]/F");
00234     HltTree->Branch("ohJetCalEta",jhcaleta,"ohJetCalEta[NohJetCal]/F");
00235     HltTree->Branch("ohJetCalE",jhcale,"ohJetCalE[NohJetCal]/F");
00236     HltTree->Branch("ohJetCalEMF",jhcalemf,"ohJetCalEMF[NohJetCal]/F");
00237     HltTree->Branch("ohJetCalN90",jhcaln90,"ohJetCalN90[NohJetCal]/F");
00238     HltTree->Branch("ohJetCalN90hits",jhcaln90hits,"ohJetCalN90hits[NohJetCal]/F");
00239 
00240     HltTree->Branch("NohJetCorCal",&nhcorjetcal,"NohJetCorCal/I");
00241     HltTree->Branch("ohJetCorCalPt",jhcorcalpt,"ohJetCorCalPt[NohJetCorCal]/F");
00242     HltTree->Branch("ohJetCorCalPhi",jhcorcalphi,"ohJetCorCalPhi[NohJetCorCal]/F");
00243     HltTree->Branch("ohJetCorCalEta",jhcorcaleta,"ohJetCorCalEta[NohJetCorCal]/F");
00244     HltTree->Branch("ohJetCorCalE",jhcorcale,"ohJetCorCalE[NohJetCorCal]/F");
00245     HltTree->Branch("ohJetCorCalEMF",jhcorcalemf,"ohJetCorCalEMF[NohJetCorCal]/F");
00246     HltTree->Branch("ohJetCorCalN90",jhcorcaln90,"ohJetCorCalN90[NohJetCorCal]/F");
00247     HltTree->Branch("ohJetCorCalN90hits",jhcorcaln90hits,"ohJetCorCalN90hits[NohJetCorCal]/F");
00248 
00249     HltTree->Branch("NohJetCorL1L2L3Cal",&nhcorL1L2L3jetcal,"NohJetCorL1L2L3Cal/I");
00250     HltTree->Branch("ohJetCorL1L2L3CalPt",jhcorL1L2L3calpt,"ohJetCorL1L2L3CalPt[NohJetCorL1L2L3Cal]/F");
00251     HltTree->Branch("ohJetCorL1L2L3CalPhi",jhcorL1L2L3calphi,"ohJetCorL1L2L3CalPhi[NohJetCorL1L2L3Cal]/F");
00252     HltTree->Branch("ohJetCorL1L2L3CalEta",jhcorL1L2L3caleta,"ohJetCorL1L2L3CalEta[NohJetCorL1L2L3Cal]/F");
00253     HltTree->Branch("ohJetCorL1L2L3CalE",jhcorL1L2L3cale,"ohJetCorL1L2L3CalE[NohJetCorL1L2L3Cal]/F");
00254     HltTree->Branch("ohJetCorL1L2L3CalEMF",jhcorL1L2L3calemf,"ohJetCorL1L2L3CalEMF[NohJetCorL1L2L3Cal]/F");
00255     HltTree->Branch("ohJetCorL1L2L3CalN90",jhcorL1L2L3caln90,"ohJetCorL1L2L3CalN90[NohJetCorL1L2L3Cal]/F");
00256     HltTree->Branch("ohJetCorL1L2L3CalN90hits",jhcorL1L2L3caln90hits,"ohJetCorL1L2L3CalN90hits[NohJetCorL1L2L3Cal]/F");
00257 
00258     //ccla GenJets
00259     HltTree->Branch("recoJetGenPt",jgenpt,"recoJetGenPt[NrecoJetGen]/F");
00260     HltTree->Branch("recoJetGenPhi",jgenphi,"recoJetGenPhi[NrecoJetGen]/F");
00261     HltTree->Branch("recoJetGenEta",jgeneta,"recoJetGenEta[NrecoJetGen]/F");
00262     HltTree->Branch("recoJetGenE",jgene,"recoJetGenE[NrecoJetGen]/F");
00263 
00264     HltTree->Branch("recoTowEt",towet,"recoTowEt[NrecoTowCal]/F");
00265     HltTree->Branch("recoTowEta",toweta,"recoTowEta[NrecoTowCal]/F");
00266     HltTree->Branch("recoTowPhi",towphi,"recoTowPhi[NrecoTowCal]/F");
00267     HltTree->Branch("recoTowE",towen,"recoTowE[NrecoTowCal]/F");
00268     HltTree->Branch("recoTowEm",towem,"recoTowEm[NrecoTowCal]/F");
00269     HltTree->Branch("recoTowHad",towhd,"recoTowHad[NrecoTowCal]/F");
00270     HltTree->Branch("recoTowOE",towoe,"recoTowOE[NrecoTowCal]/F");
00271     HltTree->Branch("recoMetCal",&mcalmet,"recoMetCal/F");
00272     HltTree->Branch("recoMetCalPhi",&mcalphi,"recoMetCalPhi/F");
00273     HltTree->Branch("recoMetCalSum",&mcalsum,"recoMetCalSum/F");
00274     HltTree->Branch("recoMetGen",&mgenmet,"recoMetGen/F");
00275     HltTree->Branch("recoMetGenPhi",&mgenphi,"recoMetGenPhi/F");
00276     HltTree->Branch("recoMetGenSum",&mgensum,"recoMetGenSum/F");
00277     HltTree->Branch("recoHTCal",&htcalet,"recoHTCal/F");
00278     HltTree->Branch("recoHTCalPhi",&htcalphi,"recoHTCalPhi/F");
00279     HltTree->Branch("recoHTCalSum",&htcalsum,"recoHTCalSum/F");
00280     HltTree->Branch("recoMetPF", &pfmet, "recoMetPF/F");
00281     HltTree->Branch("recoMetPFSum", &pfsumet, "recoMetPFSum/F");
00282     HltTree->Branch("recoMetPFPhi", &pfmetphi, "recoMetPFPhi/F");
00283 
00284     //for(int ieta=0;ieta<NETA;ieta++){std::cout << " ieta " << ieta << " eta min " << CaloTowerEtaBoundries[ieta] <<std::endl;}
00285     
00286     
00287     // Taus
00288     nohtau = 0;
00289     HltTree->Branch("NohTau",&nohtau,"NohTau/I");
00290     HltTree->Branch("ohTauEta",tauEta,"ohTauEta[NohTau]/F");
00291     HltTree->Branch("ohTauPhi",tauPhi,"ohTauPhi[NohTau]/F");
00292     HltTree->Branch("ohTauPt",tauPt,"ohTauPt[NohTau]/F");
00293     HltTree->Branch("ohTauEiso",l2tauemiso,"ohTauEiso[NohTau]/F");
00294     HltTree->Branch("ohTauL25Tpt",l25tauPt,"ohTauL25Tpt[NohTau]/F");
00295     HltTree->Branch("ohTauL3Tiso",l3tautckiso,"ohTauL3Tiso[NohTau]/I");
00296 
00297     //ohpfTaus
00298     nohPFTau = 0;
00299     HltTree->Branch("NohpfTau",&nohPFTau,"NohpfTau/I");
00300     HltTree->Branch("ohpfTauPt",ohpfTauPt,"ohpfTauPt[NohpfTau]/F");
00301     HltTree->Branch("ohpfTauProngs",ohpfTauProngs,"ohpfTauProngs[NohpfTau]/I");
00302     HltTree->Branch("ohpfTauEta",ohpfTauEta,"ohpfTauEta[NohpfTau]/F");
00303     HltTree->Branch("ohpfTauPhi",ohpfTauPhi,"ohpfTauPhi[NohpfTau]/F");
00304     HltTree->Branch("ohpfTauLeadTrackPt",ohpfTauLeadTrackPt,"ohpfTauLeadTrackPt[NohpfTau]/F");
00305     HltTree->Branch("ohpfTauLeadPionPt",ohpfTauLeadPionPt,"ohpfTauLeadPionPt[NohpfTau]/F");
00306     HltTree->Branch("ohpfTauTrkIso",ohpfTauTrkIso,"ohpfTauTrkIso[NohpfTau]/F");
00307     HltTree->Branch("ohpfTauGammaIso",ohpfTauGammaIso,"ohpfTauGammaIso[NohpfTau]/F");
00308     HltTree->Branch("ohpfTauJetPt",ohpfTauJetPt,"ohpfTauJetPt[NohpfTau]/F");    
00309 
00310     //ohpfTaus tight cone
00311     nohPFTauTightCone = 0;
00312     HltTree->Branch("NohpfTauTightCone",&nohPFTauTightCone,"NohpfTauTightCone/I");
00313     HltTree->Branch("ohpfTauTightConePt",ohpfTauTightConePt,"ohpfTauTightConePt[NohpfTauTightCone]/F");
00314     HltTree->Branch("ohpfTauTightConeProngs",ohpfTauTightConeProngs,"ohpfTauProngs[NohpfTauTightCone]/I");
00315     HltTree->Branch("ohpfTauTightConeEta",ohpfTauTightConeEta,"ohpfTauEta[NohpfTauTightCone]/F");
00316     HltTree->Branch("ohpfTauTightConePhi",ohpfTauTightConePhi,"ohpfTauPhi[NohpfTauTightCone]/F");
00317     HltTree->Branch("ohpfTauTightConeLeadTrackPt",ohpfTauTightConeLeadTrackPt,"ohpfTauTightConeLeadTrackPt[NohpfTauTightCone]/F");
00318     HltTree->Branch("ohpfTauTightConeLeadPionPt",ohpfTauTightConeLeadPionPt,"ohpfTauTightConeLeadPionPt[NohpfTauTightCone]/F");
00319     HltTree->Branch("ohpfTauTightConeTrkIso",ohpfTauTightConeTrkIso,"ohpfTauTightConeTrkIso[NohpfTauTightCone]/F");
00320     HltTree->Branch("ohpfTauTightConeGammaIso",ohpfTauTightConeGammaIso,"ohpfTauTightConeGammaIso[NohpfTauTightCone]/F");
00321     HltTree->Branch("ohpfTauTightConeJetPt",ohpfTauTightConeJetPt,"ohpfTauTightConeJetPt[NohpfTauTightCone]/F");
00322    
00323    //Reco PFTaus
00324     nRecoPFTau = 0;
00325     HltTree->Branch("NRecoPFTau",&nRecoPFTau,"NRecoPFTau/I");
00326     HltTree->Branch("recopfTauPt",recopfTauPt,"recopfTauPt[NRecoPFTau]/F");
00327     HltTree->Branch("recopfTauEta",recopfTauEta,"recopfTauEta[NRecoPFTau]/F");
00328     HltTree->Branch("recopfTauPhi",recopfTauPhi,"recopfTauPhi[NRecoPFTau]/F");
00329     HltTree->Branch("recopfTauLeadTrackPt",recopfTauLeadTrackPt,"recopfTauLeadTrackPt[NRecoPFTau]/F");
00330     HltTree->Branch("recopfTauLeadPionPt",recopfTauLeadPionPt,"recopfTauLeadPionPt[NRecoPFTau]/F");
00331     HltTree->Branch("recopfTauTrkIso",recopfTauTrkIso,"recopfTauTrkIso[NRecoPFTau]/I");
00332     HltTree->Branch("recopfTauGammaIso",recopfTauGammaIso,"recopfTauGammaIso[NRecoPFTau]/I");
00333     HltTree->Branch("recopfTauJetPt",recopfTauJetPt,"recopfTauJetPt[NRecoPFTau]/F");   
00334     HltTree->Branch("recopfTauDiscrByTancOnePercent",recopfTauDiscrByTancOnePercent,"recopfTauDiscrByTancOnePercent[NRecoPFTau]/F");   
00335     HltTree->Branch("recopfTauDiscrByTancHalfPercent",recopfTauDiscrByTancHalfPercent,"recopfTauDiscrByTancHalfPercent[NRecoPFTau]/F");   
00336     HltTree->Branch("recopfTauDiscrByTancQuarterPercent",recopfTauDiscrByTancQuarterPercent,"recopfTauDiscrByTancQuarterPercent[NRecoPFTau]/F");   
00337     HltTree->Branch("recopfTauDiscrByTancTenthPercent",recopfTauDiscrByTancTenthPercent,"recopfTauDiscrByTancTenthPercent[NRecoPFTau]/F");       
00338     HltTree->Branch("recopfTauDiscrByIso",recopfTauDiscrByIso,"recopfTauDiscrByIso[NRecoPFTau]/F");   
00339     HltTree->Branch("recopfTauDiscrAgainstMuon",recopfTauDiscrAgainstMuon,"recopfTauDiscrAgainstMuon[NRecoPFTau]/F");
00340     HltTree->Branch("recopfTauDiscrAgainstElec",recopfTauDiscrAgainstElec,"recopfTauDiscrAgainstElec[NRecoPFTau]/F");
00341   
00342     //PFJets
00343     nohPFJet = 0;
00344     HltTree->Branch("pfHT",&pfHT,"pfHT/F");
00345     HltTree->Branch("pfMHT",&pfMHT,"pfMHT/F");
00346     HltTree->Branch("NohPFJet",&nohPFJet,"NohPFJet/I");
00347     HltTree->Branch("pfJetPt",pfJetPt,"pfJetPt[NohPFJet]/F");
00348     HltTree->Branch("pfJetE",pfJetE,"pfJetE[NohPFJet]/F");
00349     HltTree->Branch("pfJetEta",pfJetEta,"pfJetEta[NohPFJet]/F");
00350     HltTree->Branch("pfJetPhi",pfJetPhi,"pfJetPhi[NohPFJet]/F");
00351 
00352     //RECO PFJets
00353     HltTree->Branch("nrpj",&nrpj,"nrpj/I");
00354     HltTree->Branch("recopfJetpt",                    jpfrecopt,                     "recopfJetpt[nrpj]/F");
00355     HltTree->Branch("recopfJete",                     jpfrecoe,                      "recopfJete[nrpj]/F");
00356     HltTree->Branch("recopfJetphi",                   jpfrecophi,                    "recopfJetphi[nrpj]/F");
00357     HltTree->Branch("recopfJeteta",                   jpfrecoeta,                    "recopfJeteta[nrpj]/F");
00358     HltTree->Branch("recopfJetneutralHadronFraction", jpfreconeutralHadronFraction,  "recopfJetneutralHadronFraction[nrpj]/F");
00359     HltTree->Branch("recopfJetneutralEMFraction",     jpfreconeutralEMFraction,      "recopfJetneutralEMFraction[nrpj]/F");
00360     HltTree->Branch("recopfJetchargedHadronFraction", jpfrecochargedHadronFraction,  "recopfJetchargedHadronFraction[nrpj]/F");
00361     HltTree->Branch("recopfJetchargedEMFraction",     jpfrecochargedEMFraction,      "recopfJetchargedEMFraction[nrpj]/F");
00362     HltTree->Branch("recopfJetneutralMultiplicity",   jpfreconeutralMultiplicity,    "recopfJetneutralMultiplicity[nrpj]/I");
00363     HltTree->Branch("recopfJetchargedMultiplicity",   jpfrecochargedMultiplicity,    "recopfJetchargedMultiplicity[nrpj]/I"); 
00364     
00365 }
00366 
00367 /* **Analyze the event** */
00368 void HLTJets::analyze(edm::Event const& iEvent,
00369                       const edm::Handle<reco::CaloJetCollection>      & ohcalojets,
00370                       const edm::Handle<reco::CaloJetCollection>      & ohcalocorjets,
00371                       const edm::Handle<reco::CaloJetCollection>      & ohcalocorL1L2L3jets,
00372                       const edm::Handle<reco::CaloJetCollection>      & rcalojets,
00373                       const edm::Handle<reco::CaloJetCollection>      & rcalocorjets,
00374                       const edm::Handle<reco::GenJetCollection>       & genjets,
00375                       const edm::Handle<reco::CaloMETCollection>      & recmets,
00376                       const edm::Handle<reco::GenMETCollection>       & genmets,
00377                       const edm::Handle<reco::METCollection>          & ht,
00378                       const edm::Handle<reco::HLTTauCollection>       & taujets,
00379                       const edm::Handle<reco::PFTauCollection>        & pfTaus,
00380                       const edm::Handle<reco::PFTauCollection>        & pfTausTightCone,
00381                       const edm::Handle<reco::PFJetCollection>        & pfJets,
00382                       const edm::Handle<reco::PFTauCollection>        & recoPfTaus,  
00383                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCOnePercent,
00384                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCHalfPercent,
00385                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCQuarterPercent,
00386                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByTanCTenthPercent,                      
00387                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrByIsolation,
00388                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrAgainstElec,
00389                       const edm::Handle<reco::PFTauDiscriminator>             & theRecoPFTauDiscrAgainstMuon,
00390                       const edm::Handle<reco::PFJetCollection>        & recoPFJets,
00391                       const edm::Handle<CaloTowerCollection>          & caloTowers,
00392                       const edm::Handle<reco::PFMETCollection>        & pfmets, 
00393                       double thresholdForSavingTowers, 
00394                       double                minPtCH,
00395                       double               minPtGamma,
00396                       TTree * HltTree) {
00397     
00398     if (_Debug) std::cout << " Beginning HLTJets " << std::endl;
00399     
00400     //initialize branch variables
00401     nhjetcal=0; nhcorjetcal=0;nhcorL1L2L3jetcal=0; njetgen=0;ntowcal=0;
00402     mcalmet=0.; mcalphi=0.;
00403     mgenmet=0.; mgenphi=0.;
00404     htcalet=0.,htcalphi=0.,htcalsum=0.;
00405 
00406     noRecoPFTausSignal = 0; noRecoPFTausIso =0;
00407     noHLTPFTausSignal = 0; noHLTPFTausIso = 0;
00408 
00409 
00410 
00411     if (rcalojets.isValid()) {
00412       reco::CaloJetCollection mycalojets;
00413       mycalojets=*rcalojets;
00414       std::sort(mycalojets.begin(),mycalojets.end(),PtGreater());
00415       typedef reco::CaloJetCollection::const_iterator cjiter;
00416       int jrcal=0;
00417       for ( cjiter i=mycalojets.begin(); i!=mycalojets.end(); i++) {
00418     
00419         if (i->pt()>_CalJetMin){
00420           jrcalpt[jrcal] = i->pt();
00421           jrcalphi[jrcal] = i->phi();
00422           jrcaleta[jrcal] = i->eta();
00423           jrcale[jrcal] = i->energy();
00424           jrcalemf[jrcal] = i->emEnergyFraction();
00425           jrcaln90[jrcal] = i->n90();
00426           jetID->calculate( iEvent, *i );
00427           jrcaln90hits[jrcal] = jetID->n90Hits();
00428           jrcal++;
00429         }
00430       }
00431       nrjetcal = jrcal;
00432     }
00433     else {nrjetcal = 0;}
00434     
00435     if (rcalocorjets.isValid()) {
00436       reco::CaloJetCollection mycalojets;
00437       mycalojets=*rcalocorjets;
00438       std::sort(mycalojets.begin(),mycalojets.end(),PtGreater());
00439       typedef reco::CaloJetCollection::const_iterator cjiter;
00440       int jrcal=0;
00441       for ( cjiter i=mycalojets.begin(); i!=mycalojets.end(); i++) {
00442     
00443         if (i->pt()>_CalJetMin){
00444           jrcorcalpt[jrcal] = i->pt();
00445           jrcorcalphi[jrcal] = i->phi();
00446           jrcorcaleta[jrcal] = i->eta();
00447           jrcorcale[jrcal] = i->energy();
00448           jrcorcalemf[jrcal] = i->emEnergyFraction();
00449           jrcorcaln90[jrcal] = i->n90();
00450           jetID->calculate( iEvent, *i );
00451           jrcorcaln90hits[jrcal] = jetID->n90Hits();
00452           jrcal++;
00453         }
00454       }
00455       nrcorjetcal = jrcal;
00456     }
00457     else {nrcorjetcal = 0;}
00458     
00459     if (ohcalojets.isValid()) {
00460         reco::CaloJetCollection mycalojets;
00461         mycalojets=*ohcalojets;
00462         std::sort(mycalojets.begin(),mycalojets.end(),PtGreater());
00463         typedef reco::CaloJetCollection::const_iterator cjiter;
00464         int jhcal=0;
00465         for ( cjiter i=mycalojets.begin(); i!=mycalojets.end(); i++) {
00466             
00467             if (i->pt()>_CalJetMin){
00468                 jhcalpt[jhcal] = i->pt();
00469                 jhcalphi[jhcal] = i->phi();
00470                 jhcaleta[jhcal] = i->eta();
00471                 jhcale[jhcal] = i->energy();
00472                 jhcalemf[jhcal] = i->emEnergyFraction();
00473                 jhcaln90[jhcal] = i->n90();
00474                 jetID->calculate( iEvent, *i );
00475                 jhcaln90hits[jhcal] = jetID->n90Hits();
00476                 jhcal++;
00477             }
00478             
00479         }
00480         nhjetcal = jhcal;
00481     }
00482     else {nhjetcal = 0;}
00483     
00484     if (ohcalocorjets.isValid()) {
00485         reco::CaloJetCollection mycalocorjets;
00486         mycalocorjets=*ohcalocorjets;
00487         std::sort(mycalocorjets.begin(),mycalocorjets.end(),PtGreater());
00488         typedef reco::CaloJetCollection::const_iterator ccorjiter;
00489         int jhcorcal=0;
00490         for ( ccorjiter i=mycalocorjets.begin(); i!=mycalocorjets.end(); i++) {
00491             
00492             if (i->pt()>_CalJetMin){
00493                 jhcorcalpt[jhcorcal] = i->pt();
00494                 jhcorcalphi[jhcorcal] = i->phi();
00495                 jhcorcaleta[jhcorcal] = i->eta();
00496                 jhcorcale[jhcorcal] = i->energy();
00497                 jhcorcalemf[jhcorcal] = i->emEnergyFraction();
00498                 jhcorcaln90[jhcorcal] = i->n90();
00499                 jetID->calculate( iEvent, *i );
00500                 jhcorcaln90hits[jhcorcal] = jetID->n90Hits();
00501                 jhcorcal++;
00502             }
00503             
00504         }
00505         nhcorjetcal = jhcorcal;
00506     }
00507     else {nhcorjetcal = 0;}
00508     
00509  if (ohcalocorL1L2L3jets.isValid()) {
00510         reco::CaloJetCollection mycalocorL1L2L3jets;
00511         mycalocorL1L2L3jets=*ohcalocorL1L2L3jets;
00512         std::sort(mycalocorL1L2L3jets.begin(),mycalocorL1L2L3jets.end(),PtGreater());
00513         typedef reco::CaloJetCollection::const_iterator ccorL1L2L3jiter;
00514         int jhcorL1L2L3cal=0;
00515         for ( ccorL1L2L3jiter i=mycalocorL1L2L3jets.begin(); i!=mycalocorL1L2L3jets.end(); i++) {
00516             
00517             if (i->pt()>_CalJetMin){
00518                 jhcorL1L2L3calpt[jhcorL1L2L3cal] = i->pt();
00519                 jhcorL1L2L3calphi[jhcorL1L2L3cal] = i->phi();
00520                 jhcorL1L2L3caleta[jhcorL1L2L3cal] = i->eta();
00521                 jhcorL1L2L3cale[jhcorL1L2L3cal] = i->energy();
00522                 jhcorL1L2L3calemf[jhcorL1L2L3cal] = i->emEnergyFraction();
00523                 jhcorL1L2L3caln90[jhcorL1L2L3cal] = i->n90();
00524                 jetID->calculate( iEvent, *i );
00525                 jhcorL1L2L3caln90hits[jhcorL1L2L3cal] = jetID->n90Hits();
00526                 jhcorL1L2L3cal++;
00527             }
00528             
00529         }
00530         nhcorL1L2L3jetcal = jhcorL1L2L3cal;
00531     }
00532     else {nhcorL1L2L3jetcal = 0;}
00533 
00534     if (caloTowers.isValid()) {
00535         //    ntowcal = caloTowers->size();
00536         int jtow = 0;
00537         for ( CaloTowerCollection::const_iterator tower=caloTowers->begin(); tower!=caloTowers->end(); tower++) {
00538             if(tower->energy() > thresholdForSavingTowers)
00539             {
00540                 towet[jtow] = tower->et();
00541                 toweta[jtow] = tower->eta();
00542                 towphi[jtow] = tower->phi();
00543                 towen[jtow] = tower->energy();
00544                 towem[jtow] = tower->emEnergy();
00545                 towhd[jtow] = tower->hadEnergy();
00546                 towoe[jtow] = tower->outerEnergy();
00547                 jtow++;
00548             }
00549         }
00550         ntowcal = jtow;
00551     }
00552     else {ntowcal = 0;}
00553     
00554     if (recmets.isValid()) {
00555         typedef reco::CaloMETCollection::const_iterator cmiter;
00556         for ( cmiter i=recmets->begin(); i!=recmets->end(); i++) {
00557             mcalmet = i->pt();
00558             mcalphi = i->phi();
00559             mcalsum = i->sumEt();
00560         }
00561     }
00562 
00563     if (pfmets.isValid()) {
00564       typedef reco::PFMETCollection::const_iterator pfmetiter;
00565       for( pfmetiter i=pfmets->begin(); i!=pfmets->end(); i++) {
00566         pfmet = i->pt();
00567         pfsumet = i->sumEt();
00568         pfmetphi = i->phi();
00569       }
00570     }
00571     
00572     if (ht.isValid()) {
00573         typedef reco::METCollection::const_iterator iter;
00574         for ( iter i=ht->begin(); i!=ht->end(); i++) {
00575             htcalet = i->pt();
00576             htcalphi = i->phi();
00577             htcalsum = i->sumEt();
00578         }
00579     }
00580     
00581     if (_Monte){
00582         
00583         if (genjets.isValid()) {
00584             reco::GenJetCollection mygenjets;
00585             mygenjets=*genjets;
00586             std::sort(mygenjets.begin(),mygenjets.end(),PtGreater());
00587             typedef reco::GenJetCollection::const_iterator gjiter;
00588             int jgen=0;
00589             for ( gjiter i=mygenjets.begin(); i!=mygenjets.end(); i++) {
00590                 
00591                 if (i->pt()>_GenJetMin){
00592                     jgenpt[jgen] = i->pt();
00593                     jgenphi[jgen] = i->phi();
00594                     jgeneta[jgen] = i->eta();
00595                     jgene[jgen] = i->energy();
00596                     jgen++;
00597                 }
00598                 
00599             }
00600             njetgen = jgen;
00601         }
00602         else {njetgen = 0;}
00603         
00604         if (genmets.isValid()) {
00605             typedef reco::GenMETCollection::const_iterator gmiter;
00606             for ( gmiter i=genmets->begin(); i!=genmets->end(); i++) {
00607                 mgenmet = i->pt();
00608                 mgenphi = i->phi();
00609                 mgensum = i->sumEt();
00610             }
00611         }
00612         
00613     }
00614     
00615     
00617     if (taujets.isValid()) {      
00618         nohtau = taujets->size();
00619         reco::HLTTauCollection mytaujets;
00620         mytaujets=*taujets;
00621         std::sort(mytaujets.begin(),mytaujets.end(),GetPtGreater());
00622         typedef reco::HLTTauCollection::const_iterator tauit;
00623         int itau=0;
00624         for(tauit i=mytaujets.begin(); i!=mytaujets.end(); i++){
00625             //Ask for Eta,Phi and Et of the tau:
00626             tauEta[itau] = i->getEta();
00627             tauPhi[itau] = i->getPhi();
00628             tauPt[itau] = i->getPt();
00629             //Ask for L2 EMIsolation cut: Nominal cut : < 5
00630             l2tauemiso[itau] = i->getEMIsolationValue();
00631             //Get L25 LeadTrackPt : Nominal cut : > 20 GeV
00632             l25tauPt[itau] = i->getL25LeadTrackPtValue();
00633             //Get TrackIsolation response (returns 0 = failed or 1= passed)
00634             l3tautckiso[itau] = i->getL3TrackIsolationResponse();
00635             //MET : > 65
00636             itau++;
00637         }      
00638     }
00639     else {nohtau = 0;}
00640 
00641     
00643     if(pfTaus.isValid()) {
00644         //float minTrkPt = minPtCH;
00645         //float minGammaPt = minPtGamma;
00646         nohPFTau  = pfTaus->size();
00647         reco::PFTauCollection taus = *pfTaus;
00648         std::sort(taus.begin(),taus.end(),GetPFPtGreater());
00649         typedef reco::PFTauCollection::const_iterator pftauit;
00650         int ipftau=0;
00651         for(pftauit i=taus.begin(); i!=taus.end(); i++){
00652             //Ask for Eta,Phi and Et of the tau:
00653             ohpfTauProngs[ipftau] = i->signalPFChargedHadrCands().size();
00654             ohpfTauEta[ipftau] = i->eta();
00655             ohpfTauPhi[ipftau] = i->phi();
00656             ohpfTauPt[ipftau] = i->pt();
00657             ohpfTauJetPt[ipftau] = i->pfTauTagInfoRef()->pfjetRef()->pt();
00658 
00659             
00660   /*
00661             if( (i->leadPFCand()).isNonnull())
00662                 pfTauLeadPionPt[ipftau] = i->leadPFCand()->pt();            
00663 */
00664             if( (i->leadPFNeutralCand()).isNonnull())
00665                 ohpfTauLeadPionPt[ipftau] = i->leadPFNeutralCand()->pt();        
00666             else 
00667                 ohpfTauLeadPionPt[ipftau] = -999.0;
00668 
00669             if((i->leadPFChargedHadrCand()).isNonnull())
00670                 ohpfTauLeadTrackPt[ipftau] = i->leadPFChargedHadrCand()->pt();
00671             else
00672                 ohpfTauLeadTrackPt[ipftau] = -999.0;
00673   
00674 
00675             float maxPtTrkIso = 0;
00676             for (unsigned int iTrk = 0; iTrk < i->isolationPFChargedHadrCands().size(); iTrk++)
00677             {
00678                 if(i->isolationPFChargedHadrCands()[iTrk]->pt() > maxPtTrkIso) maxPtTrkIso = i->isolationPFChargedHadrCands()[iTrk]->pt();
00679 
00680                 if (i->isolationPFChargedHadrCands()[iTrk]->trackRef().isNonnull()){
00681                   hltpftauIsoTrToPFTauMatch[noHLTPFTausIso]=ipftau;
00682                   HLTPFTauIsoTrDz[noHLTPFTausIso]=i->isolationPFChargedHadrCands()[iTrk]->trackRef()->dz(); // dz wrt (0,0,0), to compare offline with HLT 
00683                   HLTPFTauIsoTrPt[noHLTPFTausIso]=i->isolationPFChargedHadrCands()[iTrk]->pt();
00684                   /*
00685                   std::cout << "Adding isocand for hltpftau " << ipftau
00686                       << " pt " << HLTPFTauIsoTrPt[noHLTPFTausIso]
00687                       << " dz " << HLTPFTauIsoTrDz[noHLTPFTausIso]
00688                       << std::endl; // */
00689                   ++noHLTPFTausIso;
00690                 }
00691 
00692             }
00693                 
00694             ohpfTauTrkIso[ipftau] = maxPtTrkIso;
00695             float maxPtGammaIso = 0;
00696             for (unsigned int iGamma = 0; iGamma < i->isolationPFGammaCands().size(); iGamma++)
00697             {
00698                 if(i->isolationPFGammaCands()[iGamma]->pt() > maxPtGammaIso) maxPtGammaIso = i->isolationPFGammaCands()[iGamma]->pt();
00699             }                        
00700 
00701 
00702 
00703             for (unsigned int iTrk = 0; iTrk < i->signalPFChargedHadrCands().size(); iTrk++)
00704             {
00705               if (i->signalPFChargedHadrCands ()[iTrk]->trackRef().isNonnull()){
00706                 hltpftauSignalTrToPFTauMatch[noHLTPFTausSignal]=ipftau;
00707                 HLTPFTauSignalTrDz[noHLTPFTausSignal]=i->signalPFChargedHadrCands()[iTrk]->trackRef()->dz(); // dz wrt (0,0,0), to compare offline with HLT
00708                 HLTPFTauSignalTrPt[noHLTPFTausSignal]=i->signalPFChargedHadrCands()[iTrk]->pt();
00709                 /*
00710                   std::cout << "Adding sigcand for hltpftau " << ipftau
00711                       << " pt " << HLTPFTauSignalTrPt[noHLTPFTausSignal]
00712                       << " dz " << HLTPFTauSignalTrDz[noHLTPFTausSignal]
00713                       << std::endl; // */
00714                 ++noHLTPFTausSignal;
00715               }
00716             }
00717 
00718 
00719 
00720 
00721 
00722             ohpfTauGammaIso[ipftau] = maxPtGammaIso;
00723             ipftau++;
00724         } 
00725       
00726     }
00727 
00728     if(pfTausTightCone.isValid()) {
00729         //float minTrkPt = minPtCH;
00730         //float minGammaPt = minPtGamma;
00731         nohPFTauTightCone = pfTaus->size();
00732         reco::PFTauCollection taus = *pfTausTightCone;
00733         std::sort(taus.begin(),taus.end(),GetPFPtGreater());
00734         typedef reco::PFTauCollection::const_iterator pftauit;
00735         int ipftau=0;
00736         for(pftauit i=taus.begin(); i!=taus.end(); i++){
00737             //Ask for Eta,Phi and Et of the tau:
00738             ohpfTauTightConeProngs[ipftau] = i->signalPFChargedHadrCands().size();
00739             ohpfTauTightConeEta[ipftau] = i->eta();
00740             ohpfTauTightConePhi[ipftau] = i->phi();
00741             ohpfTauTightConePt[ipftau] = i->pt();
00742             ohpfTauTightConeJetPt[ipftau] = i->pfTauTagInfoRef()->pfjetRef()->pt();
00743 
00744     
00745             if( (i->leadPFNeutralCand()).isNonnull())
00746                 ohpfTauTightConeLeadPionPt[ipftau] = i->leadPFNeutralCand()->pt();
00747             else 
00748               ohpfTauTightConeLeadPionPt[ipftau] = -999.0;
00749 
00750 
00751             if((i->leadPFChargedHadrCand()).isNonnull())
00752                 ohpfTauTightConeLeadTrackPt[ipftau] = i->leadPFChargedHadrCand()->pt();
00753             else
00754               ohpfTauTightConeLeadTrackPt[ipftau] = -999.0;
00755 
00756             float maxPtTrkIso = 0;
00757             for (unsigned int iTrk = 0; iTrk < i->isolationPFChargedHadrCands().size(); iTrk++)
00758             {
00759                 if(i->isolationPFChargedHadrCands()[iTrk]->pt() > maxPtTrkIso) maxPtTrkIso = i->isolationPFChargedHadrCands()[iTrk]->pt();
00760             }
00761 
00762             ohpfTauTightConeTrkIso[ipftau] = maxPtTrkIso;
00763             float maxPtGammaIso = 0;
00764             for (unsigned int iGamma = 0; iGamma < i->isolationPFGammaCands().size(); iGamma++)
00765             {
00766                 if(i->isolationPFGammaCands()[iGamma]->pt() > maxPtGammaIso) maxPtGammaIso = i->isolationPFGammaCands()[iGamma]->pt();
00767             }
00768             ohpfTauTightConeGammaIso[ipftau] = maxPtGammaIso;
00769             ipftau++;
00770         }
00771 
00772     }
00773     
00775       
00776     if(recoPfTaus.isValid()) {
00777         float minTrkPt = minPtCH;
00778         float minGammaPt = minPtGamma;
00779         nRecoPFTau  = recoPfTaus->size();
00780         reco::PFTauCollection taus = *recoPfTaus;
00781 
00782         // disable sorting for proper access to discriminators
00783         //std::sort(taus.begin(),taus.end(),GetPFPtGreater());
00784         typedef reco::PFTauCollection::const_iterator pftauit;
00785         int ipftau=0;
00786         
00787         for(pftauit i=taus.begin(); i!=taus.end(); i++){
00788             //Ask for Eta,Phi and Et of the tau:
00789             recopfTauEta[ipftau] = i->eta();
00790             recopfTauPhi[ipftau] = i->phi();
00791             recopfTauPt[ipftau]  = i->pt();
00792 
00793             if( (i->leadPFNeutralCand()).isNonnull())
00794                 recopfTauLeadPionPt[ipftau] = i->leadPFNeutralCand()->pt();
00795             else
00796               recopfTauLeadPionPt[ipftau] = -999.0;  
00797 
00798 
00799             if((i->leadPFChargedHadrCand()).isNonnull())
00800                 recopfTauLeadTrackPt[ipftau] = i->leadPFChargedHadrCand()->pt();
00801             else
00802               recopfTauLeadTrackPt[ipftau]  = -999.0;
00803 
00804             int myTrks=0;
00805             for (unsigned int iTrk = 0; iTrk < i->isolationPFChargedHadrCands().size(); iTrk++)
00806             {
00807                 if(i->isolationPFChargedHadrCands()[iTrk]->pt() > minTrkPt) myTrks++;
00808                 if (i->isolationPFChargedHadrCands()[iTrk]->trackRef().isNonnull()){
00809                   isoTrToPFTauMatch[noRecoPFTausIso]=ipftau;
00810                   recoPFTauIsoTrDz[noRecoPFTausIso]=i->isolationPFChargedHadrCands()[iTrk]->trackRef()->dz(); // dz wrt (0,0,0), to compare offline with HLT
00811                   recoPFTauIsoTrPt[noRecoPFTausIso]=i->isolationPFChargedHadrCands()[iTrk]->pt();
00812                   /*
00813                   std::cout << "Adding isocand for tau " << ipftau
00814                             << " pt " << recoPFTauIsoTrPt[noRecoPFTausIso]
00815                             << " dz " << recoPFTauIsoTrDz[noRecoPFTausIso]
00816                             << std::endl;// */
00817                   ++noRecoPFTausIso;
00818                 }
00819 
00820             }
00821                
00822             recopfTauTrkIso[ipftau] = myTrks;
00823             int myGammas=0;
00824             for (unsigned int iGamma = 0; iGamma < i->isolationPFGammaCands().size(); iGamma++)
00825             {
00826                 if(i->isolationPFGammaCands()[iGamma]->pt() > minGammaPt) myGammas++;
00827             }                        
00828             recopfTauGammaIso[ipftau] = myGammas;
00829             
00830 
00831             for (unsigned int iTrk = 0; iTrk < i->signalPFChargedHadrCands().size(); iTrk++)
00832             {
00833                 if (i->signalPFChargedHadrCands ()[iTrk]->trackRef().isNonnull()){
00834                   signalTrToPFTauMatch[noRecoPFTausSignal]=ipftau;
00835                   recoPFTauSignalTrDz[noRecoPFTausSignal]=i->signalPFChargedHadrCands()[iTrk]->trackRef()->dz(); // dz wrt (0,0,0), to compare offline with HLT
00836                   recoPFTauSignalTrPt[noRecoPFTausSignal]=i->signalPFChargedHadrCands()[iTrk]->pt();
00837                   /*
00838                   std::cout << "Adding sigcand for tau " << ipftau
00839                             << " pt " << recoPFTauSignalTrPt[noRecoPFTausSignal]
00840                             << " dz " << recoPFTauSignalTrDz[noRecoPFTausSignal]
00841                             << std::endl;// */
00842                   ++noRecoPFTausSignal;
00843                 }
00844             }
00845 
00846             const reco::PFTauRef thisTauRef(recoPfTaus,ipftau);
00847             
00848             if(theRecoPFTauDiscrByTanCOnePercent.isValid()){
00849             recopfTauDiscrByTancOnePercent[ipftau] = (*theRecoPFTauDiscrByTanCOnePercent)[thisTauRef];}
00850             if(theRecoPFTauDiscrByIsolation.isValid()){ 
00851             recopfTauDiscrByIso[ipftau] = (*theRecoPFTauDiscrByIsolation)[thisTauRef];} 
00852             if(theRecoPFTauDiscrAgainstMuon.isValid()){
00853             recopfTauDiscrAgainstMuon[ipftau] = (*theRecoPFTauDiscrAgainstMuon)[thisTauRef];}
00854             if(theRecoPFTauDiscrAgainstElec.isValid()){
00855             recopfTauDiscrAgainstElec[ipftau] = (*theRecoPFTauDiscrAgainstElec)[thisTauRef];}
00856             if(theRecoPFTauDiscrByTanCHalfPercent.isValid()){
00857             recopfTauDiscrByTancHalfPercent[ipftau] = (*theRecoPFTauDiscrByTanCHalfPercent)[thisTauRef];}
00858             if(theRecoPFTauDiscrByTanCQuarterPercent.isValid()){
00859             recopfTauDiscrByTancQuarterPercent[ipftau] = (*theRecoPFTauDiscrByTanCQuarterPercent)[thisTauRef];}
00860             if(theRecoPFTauDiscrByTanCTenthPercent.isValid()){
00861             recopfTauDiscrByTancTenthPercent[ipftau] = (*theRecoPFTauDiscrByTanCTenthPercent)[thisTauRef];}
00862 
00863             ipftau++;
00864         }        
00865     }
00866    
00868     if(pfJets.isValid()) {
00869         nohPFJet  = pfJets->size();
00870         reco::PFJetCollection Jets = *pfJets;
00871         std::sort(Jets.begin(),Jets.end(),GetPFPtGreater());
00872         typedef reco::PFJetCollection::const_iterator pfJetit;
00873         int ipfJet=0;
00874         float pfMHTx = 0.;
00875         float pfMHTy = 0.;
00876         pfHT         = 0.;
00877 
00878         for(pfJetit i=Jets.begin(); i!=Jets.end(); i++){
00879             //Ask for Eta,Phi and Et of the Jet:
00880             pfJetEta[ipfJet] = i->eta();
00881             pfJetPhi[ipfJet] = i->phi();
00882             pfJetPt[ipfJet] = i->pt();           
00883             pfJetE[ipfJet] = i->energy();
00884             
00885             if (i->pt() > 40. && abs(i->eta())<3.0)
00886               pfHT  += i -> pt();
00887             if (i->pt() > 30.){
00888               pfMHTx = pfMHTx + i->px();
00889               pfMHTy = pfMHTy + i->py();
00890             }
00891             ipfJet++;   
00892         } 
00893         pfMHT = sqrt(pfMHTx*pfMHTx + pfMHTy*pfMHTy);
00894         
00895     }
00897     nrpj = 0;
00898     if(recoPFJets.isValid()){
00899             nrpj = recoPFJets->size();
00900             reco::PFJetCollection Jets = *recoPFJets;
00901             std::sort(Jets.begin(),Jets.end(),GetPFPtGreater());
00902             typedef reco::PFJetCollection::const_iterator pfJetit;
00903             int ipfJet=0;
00904             for(pfJetit i=Jets.begin(); i!=Jets.end(); i++){
00905                     //Ask for Eta,Phi and Et of the Jet:
00906                     jpfrecoeta[ipfJet] = i->eta();
00907                     jpfrecophi[ipfJet] = i->phi();
00908                     jpfrecopt[ipfJet] = i->pt();           
00909                     jpfrecoe[ipfJet] = i->energy();           
00910                     jpfreconeutralHadronFraction[ipfJet] = i->neutralHadronEnergyFraction ();
00911                     jpfrecochargedHadronFraction[ipfJet] = i->chargedHadronEnergyFraction ();
00912                     jpfreconeutralMultiplicity[ipfJet] = i->neutralMultiplicity ();
00913                     jpfrecochargedMultiplicity[ipfJet] = i->chargedMultiplicity ();
00914                     jpfreconeutralEMFraction[ipfJet] = i->neutralEmEnergyFraction ();
00915                     jpfrecochargedEMFraction[ipfJet] = i->chargedEmEnergyFraction ();
00916 
00917                     ipfJet++;   
00918             } 
00919 
00920     }
00921     
00922 }