CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
HFLightCalRand Class Reference

#include <HFLightCalRand.h>

Inheritance diagram for HFLightCalRand:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &fEvent, const edm::EventSetup &fSetup)
 
virtual void beginJob ()
 
virtual void endJob (void)
 
 HFLightCalRand (const edm::ParameterSet &fConfiguration)
 
virtual ~HFLightCalRand ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

edm::InputTag hcalCalibDigiCollectionTag_
 
edm::InputTag hfDigiCollectionTag_
 
std::string histfile
 
TH2F * hnpemapM
 
TH2F * hnpemapP
 
TH1F * hnpevar
 
TH1F * hped [26][36][2]
 
TH1F * hpedmean
 
TH1F * hpedpin [8][3]
 
TH1F * hpedrms
 
TH2F * hsignalmapM
 
TH2F * hsignalmapP
 
TH1F * hsignalmean
 
TH1F * hsignalrms
 
TH2F * hsignalRMSmapM
 
TH2F * hsignalRMSmapP
 
TH1F * hsp [26][36][2]
 
TH1F * hspe [26][36][2]
 
TH1F * hspepin [8][3]
 
TH1F * hspes
 
TH1F * hsppin [8][3]
 
TH1F * htmax
 
TH1F * htmean
 
TH1F * hts [26][36][2]
 
TH1F * htsm [26][36][2]
 
TH1F * htsmpin [8][3]
 
TH1F * htspin [8][3]
 
TFile * mFile
 
std::string prefile
 
FILE * preFile
 
std::string textfile
 
FILE * tFile
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 12 of file HFLightCalRand.h.

Constructor & Destructor Documentation

HFLightCalRand::HFLightCalRand ( const edm::ParameterSet fConfiguration)

Definition at line 43 of file HFLightCalRand.cc.

References edm::ParameterSet::getUntrackedParameter(), histfile, and textfile.

