CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFPreLightCal.cc
Go to the documentation of this file.
1 // Pre-Analysis for HFLightCal:
2 // finding time position of signal in TS
3 //
4 #include <memory>
5 #include <string>
6 #include <iostream>
7 
8 #include "TH1F.h"
9 #include "TH2F.h"
10 #include "TFile.h"
11 #include "math.h"
12 #include "TMath.h"
13 #include "TF1.h"
14 
16 
19 
24 
29 
30 using namespace std;
31 Int_t run_NN=0, event_NN=0;
32 
33 namespace {
34  //bool verbose = true;
35  bool verbose = false;
36 }
37 
39  hfDigiCollectionTag_(fConfiguration.getParameter<edm::InputTag>("hfDigiCollectionTag")),
40  hcalCalibDigiCollectionTag_(fConfiguration.getParameter<edm::InputTag>("hcalCalibDigiCollectionTag")) {
41 
42  //std::string histfile = fConfiguration.getUntrackedParameter<string>("rootFile");
43  histfile = fConfiguration.getUntrackedParameter<string>("rootPreFile");
44  textfile = fConfiguration.getUntrackedParameter<string>("textPreFile");
45 }
46 
48  //delete mFile;
49 }
50 
52 
53  char htit[64];
54  std::cout<<std::endl<<"HFPreLightCal beginJob: --> ";
55 
56  mFile = new TFile (histfile.c_str(),"RECREATE");
57  if ((tFile = fopen(textfile.c_str(),"w"))==NULL) {
58  printf("\nNo textfile open\n\n");
59  std::cout<<"Problem with output Pre-textFILE => exit"<<std::endl;
60  exit(1);
61  }
62  // Histos
63  htsmax = new TH1F("htsmax","Max TS",100,0,10);
64  htspinmax = new TH1F("htspinmax","Max TS PIN",100,0,10);
65  // Channel-by-channel histos
66  for (int i=0;i<13;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
67  if (i>10 && j%2==0) continue;
68  sprintf(htit,"tspre_+%d_%d_%d",i+29,j*2+1,k+1);
69  hts[i][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape)
70  sprintf(htit,"tspre_-%d_%d_%d",i+29,j*2+1,k+1);
71  hts[i+13][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape)
72  }
73  // PIN-diodes histos
74  for (int i=0;i<4;i++) for (int j=0;j<3;j++) {
75  sprintf(htit,"tspre_PIN%d_+Q%d",j+1,i+1);
76  htspin[i][j] = new TH1F(htit,htit,10,-0.5,9.5);
77  sprintf(htit,"tspre_PIN%d_-Q%d",j+1,i+1);
78  htspin[i+4][j] = new TH1F(htit,htit,10,-0.5,9.5);
79  }
80  std::cout<<"histfile="<<histfile.c_str()<<" textfile="<<textfile.c_str()<<std::endl;
81  return;
82 }
83 
85 {
86  Double_t sum,cont;
87  Int_t tsmax;
88 
89  std::cout<<std::endl<<"HFPreLightCal endJob --> ";
90 
91  for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
92  if (i>10 && i<13 && j%2==0) continue;
93  if (i>23 && j%2==0) continue;
94  sum=tsmax=0;
95  for (int ii=1; ii<=10; ii++) {
96  cont = hts[i][j][k]->GetBinContent(ii);
97  if (ii<3) cont=cont-(hts[i][j][k]->GetBinContent(ii+4)+hts[i][j][k]->GetBinContent(ii+8))/2;
98  else if (ii<5) cont=cont-hts[i][j][k]->GetBinContent(ii+4);
99  else if (ii<7) cont=cont-(hts[i][j][k]->GetBinContent(ii-4)+hts[i][j][k]->GetBinContent(ii+4))/2;
100  else if (ii<9) cont=cont-hts[i][j][k]->GetBinContent(ii-4);
101  else cont=cont-(hts[i][j][k]->GetBinContent(ii-4)+hts[i][j][k]->GetBinContent(ii-8))/2;
102  if (cont>sum) {
103  sum = cont;
104  tsmax=ii;
105  }
106  }
107  htsmax->Fill(tsmax);
108  if (i<13) fprintf(tFile," %d %d %d %d\n",i+29,j*2+1,k+1,tsmax);
109  else fprintf(tFile," %d %d %d %d\n",13-i-29,j*2+1,k+1,tsmax);
110  }
111 
112  for (int i=0;i<8;i++) for (int j=0;j<3;j++) {
113  sum=tsmax=0;
114  tsmax = htspin[i][j]->GetMaximumBin();
115  htspinmax->Fill(tsmax);
116  if (i<4) fprintf(tFile,"%d %d %d\n",i+1,j+1,tsmax);
117  else fprintf(tFile,"%d %d %d\n",-i+3,j+1,tsmax);
118  }
119 
120  mFile->Write();
121  mFile->Close();
122  fclose(tFile);
123  std::cout<<" Nevents = "<<event_NN<<std::endl;
124  return;
125 }
126 
127 void HFPreLightCal::analyze(const edm::Event& fEvent, const edm::EventSetup& fSetup) {
128 
129  // event ID
130  edm::EventID eventId = fEvent.id();
131  int runNumber = eventId.run ();
132  int eventNumber = eventId.event ();
133  if (run_NN==0) run_NN=runNumber;
134  event_NN++;
135  if (verbose) std::cout << "========================================="<<std::endl
136  << "run/event: "<<runNumber<<'/'<<eventNumber<<std::endl;
137 
138  // HF PIN-diodes
141  if (verbose) std::cout<<"Analysis-> total CAL digis= "<<calib->size()<<std::endl;
142  /* COMMENTED OUT by J. Mans (7-28-2008) as major changes needed with new Calib DetId
143 
144  for (unsigned j = 0; j < calib->size (); ++j) {
145  const HcalCalibDataFrame digi = (*calib)[j];
146  HcalElectronicsId elecId = digi.elecId();
147  HcalCalibDetId calibId = digi.id();
148  if (verbose) std::cout<<calibId.sectorString().c_str()<<" "<<calibId.rbx()<<" "<<elecId.fiberChanId()<<std::endl;
149  int isector = calibId.rbx()-1;
150  int ipin = elecId.fiberChanId();
151  int iside = -1;
152  if (calibId.sectorString() == "HFP") iside = 0;
153  else if (calibId.sectorString() == "HFM") iside = 4;
154 
155  if (iside != -1) {
156  for (int isample = 0; isample < digi.size(); ++isample) {
157  int adc = digi[isample].adc();
158  int capid = digi[isample].capid ();
159  double linear_ADC = digi[isample].nominal_fC();
160  if (verbose) std::cout<<"PIN linear_ADC = "<<linear_ADC<<std::endl;
161  htspin[isector+iside][ipin]->Fill(isample,linear_ADC);
162  }
163  }
164  }
165  */
166  // HF
168  fEvent.getByLabel(hfDigiCollectionTag_, hf_digi);
169  if (verbose) std::cout<<"Analysis-> total HF digis= "<<hf_digi->size()<<std::endl;
170 
171  for (unsigned ihit = 0; ihit < hf_digi->size (); ++ihit) {
172  const HFDataFrame& frame = (*hf_digi)[ihit];
173  HcalDetId detId = frame.id();
174  int ieta = detId.ieta();
175  int iphi = detId.iphi();
176  int depth = detId.depth();
177  if (verbose) std::cout <<"HF digi # " <<ihit<<": eta/phi/depth: "
178  <<ieta<<'/'<<iphi<<'/'<< depth << std::endl;
179 
180  if (ieta>0) ieta = ieta-29;
181  else ieta = 13-ieta-29;
182 
183  for (int isample = 0; isample < frame.size(); ++isample) {
184  int adc = frame[isample].adc();
185  int capid = frame[isample].capid ();
186  double linear_ADC = frame[isample].nominal_fC();
187  double nominal_fC = detId.subdet () == HcalForward ? 2.6 * linear_ADC : linear_ADC;
188 
189  if (verbose) std::cout << "Analysis-> HF sample # " << isample
190  << ", capid=" << capid
191  << ": ADC=" << adc
192  << ", linearized ADC=" << linear_ADC
193  << ", nominal fC=" << nominal_fC << std::endl;
194 
195  hts[ieta][(iphi-1)/2][depth-1]->Fill(isample,linear_ADC);
196  }
197  }
198  return;
199 }
200 
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
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual void analyze(const edm::Event &fEvent, const edm::EventSetup &fSetup)
edm::InputTag hcalCalibDigiCollectionTag_
Definition: HFPreLightCal.h:36
Int_t run_NN
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:32
std::string histfile
Definition: HFPreLightCal.h:27
Int_t event_NN
virtual void endJob(void)
#define NULL
Definition: scimark2.h:8
int ii
Definition: cuy.py:588
TH1F * hts[26][36][2]
Definition: HFPreLightCal.h:31
int depth() const
get the tower depth
Definition: HcalDetId.h:42
MVATrainerComputer * calib
Definition: MVATrainer.cc:64
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
int ieta() const
get the cell ieta
Definition: HcalDetId.h:38
int j
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
TH1F * htspin[8][3]
Definition: HFPreLightCal.h:33
int k[5][pyjets_maxn]
int iphi() const
get the cell iphi
Definition: HcalDetId.h:40
std::string textfile
Definition: HFPreLightCal.h:28
int size() const
total number of samples in the digi
Definition: HFDataFrame.h:26
edm::InputTag hfDigiCollectionTag_
Definition: HFPreLightCal.h:35
int cont
virtual void beginJob()
edm::EventID id() const
Definition: EventBase.h:56
tuple cout
Definition: gather_cfg.py:121
const HcalDetId & id() const
Definition: HFDataFrame.h:22
virtual ~HFPreLightCal()
TH1F * htspinmax
Definition: HFPreLightCal.h:32
HFPreLightCal(const edm::ParameterSet &fConfiguration)