00001 #include "HLTriggerOffline/Tau/interface/L2TauAnalyzer.h"
00002 #include "Math/GenVector/VectorUtil.h"
00003 #include <iostream>
00004 #include <iomanip>
00005 #include <fstream>
00006
00007 L2TauAnalyzer::L2TauAnalyzer(const edm::ParameterSet& iConfig):
00008 l2TauInfoAssoc_(iConfig.getParameter<edm::InputTag>("L2InfoAssociationInput")),
00009 l1Taus_(iConfig.getParameter<edm::InputTag>("L1TauCollection")),
00010 l1Jets_(iConfig.getParameter<edm::InputTag>("L1JetCollection")),
00011 rootFile_(iConfig.getParameter<std::string>("outputFileName")),
00012 IsSignal_(iConfig.getParameter<bool>("IsSignal")),
00013 mcColl_(iConfig.getParameter<edm::InputTag>("MatchedCollection"))
00014 {
00015
00016 l2file = new TFile(rootFile_.c_str(),"recreate");
00017
00018 l2tree = new TTree("l2tree","Level 2 Tau Tree");
00019
00020
00021
00022 ecalIsol_Et=0.;
00023 towerIsol_Et=0.;
00024 cl_etaRMS=0.;
00025 cl_phiRMS=0.;
00026 cl_drRMS=0.;
00027 MCeta=0.;
00028 MCet=0.;
00029 cl_Nclusters=0;
00030 seedTowerEt = 0.;
00031 JetEt=0.;
00032 JetEta=0.;
00033 L1et=0.;
00034 L1eta=0.;
00035
00036
00037 l2tree->Branch("ecalIsolEt",&ecalIsol_Et,"ecalIsolEt/F");
00038 l2tree->Branch("towerIsolEt",&towerIsol_Et,"towerIsolEt/F");
00039 l2tree->Branch("clEtaRMS",&cl_etaRMS,"clEtaRMS/F");
00040 l2tree->Branch("clPhiRMS",&cl_phiRMS,"clPhiRMS/F");
00041 l2tree->Branch("clDrRMS",&cl_drRMS,"clDrRMS/F");
00042 l2tree->Branch("mcEta",&MCeta,"mcEta/F");
00043 l2tree->Branch("mcEt",&MCet,"mcEt/F");
00044 l2tree->Branch("clNclusters",&cl_Nclusters,"clNclusters/I");
00045 l2tree->Branch("seedTowerEt",&seedTowerEt,"seedTowerEt/F");
00046 l2tree->Branch("jetEt",&JetEt,"jetEt/F");
00047 l2tree->Branch("jetEta",&JetEta,"jetEta/F");
00048 l2tree->Branch("L1Et",&L1et,"L1Et/F");
00049 l2tree->Branch("L1Eta",&L1eta,"L1Eta/F");
00050
00051 }
00052
00053
00054 L2TauAnalyzer::~L2TauAnalyzer()
00055 {
00056 }
00057
00058
00059
00060 void
00061 L2TauAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00062 {
00063 using namespace edm;
00064 using namespace reco;
00065
00066 Handle<L2TauInfoAssociation> l2TauInfoAssoc;
00067 Handle<LVColl> McInfo;
00068 Handle<l1extra::L1JetParticleCollection> L1Taus;
00069 Handle<l1extra::L1JetParticleCollection> L1Jets;
00070
00071 if(iEvent.getByLabel(l2TauInfoAssoc_,l2TauInfoAssoc))
00072 {
00073 if(l2TauInfoAssoc->size()>0)
00074 for(L2TauInfoAssociation::const_iterator p = l2TauInfoAssoc->begin();p!=l2TauInfoAssoc->end();++p)
00075 {
00076 const L2TauIsolationInfo l2info = p->val;
00077 const CaloJet& jet =*(p->key);
00078
00079 MatchElementL2 mcMatch;
00080 mcMatch.matched=false;
00081 mcMatch.mcEt=0;
00082 mcMatch.mcEta=0;
00083 mcMatch.deltar=0;
00084
00085 if(IsSignal_)
00086 {
00087 if(iEvent.getByLabel(mcColl_,McInfo))
00088 mcMatch=match(jet,*McInfo);
00089 }
00090
00091 if((mcMatch.matched&&IsSignal_)||(!IsSignal_))
00092 {
00093
00094 ecalIsol_Et=l2info.ECALIsolConeCut;
00095 towerIsol_Et=l2info.TowerIsolConeCut;
00096 cl_Nclusters=l2info.ECALClusterNClusters;
00097 cl_etaRMS=l2info.ECALClusterEtaRMS;
00098 cl_phiRMS=l2info.ECALClusterPhiRMS;
00099 cl_drRMS=l2info.ECALClusterDRRMS;
00100 seedTowerEt = l2info.SeedTowerEt;
00101 MCeta =mcMatch.mcEta;
00102 MCet=mcMatch.mcEt;
00103 JetEt = jet.et();
00104 JetEta = jet.eta();
00105
00106
00107 L1et=0;
00108 L1eta=0;
00109 if(iEvent.getByLabel(l1Taus_,L1Taus))
00110 {
00111 MatchElementL2 l1Match;
00112 l1Match.matched=false;
00113 l1Match.mcEt=0;
00114 l1Match.mcEta=0;
00115 l1Match.deltar=0;
00116 l1Match=match(jet,*L1Taus);
00117 if(l1Match.matched)
00118 {
00119 L1et=l1Match.mcEt;
00120 L1eta=l1Match.mcEta;
00121 }
00122
00123 else
00124 {
00125 if(iEvent.getByLabel(l1Jets_,L1Jets))
00126 {
00127 l1Match=match(jet,*L1Taus);
00128 if(l1Match.matched)
00129 {
00130 L1et=l1Match.mcEt;
00131 L1eta=l1Match.mcEta;
00132 }
00133
00134 }
00135 }
00136
00137 }
00138
00139 l2tree->Fill();
00140 }
00141
00142 }
00143 }
00144 }
00145
00146
00147
00148 void
00149 L2TauAnalyzer::beginJob(const edm::EventSetup&)
00150 {
00151
00152 }
00153
00154
00155 void
00156 L2TauAnalyzer::endJob() {
00157 l2file->Write();
00158
00159 }
00160
00161 MatchElementL2
00162 L2TauAnalyzer::match(const reco::Jet& jet,const LVColl& McInfo)
00163 {
00164
00165
00166
00167
00168 bool matched=false;
00169 double delta_min=100.;
00170 double mceta=0;
00171 double mcet=0;
00172
00173 double matchingDR=0.3;
00174
00175
00176
00177
00178 if(McInfo.size()>0)
00179 for(std::vector<LV>::const_iterator it = McInfo.begin();it!=McInfo.end();++it)
00180 {
00181 double delta = ROOT::Math::VectorUtil::DeltaR(jet.p4().Vect(),*it);
00182 if(delta<matchingDR)
00183 {
00184 matched=true;
00185 if(delta<delta_min)
00186 {
00187 delta_min=delta;
00188 mceta=it->eta();
00189 mcet=it->Et();
00190 }
00191 }
00192 }
00193
00194
00195 MatchElementL2 match;
00196 match.matched=matched;
00197 match.deltar=delta_min;
00198 match.mcEta = mceta;
00199 match.mcEt = mcet;
00200
00201
00202 return match;
00203 }
00204
00205 MatchElementL2
00206 L2TauAnalyzer::match(const reco::Jet& jet,const l1extra::L1JetParticleCollection& McInfo)
00207 {
00208
00209
00210
00211
00212 bool matched=false;
00213 double delta_min=100.;
00214 double mceta=0;
00215 double mcet=0;
00216
00217 double matchingDR=0.5;
00218
00219
00220
00221
00222 if(McInfo.size()>0)
00223 for(l1extra::L1JetParticleCollection::const_iterator it = McInfo.begin();it!=McInfo.end();++it)
00224 {
00225 double delta = ROOT::Math::VectorUtil::DeltaR(jet.p4().Vect(),it->p4().Vect());
00226 if(delta<matchingDR)
00227 {
00228 matched=true;
00229 if(delta<delta_min)
00230 {
00231 delta_min=delta;
00232 mceta=it->eta();
00233 mcet=it->et();
00234 }
00235 }
00236 }
00237
00238
00239 MatchElementL2 match;
00240 match.matched=matched;
00241 match.deltar=delta_min;
00242 match.mcEta = mceta;
00243 match.mcEt = mcet;
00244
00245
00246 return match;
00247 }
00248
00249
00250
00251
00252