CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalABAnalyzer.cc
Go to the documentation of this file.
1 /*
2  * \class EcalABAnalyzer
3  *
4  * $Date: 2012/02/09 10:07:36 $
5  * primary author: Julie Malcles - CEA/Saclay
6  * author: Gautier Hamel De Monchenault - CEA/Saclay
7  */
8 
9 #include <TAxis.h>
10 #include <TH1.h>
11 #include <TProfile.h>
12 #include <TTree.h>
13 #include <TChain.h>
14 #include <TFile.h>
15 #include <TMath.h>
16 
18 
19 #include <sstream>
20 #include <iostream>
21 #include <fstream>
22 #include <iomanip>
23 #include <vector>
24 
27 
30 
33 
41 
49 
50 using namespace std;
51 
52 
53 //========================================================================
55 //========================================================================
56  :
57 iEvent(0),
58 
59 // framework parameters with default values
60 _nsamples( iConfig.getUntrackedParameter< unsigned int >( "nSamples", 10 ) ),
61 _presample( iConfig.getUntrackedParameter< unsigned int >( "nPresamples", 2 ) ),
62 _firstsample( iConfig.getUntrackedParameter< unsigned int >( "firstSample", 1 ) ),
63 _lastsample( iConfig.getUntrackedParameter< unsigned int >( "lastSample", 2 ) ),
64 _timingcutlow( iConfig.getUntrackedParameter< unsigned int >( "timingCutLow", 2 ) ),
65 _timingcuthigh( iConfig.getUntrackedParameter< unsigned int >( "timingCutHigh", 9 ) ),
66 _timingquallow( iConfig.getUntrackedParameter< unsigned int >( "timingQualLow", 3 ) ),
67 _timingqualhigh( iConfig.getUntrackedParameter< unsigned int >( "timingQualHigh", 8 ) ),
68 _ratiomincutlow( iConfig.getUntrackedParameter< double >( "ratioMinCutLow", 0.4 ) ),
69 _ratiomincuthigh( iConfig.getUntrackedParameter< double >( "ratioMinCutHigh", 0.95 ) ),
70 _ratiomaxcutlow( iConfig.getUntrackedParameter< double >( "ratioMaxCutLow", 0.8 ) ),
71 _presamplecut( iConfig.getUntrackedParameter< double >( "presampleCut", 5.0 ) ),
72 _niter( iConfig.getUntrackedParameter< unsigned int >( "nIter", 3 ) ),
73 _alpha( iConfig.getUntrackedParameter< double >( "alpha", 1.5076494 ) ),
74 _beta( iConfig.getUntrackedParameter< double >( "beta", 1.5136036 ) ),
75 _nevtmax( iConfig.getUntrackedParameter< unsigned int >( "nEvtMax", 200 ) ),
76 _noise( iConfig.getUntrackedParameter< double >( "noise", 2.0 ) ),
77 _chi2cut( iConfig.getUntrackedParameter< double >( "chi2cut", 100.0 ) ),
78 _ecalPart( iConfig.getUntrackedParameter< std::string >( "ecalPart", "EB" ) ),
79 _qualpercent( iConfig.getUntrackedParameter< double >( "qualPercent", 0.2 ) ),
80 _debug( iConfig.getUntrackedParameter< int >( "debug", 0 ) ),
81 nCrys( NCRYSEB),
82 runType(-1), runNum(0), fedID(-1), dccID(-1), side(2), lightside(2), iZ(1),
83 phi(-1), eta(-1), event(0), color(-1), channelIteratorEE(0)
84 
85  //========================================================================
86 
87 {
88 
89 
90  // Initialization from cfg file
91 
92  resdir_ = iConfig.getUntrackedParameter<std::string>("resDir");
93 
94  digiCollection_ = iConfig.getParameter<std::string>("digiCollection");
95  digiProducer_ = iConfig.getParameter<std::string>("digiProducer");
96 
97  eventHeaderCollection_ = iConfig.getParameter<std::string>("eventHeaderCollection");
98  eventHeaderProducer_ = iConfig.getParameter<std::string>("eventHeaderProducer");
99 
100 
101 
102  // Geometrical constants initialization
103 
104 
105  if (_ecalPart == "EB") {
106  nCrys = NCRYSEB;
107  } else {
108  nCrys = NCRYSEE;
109  }
110  iZ = 1;
111  if(_fedid <= 609 )
112  iZ = -1;
113 
114  for(unsigned int j=0;j<nCrys;j++){
115  iEta[j]=-1;
116  iPhi[j]=-1;
117  iTowerID[j]=-1;
118  iChannelID[j]=-1;
119  idccID[j]=-1;
120  iside[j]=-1;
121  wasTimingOK[j]=true;
122  wasGainOK[j]=true;
123  nevtAB[j]=0 ;
124  }
125 
126  // Quality check flags
127 
128  isGainOK=true;
129  isTimingOK=true;
130 
131  // Objects dealing with pulses
132 
136 
137 
138  // Objects needed for npresample calculation
139 
140  Delta01=new TMom();
141  Delta12=new TMom();
142  _fitab=true;
143 
144 }
145 
146 //========================================================================
148  //========================================================================
149 
150 
151  // do anything here that needs to be done at desctruction time
152  // (e.g. close files, deallocate resources etc.)
153 
154 }
155 
156 
157 
158 //========================================================================
160  //========================================================================
161 
162 
163  //Calculate alpha and beta
164 
165 
166  // Define output results filenames and shape analyzer object (alpha,beta)
167  //=====================================================================
168 
169  // 1) AlphaBeta files
170 
171  doesABTreeExist=true;
172 
173  std::stringstream nameabinitfile;
174  nameabinitfile << resdir_ <<"/ABInit.root";
175  alphainitfile=nameabinitfile.str();
176 
177  std::stringstream nameabfile;
178  std::stringstream link;
179  nameabfile << resdir_ <<"/AB.root";
180 
181  FILE *test;
182  test = fopen(nameabinitfile.str().c_str(),"r");
183  if(test == NULL) {
184  doesABTreeExist=false;
185  _fitab=true;
186  };
187  delete test;
188 
189 
190  TFile *fAB=0; TTree *ABInit=0;
191  if(doesABTreeExist){
192  fAB=new TFile(nameabinitfile.str().c_str());
193  }
194  if(doesABTreeExist && fAB){
195  ABInit = (TTree*) fAB->Get("ABCol0");
196  }
197 
198  // 2) Shape analyzer
199 
200  if(doesABTreeExist && fAB && ABInit && ABInit->GetEntries()!=0){
201  shapana= new TShapeAnalysis(ABInit, _alpha, _beta, 5.5, 1.0);
202  doesABTreeExist=true;
203  }else{
204  shapana= new TShapeAnalysis(_alpha, _beta, 5.5, 1.0);
205  doesABTreeExist=false;
206  _fitab=true;
207  }
210 
211  if(doesABTreeExist && fAB ) fAB->Close();
212 
213  if(_fitab){
214  alphafile=nameabfile.str();
215  }else{
217  link<< "ln -s "<<resdir_<<"/ABInit.root "<< resdir_<<"/AB.root";
218  system(link.str().c_str());
219  }
220 
221  // Define output results files' names
222 
223  std::stringstream namefile;
224  namefile << resdir_ <<"/AB.root";
225  alphafile=namefile.str();
226 
227 }
228 
229 //========================================================================
231 //========================================================================
232 
233  ++iEvent;
234 
235  // retrieving DCC header
237  const EcalRawDataCollection* DCCHeader=0;
238  try {
240  DCCHeader=pDCCHeader.product();
241  }catch ( std::exception& ex ) {
242  std::cerr << "Error! can't get the product retrieving DCC header" << eventHeaderCollection_.c_str() <<" "<< eventHeaderProducer_.c_str() << std::endl;
243  }
244 
245  //retrieving crystal data from Event
247  const EBDigiCollection* EBDigi=0;
249  const EEDigiCollection* EEDigi=0;
250  if (_ecalPart == "EB") {
251  try {
253  EBDigi=pEBDigi.product();
254  }catch ( std::exception& ex ) {
255  std::cerr << "Error! can't get the product retrieving EB crystal data " << digiCollection_.c_str() << std::endl;
256  }
257  } else if (_ecalPart == "EE") {
258  try {
260  EEDigi=pEEDigi.product();
261  }catch ( std::exception& ex ) {
262  std::cerr << "Error! can't get the product retrieving EE crystal data " << digiCollection_.c_str() << std::endl;
263  }
264  } else {
265  std::cout <<" Wrong ecalPart in cfg file " << std::endl;
266  return;
267  }
268 
269  // retrieving electronics mapping
271  const EcalElectronicsMapping* TheMapping=0;
272  try{
273  c.get< EcalMappingRcd >().get(ecalmapping);
274  TheMapping = ecalmapping.product();
275  }catch ( std::exception& ex ) {
276  std::cerr << "Error! can't get the product EcalMappingRcd"<< std::endl;
277  }
278 
279  // =============================
280  // Decode DCCHeader Information
281  // =============================
282 
283 
284  for ( EcalRawDataCollection::const_iterator headerItr= DCCHeader->begin();headerItr != DCCHeader->end();
285  ++headerItr ) {
286 
287  // Get run type and run number
288 
289  int fed = headerItr->fedId();
290  if(fed!=_fedid && _fedid!=-999) continue;
291 
292  runType=headerItr->getRunType();
293  runNum=headerItr->getRunNumber();
294  event=headerItr->getLV1();
295 
296  dccID=headerItr->getDccInTCCCommand();
297  fedID=headerItr->fedId();
298  lightside=headerItr->getRtHalf();
299 
300  // Check fed corresponds to the DCC in TCC
301 
302  if( 600+dccID != fedID ) continue;
303 
304  // Cut on runType
305 
310 
311  // Retrieve laser color and event number
312 
313  EcalDCCHeaderBlock::EcalDCCEventSettings settings = headerItr->getEventSettings();
314  color = settings.wavelength;
315  if( color <0 ) return;
316 
317  std::vector<int>::iterator iter = find( colors.begin(), colors.end(), color );
318  if( iter==colors.end() ){
319  colors.push_back( color );
320  }
321  }
322 
323 
324  // Cut on fedID
325 
326  if(fedID!=_fedid && _fedid!=-999) return;
327 
328 
329  // ===========================
330  // Decode EBDigis Information
331  // ===========================
332 
333  int adcGain=0;
334 
335  if (EBDigi){
336 
337  // Loop on crystals
338  //===================
339 
340  for ( EBDigiCollection::const_iterator digiItr= EBDigi->begin(); digiItr != EBDigi->end(); ++digiItr ) { // Loop on crystals
341 
342  // Retrieve geometry
343  //===================
344 
345  EBDetId id_crystal(digiItr->id()) ;
346  EBDataFrame df( *digiItr );
347 
348  int etaG = id_crystal.ieta() ; // global
349  int phiG = id_crystal.iphi() ; // global
350 
351  int etaL ; // local
352  int phiL ; // local
353  std::pair<int, int> LocalCoord=MEEBGeom::localCoord( etaG , phiG );
354 
355  etaL=LocalCoord.first ;
356  phiL=LocalCoord.second ;
357 
358  eta = etaG;
359  phi = phiG;
360 
361  side=MEEBGeom::side(etaG,phiG);
362 
363  // Recover the TT id and the electronic crystal numbering from EcalElectronicsMapping
364 
365  EcalElectronicsId elecid_crystal = TheMapping->getElectronicsId(id_crystal);
366 
367  int towerID=elecid_crystal.towerId();
368  int strip=elecid_crystal.stripId();
369  int xtal=elecid_crystal.xtalId();
370  int channelID= 5*(strip-1) + xtal-1;
371 
372  unsigned int channel=MEEBGeom::electronic_channel( etaL, phiL );
373 
374  assert( channel < nCrys );
375 
376  iEta[channel]=eta;
377  iPhi[channel]=phi;
378  iTowerID[channel]=towerID;
379  iChannelID[channel]=channelID;
380  idccID[channel]=dccID;
381  iside[channel]=side;
382 
383  // APD Pulse
384  //===========
385 
386  // Loop on adc samples
387 
388  for (unsigned int i=0; i< (*digiItr).size() ; ++i ) {
389 
390  EcalMGPASample samp_crystal(df.sample(i));
391  adc[i]=samp_crystal.adc() ;
392  adcG[i]=samp_crystal.gainId();
393  adc[i]*=adcG[i];
394  if (i==0) adcGain=adcG[i];
395  if (i>0) adcGain=TMath::Max(adcG[i],adcGain);
396  }
397 
399 
400 
401  // Quality checks
402  //================
403 
404  if(adcGain!=1) nEvtBadGain[channel]++;
405  if(!APDPulse->isTimingQualOK()) nEvtBadTiming[channel]++;
406  nEvtTot[channel]++;
407 
408 
409  // Fill if Pulse is fine
410  //=======================
411 
412  if( APDPulse->isPulseOK() && lightside==side){
413 
416 
417  if( nevtAB[channel] < _nevtmax && _fitab ){
418  if(doesABTreeExist) shapana -> putAllVals(channel, adc, eta, phi);
419  else shapana -> putAllVals(channel, adc, eta, phi, dccID, side, towerID, channelID);
420  nevtAB[channel]++ ;
421  }
422  }
423  }
424 
425  } else if (EEDigi) {
426 
427  // Loop on crystals
428  //===================
429 
430  for ( EEDigiCollection::const_iterator digiItr= EEDigi->begin(); digiItr != EEDigi->end(); ++digiItr ) { // Loop on crystals
431 
432  // Retrieve geometry
433  //===================
434 
435  EEDetId id_crystal(digiItr->id()) ;
436  EEDataFrame df( *digiItr );
437 
438  phi = id_crystal.ix() ;
439  eta = id_crystal.iy() ;
440 
441  int iX = (phi-1)/5+1;
442  int iY = (eta-1)/5+1;
443 
444  side=MEEEGeom::side( iX, iY ,iZ);
445  EcalElectronicsId elecid_crystal = TheMapping->getElectronicsId(id_crystal);
446 
447  int towerID=elecid_crystal.towerId();
448  int channelID=elecid_crystal.channelId()-1;
449 
450  int hashedIndex=100000*eta+phi;
451 
452  if( channelMapEE.count(hashedIndex) == 0 ){
455  }
456 
457  unsigned int channel=channelMapEE[hashedIndex];
458 
459  assert ( channel < nCrys );
460 
461  iEta[channel]=eta;
462  iPhi[channel]=phi;
463  iTowerID[channel]=towerID;
464  iChannelID[channel]=channelID;
465  idccID[channel]=dccID;
466  iside[channel]=side;
467 
468  // APD Pulse
469  //===========
470 
471  if( (*digiItr).size()>10) std::cout <<"SAMPLES SIZE > 10!" << (*digiItr).size()<< std::endl;
472 
473  // Loop on adc samples
474 
475  for (unsigned int i=0; i< (*digiItr).size() ; ++i ) {
476 
477  EcalMGPASample samp_crystal(df.sample(i));
478  adc[i]=samp_crystal.adc() ;
479  adcG[i]=samp_crystal.gainId();
480  adc[i]*=adcG[i];
481 
482  if (i==0) adcGain=adcG[i];
483  if (i>0) adcGain=TMath::Max(adcG[i],adcGain);
484  }
485 
487 
488  // Quality checks
489  //================
490 
491  if(adcGain!=1) nEvtBadGain[channel]++;
492  if(!APDPulse->isTimingQualOK()) nEvtBadTiming[channel]++;
493  nEvtTot[channel]++;
494 
495  // Fill if Pulse is fine
496  //=======================
497 
498  if( APDPulse->isPulseOK() && lightside==side){
499 
502 
503  if( nevtAB[channel] < _nevtmax && _fitab ){
504  if(doesABTreeExist) shapana -> putAllVals(channel, adc, eta, phi);
505  else shapana -> putAllVals(channel, adc, eta, phi, dccID, side, towerID, channelID);
506  nevtAB[channel]++ ;
507  }
508  }
509  }
510  }
511 }// analyze
512 
513 
514 //========================================================================
516 //========================================================================
517 
518  std::cout << "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << std::endl;
519  std::cout << "\t+=+ Analyzing data: getting (alpha, beta) +=+" << std::endl;
520 
521  // Adjust channel numbers for EE
522  //===============================
523 
524  if( _ecalPart == "EE" ) {
525  nCrys=channelMapEE.size();
527  }
528 
529 
530  // Set presamples number
531  //======================
532 
533  double delta01=Delta01->getMean();
534  double delta12=Delta12->getMean();
535  if(delta12>_presamplecut) {
536  _presample=2;
537  if(delta01>_presamplecut) _presample=1;
538  }
539 
542 
543  // Get alpha and beta
544  //======================
545 
546  if(_fitab){
547  std::cout << "\n\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << std::endl;
548  std::cout << "\t+=+ Analyzing data: getting (alpha, beta) +=+" << std::endl;
549  TFile *fAB=0; TTree *ABInit=0;
550  if(doesABTreeExist){
551  fAB=new TFile(alphainitfile.c_str());
552  }
553  if(doesABTreeExist && fAB){
554  ABInit = (TTree*) fAB->Get("ABCol0");
555  }
556  shapana->computeShape(alphafile, ABInit);
557 
558  // Set quality flags for gains and timing
559 
560  double BadGainEvtPercentage=0.0;
561  double BadTimingEvtPercentage=0.0;
562 
563  int nChanBadGain=0;
564  int nChanBadTiming=0;
565 
566  for (unsigned int i=0;i<nCrys;i++){
567  if(nEvtTot[i]!=0){
568  BadGainEvtPercentage=double(nEvtBadGain[i])/double(nEvtTot[i]);
569  BadTimingEvtPercentage=double(nEvtBadTiming[i])/double(nEvtTot[i]);
570  }
571  if(BadGainEvtPercentage>_qualpercent) {
572  wasGainOK[i]=false;
573  nChanBadGain++;
574  }
575  if(BadTimingEvtPercentage>_qualpercent){
576  wasTimingOK[i]=false;
577  nChanBadTiming++;
578  }
579  }
580 
581  double BadGainChanPercentage=double(nChanBadGain)/double(nCrys);
582  double BadTimingChanPercentage=double(nChanBadTiming)/double(nCrys);
583 
584  if(BadGainChanPercentage>_qualpercent) isGainOK = false;
585  if(BadTimingChanPercentage>_qualpercent) isTimingOK = false;
586 
587 
588  if( !isGainOK )
589  std::cout << "\t+=+ ............................ WARNING! APD GAIN WAS NOT 1 +=+" << std::endl;
590  if( !isTimingOK )
591  std::cout << "\t+=+ ............................ WARNING! TIMING WAS BAD +=+" << std::endl;
592 
593  std::cout << "\t+=+ .................................... done +=+" << std::endl;
594  std::cout << "\t+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+=+" << std::endl;
595  }
596 
597 
598 }
599 
601 
void addEntry(double val)
Definition: TMom.cc:111
static XYCoord localCoord(int icr)
Definition: MEEBGeom.cc:153
T getParameter(std::string const &) const
unsigned int nCrys
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::string digiProducer_
std::string alphainitfile
int xtalId() const
get the channel id
std::string eventHeaderProducer_
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
std::string _ecalPart
int stripId() const
get the tower id
double getDelta(int, int)
Definition: TAPDPulse.cc:96
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Ecal readout channel identification [32:20] Unused (so far) [19:13] DCC id [12:6] tower [5:3] strip [...
double _ratiomaxcutlow
#define NCRYSEE
bool setPulse(double *)
Definition: TAPDPulse.cc:65
static int side(SuperCrysCoord iX, SuperCrysCoord iY, int iz)
Definition: MEEEGeom.cc:939
double _ratiomincuthigh
unsigned int nevtAB[NCRYSEB]
std::vector< T >::const_iterator const_iterator
unsigned int _nsamples
int nEvtBadGain[NCRYSEB]
#define NULL
Definition: scimark2.h:8
void computeShape(std::string namefile, TTree *)
int towerId() const
get the tower id
unsigned int _firstsample
int nEvtBadTiming[NCRYSEB]
const_iterator begin() const
T eta() const
Definition: TMom.h:7
virtual void beginJob()
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
unsigned int _nevtmax
double adc[10]
std::vector< int > colors
std::string alphafile
std::string link(std::string &nm, std::string &ns)
Definition: hierarchy.cc:24
int hashedIndex(int ieta, int iphi)
Definition: EcalPyUtils.cc:42
EcalElectronicsId getElectronicsId(const DetId &id) const
Get the electronics id for this det id.
double _presamplecut
int iEvent
Definition: GenABIO.cc:243
static int electronic_channel(EBLocalCoord ix, EBLocalCoord iy)
Definition: MEEBGeom.cc:341
#define NCRYSEB
int nEvtTot[NCRYSEB]
int j
Definition: DBlmapReader.cc:9
EcalABAnalyzer(const edm::ParameterSet &iConfig)
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
TShapeAnalysis * shapana
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
bool wasTimingOK[NCRYSEB]
int iside[NCRYSEB]
unsigned int _timingcuthigh
std::string resdir_
const_iterator end() const
static int side(EBGlobalCoord ieta, EBGlobalCoord iphi)
Definition: MEEBGeom.cc:114
bool isTimingQualOK()
Definition: TAPDPulse.cc:124
int idccID[NCRYSEB]
unsigned int _lastsample
int iEta[NCRYSEB]
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::map< int, int > channelMapEE
T const * product() const
Definition: Handle.h:74
TAPDPulse * APDPulse
const_iterator end() const
std::string digiCollection_
virtual void endJob()
unsigned int _timingqualhigh
void setPresamples(int)
Definition: TAPDPulse.cc:223
bool isPulseOK()
Definition: TAPDPulse.cc:140
unsigned int _presample
std::string eventHeaderCollection_
double _ratiomincutlow
virtual void analyze(const edm::Event &e, const edm::EventSetup &c)
void set_presample(int)
unsigned int _timingcutlow
int iTowerID[NCRYSEB]
tuple cout
Definition: gather_cfg.py:121
double getMean()
Definition: TMom.cc:148
int iPhi[NCRYSEB]
unsigned int _timingquallow
int channelId() const
so far for EndCap only :
const_iterator begin() const
bool wasGainOK[NCRYSEB]
Definition: DDAxes.h:10
int iChannelID[NCRYSEB]