File indexing completed on 2024-09-07 04:38:06
0001
0002 #include <memory>
0003 #include <utility>
0004
0005
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008
0009 #include "Geometry/CaloTopology/interface/EcalTrigTowerConstituentsMap.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013
0014 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
0015
0016 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
0017 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0018
0019 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0020
0021 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
0022 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0023
0024
0025
0026 #include "EcalEBTrigPrimAnalyzer.h"
0027
0028 #include <TMath.h>
0029
0030 using namespace edm;
0031 class CaloSubdetectorGeometry;
0032
0033 EcalEBTrigPrimAnalyzer::EcalEBTrigPrimAnalyzer(const edm::ParameterSet& iConfig)
0034 : tokens_(consumesCollector())
0035
0036 {
0037 ecal_parts_.push_back("Barrel");
0038
0039 histfile_ = new TFile("histos.root", "RECREATE");
0040 tree_ = new TTree("TPGtree", "TPGtree");
0041 tree_->Branch("tpIphi", &tpIphi_, "tpIphi/I");
0042 tree_->Branch("tpIeta", &tpIeta_, "tpIeta/I");
0043 tree_->Branch("rhIphi", &rhIphi_, "rhIphi/I");
0044 tree_->Branch("rhIeta", &rhIeta_, "rhIeta/I");
0045 tree_->Branch("eRec", &eRec_, "eRec/F");
0046 tree_->Branch("tpgADC", &tpgADC_, "tpgADC/I");
0047 tree_->Branch("tpgGeV", &tpgGeV_, "tpgGeV/F");
0048 tree_->Branch("ttf", &ttf_, "ttf/I");
0049 tree_->Branch("fg", &fg_, "fg/I");
0050 for (unsigned int i = 0; i < ecal_parts_.size(); ++i) {
0051 char title[30];
0052 sprintf(title, "%s_Et", ecal_parts_[i].c_str());
0053 ecal_et_[i] = new TH1I(title, "Et", 255, 0, 255);
0054 sprintf(title, "%s_ttf", ecal_parts_[i].c_str());
0055 ecal_tt_[i] = new TH1I(title, "TTF", 10, 0, 10);
0056 sprintf(title, "%s_fgvb", ecal_parts_[i].c_str());
0057 ecal_fgvb_[i] = new TH1I(title, "FGVB", 10, 0, 10);
0058 }
0059 eTTmapToken_ = esConsumes<EcalTrigTowerConstituentsMap, IdealGeometryRecord>();
0060 recHits_ = iConfig.getParameter<bool>("AnalyzeRecHits");
0061 debug_ = iConfig.getParameter<bool>("Debug");
0062 rechits_labelEB_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("inputRecHitsEB"));
0063 primToken_ = consumes<EcalEBTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("inputTP"));
0064 barrelGeomToken_ = esConsumes<CaloSubdetectorGeometry, EcalBarrelGeometryRecord>(edm::ESInputTag("", "EcalBarrel"));
0065 tokenEBdigi_ = consumes<EBDigiCollection>(iConfig.getParameter<edm::InputTag>("barrelEcalDigis"));
0066 nEvents_ = 0;
0067
0068 hTPvsTow_eta_ = new TH2F("TP_vs_Tow_eta", "TP vs Tow eta ; #eta(tow); #eta(tp)", 50, -2.5, 2.5, 50, -2.5, 2.5);
0069 hAllTPperEvt_ = new TH1F("AllTPperEvt", "TP per Event; N_{TP}; ", 100, 0., 20000.);
0070 hTPperEvt_ = new TH1F("TPperEvt", "N_{TP} per Event; N_{TP}; ", 100, 0., 500.);
0071 hTP_iphiVsieta_ = new TH2F("TP_iphiVsieta", "TP i#phi vs i#eta ; i#eta(tp); i#phi(tp)", 10, 70, 80, 10, 340, 350);
0072 hTP_iphiVsieta_fullrange_ =
0073 new TH2F("TP_iphiVsieta_fullrange", "TP i#phi vs i#eta ; i#eta(tp); i#phi(tp)", 200, -100, 100, 350, 0, 350);
0074
0075 if (recHits_) {
0076 hTPvsTow_ieta_ =
0077 new TH2F("TP_vs_Tow_ieta", "TP vs Tow ieta ; i#eta(tow); i#eta(tp)", 200, -100, 100, 200, -100, 100);
0078
0079 hTPvsRechit_ = new TH2F("TP_vs_RecHit", "TP vs rechit Et;E_{T}(rh) (GeV);E_{T}(tp) (GeV)", 100, 0, 50, 100, 0, 50);
0080 hDeltaEt_ = new TH1F("DeltaEt", "[Et(rh)-Et(TP)]/Et(rh); [E_{T}(rh)-E_{T}(tp)]/E_{T}(rh); Counts", 200, -1, 1);
0081 hTPoverRechit_ = new TH1F("TP_over_RecHit", "Et(TP/rechit); E_{T}(tp)/E_{T}(rh); Counts", 200, 0, 2);
0082 hRechitEt_ = new TH1F("RecHitEt", "E_{T};E_{T}(rh) (GeV);Counts", 100, 0, 50);
0083 hTPEt_ = new TH1F("TPEt", "E_{T}{tp);E_{T}(rh) (GeV);Count", 100, 0, 50);
0084 hRatioEt_ = new TH1F("RatioTPoverRH", "Et", 100, 0, 50);
0085 hAllRechitEt_ = new TH1F("AllRecHit", "Et", 100, 0, 50);
0086
0087 hRH_iphiVsieta_ = new TH2F("RH_iphiVsieta", "RH i#phi vs i#eta ; i#eta(rh); i#phi(rh)", 10, 70, 80, 10, 340, 350);
0088 hRH_iphiVsieta_fullrange_ =
0089 new TH2F("RH_iphiVsieta_fullrange", "RH i#phi vs i#eta ; i#eta(rh); i#phi(rh)", 200, -100, 100, 350, 0, 350);
0090 }
0091 }
0092
0093 EcalEBTrigPrimAnalyzer::~EcalEBTrigPrimAnalyzer() {
0094
0095
0096
0097 histfile_->Write();
0098 histfile_->Close();
0099 }
0100
0101 void EcalEBTrigPrimAnalyzer::init(const edm::EventSetup& iSetup) { eTTmap_ = iSetup.getHandle(eTTmapToken_); }
0102
0103
0104
0105
0106
0107
0108 void EcalEBTrigPrimAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0109 using namespace edm;
0110 using namespace std;
0111 nEvents_++;
0112
0113 if (nEvents_ == 1)
0114 this->init(iSetup);
0115
0116
0117 edm::Handle<EcalEBTrigPrimDigiCollection> tp;
0118 iEvent.getByToken(primToken_, tp);
0119
0120
0121
0122
0123
0124
0125
0126
0127 for (unsigned int i = 0; i < tp.product()->size(); i++) {
0128 EcalEBTriggerPrimitiveDigi d = (*(tp.product()))[i];
0129 int subdet = 0;
0130 if (subdet == 0) {
0131 ecal_et_[subdet]->Fill(d.encodedEt());
0132 }
0133 }
0134
0135
0136
0137 edm::Handle<EcalRecHitCollection> rechit_EB_col;
0138 if (recHits_) {
0139
0140 iEvent.getByToken(rechits_labelEB_, rechit_EB_col);
0141 }
0142
0143 edm::ESHandle<CaloSubdetectorGeometry> theBarrelGeometry_handle = iSetup.getHandle(barrelGeomToken_);
0144 const CaloSubdetectorGeometry* theBarrelGeometry = theBarrelGeometry_handle.product();
0145
0146 if (debug_) {
0147 std::cout << " TP analyzer =================> Treating event " << iEvent.id() << " Number of TPs "
0148 << tp.product()->size() << std::endl;
0149 if (recHits_)
0150 std::cout << " Number of EB rechits " << rechit_EB_col.product()->size() << std::endl;
0151 }
0152 hAllTPperEvt_->Fill(float(tp.product()->size()));
0153
0154
0155
0156
0157 EcalTPGScale ecalScale(tokens_, iSetup);
0158
0159
0160
0161
0162
0163 int nTP = 0;
0164 for (unsigned int i = 0; i < tp.product()->size(); i++) {
0165 EcalEBTriggerPrimitiveDigi d = (*(tp.product()))[i];
0166 const EBDetId TPid = d.id();
0167
0168
0169
0170
0171
0172
0173
0174
0175 float Et = d.encodedEt() / 8.;
0176 if (Et <= 5)
0177 continue;
0178
0179 nTP++;
0180
0181 std::cout << " TP digi size " << d.size() << std::endl;
0182 for (int iBx = 0; iBx < d.size(); iBx++) {
0183 std::cout << " TP samples " << d.sample(iBx) << std::endl;
0184 }
0185
0186
0187
0188
0189
0190
0191
0192 tpIphi_ = TPid.iphi();
0193 tpIeta_ = TPid.ieta();
0194 tpgADC_ = d.encodedEt();
0195 tpgGeV_ = Et;
0196
0197 hTP_iphiVsieta_->Fill(TPid.ieta(), TPid.iphi(), Et);
0198 hTP_iphiVsieta_fullrange_->Fill(TPid.ieta(), TPid.iphi(), Et);
0199
0200 if (recHits_) {
0201 for (unsigned int j = 0; j < rechit_EB_col.product()->size(); j++) {
0202 const EBDetId& myid1 = (*rechit_EB_col.product())[j].id();
0203 float theta = theBarrelGeometry->getGeometry(myid1)->getPosition().theta();
0204 float rhEt = ((*rechit_EB_col.product())[j].energy()) * sin(theta);
0205 if (myid1 == TPid) {
0206 if (debug_)
0207 std::cout << " Analyzer same cristal " << myid1 << " " << TPid << std::endl;
0208
0209
0210
0211
0212
0213
0214
0215
0216
0217
0218
0219 rhIphi_ = myid1.iphi();
0220 rhIeta_ = myid1.ieta();
0221 hRH_iphiVsieta_->Fill(myid1.ieta(), myid1.iphi(), rhEt);
0222 hRH_iphiVsieta_fullrange_->Fill(myid1.ieta(), myid1.iphi(), rhEt);
0223
0224 hTPvsRechit_->Fill(rhEt, Et);
0225 hTPoverRechit_->Fill(Et / rhEt);
0226 hDeltaEt_->Fill((rhEt - Et) / rhEt);
0227 if (debug_)
0228 std::cout << " TP compressed et " << d.encodedEt() << " Et in GeV " << Et << " RH Et " << rhEt
0229 << " Et/rhEt " << Et / rhEt << std::endl;
0230 hRechitEt_->Fill(rhEt);
0231 hTPEt_->Fill(Et);
0232 if (rhEt < 1000000)
0233 eRec_ = rhEt;
0234 tree_->Fill();
0235 }
0236
0237 }
0238 }
0239
0240 }
0241
0242
0243
0244 hTPperEvt_->Fill(float(nTP));
0245
0246 if (recHits_) {
0247 hRatioEt_->Divide(hTPEt_, hRechitEt_);
0248 for (unsigned int j = 0; j < rechit_EB_col.product()->size(); j++) {
0249 const EBDetId& myid1 = (*rechit_EB_col.product())[j].id();
0250 float theta = theBarrelGeometry->getGeometry(myid1)->getPosition().theta();
0251 float rhEt = ((*rechit_EB_col.product())[j].energy()) * sin(theta);
0252 if (rhEt > 0)
0253 hAllRechitEt_->Fill(rhEt);
0254 }
0255 }
0256 }
0257
0258 void EcalEBTrigPrimAnalyzer::endJob() {
0259 for (unsigned int i = 0; i < ecal_parts_.size(); ++i) {
0260 ecal_et_[i]->Write();
0261 ecal_tt_[i]->Write();
0262 ecal_fgvb_[i]->Write();
0263 }
0264
0265 hAllTPperEvt_->Write();
0266 hTPperEvt_->Write();
0267
0268 if (recHits_) {
0269 hTPvsTow_ieta_->Write();
0270 hTPvsTow_eta_->Write();
0271 hTPvsRechit_->Write();
0272 hTPoverRechit_->Write();
0273 hAllRechitEt_->Write();
0274 hRechitEt_->Write();
0275 hDeltaEt_->Write();
0276 hTPEt_->Write();
0277 hRatioEt_->Write();
0278 hTP_iphiVsieta_->Write();
0279 hRH_iphiVsieta_->Write();
0280 hTP_iphiVsieta_fullrange_->Write();
0281 hRH_iphiVsieta_fullrange_->Write();
0282 }
0283 }