43  :
44  hfDigiCollectionTag_(fConfiguration.getParameter<edm::InputTag>("hfDigiCollectionTag")),
45  hcalCalibDigiCollectionTag_(fConfiguration.getParameter<edm::InputTag>("hcalCalibDigiCollectionTag")) {
46 
47  //std::string histfile = fConfiguration.getUntrackedParameter<string>("rootFile");
48  histfile = fConfiguration.getUntrackedParameter<string>("rootFile");
49  textfile = fConfiguration.getUntrackedParameter<string>("textFile");
50 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::string textfile
edm::InputTag hfDigiCollectionTag_
std::string histfile
edm::InputTag hcalCalibDigiCollectionTag_
HFLightCalRand::~HFLightCalRand ( )
virtual

Definition at line 52 of file HFLightCalRand.cc.

52  {
53  //delete mFile;
54 }

Member Function Documentation

void HFLightCalRand::analyze ( const edm::Event fEvent,
const edm::EventSetup fSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 373 of file HFLightCalRand.cc.

References ecalMGPA::adc(), calib, gather_cfg::cout, HcalDetId::depth(), edm::EventID::event(), EventN, HcalObjRepresent::Fill(), edm::Event::getByLabel(), hcalCalibDigiCollectionTag_, HcalForward, hfDigiCollectionTag_, hped, hsp, hspe, htmax, htmean, hts, htsm, HFDataFrame::id(), edm::EventBase::id(), HcalDetId::ieta(), cuy::ii, HcalDetId::iphi(), edm::EventID::run(), runNumb, convertSQLiteXML::runNumber, HFDataFrame::size(), and HcalDetId::subdet().

373  {
374 
375  // event ID
376  edm::EventID eventId = fEvent.id();
377  int runNumber = eventId.run ();
378  int eventNumber = eventId.event ();
379  if (runNumb==0) runNumb=runNumber;
380  EventN++;
381  if (verbose) std::cout << "========================================="<<std::endl
382  << "run/event: "<<runNumber<<'/'<<eventNumber<<std::endl;
383 
384  Double_t buf[20];
385  Double_t maxADC,signal,ped=0,meant;
386  Int_t maxisample=0,i1=3,i2=6;
387 
388  // HF PIN-diodes
391  if (verbose) std::cout<<"Analysis-> total CAL digis= "<<calib->size()<<std::endl;
392  /* COMMENTED OUT by J. Mans (7-28-2008) as major changes needed with new Calib DetId
393  re-commented out by R.Ofierzynski (11.Nov.2008) - to be able to provide a consistent code for CMSSW_3_0_0_pre3:
394  major changes are needed for the new Calib DetId which does not have the old methods any more
395 
396  for (unsigned j = 0; j < calib->size (); ++j) {
397  const HcalCalibDataFrame digi = (*calib)[j];
398  HcalElectronicsId elecId = digi.elecId();
399  HcalCalibDetId calibId = digi.id();
400  if (verbose) std::cout<<calibId.sectorString().c_str()<<" "<<calibId.rbx()<<" "<<elecId.fiberChanId()<<std::endl;
401  int isector = calibId.rbx()-1;
402  int ipin = elecId.fiberChanId();
403  int iside = -1;
404  if (calibId.sectorString() == "HFP") iside = 0;
405  else if (calibId.sectorString() == "HFM") iside = 4;
406 
407  if (iside != -1) {
408  maxADC=-99;
409  for (int isample = 0; isample < digi.size(); ++isample) {
410  int adc = digi[isample].adc();
411  int capid = digi[isample].capid ();
412  double linear_ADC = digi[isample].nominal_fC();
413  if (verbose) std::cout<<"PIN linear_ADC = "<<linear_ADC<<std::endl;
414  htspin[isector+iside][ipin]->Fill(isample,linear_ADC);
415  buf[isample]=linear_ADC;
416  if (maxADC<linear_ADC) {
417  maxADC=linear_ADC;
418  maxisample=isample;
419  }
420  }
421  i1=maxisample-1;
422  i2=maxisample+2;
423  if (i1<0) {i1=0;i2=3;}
424  else if (i2>9) {i1=6;i2=9;}
425  if (i1==0) ped=buf[8]+buf[9]+buf[6]+buf[7];
426  else if (i1==1) ped=buf[8]+buf[9]+buf[6]+buf[7];
427  else if (i1==2) ped=buf[0]+buf[1]+buf[6]+buf[7];
428  else if (i1==3) ped=buf[0]+buf[1]+buf[2]+buf[7];
429  else if (i1>=4) ped=buf[0]+buf[1]+buf[2]+buf[3];
430  signal=0;
431  for (ii=0;ii<4;ii++) signal+=TMath::Max(buf[ii+i1],ped/4.0);
432  hsppin[isector+iside][ipin]->Fill(signal);
433  hspepin[isector+iside][ipin]->Fill(signal);
434  hpedpin[isector+iside][ipin]->Fill(ped);
435 
436  // Mean PIN signal time estimation
437  ped=ped/4;
438  meant=0;
439  for (ii=0;ii<4;ii++) meant+=(TMath::Max(buf[ii+i1],ped)-ped)*(ii+i1);
440  if (signal-ped*4>0) meant/=(signal-ped*4);
441  else meant=i1+1;
442  htsmpin[isector+iside][ipin]->Fill(meant);
443  }
444  }
445  */
446  // HF
448  fEvent.getByLabel(hfDigiCollectionTag_, hf_digi);
449  if (verbose) std::cout<<"Analysis-> total HF digis= "<<hf_digi->size()<<std::endl;
450 
451  for (unsigned ihit = 0; ihit < hf_digi->size (); ++ihit) {
452  const HFDataFrame& frame = (*hf_digi)[ihit];
453  HcalDetId detId = frame.id();
454  int ieta = detId.ieta();
455  int iphi = detId.iphi();
456  int depth = detId.depth();
457  if (verbose) std::cout <<"HF digi # " <<ihit<<": eta/phi/depth: "
458  <<ieta<<'/'<<iphi<<'/'<< depth << std::endl;
459 
460  if (ieta>0) ieta = ieta-29;
461  else ieta = 13-ieta-29;
462 
463  for (int ii=0; ii<10; ii++) buf[ii]=0;
464  maxADC=-99;
465  for (int isample = 0; isample < frame.size(); ++isample) {
466  int adc = frame[isample].adc();
467  int capid = frame[isample].capid ();
468  double linear_ADC = frame[isample].nominal_fC();
469  double nominal_fC = detId.subdet () == HcalForward ? 2.6 * linear_ADC : linear_ADC;
470 
471  if (verbose) std::cout << "Analysis-> HF sample # " << isample
472  << ", capid=" << capid
473  << ": ADC=" << adc
474  << ", linearized ADC=" << linear_ADC
475  << ", nominal fC=" << nominal_fC << std::endl;
476 
477  hts[ieta][(iphi-1)/2][depth-1]->Fill(isample,linear_ADC);
478  buf[isample]=linear_ADC;
479  if (maxADC<linear_ADC) {
480  maxADC=linear_ADC;
481  maxisample=isample;
482  }
483  }
484 
485  // Signal is four capIDs around maxTS, Pedestal is four capID off the signal
486  maxADC=-99;
487  for (int ibf=0; ibf<frame.size()-1; ibf++) {
488  Double_t sumbuf=0;
489  for (int jbf=0; jbf<2; jbf++) {
490  sumbuf += buf[ibf+jbf];
491  if (ibf+jbf<2) sumbuf -= (buf[ibf+jbf+4]+buf[ibf+jbf+8])/2.0;
492  else if (ibf+jbf<4) sumbuf -= buf[ibf+jbf+4];
493  else if (ibf+jbf<6) sumbuf -= (buf[ibf+jbf-4]+buf[ibf+jbf+4])/2.0;
494  else if (ibf+jbf<8) sumbuf -= buf[ibf+jbf-4];
495  else if (ibf+jbf<10) sumbuf -= (buf[ibf+jbf-4]+buf[ibf+jbf-8])/2.0;
496  }
497  if (sumbuf>maxADC) {
498  maxADC=sumbuf;
499  if (buf[ibf]>buf[ibf+1]) maxisample=ibf;
500  else maxisample=ibf+1;
501  //maxisample=ibf;
502  }
503  }
504  htmax->Fill(maxisample);
505  i1=maxisample-1;
506  i2=maxisample+2;
507  if (i1<0) {i1=0;i2=3;}
508  else if (i2>9) {i1=6;i2=9;}
509  else if (i2<9) {
510  if (buf[i1]<buf[i2+1]) {i1=i1+1;i2=i2+1;}
511  }
512  signal=buf[i1]+buf[i1+1]+buf[i1+2]+buf[i1+3];
513  hsp[ieta][(iphi-1)/2][depth-1]->Fill(signal);
514  hspe[ieta][(iphi-1)/2][depth-1]->Fill(signal);
515 
516  if (i1==0) ped=(buf[4]+buf[8])/2.0+(buf[5]+buf[9])/2.0+buf[6]+buf[7];
517  else if (i1==1) ped=(buf[0]+buf[8])/2.0+(buf[5]+buf[9])/2.0+buf[6]+buf[7];
518  else if (i1==2) ped=(buf[0]+buf[8])/2.0+(buf[1]+buf[9])/2.0+buf[6]+buf[7];
519  else if (i1==3) ped=(buf[0]+buf[8])/2.0+(buf[1]+buf[9])/2.0+buf[2]+buf[7];
520  else if (i1==4) ped=(buf[0]+buf[8])/2.0+(buf[1]+buf[9])/2.0+buf[2]+buf[3];
521  else if (i1==5) ped=(buf[0]+buf[4])/2.0+(buf[1]+buf[9])/2.0+buf[2]+buf[3];
522  else if (i1==6) ped=(buf[0]+buf[4])/2.0+(buf[1]+buf[5])/2.0+buf[2]+buf[3];
523  /*
524  if (i1==0) ped=buf[8]+buf[9]+buf[6]+buf[7];
525  else if (i1==1) ped=buf[0]+buf[9]+buf[6]+buf[7];
526  else if (i1==2) ped=buf[0]+buf[1]+buf[6]+buf[7];
527  else if (i1==3) ped=buf[0]+buf[1]+buf[2]+buf[7];
528  else if (i1>=4) ped=buf[0]+buf[1]+buf[2]+buf[3];
529  */
530  hped[ieta][(iphi-1)/2][depth-1]->Fill(ped);
531 
532  // Mean signal time estimation
533  ped=ped/4;
534  meant=(buf[i1]-ped)*i1+(buf[i1+1]-ped)*(i1+1)+(buf[i1+2]-ped)*(i1+2)+(buf[i1+3]-ped)*(i1+3);
535  meant /= (buf[i1]-ped)+(buf[i1+1]-ped)+(buf[i1+2]-ped)+(buf[i1+3]-ped);
536  htmean->Fill(meant);
537  htsm[ieta][(iphi-1)/2][depth-1]->Fill(meant);
538  }
539 }
int adc(sample_type sample)
get the ADC sample (12 bits)
RunNumber_t run() const
Definition: EventID.h:42
EventNumber_t event() const
Definition: EventID.h:44
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:32
int ii
Definition: cuy.py:588
Int_t runNumb
int depth() const
get the tower depth
Definition: HcalDetId.h:42
MVATrainerComputer * calib
Definition: MVATrainer.cc:64
Int_t EventN
TH1F * hped[26][36][2]
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
TH1F * hspe[26][36][2]
int ieta() const
get the cell ieta
Definition: HcalDetId.h:38
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
edm::InputTag hfDigiCollectionTag_
int iphi() const
get the cell iphi
Definition: HcalDetId.h:40
int size() const
total number of samples in the digi
Definition: HFDataFrame.h:26
TH1F * hsp[26][36][2]
edm::EventID id() const
Definition: EventBase.h:56
TH1F * htsm[26][36][2]
tuple cout
Definition: gather_cfg.py:121
const HcalDetId & id() const
Definition: HFDataFrame.h:22
edm::InputTag hcalCalibDigiCollectionTag_
TH1F * hts[26][36][2]
void HFLightCalRand::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 56 of file HFLightCalRand.cc.

References gather_cfg::cout, histfile, hnpemapM, hnpemapP, hnpevar, hped, hpedmean, hpedpin, hpedrms, hsignalmapM, hsignalmapP, hsignalmean, hsignalrms, hsignalRMSmapM, hsignalRMSmapP, hsp, hspe, hspepin, hspes, hsppin, htmax, htmean, hts, htsm, htsmpin, htspin, i, j, gen::k, mFile, NULL, textfile, and tFile.

56  {
57 
58  char htit[64];
59  mFile = new TFile (histfile.c_str(),"RECREATE");
60  if ((tFile = fopen(textfile.c_str(),"w"))==NULL)printf("\nNo textfile open\n\n");
61  // General Histos
62  htmax = new TH1F("htmax","Max TS",10,-0.5,9.5);
63  htmean = new TH1F("htmean","Mean signal TS",100,0,10);
64  hsignalmean = new TH1F("hsignalmean","Mean ADC 4maxTS",1201,-25,30000);
65  hsignalrms = new TH1F("hsignalrms","RMS ADC 4maxTS",500,0,500);
66  hpedmean = new TH1F("hpedmean","Mean ADC 4lowTS",200,-10,90);
67  hpedrms = new TH1F("hpedrms","RMS ADC 4lowTS",200,0,100);
68  hspes = new TH1F("hspes","SPE if measured",200,0,40);
69  hnpevar = new TH1F("hnpevar","~N PE input",500,0,500);
70  hsignalmapP = new TH2F("hsignalmapP","Mean(Response) - Mean(Pedestal) HFP;#eta;#phi",26,28.5,41.5,36,0,72);
71  hsignalRMSmapP = new TH2F("hsignalRMSmapP","RMS Response HFP;#eta;#phi",26,28.5,41.5,36,0,72);
72  hnpemapP = new TH2F("hnpemapP","~N PE input HFP;#eta;#phi",26,28.5,41.5,36,0,72);
73  hnpemapP->SetOption("COLZ");hsignalmapP->SetOption("COLZ");hsignalRMSmapP->SetOption("COLZ");
74  hsignalmapM = new TH2F("hsignalmapM","Mean(Response) - Mean(Pedestal) HFM;#eta;#phi",26,-41.5,-28.5,36,0,72);
75  hsignalRMSmapM = new TH2F("hsignalRMSmapM","RMS Response HFM;#eta;#phi",26,-41.5,-28.5,36,0,72);
76  hnpemapM = new TH2F("hnpemapM","~N PE input HFM;#eta;#phi",26,-41.5,-28.5,36,0,72);
77  hnpemapM->SetOption("COLZ");hsignalmapM->SetOption("COLZ");hsignalRMSmapM->SetOption("COLZ");
78  // Channel-by-channel histos
79  for (int i=0;i<13;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
80  if (i>10 && j%2==0) continue;
81  sprintf(htit,"ts_+%d_%d_%d",i+29,j*2+1,k+1);
82  hts[i][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape)
83  sprintf(htit,"tsmean_+%d_%d_%d",i+29,j*2+1,k+1);
84  htsm[i][j][k] = new TH1F(htit,htit,100,0,10); // Mean signal time estimated from TS
85  sprintf(htit,"sp_+%d_%d_%d",i+29,j*2+1,k+1);
86  hsp[i][j][k] = new TH1F(htit,htit,1201,-25,30000); // Big-scale spectrum (linear ADC)
87  sprintf(htit,"spe_+%d_%d_%d",i+29,j*2+1,k+1);
88  hspe[i][j][k] = new TH1F(htit,htit,200,-9.5,190.5); // Small-scale spectrum (linear ADC)
89  sprintf(htit,"ped_+%d_%d_%d",i+29,j*2+1,k+1);
90  hped[i][j][k] = new TH1F(htit,htit,200,-9.5,190.5); // Pedestal spectrum
91  sprintf(htit,"ts_-%d_%d_%d",i+29,j*2+1,k+1);
92  hts[i+13][j][k] = new TH1F(htit,htit,10,-0.5,9.5);
93  sprintf(htit,"tsmean_-%d_%d_%d",i+29,j*2+1,k+1);
94  htsm[i+13][j][k] = new TH1F(htit,htit,100,0,10);
95  sprintf(htit,"sp_-%d_%d_%d",i+29,j*2+1,k+1);
96  hsp[i+13][j][k] = new TH1F(htit,htit,1201,-25,30000);
97  sprintf(htit,"spe_-%d_%d_%d",i+29,j*2+1,k+1);
98  hspe[i+13][j][k] = new TH1F(htit,htit,200,-9.5,190.5);
99  sprintf(htit,"ped_-%d_%d_%d",i+29,j*2+1,k+1);
100  hped[i+13][j][k] = new TH1F(htit,htit,200,-9.5,190.5);
101  }
102  // PIN-diodes histos
103  for (int i=0;i<4;i++) for (int j=0;j<3;j++) {
104  sprintf(htit,"ts_PIN%d_+Q%d",j+1,i+1);
105  htspin[i][j] = new TH1F(htit,htit,10,-0.5,9.5);
106  sprintf(htit,"sp_PIN%d_+Q%d",j+1,i+1);
107  hsppin[i][j] = new TH1F(htit,htit,1601,-25,40000);
108  sprintf(htit,"spe_PIN%d_+Q%d",j+1,i+1);
109  hspepin[i][j] = new TH1F(htit,htit,200,-9.5,190.5);
110  sprintf(htit,"ped_PIN%d_+Q%d",j+1,i+1);
111  hpedpin[i][j] = new TH1F(htit,htit,200,-9.5,190.5);
112  sprintf(htit,"tsmean_PIN%d_+Q%d",j+1,i+1);
113  htsmpin[i][j] = new TH1F(htit,htit,100,0,10);
114  sprintf(htit,"ts_PIN%d_-Q%d",j+1,i+1);
115  htspin[i+4][j] = new TH1F(htit,htit,10,-0.5,9.5);
116  sprintf(htit,"sp_PIN%d_-Q%d",j+1,i+1);
117  hsppin[i+4][j] = new TH1F(htit,htit,1601,-25,40000);
118  sprintf(htit,"spe_PIN%d_-Q%d",j+1,i+1);
119  hspepin[i+4][j] = new TH1F(htit,htit,200,-9.5,190.5);
120  sprintf(htit,"ped_PIN%d_-Q%d",j+1,i+1);
121  hpedpin[i+4][j] = new TH1F(htit,htit,200,-9.5,190.5);
122  sprintf(htit,"tsmean_PIN%d_-Q%d",j+1,i+1);
123  htsmpin[i+4][j] = new TH1F(htit,htit,100,0,10);
124  }
125  std::cout<<std::endl<<"beginJob: histfile="<<histfile.c_str()<<" textfile="<<textfile.c_str()<<std::endl;
126  return;
127 }
int i
Definition: DBlmapReader.cc:9
TH1F * hpedpin[8][3]
#define NULL
Definition: scimark2.h:8
TH1F * hsppin[8][3]
TH1F * hped[26][36][2]
TH1F * hspe[26][36][2]
int j
Definition: DBlmapReader.cc:9
std::string textfile
std::string histfile
int k[5][pyjets_maxn]
TH1F * htsmpin[8][3]
TH2F * hsignalRMSmapM
TH1F * hsp[26][36][2]
TH2F * hsignalRMSmapP
TH1F * htsm[26][36][2]
tuple cout
Definition: gather_cfg.py:121
TH1F * htspin[8][3]
TH1F * hspepin[8][3]
TH1F * hts[26][36][2]
void HFLightCalRand::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 193 of file HFLightCalRand.cc.

References gather_cfg::cout, alignCSCRings::e, EventN, Fit3Peak(), HistSpec(), hnpemapM, hnpemapP, hnpevar, hped, hpedmean, hpedpin, hpedrms, hsignalmapM, hsignalmapP, hsignalmean, hsignalrms, hsignalRMSmapM, hsignalRMSmapP, hsp, hspe, hspepin, hspes, hsppin, hts, htspin, i, j, gen::k, timingPdfMaker::mean, mFile, NEvents, plotscripts::rms(), runNumb, and tFile.

194 {
195  Double_t mean,rms,meanped,rmsped,npevar;
196  Double_t par[5],parm[5],dspe=0,dnpe;
197  Int_t tsmax,intspe;
198  fprintf(tFile,"#RunN %d Events processed %d",runNumb,EventN);
199  std::cout<<"endJob: histos processing..."<<std::endl;
200  std::cout<<"RunN= "<<runNumb<<" Events processed= "<<EventN<<std::endl;
201 
202  for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
203  if (i>10 && i<13 && j%2==0) continue;
204  if (i>23 && j%2==0) continue;
205  meanped=rmsped=mean=rms=0;
206  if (hsp[i][j][k]->Integral()>0) {
207  HistSpec(hped[i][j][k],meanped,rmsped);
208  HistSpec(hsp[i][j][k],mean,rms);
209  if (hspe[i][j][k]->Integral()>hsp[i][j][k]->Integral()*0.9 || mean<100) {
210  HistSpec(hspe[i][j][k],mean,rms);
211  }
212  hsignalmean->Fill(mean);
213  hsignalrms->Fill(rms);
214  hpedmean->Fill(meanped);
215  hpedrms->Fill(rmsped);
216  }
217  }
218 
219  meanped=hpedmean->GetMean();
220  rmsped=hpedrms->GetMean();
221  mean=hsignalmean->GetMean();
222  rms=hsignalrms->GetMean();
223  fprintf(tFile," MeanInput=<%.2f> [linADCcount] RMS=<%.2f>\n",mean,rms);
224  fprintf(tFile,"#eta/phi/depth sum4maxTS RMS ~N_PE sum4lowTS RMS maxTS SPE +/- Err Comment\n");
225  TF1* fPed = new TF1("fPed","gaus",0,120);
226  fPed->SetNpx(200);
227  TF1 *fTot = new TF1("fTot",Fit3Peak ,0,200,5);
228  fTot->SetNpx(800);
229  for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
230  if (i>10 && i<13 && j%2==0) continue;
231  if (i>23 && j%2==0) continue;
232  HistSpec(hped[i][j][k],meanped,rmsped);
233  HistSpec(hsp[i][j][k],mean,rms);
234  par[3]=0;
235  if (hspe[i][j][k]->Integral()>hsp[i][j][k]->Integral()*0.9 || mean<100) {
236  HistSpec(hspe[i][j][k],mean,rms);
237  if (hspe[i][j][k]->Integral(1,(int) (meanped+3*rmsped+12))/NEvents>0.1) {
238  //if (hspe[i][j][k]->Integral()>100 && mean-meanped<100) {
239  if (mean+rms*3-meanped-rmsped*3>1 && rmsped>0) { // SPE fit if low intensity>0
240  par[1] = meanped;
241  par[2] = rmsped;
242  par[0] = hped[i][j][k]->GetMaximum();
243  fPed->SetParameters(par);
244  hped[i][j][k]->Fit(fPed,"BQ0");
245  fPed->GetParameters(&par[0]);
246  hped[i][j][k]->Fit(fPed,"B0Q","",par[1]-par[2]*3,par[1]+par[2]*3);
247  fPed->GetParameters(par);
248  hped[i][j][k]->Fit(fPed,"BLIQ","",par[1]-par[2]*3,par[1]+par[2]*3);
249  fPed->GetParameters(&par[0]);
250  fPed->GetParameters(&parm[0]);
251  NEvents = (int) hspe[i][j][k]->Integral();
252  par[0]=0.1;
253  par[3]=10;
254  par[4]=6;
255  par[2]=par[2]/0.93;
256  fTot->SetParameters(par);
257  fTot->SetParLimits(0,0,2);
258  //fTot->FixParameter(1,par[1]);
259  //fTot->FixParameter(2,par[2]);
260  fTot->SetParLimits(1,par[1]-4,par[1]+4);
261  fTot->SetParLimits(2,par[2]*0.9,par[2]);
262  fTot->SetParLimits(3,1.2,50);
263  //fTot->SetParLimits(4,par[2]*1.1,100);
264  par[4]=0;
265  fTot->SetParLimits(4,-1.64,1.64);
266  hspe[i][j][k]->Fit(fTot,"BLEQ","");
267  fTot->GetParameters(par);
268  dspe=fTot->GetParError(3);
269  dnpe=fTot->GetParError(0);
270  if (par[3]<1.21 || dnpe>par[0]) par[3]=-1;
271  else if (par[0]>1.96 || par[3]>49) par[3]=0;
272  else {
273  hspes->Fill(par[3]);
274  }
275  }
276  }
277  }
278 
279  // NPE
280  npevar=0;
281  if (par[3]>0) npevar=par[0]; // NPE from SPE fit
282  else { // NPE from high intensity signal
283  if (hspe[i][j][k]->Integral()>hsp[i][j][k]->Integral()*0.98) {
284  HistSpec(hspe[i][j][k],mean,rms,3);
285  }
286  else {
287  HistSpec(hsp[i][j][k],mean,rms,3);
288  }
289  if (rmsped>0) {
290  if ((rms*rms-rmsped*rmsped/0.8836)>1 && mean>meanped+2) {
291  npevar=(mean-meanped-2)*(mean-meanped-2)/(rms*rms-rmsped*rmsped/0.8836);
292  }
293  else if (mean<100) {
294  intspe=int(hspe[i][j][k]->Integral());
295  hspe[i][j][k]->SetAxisRange(meanped+2+rmsped*4,300);
296  npevar=hspe[i][j][k]->Integral()/intspe;
297  if (npevar>0.01) npevar=-1;
298  else npevar=0;
299  hspe[i][j][k]->SetAxisRange(-20,300);
300  }
301  }
302  }
303  if (npevar>5.0e-5) hnpevar->Fill(npevar);
304 
305  if (i<13) {
306  hsignalmapP->Fill(i+28.6+k/2.0,j*2+1,mean-meanped-1.8);
307  hsignalRMSmapP->Fill(i+28.6+k/2.0,j*2+1,rms);
308  if (npevar>0) hnpemapP->Fill(i+28.6+k/2.0,j*2+1,npevar);
309  fprintf(tFile,"%3d%4d%5d %11.2f%8.2f",i+29,j*2+1,k+1,mean,rms);
310  }
311  else {
312  fprintf(tFile,"%3d%4d%5d %11.2f%8.2f",13-i-29,j*2+1,k+1,mean,rms);
313  hsignalmapM->Fill(13-i-28.6-k/2.0,j*2+1,mean-meanped-1.8);
314  hsignalRMSmapM->Fill(13-i-28.6-k/2.0,j*2+1,rms);
315  if (npevar>0) hnpemapM->Fill(13-i-28.6-k/2.0,j*2+1,npevar);
316  }
317  if (npevar>0) fprintf(tFile," %9.4f",npevar);
318  else fprintf(tFile," 0 ");
319  fprintf(tFile," %8.2f%8.2f",meanped,rmsped);
320  tsmax=hts[i][j][k]->GetMaximumBin()-1;
321  fprintf(tFile," %4d",tsmax);
322  if (par[3]>0 && par[3]<99) fprintf(tFile,"%8.2f%7.2f",par[3]+1,dspe);
323  else if (npevar>0) fprintf(tFile,"%8.2f 0 ",(mean-meanped-1)/npevar);
324  else fprintf(tFile," 0 0 ");
325 
326  // Diagnostics
327  fprintf(tFile," ");
328  if (hsp[i][j][k]->GetEntries()<=0) fprintf(tFile,"NoSignal\n");
329  else if (hsp[i][j][k]->GetEntries()<=10) fprintf(tFile,"Nev<10\n");
330  else {
331  if (hsp[i][j][k]->Integral()<=10 || mean>12000) fprintf(tFile,"SignalOffRange\n");
332  else {
333  if (hsp[i][j][k]->Integral()<100) fprintf(tFile,"Nev<100/");
334  if (npevar>0 && par[3]>0 && (npevar*NEvents<10 || npevar<0.001))
335  fprintf(tFile,"LowSignal/");
336  else if (npevar==0 && fabs(mean-meanped)<3) fprintf(tFile,"LowSignal/");
337  if (par[3]<0) fprintf(tFile,"BadFit/");
338  else if (par[3]==0) fprintf(tFile,"NoSPEFit/");
339  else if (par[3]>0 && npevar>1) fprintf(tFile,"NPE>1/");
340  if (npevar<0) fprintf(tFile,"Problem/");
341  if (mean<2) fprintf(tFile,"LowMean/");
342  if (rms<0.5) fprintf(tFile,"LowRMS/");
343  if (meanped<-1) fprintf(tFile,"LowPed/");
344  else if (meanped>25) fprintf(tFile,"HighPed/");
345  if (rmsped<0.5 && rmsped>0) fprintf(tFile,"NarrowPed/");
346  else if (rmsped>10) fprintf(tFile,"WidePed/");
347  if (hped[i][j][k]->GetBinContent(201)>10) fprintf(tFile,"PedOffRange");
348  fprintf(tFile,"-\n");
349  }
350  }
351  }
352 
353  for (int i=0;i<8;i++) for (int j=0;j<3;j++) {
354  HistSpec(hpedpin[i][j],meanped,rmsped);
355  HistSpec(hsppin[i][j],mean,rms);
356  if (hspepin[i][j]->Integral()>hsppin[i][j]->Integral()*0.9 || mean<100) {
357  HistSpec(hspepin[i][j],mean,rms);
358  }
359  if (i<4) fprintf(tFile," PIN%d +Q%d %12.2f %6.2f",j+1,i+1,mean,rms);
360  else fprintf(tFile," PIN%d -Q%d %12.2f %6.2f",j+1,i-3,mean,rms);
361  fprintf(tFile," %15.2f %6.2f",meanped,rmsped);
362  tsmax=htspin[i][j]->GetMaximumBin()-1;
363  fprintf(tFile," %4d\n",tsmax);
364  }
365 
366  mFile->Write();
367  mFile->Close();
368  fclose(tFile);
369  std::cout<<std::endl<<" --endJob-- done"<<std::endl;
370  return;
371 }
int i
Definition: DBlmapReader.cc:9
Double_t Fit3Peak(Double_t *x, Double_t *par)
TH1F * hpedpin[8][3]
TH1F * hsppin[8][3]
Int_t runNumb
Int_t EventN
TH1F * hped[26][36][2]
TH1F * hspe[26][36][2]
int j
Definition: DBlmapReader.cc:9
int k[5][pyjets_maxn]
TH2F * hsignalRMSmapM
TH1F * hsp[26][36][2]
TH2F * hsignalRMSmapP
Int_t NEvents
tuple cout
Definition: gather_cfg.py:121
TH1F * htspin[8][3]
TH1F * hspepin[8][3]
TH1F * hts[26][36][2]
void HistSpec(TH1F *hist, Double_t &mean, Double_t &rms, Double_t range=4)

