CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  //rootfile_ = config.getUntrackedParameter<std::string>("rootfile","EcalTBValidation.root");
34 
35  // verbosity...
36  verbose_ = config.getUntrackedParameter<bool>("verbose", false);
37 
39  if( dbe_ ) {
40  if( verbose_ ) {
41  dbe_->setVerbose(1);
42  dbe_->showDirStructure();
43  }
44  else {
45  dbe_->setVerbose(0);
46  }
47  }
48 
49  meETBxib_ = 0;
50  meETBampltdc_ = 0;
51  meETBShape_ = 0;
52  meETBhodoX_ = 0;
53  meETBhodoY_ = 0;
54  meETBe1x1_ = 0;
55  meETBe3x3_ = 0;
56  meETBe5x5_ = 0;
57  meETBe1e9_ = 0;
58  meETBe1e25_ = 0;
59  meETBe9e25_ = 0;
63  meETBe1vsX_ = 0;
64  meETBe1vsY_ = 0;
65  meETBe1e9vsX_ = 0;
66  meETBe1e9vsY_ = 0;
67  meETBe1e25vsX_ = 0;
68  meETBe1e25vsY_ = 0;
69  meETBe9e25vsX_ = 0;
70  meETBe9e25vsY_ = 0;
71 
72  if( dbe_ ) {
73 
74  std::string hname;
75  dbe_->setCurrentFolder( "EcalRecHitsV/EcalTBValidationTask" );
76 
77  hname = "xtal in beam position";
78  meETBxib_ = dbe_->book2D( hname, hname, 85, 0., 85., 20,0., 20. );
79  hname = "Max Amplitude vs TDC offset";
80  meETBampltdc_ = dbe_->book2D( hname, hname, 100, 0., 1., 1000, 0., 4000. );
81  hname = "Beam Profile X";
82  meETBhodoX_ = dbe_->book1D( hname, hname, 100, -20., 20. );
83  hname = "Beam Profile Y";
84  meETBhodoY_ = dbe_->book1D( hname, hname, 100, -20., 20. );
85  hname = "E1x1 energy";
86  meETBe1x1_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
87  hname = "E3x3 energy";
88  meETBe3x3_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
89  hname = "E5x5 energy";
90  meETBe5x5_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
91  hname = "E1x1 energy center";
92  meETBe1x1_center_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
93  hname = "E3x3 energy center";
94  meETBe3x3_center_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
95  hname = "E5x5 energy center";
96  meETBe5x5_center_ = dbe_->book1D( hname, hname, 1000, 0., 4000. );
97  hname = "E1 over E9 ratio";
98  meETBe1e9_ = dbe_->book1D( hname, hname, 600, 0., 1.2 );
99  hname = "E1 over E25 ratio";
100  meETBe1e25_ = dbe_->book1D( hname, hname, 600, 0., 1.2 );
101  hname = "E9 over E25 ratio";
102  meETBe9e25_ = dbe_->book1D( hname, hname, 600, 0., 1.2 );
103  hname = "E1 vs X";
104  meETBe1vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 1000, 0., 4000. );
105  hname = "E1 vs Y";
106  meETBe1vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 1000, 0., 4000. );
107  hname = "E1 over E9 vs X";
108  meETBe1e9vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
109  hname = "E1 over E9 vs Y";
110  meETBe1e9vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
111  hname = "E1 over E25 vs X";
112  meETBe1e25vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
113  hname = "E1 over E25 vs Y";
114  meETBe1e25vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
115  hname = "E9 over E25 vs X";
116  meETBe9e25vsX_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
117  hname = "E9 over E25 vs Y";
118  meETBe9e25vsY_ = dbe_->book2D( hname, hname, 80, -20, 20, 600, 0., 1.2 );
119  hname = "Xtal in Beam Shape";
120  meETBShape_ = dbe_->book2D( hname, hname, 250, 0, 10, 350, 0, 3500 );
121  }
122 
123 }
124 
125 
127 
129 
131 
133 
134  using namespace edm;
135  using namespace cms;
136 
137  // digis
138  const EBDigiCollection* theDigis=0;
140  event.getByLabel(digiProducer_, digiCollection_, pdigis);
141  if(pdigis.isValid()){
142  theDigis = pdigis.product();
143  }
144  else {
145  std::cerr << "Error! can't get the product " << digiCollection_.c_str() << std::endl;
146  return;
147  }
148 
149  // rechits
150  const EBUncalibratedRecHitCollection* theHits=0;
152  event.getByLabel(hitProducer_, hitCollection_, phits);
153  if(phits.isValid()){
154  theHits = phits.product();
155  }
156  else {
157  std::cerr << "Error! can't get the product " << hitCollection_.c_str() << std::endl;
158  return;
159  }
160 
161  // hodoscopes
162  const EcalTBHodoscopeRecInfo* theHodo=0;
164  event.getByLabel(hodoRecInfoProducer_, hodoRecInfoCollection_, pHodo);
165  if(pHodo.isValid()){
166  theHodo = pHodo.product();
167  }
168  else{
169  std::cerr << "Error! can't get the product " << hodoRecInfoCollection_.c_str() << std::endl;
170  return;
171  }
172 
173  // tdc
174  const EcalTBTDCRecInfo* theTDC=0;
176  event.getByLabel(tdcRecInfoProducer_, tdcRecInfoCollection_, pTDC);
177  if(pTDC.isValid()){
178  theTDC = pTDC.product();
179  }
180  else{
181  std::cerr << "Error! can't get the product " << tdcRecInfoCollection_.c_str() << std::endl;
182  return;
183  }
184 
185  // event header
186  const EcalTBEventHeader* evtHeader=0;
187  Handle<EcalTBEventHeader> pEventHeader;
188  event.getByLabel(eventHeaderProducer_ , pEventHeader);
189  if(pEventHeader.isValid()){
190  evtHeader = pEventHeader.product();
191  }
192  else{
193  std::cerr << "Error! can't get the product " << eventHeaderProducer_.c_str() << std::endl;
194  return;
195  }
196 
197 
198  // -----------------------------------------------------------------------
199  // xtal-in-beam
200  EBDetId xtalInBeamId(1,xtalInBeam_, EBDetId::SMCRYSTALMODE);
201  if (xtalInBeamId==EBDetId(0)){ return; }
202  int xibEta = xtalInBeamId.ieta();
203  int xibPhi = xtalInBeamId.iphi();
204 
205  // skipping events with moving table (in data)
206  if (data_ && (evtHeader->tableIsMoving())) return;
207 
208  // amplitudes
209  EBDetId Xtals5x5[25];
210  for (unsigned int icry=0;icry<25;icry++){
211  unsigned int row = icry/5;
212  unsigned int column = icry%5;
213  int ieta = xtalInBeamId.ieta()+column-2;
214  int iphi = xtalInBeamId.iphi()+row-2;
215  if(EBDetId::validDetId(ieta, iphi)){
216  EBDetId tempId(ieta, iphi,EBDetId::ETAPHIMODE);
217  if (tempId.ism()==1)
218  Xtals5x5[icry] = tempId;
219  else
220  Xtals5x5[icry] = EBDetId(0);
221  } else {
222  Xtals5x5[icry] = EBDetId(0);
223  }
224  }
225 
226  // matrices
227  double ampl1x1 = 0.;
228  double ampl3x3 = 0.;
229  double ampl5x5 = 0.;
230  for (unsigned int icry=0;icry<25;icry++) {
231  if (!Xtals5x5[icry].null()){
232  double theAmpl = (theHits->find(Xtals5x5[icry]))->amplitude();
233  ampl5x5 += theAmpl;
234  if (icry==12){ampl1x1 = theAmpl;}
235  if (icry==6 || icry==7 || icry==8 || icry==11 || icry==12 || icry==13 || icry==16 || icry==17 || icry==18){ampl3x3 += theAmpl;}
236  }}
237 
238 
239  // pulse shape
240  double sampleSave[10];
241  for(int ii=0; ii < 10; ++ii){ sampleSave[ii] = 0.0; }
242  EBDigiCollection::const_iterator thisDigi = theDigis->find(xtalInBeamId);
243  // int sMax = -1; // UNUSED
244  double eMax = 0.;
245  if (thisDigi != theDigis->end()){
246  EBDataFrame myDigi = (*thisDigi);
247  for (int sample=0; sample < myDigi.size(); ++sample){
248  double analogSample = myDigi.sample(sample).adc();
249  sampleSave[sample] = analogSample;
250  if ( eMax < analogSample ) {
251  eMax = analogSample;
252  // sMax = sample; // UNUSED
253  }
254  }
255  }
256 
257  // beam profile
258  double xBeam = theHodo->posX();
259  double yBeam = theHodo->posY();
260 
261  // filling histos
262  meETBxib_ -> Fill(xibEta, xibPhi);
263  meETBhodoX_ -> Fill(xBeam);
264  meETBhodoY_ -> Fill(yBeam);
265  meETBampltdc_ -> Fill(theTDC->offset(),ampl1x1);
266  meETBe1x1_ -> Fill(ampl1x1);
267  meETBe3x3_ -> Fill(ampl3x3);
268  meETBe5x5_ -> Fill(ampl5x5);
269  meETBe1e9_ -> Fill(ampl1x1/ampl3x3);
270  meETBe1e25_ -> Fill(ampl1x1/ampl5x5);
271  meETBe9e25_ -> Fill(ampl3x3/ampl5x5);
272  meETBe1vsX_ -> Fill(xBeam,ampl1x1);
273  meETBe1vsY_ -> Fill(yBeam,ampl1x1);
274  meETBe1e9vsX_ -> Fill(xBeam,ampl1x1/ampl3x3);
275  meETBe1e9vsY_ -> Fill(yBeam,ampl1x1/ampl3x3);
276  meETBe1e25vsX_ -> Fill(xBeam,ampl1x1/ampl5x5);
277  meETBe1e25vsY_ -> Fill(yBeam,ampl1x1/ampl5x5);
278  meETBe9e25vsX_ -> Fill(xBeam,ampl3x3/ampl5x5);
279  meETBe9e25vsY_ -> Fill(yBeam,ampl3x3/ampl5x5);
280 
281  for(int ii=0; ii < 10; ++ii){ meETBShape_->Fill(double(ii)+theTDC->offset(),sampleSave[ii]); }
282 
283  if ( (fabs(xBeam)<2.5) && (fabs(yBeam)<2.5) ){
284  meETBe1x1_center_ -> Fill(ampl1x1);
285  meETBe3x3_center_ -> Fill(ampl3x3);
286  meETBe5x5_center_ -> Fill(ampl5x5);
287  }
288 }
289 
290 
MonitorElement * meETBShape_
static bool validDetId(int i, int j)
check if a valid index combination
Definition: EBDetId.cc:59
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_
MonitorElement * meETBe9e25vsX_
EcalMGPASample sample(int i) const
Definition: EcalDataFrame.h:30
MonitorElement * meETBe1e25vsX_
std::string eventHeaderProducer_
int ism() const
get the ECAL/SM id
Definition: EBDetId.cc:79
int size() const
Definition: EcalDataFrame.h:27
MonitorElement * meETBxib_
void Fill(long long x)
int iphi() const
get the crystal iphi
Definition: EBDetId.h:46
EcalTBValidation(const edm::ParameterSet &)
MonitorElement * meETBe9e25_
MonitorElement * meETBe9e25vsY_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
std::string eventHeaderCollection_
MonitorElement * meETBe1vsX_
std::string hitProducer_
std::string hodoRecInfoProducer_
std::string digiProducer_
int ieta() const
get the crystal ieta
Definition: EBDetId.h:44
MonitorElement * meETBe5x5_center_
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
static const int ETAPHIMODE
Definition: EBDetId.h:145
std::string hitCollection_
MonitorElement * meETBe1e9vsX_
T const * product() const
Definition: Handle.h:74
MonitorElement * meETBe1e25vsY_
MonitorElement * meETBe1x1_
std::string tdcRecInfoCollection_
std::string tdcRecInfoProducer_
MonitorElement * meETBe1e9_
virtual void beginJob()
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * meETBe3x3_center_
static const int SMCRYSTALMODE
Definition: EBDetId.h:146
MonitorElement * meETBe1e9vsY_
std::string digiCollection_
MonitorElement * meETBhodoY_
MonitorElement * meETBe1vsY_
std::string hodoRecInfoCollection_
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
MonitorElement * meETBe1x1_center_
MonitorElement * meETBhodoX_
int adc() const
get the ADC sample (12 bits)
virtual void endJob()