CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ESTimingTask.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <fstream>
3 #include <iostream>
4 
18 
19 #include "TMath.h"
20 #include "TGraph.h"
21 
22 using namespace cms;
23 using namespace edm;
24 using namespace std;
25 
26 // fit function
27 double fitf(double *x, double *par) {
28 
29  double wc = par[2];
30  double n = par[3]; // n-1 (in fact)
31  double v1 = pow(wc/n*(x[0]-par[1]), n);
32  double v2 = TMath::Exp(n-wc*(x[0]-par[1]));
33  double v = par[0]*v1*v2;
34 
35  if (x[0] < par[1]) v = 0;
36 
37  return v;
38 }
39 
41 
42  digilabel_ = consumes<ESDigiCollection>(ps.getParameter<InputTag>("DigiLabel"));
43  prefixME_ = ps.getUntrackedParameter<string>("prefixME", "EcalPreshower");
44 
45  eCount_ = 0;
46 
47  fit_ = new TF1("fitShape", fitf, -200, 200, 4);
48  fit_->SetParameters(50, 10, 0, 0);
49 
50  //Histogram init
51  for (int i = 0; i < 2; ++i)
52  for (int j = 0; j < 2; ++j)
53  hTiming_[i][j] = 0;
54 
55  htESP_ = new TH1F("htESP", "Timing ES+", 81, -20.5, 20.5);
56  htESM_ = new TH1F("htESM", "Timing ES-", 81, -20.5, 20.5);
57 }
58 
59 void
61 {
62  iBooker.setCurrentFolder(prefixME_ + "/ESTimingTask");
63 
64  //Booking Histograms
65  //Notice: Change ESRenderPlugin under DQM/RenderPlugins/src if you change this histogram name.
66  char histo[200];
67  for (int i=0 ; i<2; ++i)
68  for (int j=0 ; j<2; ++j) {
69  int iz = (i==0)? 1:-1;
70  sprintf(histo, "ES Timing Z %d P %d", iz, j+1);
71  hTiming_[i][j] = iBooker.book1D(histo, histo, 81, -20.5, 20.5);
72  hTiming_[i][j]->setAxisTitle("ES Timing (ns)", 1);
73  }
74 
75  sprintf(histo, "ES 2D Timing");
76  h2DTiming_ = iBooker.book2D(histo, histo, 81, -20.5, 20.5, 81, -20.5, 20.5);
77  h2DTiming_->setAxisTitle("ES- Timing (ns)", 1);
78  h2DTiming_->setAxisTitle("ES+ Timing (ns)", 2);
79 }
80 
82  delete fit_;
83  delete htESP_;
84  delete htESM_;
85 }
86 
87 void ESTimingTask::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
88 
89  set(iSetup);
90 
91  runNum_ = e.id().run();
92  eCount_++;
93 
94  htESP_->Reset();
95  htESM_->Reset();
96 
97  //Digis
98  int zside, plane, ix, iy, is;
99  double adc[3];
100  // double para[10];
101  //double tx[3] = {-5., 20., 45.};
103  if ( e.getByToken(digilabel_, digis) ) {
104 
105  for (ESDigiCollection::const_iterator digiItr = digis->begin(); digiItr != digis->end(); ++digiItr) {
106 
107  ESDataFrame dataframe = (*digiItr);
108  ESDetId id = dataframe.id();
109 
110  zside = id.zside();
111  plane = id.plane();
112  ix = id.six();
113  iy = id.siy();
114  is = id.strip();
115 
116  //if (zside==1 && plane==1 && ix==15 && iy==6) continue;
117  if (zside==1 && plane==1 && ix==7 && iy==28) continue;
118  if (zside==1 && plane==1 && ix==24 && iy==9 && is==21) continue;
119  if (zside==-1 && plane==2 && ix==35 && iy==17 && is==23) continue;
120 
121  int i = (zside==1)? 0:1;
122  int j = plane-1;
123 
124  for (int k=0; k<dataframe.size(); ++k)
125  adc[k] = dataframe.sample(k).adc();
126 
127  double status = 0;
128  if (adc[1] < 200) status = 1;
129  if (fabs(adc[0]) > 10) status = 1;
130  if (adc[1] < 0 || adc[2] < 0) status = 1;
131  if (adc[0] > adc[1] || adc[0] > adc[2]) status = 1;
132  if (adc[2] > adc[1]) status = 1;
133 
134  if (int(status) == 0) {
135 
136  double A1 = adc[1];
137  double A2 = adc[2];
138  double DeltaT = 25.;
139  double aaa = (A2 > 0 && A1 > 0) ? log(A2/A1)/n_ : 20.; // if A1=0, t0=20
140  double bbb = wc_/n_*DeltaT;
141  double ccc= exp(aaa+bbb);
142 
143  double t0 = (2.-ccc)/(1.-ccc) * DeltaT - 5;
144  hTiming_[i][j]->Fill(t0);
145  //cout<<"t0 : "<<t0<<endl;
146  /*
147  TGraph *gr = new TGraph(3, tx, adc);
148  fit_->SetParameters(50, 10, wc_, n_);
149  fit_->FixParameter(2, wc_);
150  fit_->FixParameter(3, n_);
151  fit_->Print();
152  gr->Fit("fitShape", "MQ");
153  fit_->GetParameters(para);
154  delete gr;
155  //hTiming_[i][j]->Fill(para[1]);
156  */
157  //cout<<"ADC : "<<zside<<" "<<plane<<" "<<ix<<" "<<iy<<" "<<is<<" "<<adc[0]<<" "<<adc[1]<<" "<<adc[2]<<" "<<para[1]<<" "<<wc_<<" "<<n_<<endl;
158 
159  if (zside == 1) htESP_->Fill(t0);
160  else if (zside == -1) htESM_->Fill(t0);
161  }
162 
163  }
164  } else {
165  LogWarning("ESTimingTask") << "DigiCollection not available";
166  }
167 
168  if (htESP_->GetEntries() > 0 && htESM_->GetEntries() > 0)
169  h2DTiming_->Fill(htESM_->GetMean(), htESP_->GetMean());
170 
171 }
172 
174 
175 
176  es.get<ESGainRcd>().get(esgain_);
177  const ESGain *gain = esgain_.product();
178 
179  int ESGain = (int) gain->getESGain();
180 
181  if (ESGain == 1) { // LG
182  wc_ = 0.0837264;
183  n_ = 2.016;
184  } else { // HG
185  wc_ = 0.07291;
186  n_ = 1.798;
187  }
188 
189  //cout<<"gain : "<<ESGain<<endl;
190  //cout<<wc_<<" "<<n_<<endl;
191 
192 }
193 
194 //define this as a plug-in
int adc(sample_type sample)
get the ADC sample (12 bits)
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ESTimingTask.cc:60
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
double fitf(double *x, double *par)
Definition: ESTimingTask.cc:27
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
Definition: ESGain.h:7
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const ESDetId & id() const
Definition: ESDataFrame.h:21
void set(const edm::EventSetup &es)
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: ESTimingTask.cc:87
int zside(DetId const &)
int size() const
Definition: ESDataFrame.h:23
virtual ~ESTimingTask()
Definition: ESTimingTask.cc:81
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
const ESSample & sample(int i) const
Definition: ESDataFrame.h:26
int j
Definition: DBlmapReader.cc:9
int zside() const
Definition: ESDetId.h:44
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
float getESGain() const
Definition: ESGain.h:13
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
const T & get() const
Definition: EventSetup.h:55
edm::EventID id() const
Definition: EventBase.h:60
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:18
Definition: DDAxes.h:10
tuple status
Definition: ntuplemaker.py:245
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
ESTimingTask(const edm::ParameterSet &ps)
Definition: ESTimingTask.cc:40
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Definition: Run.h:41