Member Data Documentation

edm::InputTag HFLightCalRand::hcalCalibDigiCollectionTag_
private

Definition at line 47 of file HFLightCalRand.h.

Referenced by analyze().

edm::InputTag HFLightCalRand::hfDigiCollectionTag_
private

Definition at line 46 of file HFLightCalRand.h.

Referenced by analyze().

std::string HFLightCalRand::histfile
private

Definition at line 27 of file HFLightCalRand.h.

Referenced by beginJob(), and HFLightCalRand().

TH2F * HFLightCalRand::hnpemapM
private

Definition at line 38 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH2F* HFLightCalRand::hnpemapP
private

Definition at line 38 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F * HFLightCalRand::hnpevar
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F* HFLightCalRand::hped[26][36][2]
private

Definition at line 37 of file HFLightCalRand.h.

Referenced by analyze(), beginJob(), and endJob().

TH1F * HFLightCalRand::hpedmean
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F* HFLightCalRand::hpedpin[8][3]
private

Definition at line 43 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F * HFLightCalRand::hpedrms
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH2F * HFLightCalRand::hsignalmapM
private

Definition at line 38 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH2F * HFLightCalRand::hsignalmapP
private

Definition at line 38 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F* HFLightCalRand::hsignalmean
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F * HFLightCalRand::hsignalrms
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH2F * HFLightCalRand::hsignalRMSmapM
private

