CMS 3D CMS Logo

EcalTBValidation.cc
Go to the documentation of this file.
1 // Ecal H4 tesbeam analysis
11 
12 
13 #include <iostream>
14 #include <string>
15 
16 
18 
19  data_ = config.getUntrackedParameter<int>("data",-1000);
20  xtalInBeam_ = config.getUntrackedParameter<int>("xtalInBeam",-1000);
21 
22  digiCollection_ = config.getParameter<std::string>("digiCollection");
23  digiProducer_ = config.getParameter<std::string>("digiProducer");
24  hitCollection_ = config.getParameter<std::string>("hitCollection");
25  hitProducer_ = config.getParameter<std::string>("hitProducer");
26  hodoRecInfoCollection_ = config.getParameter<std::string>("hodoRecInfoCollection");
27  hodoRecInfoProducer_ = config.getParameter<std::string>("hodoRecInfoProducer");
28  tdcRecInfoCollection_ = config.getParameter<std::string>("tdcRecInfoCollection");
29  tdcRecInfoProducer_ = config.getParameter<std::string>("tdcRecInfoProducer");
30  eventHeaderCollection_ = config.getParameter<std::string>("eventHeaderCollection");
31  eventHeaderProducer_ = config.getParameter<std::string>("eventHeaderProducer");
32 
33  digi_Token_ = consumes<EBDigiCollection>(edm::InputTag(digiProducer_, digiCollection_));
34  hit_Token_ = consumes<EBUncalibratedRecHitCollection>(edm::InputTag(hitProducer_, hitCollection_));
35  hodoRec_Token_ = consumes<EcalTBHodoscopeRecInfo>(edm::InputTag(hodoRecInfoProducer_, hodoRecInfoCollection_));
36  tdcRec_Token_ = consumes<EcalTBTDCRecInfo>(edm::InputTag(tdcRecInfoProducer_, tdcRecInfoCollection_));
37  eventHeader_Token_ = consumes<EcalTBEventHeader>(edm::InputTag(eventHeaderProducer_));
38  //rootfile_ = config.getUntrackedParameter<std::string>("rootfile","EcalTBValidation.root");
39 
40  // verbosity...
41  verbose_ = config.getUntrackedParameter<bool>("verbose", false);
42 
43  meETBxib_ = nullptr;
44  meETBampltdc_ = nullptr;
45  meETBShape_ = nullptr;
46  meETBhodoX_ = nullptr;
47  meETBhodoY_ = nullptr;
48  meETBe1x1_ = nullptr;
49  meETBe3x3_ = nullptr;
50  meETBe5x5_ = nullptr;
51  meETBe1e9_ = nullptr;
52  meETBe1e25_ = nullptr;
53  meETBe9e25_ = nullptr;
54  meETBe1x1_center_ = nullptr;
55  meETBe3x3_center_ = nullptr;
56  meETBe5x5_center_ = nullptr;
57  meETBe1vsX_ = nullptr;
58  meETBe1vsY_ = nullptr;
59  meETBe1e9vsX_ = nullptr;
60  meETBe1e9vsY_ = nullptr;
61  meETBe1e25vsX_ = nullptr;
62  meETBe1e25vsY_ = nullptr;
63  meETBe9e25vsX_ = nullptr;
64  meETBe9e25vsY_ = nullptr;
65 
66 }
67 
68 
70 
72 
73  std::string hname;
74  ibooker.setCurrentFolder( "EcalRecHitsV/EcalTBValidationTask" );
75 
76  hname = "xtal in beam position";
77  meETBxib_ = ibooker.book2D( hname, hname, 85, 0., 85., 20,0., 20. );
78  hname = "Max Amplitude vs TDC offset";
79  meETBampltdc_ = ibooker.book2D( hname, hname, 100, 0., 1., 1000, 0., 4000. );
80  hname = "Beam Profile X";
81  meETBhodoX_ = ibooker.book1D( hname, hname, 100, -20., 20. );
82  hname = "Beam Profile Y";
83  meETBhodoY_ = ibooker.book1D( hname, hname, 100, -20., 20. );
84  hname = "E1x1 energy";
85  meETBe1x1_ = ibooker.book1D( hname, hname, 1000, 0., 4000. );
86  hname = "E3x3 energy";
87  meETBe3x3_ = ibooker.book1D( hname, hname, 1000, 0., 4000. );
88  hname = "E5x5 energy";
89  meETBe5x5_ = ibooker.book1D( hname, hname, 1000, 0., 4000. );
90  hname = "E1x1 energy center";
91  meETBe1x1_center_ = ibooker.book1D( hname, hname, 1000, 0., 4000. );
92  hname = "E3x3 energy center";
93  meETBe3x3_center_ = ibooker.book1D( hname, hname, 1000, 0., 4000. );
94  hname = "E5x5 energy center";
95  meETBe5x5_center_ = ibooker.book1D( hname, hname, 1000, 0., 4000. );
96  hname = "E1 over E9 ratio";
97  meETBe1e9_ = ibooker.book1D( hname, hname, 600, 0., 1.2 );
98  hname = "E1 over E25 ratio";
99  meETBe1e25_ = ibooker.book1D( hname, hname, 600, 0., 1.2 );
100  hname = "E9 over E25 ratio";
101  meETBe9e25_ = ibooker.book1D( hname, hname, 600, 0., 1.2 );
102  hname = "E1 vs X";
103  meETBe1vsX_ = ibooker.book2D( hname, hname, 80, -20, 20, 1000, 0., 4000. );
104  hname = "E1 vs Y";
105  meETBe1vsY_ = ibooker.book2D( hname, hname, 80, -20, 20, 1000, 0., 4000. );
106  hname = "E1 over E9 vs X";
107  meETBe1e9vsX_ = ibooker.book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
108  hname = "E1 over E9 vs Y";
109  meETBe1e9vsY_ = ibooker.book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
110  hname = "E1 over E25 vs X";
111  meETBe1e25vsX_ = ibooker.book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
112  hname = "E1 over E25 vs Y";
113  meETBe1e25vsY_ = ibooker.book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
114  hname = "E9 over E25 vs X";
115  meETBe9e25vsX_ = ibooker.book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
116  hname = "E9 over E25 vs Y";
117  meETBe9e25vsY_ = ibooker.book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
118  hname = "Xtal in Beam Shape";
119  meETBShape_ = ibooker.book2D( hname, hname, 250, 0, 10, 350, 0, 3500 );
120 
121 }
122 
124 
125  using namespace edm;
126  using namespace cms;
127 
128  // digis
129  const EBDigiCollection* theDigis=nullptr;
131  event.getByToken(digi_Token_, pdigis);
132  if(pdigis.isValid()){
133  theDigis = pdigis.product();
134  }
135  else {
136  std::cerr << "Error! can't get the product " << digiCollection_.c_str() << std::endl;
137  return;
138  }
139 
140  // rechits
141  const EBUncalibratedRecHitCollection* theHits=nullptr;
143  event.getByToken(hit_Token_, phits);
144  if(phits.isValid()){
145  theHits = phits.product();
146  }
147  else {
148  std::cerr << "Error! can't get the product " << hitCollection_.c_str() << std::endl;
149  return;
150  }
151 
152  // hodoscopes
153  const EcalTBHodoscopeRecInfo* theHodo=nullptr;
155  event.getByToken(hodoRec_Token_, pHodo);
156  if(pHodo.isValid()){
157  theHodo = pHodo.product();
158  }
159  else{
160  std::cerr << "Error! can't get the product " << hodoRecInfoCollection_.c_str() << std::endl;
161  return;
162  }
163 
164  // tdc
165  const EcalTBTDCRecInfo* theTDC=nullptr;
167  event.getByToken(tdcRec_Token_, pTDC);
168  if(pTDC.isValid()){
169  theTDC = pTDC.product();
170  }
171  else{
172  std::cerr << "Error! can't get the product " << tdcRecInfoCollection_.c_str() << std::endl;
173  return;
174  }
175 
176  // event header
177  const EcalTBEventHeader* evtHeader=nullptr;
178  Handle<EcalTBEventHeader> pEventHeader;
179  event.getByToken(eventHeader_Token_ , pEventHeader);
180  if(pEventHeader.isValid()){
181  evtHeader = pEventHeader.product();
182  }
183  else{
184  std::cerr << "Error! can't get the product " << eventHeaderProducer_.c_str() << std::endl;
185  return;
186  }
187 
188 
189  // -----------------------------------------------------------------------
190  // xtal-in-beam
191  EBDetId xtalInBeamId(1,xtalInBeam_, EBDetId::SMCRYSTALMODE);
192  if (xtalInBeamId==EBDetId(0)){ return; }
193  int xibEta = xtalInBeamId.ieta();
194  int xibPhi = xtalInBeamId.iphi();
195 
196  // skipping events with moving table (in data)
197  if (data_ && (evtHeader->tableIsMoving())) return;
198 
199  // amplitudes
200  EBDetId Xtals5x5[25];
201  for (unsigned int icry=0;icry<25;icry++){
202  unsigned int row = icry/5;
203  unsigned int column = icry%5;
204  int ieta = xtalInBeamId.ieta()+column-2;
205  int iphi = xtalInBeamId.iphi()+row-2;
206  if(EBDetId::validDetId(ieta, iphi)){
207  EBDetId tempId(ieta, iphi,EBDetId::ETAPHIMODE);
208  if (tempId.ism()==1)
209  Xtals5x5[icry] = tempId;
210  else
211  Xtals5x5[icry] = EBDetId(0);
212  } else {
213  Xtals5x5[icry] = EBDetId(0);
214  }
215  }
216 
217  // matrices
218  double ampl1x1 = 0.;
219  double ampl3x3 = 0.;
220  double ampl5x5 = 0.;
221  for (unsigned int icry=0;icry<25;icry++) {
222  if (!Xtals5x5[icry].null()){
223  double theAmpl = (theHits->find(Xtals5x5[icry]))->amplitude();
224  ampl5x5 += theAmpl;
225  if (icry==12){ampl1x1 = theAmpl;}
226  if (icry==6 || icry==7 || icry==8 || icry==11 || icry==12 || icry==13 || icry==16 || icry==17 || icry==18){ampl3x3 += theAmpl;}
227  }}
228 
229 
230  // pulse shape
231  double sampleSave[10];
232  for(int ii=0; ii < 10; ++ii){ sampleSave[ii] = 0.0; }
233  EBDigiCollection::const_iterator thisDigi = theDigis->find(xtalInBeamId);
234  // int sMax = -1; // UNUSED
235  double eMax = 0.;
236  if (thisDigi != theDigis->end()){
237  EBDataFrame myDigi = (*thisDigi);
238  for (int sample=0; sample < myDigi.size(); ++sample){
239  double analogSample = myDigi.sample(sample).adc();
240  sampleSave[sample] = analogSample;
241  if ( eMax < analogSample ) {
242  eMax = analogSample;
243  // sMax = sample; // UNUSED
244  }
245  }
246  }
247 
248  // beam profile
249  double xBeam = theHodo->posX();
250  double yBeam = theHodo->posY();
251 
252  // filling histos
253  meETBxib_ -> Fill(xibEta, xibPhi);
254  meETBhodoX_ -> Fill(xBeam);
255  meETBhodoY_ -> Fill(yBeam);
256  meETBampltdc_ -> Fill(theTDC->offset(),ampl1x1);
257  meETBe1x1_ -> Fill(ampl1x1);
258  meETBe3x3_ -> Fill(ampl3x3);
259  meETBe5x5_ -> Fill(ampl5x5);
260  meETBe1e9_ -> Fill(ampl1x1/ampl3x3);
261  meETBe1e25_ -> Fill(ampl1x1/ampl5x5);
262  meETBe9e25_ -> Fill(ampl3x3/ampl5x5);
263  meETBe1vsX_ -> Fill(xBeam,ampl1x1);
264  meETBe1vsY_ -> Fill(yBeam,ampl1x1);
265  meETBe1e9vsX_ -> Fill(xBeam,ampl1x1/ampl3x3);
266  meETBe1e9vsY_ -> Fill(yBeam,ampl1x1/ampl3x3);
267  meETBe1e25vsX_ -> Fill(xBeam,ampl1x1/ampl5x5);
268  meETBe1e25vsY_ -> Fill(yBeam,ampl1x1/ampl5x5);
269  meETBe9e25vsX_ -> Fill(xBeam,ampl3x3/ampl5x5);
270  meETBe9e25vsY_ -> Fill(yBeam,ampl3x3/ampl5x5);
271 
272  for(int ii=0; ii < 10; ++ii){ meETBShape_->Fill(double(ii)+theTDC->offset(),sampleSave[ii]); }
273 
274  if ( (fabs(xBeam)<2.5) && (fabs(yBeam)<2.5) ){
275  meETBe1x1_center_ -> Fill(ampl1x1);
276  meETBe3x3_center_ -> Fill(ampl3x3);
277  meETBe5x5_center_ -> Fill(ampl5x5);
278  }
279 }
280 
281 
MonitorElement * meETBShape_
T getParameter(std::string const &) const
MonitorElement * meETBampltdc_
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * meETBe5x5_
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
MonitorElement * meETBe1e25_
MonitorElement * meETBe3x3_
bool tableIsMoving() const
Tell if the table is Moving.
edm::EDGetTokenT< EcalTBTDCRecInfo > tdcRec_Token_
edm::EDGetTokenT< EBDigiCollection > digi_Token_
MonitorElement * meETBe9e25vsX_
edm::EDGetTokenT< EcalTBHodoscopeRecInfo > hodoRec_Token_
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
EcalMGPASample sample(int i) const
Definition: EcalDataFrame.h:29
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
MonitorElement * meETBe1e25vsX_
std::string eventHeaderProducer_
int ism() const
get the ECAL/SM id
Definition: EBDetId.h:61
int size() const
Definition: EcalDataFrame.h:26
static bool validDetId(int i, int j)
check if a valid index combination
Definition: EBDetId.h:124
MonitorElement * meETBxib_
void Fill(long long x)
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
EcalTBValidation(const edm::ParameterSet &)
edm::EDGetTokenT< EBUncalibratedRecHitCollection > hit_Token_
const_iterator find(id_type i) const
MonitorElement * meETBe9e25_
MonitorElement * meETBe9e25vsY_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
std::string eventHeaderCollection_
MonitorElement * meETBe1vsX_
std::string hitProducer_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
std::string hodoRecInfoProducer_
std::string digiProducer_
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
MonitorElement * meETBe5x5_center_
bool isValid() const
Definition: HandleBase.h:74
static const int ETAPHIMODE
Definition: EBDetId.h:166
ii
Definition: cuy.py:589
Namespace of DDCMS conversion namespace.
std::string hitCollection_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
T const * product() const
Definition: Handle.h:81
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
MonitorElement * meETBe1e9vsX_
MonitorElement * meETBe1e25vsY_
MonitorElement * meETBe1x1_
const_iterator end() const
std::string tdcRecInfoCollection_
~EcalTBValidation() override
iterator find(key_type k)
HLT enums.
std::string tdcRecInfoProducer_
edm::EDGetTokenT< EcalTBEventHeader > eventHeader_Token_
MonitorElement * meETBe1e9_
MonitorElement * meETBe3x3_center_
static const int SMCRYSTALMODE
Definition: EBDetId.h:167
MonitorElement * meETBe1e9vsY_
std::string digiCollection_
MonitorElement * meETBhodoY_
MonitorElement * meETBe1vsY_
std::string hodoRecInfoCollection_
MonitorElement * meETBe1x1_center_
float offset() const
Definition: event.py:1
Definition: Run.h:44
MonitorElement * meETBhodoX_
void analyze(const edm::Event &, const edm::EventSetup &) override
int adc() const
get the ADC sample (12 bits)