CMS 3D CMS Logo

ESTimingTask.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <fstream>
3 #include <iostream>
4 
17 
18 #include "TMath.h"
19 #include "TGraph.h"
20 
21 using namespace cms;
22 using namespace edm;
23 using namespace std;
24 
25 // fit function
26 double fitf(double* x, double* par) {
27  double wc = par[2];
28  double n = par[3]; // n-1 (in fact)
29  double v1 = pow(wc / n * (x[0] - par[1]), n);
30  double v2 = TMath::Exp(n - wc * (x[0] - par[1]));
31  double v = par[0] * v1 * v2;
32 
33  if (x[0] < par[1])
34  v = 0;
35 
36  return v;
37 }
38 
40  digilabel_ = consumes<ESDigiCollection>(ps.getParameter<InputTag>("DigiLabel"));
41  prefixME_ = ps.getUntrackedParameter<string>("prefixME", "EcalPreshower");
42  esgainToken_ = esConsumes();
43  eCount_ = 0;
44 
45  fit_ = new TF1("fitShape", fitf, -200, 200, 4);
46  fit_->SetParameters(50, 10, 0, 0);
47 
48  //Histogram init
49  for (int i = 0; i < 2; ++i)
50  for (int j = 0; j < 2; ++j)
51  hTiming_[i][j] = nullptr;
52 
53  htESP_ = new TH1F("htESP", "Timing ES+", 81, -20.5, 20.5);
54  htESM_ = new TH1F("htESM", "Timing ES-", 81, -20.5, 20.5);
55 }
56 
58  iBooker.setCurrentFolder(prefixME_ + "/ESTimingTask");
59 
60  //Booking Histograms
61  //Notice: Change ESRenderPlugin under DQM/RenderPlugins/src if you change this histogram name.
62  char histo[200];
63  for (int i = 0; i < 2; ++i)
64  for (int j = 0; j < 2; ++j) {
65  int iz = (i == 0) ? 1 : -1;
66  sprintf(histo, "ES Timing Z %d P %d", iz, j + 1);
67  hTiming_[i][j] = iBooker.book1D(histo, histo, 81, -20.5, 20.5);
68  hTiming_[i][j]->setAxisTitle("ES Timing (ns)", 1);
69  }
70 
71  sprintf(histo, "ES 2D Timing");
72  h2DTiming_ = iBooker.book2D(histo, histo, 81, -20.5, 20.5, 81, -20.5, 20.5);
73  h2DTiming_->setAxisTitle("ES- Timing (ns)", 1);
74  h2DTiming_->setAxisTitle("ES+ Timing (ns)", 2);
75 }
76 
78  delete fit_;
79  delete htESP_;
80  delete htESM_;
81 }
82 
83 void ESTimingTask::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
84  set(iSetup);
85 
86  runNum_ = e.id().run();
87  eCount_++;
88 
89  htESP_->Reset();
90  htESM_->Reset();
91 
92  //Digis
93  int zside, plane, ix, iy, is;
94  double adc[3];
95  // double para[10];
96  //double tx[3] = {-5., 20., 45.};
98  if (e.getByToken(digilabel_, digis)) {
99  for (ESDigiCollection::const_iterator digiItr = digis->begin(); digiItr != digis->end(); ++digiItr) {
100  ESDataFrame dataframe = (*digiItr);
101  ESDetId id = dataframe.id();
102 
103  zside = id.zside();
104  plane = id.plane();
105  ix = id.six();
106  iy = id.siy();
107  is = id.strip();
108 
109  //if (zside==1 && plane==1 && ix==15 && iy==6) continue;
110  if (zside == 1 && plane == 1 && ix == 7 && iy == 28)
111  continue;
112  if (zside == 1 && plane == 1 && ix == 24 && iy == 9 && is == 21)
113  continue;
114  if (zside == -1 && plane == 2 && ix == 35 && iy == 17 && is == 23)
115  continue;
116 
117  int i = (zside == 1) ? 0 : 1;
118  int j = plane - 1;
119 
120  for (int k = 0; k < dataframe.size(); ++k)
121  adc[k] = dataframe.sample(k).adc();
122 
123  double status = 0;
124  if (adc[1] < 200)
125  status = 1;
126  if (fabs(adc[0]) > 10)
127  status = 1;
128  if (adc[1] < 0 || adc[2] < 0)
129  status = 1;
130  if (adc[0] > adc[1] || adc[0] > adc[2])
131  status = 1;
132  if (adc[2] > adc[1])
133  status = 1;
134 
135  if (int(status) == 0) {
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)
160  htESP_->Fill(t0);
161  else if (zside == -1)
162  htESM_->Fill(t0);
163  }
164  }
165  } else {
166  LogWarning("ESTimingTask") << "DigiCollection not available";
167  }
168 
169  if (htESP_->GetEntries() > 0 && htESM_->GetEntries() > 0)
170  h2DTiming_->Fill(htESM_->GetMean(), htESP_->GetMean());
171 }
172 
174  const ESGain* gain = &es.getData(esgainToken_);
175 
176  int ESGain = (int)gain->getESGain();
177 
178  if (ESGain == 1) { // LG
179  wc_ = 0.0837264;
180  n_ = 2.016;
181  } else { // HG
182  wc_ = 0.07291;
183  n_ = 1.798;
184  }
185 
186  //cout<<"gain : "<<ESGain<<endl;
187  //cout<<wc_<<" "<<n_<<endl;
188 }
189 
190 //define this as a plug-in
const ESDetId & id() const
Definition: ESDataFrame.h:19
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: ESTimingTask.cc:57
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
int size() const
Definition: ESDataFrame.h:21
Definition: ESGain.h:7
~ESTimingTask() override
Definition: ESTimingTask.cc:77
void set(const edm::EventSetup &es)
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: ESTimingTask.cc:83
int zside(DetId const &)
T getUntrackedParameter(std::string const &, T const &) const
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:16
const ESSample & sample(int i) const
Definition: ESDataFrame.h:24
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double fitf(double *x, double *par)
Definition: ESTimingTask.cc:26
const_iterator end() const
Namespace of DDCMS conversion namespace.
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t ix(uint32_t id)
const_iterator begin() const
The iterator returned can not safely be used across threads.
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
HLT enums.
float x
ALPAKA_FN_ACC ALPAKA_FN_INLINE uint32_t iy(uint32_t id)
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
ESTimingTask(const edm::ParameterSet &ps)
Definition: ESTimingTask.cc:39
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Definition: Run.h:45
uint16_t *__restrict__ uint16_t const *__restrict__ adc
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)