Definition at line 38 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH2F * HFLightCalRand::hsignalRMSmapP
private

Definition at line 38 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F* HFLightCalRand::hsp[26][36][2]
private

Definition at line 35 of file HFLightCalRand.h.

Referenced by analyze(), beginJob(), and endJob().

TH1F* HFLightCalRand::hspe[26][36][2]
private

Definition at line 36 of file HFLightCalRand.h.

Referenced by analyze(), beginJob(), and endJob().

TH1F* HFLightCalRand::hspepin[8][3]
private

Definition at line 42 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F * HFLightCalRand::hspes
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F* HFLightCalRand::hsppin[8][3]
private

Definition at line 41 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TH1F * HFLightCalRand::htmax
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by analyze(), and beginJob().

TH1F * HFLightCalRand::htmean
private

Definition at line 39 of file HFLightCalRand.h.

Referenced by analyze(), and beginJob().

TH1F* HFLightCalRand::hts[26][36][2]
private

Definition at line 33 of file HFLightCalRand.h.

Referenced by analyze(), beginJob(), and endJob().

TH1F* HFLightCalRand::htsm[26][36][2]
private

Definition at line 34 of file HFLightCalRand.h.

Referenced by analyze(), and beginJob().

TH1F* HFLightCalRand::htsmpin[8][3]
private

Definition at line 44 of file HFLightCalRand.h.

Referenced by beginJob().

TH1F* HFLightCalRand::htspin[8][3]
private

Definition at line 40 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

TFile* HFLightCalRand::mFile
private

Definition at line 30 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().

std::string HFLightCalRand::prefile
private

Definition at line 29 of file HFLightCalRand.h.

FILE* HFLightCalRand::preFile
private

Definition at line 32 of file HFLightCalRand.h.

std::string HFLightCalRand::textfile
private

Definition at line 28 of file HFLightCalRand.h.

Referenced by beginJob(), and HFLightCalRand().

FILE* HFLightCalRand::tFile
private

Definition at line 31 of file HFLightCalRand.h.

Referenced by beginJob(), and endJob().