CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalTrigPrimAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Class: EcalTrigPrimAnalyzer
4 //
10 //
11 //
12 // Original Author: Ursula Berthon, Stephanie Baffioni, Pascal Paganini
13 // Created: Thu Jul 4 11:38:38 CEST 2005
14 //
15 //
16 
17 
18 // system include files
19 #include <memory>
20 #include <utility>
21 
22 // user include files
26 
31 
36 
41 
43 
44 #include "EcalTrigPrimAnalyzer.h"
45 
46 #include <TMath.h>
47 
48 using namespace edm;
50 
52 
53 {
54  ecal_parts_.push_back("Barrel");
55  ecal_parts_.push_back("Endcap");
56 
57  histfile_=new TFile("histos.root","RECREATE");
58  tree_ = new TTree("TPGtree","TPGtree");
59  tree_->Branch("iphi",&iphi_,"iphi/I");
60  tree_->Branch("ieta",&ieta_,"ieta/I");
61  tree_->Branch("eRec",&eRec_,"eRec/F");
62  tree_->Branch("tpgADC",&tpgADC_,"tpgADC/I");
63  tree_->Branch("tpgGeV",&tpgGeV_,"tpgGeV/F");
64  tree_->Branch("ttf",&ttf_,"ttf/I");
65  tree_->Branch("fg",&fg_,"fg/I");
66  for (unsigned int i=0;i<2;++i) {
67  ecal_et_[i]=new TH1I(ecal_parts_[i].c_str(),"Et",255,0,255);
68  char title[30];
69  sprintf(title,"%s_ttf",ecal_parts_[i].c_str());
70  ecal_tt_[i]=new TH1I(title,"TTF",10,0,10);
71  sprintf(title,"%s_fgvb",ecal_parts_[i].c_str());
72  ecal_fgvb_[i]=new TH1I(title,"FGVB",10,0,10);
73  }
74 
75  recHits_= iConfig.getParameter<bool>("AnalyzeRecHits");
76  label_=iConfig.getParameter<edm::InputTag>("inputTP");
77  if (recHits_) {
78  hTPvsRechit_= new TH2F("TP_vs_RecHit","TP vs rechit",256,-1,255,255,0,255);
79  hTPoverRechit_= new TH1F("TP_over_RecHit","TP over rechit",500,0,4);
80  rechits_labelEB_=iConfig.getParameter<edm::InputTag>("inputRecHitsEB");
81  rechits_labelEE_=iConfig.getParameter<edm::InputTag>("inputRecHitsEE");
82  }
83 }
84 
85 
87 {
88 
89  // do anything here that needs to be done at desctruction time
90  // (e.g. close files, deallocate resources etc.)
91 
92  histfile_->Write();
93  histfile_->Close();
94 
95 }
96 
97 
98 //
99 // member functions
100 //
101 
102 // ------------ method called to analyze the data ------------
103 void
105 {
106  using namespace edm;
107  using namespace std;
108 
109  // Get input
111  iEvent.getByLabel(label_,tp);
112  for (unsigned int i=0;i<tp.product()->size();i++) {
113  EcalTriggerPrimitiveDigi d=(*(tp.product()))[i];
114  int subdet=d.id().subDet()-1;
115  if (subdet==0) {
116  ecal_et_[subdet]->Fill(d.compressedEt());
117  }
118  else {
119  if (d.id().ietaAbs()==27 || d.id().ietaAbs()==28) {
120  if (i%2) ecal_et_[subdet]->Fill(d.compressedEt()*2.);
121  }
122  else ecal_et_[subdet]->Fill(d.compressedEt());
123  }
124  ecal_tt_[subdet]->Fill(d.ttFlag());
125  ecal_fgvb_[subdet]->Fill(d.fineGrain());
126 
127  }
128  if (!recHits_) return;
129 
130  // comparison with RecHits
131  edm::Handle<EcalRecHitCollection> rechit_EB_col;
132  iEvent.getByLabel(rechits_labelEB_,rechit_EB_col);
133 
134  edm::Handle<EcalRecHitCollection> rechit_EE_col;
135  iEvent.getByLabel(rechits_labelEE_,rechit_EE_col);
136 
137 
138  edm::ESHandle<CaloGeometry> theGeometry;
139  edm::ESHandle<CaloSubdetectorGeometry> theBarrelGeometry_handle;
140  edm::ESHandle<CaloSubdetectorGeometry> theEndcapGeometry_handle;
141  iSetup.get<CaloGeometryRecord>().get( theGeometry );
142  iSetup.get<EcalEndcapGeometryRecord>().get("EcalEndcap",theEndcapGeometry_handle);
143  iSetup.get<EcalBarrelGeometryRecord>().get("EcalBarrel",theBarrelGeometry_handle);
144 
145  const CaloSubdetectorGeometry *theEndcapGeometry,*theBarrelGeometry;
146  theEndcapGeometry = &(*theEndcapGeometry_handle);
147  theBarrelGeometry = &(*theBarrelGeometry_handle);
149  iSetup.get<IdealGeometryRecord>().get(eTTmap_);
150 
151  map<EcalTrigTowerDetId, float> mapTow_Et;
152 
153 
154  for (unsigned int i=0;i<rechit_EB_col.product()->size();i++) {
155  const EBDetId & myid1=(*rechit_EB_col.product())[i].id();
156  EcalTrigTowerDetId towid1= myid1.tower();
157  float theta = theBarrelGeometry->getGeometry(myid1)->getPosition().theta();
158  float Etsum=((*rechit_EB_col.product())[i].energy())*sin(theta);
159  bool test_alreadyin= false;
160  map<EcalTrigTowerDetId, float>::iterator ittest= mapTow_Et.find(towid1);
161  if (ittest!= mapTow_Et.end()) test_alreadyin=true;
162  if (test_alreadyin) continue;
163  unsigned int j=i+1;
164  bool loopend=false;
165  unsigned int count=0;
166  while( j<rechit_EB_col.product()->size() && !loopend){
167  count++;
168  const EBDetId & myid2=(*rechit_EB_col.product())[j].id();
169  EcalTrigTowerDetId towid2= myid2.tower();
170  if( towid1==towid2 ) {
171  float theta=theBarrelGeometry->getGeometry(myid2)->getPosition().theta();
172  Etsum += (*rechit_EB_col.product())[j].energy()*sin(theta);
173  }
174  j++;
175  if (count>1800) loopend=true;
176  }
177  mapTow_Et.insert(pair<EcalTrigTowerDetId,float>(towid1, Etsum));
178  }
179 
180 
181  for (unsigned int i=0;i<rechit_EE_col.product()->size();i++) {
182  const EEDetId & myid1=(*rechit_EE_col.product())[i].id();
183  EcalTrigTowerDetId towid1= (*eTTmap_).towerOf(myid1);
184  float theta=theEndcapGeometry->getGeometry(myid1)->getPosition().theta();
185  float Etsum=(*rechit_EE_col.product())[i].energy()*sin(theta);
186  bool test_alreadyin= false;
187  map<EcalTrigTowerDetId, float>::iterator ittest= mapTow_Et.find(towid1);
188  if (ittest!= mapTow_Et.end()) test_alreadyin=true;
189  if (test_alreadyin) continue;
190  unsigned int j=i+1;
191  bool loopend=false;
192  unsigned int count=0;
193  while( j<rechit_EE_col.product()->size() && !loopend){
194  const EEDetId & myid2=(*rechit_EE_col.product())[j].id();
195  EcalTrigTowerDetId towid2= (*eTTmap_).towerOf(myid2);
196  if( towid1==towid2 ) {
197  float theta=theEndcapGeometry->getGeometry(myid2)->getPosition().theta();
198  Etsum += (*rechit_EE_col.product())[j].energy()*sin(theta);
199  }
200  // else loopend=true;
201  j++;
202  if (count>500) loopend=true;
203  }
204  // alreadyin_EE.push_back(towid1);
205  mapTow_Et.insert(pair<EcalTrigTowerDetId,float>(towid1, Etsum));
206  }
207 
208 
209  EcalTPGScale ecalScale ;
210  ecalScale.setEventSetup(iSetup) ;
211  for (unsigned int i=0;i<tp.product()->size();i++) {
212  EcalTriggerPrimitiveDigi d=(*(tp.product()))[i];
213  const EcalTrigTowerDetId TPtowid= d.id();
214  map<EcalTrigTowerDetId, float>::iterator it= mapTow_Et.find(TPtowid);
215  float Et = ecalScale.getTPGInGeV(d.compressedEt(), TPtowid) ;
216  if (d.id().ietaAbs()==27 || d.id().ietaAbs()==28) Et*=2;
217  iphi_ = TPtowid.iphi() ;
218  ieta_ = TPtowid.ieta() ;
219  tpgADC_ = d.compressedEt() ;
220  tpgGeV_ = Et ;
221  ttf_ = d.ttFlag() ;
222  fg_ = d.fineGrain() ;
223  if (it!= mapTow_Et.end()) {
224  hTPvsRechit_->Fill(it->second,Et);
225  hTPoverRechit_->Fill(Et/it->second);
226  eRec_ = it->second ;
227  }
228  tree_->Fill() ;
229  }
230 
231 }
232 
233 void
235  for (unsigned int i=0;i<2;++i) {
236  ecal_et_[i]->Write();
237  ecal_tt_[i]->Write();
238  ecal_fgvb_[i]->Write();
239  }
240  if (recHits_) {
241  hTPvsRechit_->Write();
242  hTPoverRechit_->Write();
243  }
244 }
245 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void setEventSetup(const edm::EventSetup &evtSetup)
Definition: EcalTPGScale.cc:19
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
int ieta() const
get the tower ieta
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
int compressedEt() const
get the encoded/compressed Et of interesting sample
int iEvent
Definition: GenABIO.cc:230
int ietaAbs() const
get the absolute value of the tower ieta
int j
Definition: DBlmapReader.cc:9
EcalTrigTowerDetId tower() const
get the HCAL/trigger iphi of this crystal
Definition: EBDetId.h:59
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
EcalTrigPrimAnalyzer(const edm::ParameterSet &)
int iphi() const
get the tower iphi
const EcalTrigTowerDetId & id() const
const T & get() const
Definition: EventSetup.h:55
EcalSubdetector subDet() const
get the subDetector associated to the Trigger Tower
bool fineGrain() const
get the fine-grain bit of interesting sample
int ttFlag() const
get the Trigger tower Flag of interesting sample