CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
ESTimingTask Class Reference

#include <ESTimingTask.h>

Inheritance diagram for ESTimingTask:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 ESTimingTask (const edm::ParameterSet &ps)
 
 ~ESTimingTask () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void set (const edm::EventSetup &es)
 

Private Attributes

edm::EDGetTokenT< ESDigiCollectiondigilabel_
 
int eCount_
 
edm::ESHandle< ESGainesgain_
 
TF1 * fit_
 
MonitorElementh2DTiming_
 
TH1F * htESM_
 
TH1F * htESP_
 
MonitorElementhTiming_ [2][2]
 
Double_t n_
 
std::string prefixME_
 
int runNum_
 
Double_t wc_
 

Detailed Description

Definition at line 22 of file ESTimingTask.h.

Constructor & Destructor Documentation

ESTimingTask::ESTimingTask ( const edm::ParameterSet ps)

Definition at line 40 of file ESTimingTask.cc.

References fitf(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and mps_fire::i.

40  {
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] = nullptr;
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 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
TH1F * htESM_
Definition: ESTimingTask.h:46
std::string prefixME_
Definition: ESTimingTask.h:37
TH1F * htESP_
Definition: ESTimingTask.h:45
edm::EDGetTokenT< ESDigiCollection > digilabel_
Definition: ESTimingTask.h:36
double fitf(double *x, double *par)
Definition: ESTimingTask.cc:27
MonitorElement * hTiming_[2][2]
Definition: ESTimingTask.h:39
ESTimingTask::~ESTimingTask ( )
override

Definition at line 81 of file ESTimingTask.cc.

81  {
82  delete fit_;
83  delete htESP_;
84  delete htESM_;
85 }
TH1F * htESM_
Definition: ESTimingTask.h:46
TH1F * htESP_
Definition: ESTimingTask.h:45

Member Function Documentation

void ESTimingTask::analyze ( const edm::Event e,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 87 of file ESTimingTask.cc.

References ecalMGPA::adc(), ESSample::adc(), edm::DataFrameContainer::begin(), edm::DataFrameContainer::end(), JetChargeProducer_cfi::exp, edm::Event::getByToken(), mps_fire::i, ESDataFrame::id(), edm::EventBase::id(), gen::k, cmsBatch::log, edm::EventID::run(), ESDataFrame::sample(), ESDataFrame::size(), mps_update::status, genVertex_cff::t0, ESDetId::zside(), and ecaldqm::zside().

87  {
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 }
RunNumber_t run() const
Definition: EventID.h:39
TH1F * htESM_
Definition: ESTimingTask.h:46
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const ESDetId & id() const
Definition: ESDataFrame.h:21
const_iterator begin() const
int zside(DetId const &)
int size() const
Definition: ESDataFrame.h:23
TH1F * htESP_
Definition: ESTimingTask.h:45
void Fill(long long x)
const ESSample & sample(int i) const
Definition: ESDataFrame.h:26
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
int zside() const
Definition: ESDetId.h:39
edm::EDGetTokenT< ESDigiCollection > digilabel_
Definition: ESTimingTask.h:36
int k[5][pyjets_maxn]
MonitorElement * h2DTiming_
Definition: ESTimingTask.h:40
const_iterator end() const
edm::EventID id() const
Definition: EventBase.h:59
Double_t wc_
Definition: ESTimingTask.h:49
MonitorElement * hTiming_[2][2]
Definition: ESTimingTask.h:39
Double_t n_
Definition: ESTimingTask.h:49
int adc() const
get the ADC sample (singed 16 bits)
Definition: ESSample.h:18
void ESTimingTask::bookHistograms ( DQMStore::IBooker iBooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprivate

Definition at line 60 of file ESTimingTask.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), trackerHits::histo, mps_fire::i, MonitorElement::setAxisTitle(), and DQMStore::IBooker::setCurrentFolder().

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 }
std::string prefixME_
Definition: ESTimingTask.h:37
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * h2DTiming_
Definition: ESTimingTask.h:40
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
MonitorElement * hTiming_[2][2]
Definition: ESTimingTask.h:39
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void ESTimingTask::set ( const edm::EventSetup es)
private

Definition at line 173 of file ESTimingTask.cc.

References DEFINE_FWK_MODULE, muonCSCDigis_cfi::gain, edm::EventSetup::get(), ESGain::getESGain(), and createfilelist::int.

173  {
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 }
edm::ESHandle< ESGain > esgain_
Definition: ESTimingTask.h:42
Definition: ESGain.h:7
float getESGain() const
Definition: ESGain.h:13
Double_t wc_
Definition: ESTimingTask.h:49
Double_t n_
Definition: ESTimingTask.h:49
T get() const
Definition: EventSetup.h:71
T const * product() const
Definition: ESHandle.h:86

Member Data Documentation

edm::EDGetTokenT<ESDigiCollection> ESTimingTask::digilabel_
private

Definition at line 36 of file ESTimingTask.h.

int ESTimingTask::eCount_
private

Definition at line 48 of file ESTimingTask.h.

edm::ESHandle<ESGain> ESTimingTask::esgain_
private

Definition at line 42 of file ESTimingTask.h.

TF1* ESTimingTask::fit_
private

Definition at line 44 of file ESTimingTask.h.

MonitorElement* ESTimingTask::h2DTiming_
private

Definition at line 40 of file ESTimingTask.h.

TH1F* ESTimingTask::htESM_
private

Definition at line 46 of file ESTimingTask.h.

TH1F* ESTimingTask::htESP_
private

Definition at line 45 of file ESTimingTask.h.

MonitorElement* ESTimingTask::hTiming_[2][2]
private

Definition at line 39 of file ESTimingTask.h.

Double_t ESTimingTask::n_
private

Definition at line 49 of file ESTimingTask.h.

std::string ESTimingTask::prefixME_
private

Definition at line 37 of file ESTimingTask.h.

int ESTimingTask::runNum_
private

Definition at line 48 of file ESTimingTask.h.

Double_t ESTimingTask::wc_
private

Definition at line 49 of file ESTimingTask.h.