CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EcalPulseShapeGrapher.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EcalPulseShapeGrapher
4 // Class: EcalPulseShapeGrapher
5 //
13 //
14 // Original Author: Seth Cooper
15 // Created: Tue Feb 5 11:35:45 CST 2008
16 //
17 //
18 
20 
21 //
22 // constants, enums and typedefs
23 //
24 
25 //
26 // static data member definitions
27 //
28 
29 //
30 // constructors and destructor
31 //
33  : EBUncalibratedRecHitCollection_(iConfig.getParameter<edm::InputTag>("EBUncalibratedRecHitCollection")),
34  EBDigis_(iConfig.getParameter<edm::InputTag>("EBDigiCollection")),
35  EEUncalibratedRecHitCollection_(iConfig.getParameter<edm::InputTag>("EEUncalibratedRecHitCollection")),
36  EEDigis_(iConfig.getParameter<edm::InputTag>("EEDigiCollection")),
37  ampCut_(iConfig.getUntrackedParameter<int>("AmplitudeCutADC", 13)),
38  rootFilename_(iConfig.getUntrackedParameter<std::string>("rootFilename", "pulseShapeGrapher")) {
39  //now do what ever initialization is needed
40 
41  std::vector<int> listDefaults;
42  listDefaults.push_back(-1);
43  listChannels_ = iConfig.getUntrackedParameter<std::vector<int> >("listChannels", listDefaults);
44 
45  for (std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr) {
46  std::string title = "Amplitude of cry " + intToString(*itr);
47  std::string name = "ampOfCry" + intToString(*itr);
48  ampHistMap_[*itr] = new TH1F(name.c_str(), title.c_str(), 100, 0, 100);
49  ampHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
50 
51  title = "Amplitude (over 13 ADC) of cry " + intToString(*itr);
52  name = "cutAmpOfCry" + intToString(*itr);
53  cutAmpHistMap_[*itr] = new TH1F(name.c_str(), title.c_str(), 100, 0, 100);
54  cutAmpHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
55 
56  title = "Pulse shape of cry " + intToString(*itr);
57  name = "PulseShapeCry" + intToString(*itr);
58  pulseShapeHistMap_[*itr] = new TH2F(name.c_str(), title.c_str(), 10, 0, 10, 220, -20, 2);
59  pulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
60 
61  title = "Raw Pulse shape of cry " + intToString(*itr);
62  name = "RawPulseShapeCry" + intToString(*itr);
63  rawPulseShapeHistMap_[*itr] = new TH2F(name.c_str(), title.c_str(), 10, 0, 10, 500, 0, 500);
64  rawPulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
65 
66  title = "Amplitude of first sample, cry " + intToString(*itr);
67  name = "AmpOfFirstSampleCry" + intToString(*itr);
68  firstSampleHistMap_[*itr] = new TH1F(name.c_str(), title.c_str(), 300, 100, 400);
69  firstSampleHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
70  }
71 
72  fedMap_ = new EcalFedMap();
73 
74  for (int i = 0; i < 10; i++)
75  abscissa[i] = i;
76 }
77 
79  // do anything here that needs to be done at desctruction time
80  // (e.g. close files, deallocate resources etc.)
81 }
82 
83 //
84 // member functions
85 //
86 
87 // ------------ method called to for each event ------------
89  using namespace edm;
90  using namespace std;
91 
92  int numHitsWithActivity = 0;
93 
94  //int eventNum = iEvent.id().event();
95  //vector<EcalUncalibratedRecHit> sampleHitsPastCut;
96 
101  //cout << "event: " << eventNum << " sample hits collection size: " << sampleHits->size() << endl;
102  //Handle<EcalUncalibratedRecHitCollection> fittedHits;
103  //iEvent.getByLabel(FittedUncalibratedRecHitCollection_, fittedHits);
104  Handle<EBDigiCollection> EBdigis;
105  iEvent.getByLabel(EBDigis_, EBdigis);
106  Handle<EEDigiCollection> EEdigis;
107  iEvent.getByLabel(EEDigis_, EEdigis);
108  //cout << "event: " << eventNum << " digi collection size: " << digis->size() << endl;
109 
110  unique_ptr<EcalElectronicsMapping> ecalElectronicsMap(new EcalElectronicsMapping);
111 
112  // Loop over the hits
113  for (EcalUncalibratedRecHitCollection::const_iterator hitItr = EBHits->begin(); hitItr != EBHits->end(); ++hitItr) {
114  EcalUncalibratedRecHit hit = (*hitItr);
115  float amplitude = hit.amplitude();
116  EBDetId hitDetId = hit.id();
117 
118  // Get the Fedid
119  EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
120  int FEDid = 600 + elecId.dccId();
121  string SMname = fedMap_->getSliceFromFed(FEDid);
122 
123  vector<int>::const_iterator itr = listChannels_.begin();
124  while (itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex()) {
125  itr++;
126  }
127  if (itr == listChannels_.end())
128  continue;
129 
130  ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
131  //cout << "Cry hash:" << hitDetId.hashedIndex() << " amplitude: " << amplitude << endl;
132  if (amplitude < ampCut_)
133  continue;
134 
135  cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
136  numHitsWithActivity++;
137  EBDigiCollection::const_iterator digiItr = EBdigis->begin();
138  while (digiItr != EBdigis->end() && digiItr->id() != hitItr->id()) {
139  digiItr++;
140  }
141  if (digiItr == EBdigis->end())
142  continue;
143 
144  double sampleADC[10];
145  EBDataFrame df(*digiItr);
146  double pedestal = 200;
147 
148  if (df.sample(0).gainId() != 1 || df.sample(1).gainId() != 1)
149  continue; //goes to the next digi
150  else {
151  sampleADC[0] = df.sample(0).adc();
152  sampleADC[1] = df.sample(1).adc();
153  pedestal = (double)(sampleADC[0] + sampleADC[1]) / (double)2;
154  }
155 
156  for (int i = 0; (unsigned int)i < digiItr->size(); ++i) {
157  EBDataFrame df(*digiItr);
158  double gain = 12.;
159  if (df.sample(i).gainId() == 1)
160  gain = 1.;
161  else if (df.sample(i).gainId() == 2)
162  gain = 2.;
163  sampleADC[i] = pedestal + (df.sample(i).adc() - pedestal) * gain;
164  }
165 
166  //cout << "1) maxsample amp:" << maxSampleAmp << " maxSampleIndex:" << maxSampleIndex << endl;
167  for (int i = 0; i < 10; ++i) {
168  //cout << "Filling hist for:" << hitDetId.hashedIndex() << " with sample:" << i
169  //<< " amp:" <<(float)(df.sample(i).adc()-baseline)/(maxSampleAmp-baseline) << endl;
170  //cout << "ADC of sample:" << df.sample(i).adc() << " baseline:" << baseline << " maxSampleAmp:"
171  //<< maxSampleAmp << endl << endl;
172  pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i] - pedestal) / amplitude);
173  rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
174  }
175  firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
176  }
177 
178  // Now do the same for the EE hits
179  for (EcalUncalibratedRecHitCollection::const_iterator hitItr = EEHits->begin(); hitItr != EEHits->end(); ++hitItr) {
180  EcalUncalibratedRecHit hit = (*hitItr);
181  float amplitude = hit.amplitude();
182  EEDetId hitDetId = hit.id();
183 
184  // Get the Fedid
185  EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
186  int FEDid = 600 + elecId.dccId();
187  string SMname = fedMap_->getSliceFromFed(FEDid);
188 
189  vector<int>::const_iterator itr = listChannels_.begin();
190  while (itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex()) {
191  itr++;
192  }
193  if (itr == listChannels_.end())
194  continue;
195 
196  ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
197  //cout << "Cry hash:" << hitDetId.hashedIndex() << " amplitude: " << amplitude << endl;
198  if (amplitude < ampCut_)
199  continue;
200 
201  cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
202  numHitsWithActivity++;
203  EEDigiCollection::const_iterator digiItr = EEdigis->begin();
204  while (digiItr != EEdigis->end() && digiItr->id() != hitItr->id()) {
205  digiItr++;
206  }
207  if (digiItr == EEdigis->end())
208  continue;
209 
210  double sampleADC[10];
211  EEDataFrame df(*digiItr);
212  double pedestal = 200;
213 
214  if (df.sample(0).gainId() != 1 || df.sample(1).gainId() != 1)
215  continue; //goes to the next digi
216  else {
217  sampleADC[0] = df.sample(0).adc();
218  sampleADC[1] = df.sample(1).adc();
219  pedestal = (double)(sampleADC[0] + sampleADC[1]) / (double)2;
220  }
221 
222  for (int i = 0; (unsigned int)i < digiItr->size(); ++i) {
223  EEDataFrame df(*digiItr);
224  double gain = 12.;
225  if (df.sample(i).gainId() == 1)
226  gain = 1.;
227  else if (df.sample(i).gainId() == 2)
228  gain = 2.;
229  sampleADC[i] = pedestal + (df.sample(i).adc() - pedestal) * gain;
230  }
231 
232  //cout << "1) maxsample amp:" << maxSampleAmp << " maxSampleIndex:" << maxSampleIndex << endl;
233  for (int i = 0; i < 10; ++i) {
234  //cout << "Filling hist for:" << hitDetId.hashedIndex() << " with sample:" << i
235  //<< " amp:" <<(float)(df.sample(i).adc()-baseline)/(maxSampleAmp-baseline) << endl;
236  //cout << "ADC of sample:" << df.sample(i).adc() << " baseline:" << baseline << " maxSampleAmp:"
237  //<< maxSampleAmp << endl << endl;
238  pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i] - pedestal) / amplitude);
239  rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
240  }
241  firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
242  }
243 }
244 
245 // ------------ method called once each job just after ending the event loop ------------
247  rootFilename_ += ".root";
248  file_ = new TFile(rootFilename_.c_str(), "RECREATE");
249  TH1::AddDirectory(false);
250 
251  for (std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr) {
252  ampHistMap_[*itr]->Write();
253  cutAmpHistMap_[*itr]->Write();
254  firstSampleHistMap_[*itr]->Write();
255 
256  rawPulseShapeHistMap_[*itr]->Write();
257  TProfile* t2 = (TProfile*)(rawPulseShapeHistMap_[*itr]->ProfileX());
258  t2->Write();
259  //TODO: fix the normalization so these are correct
260  //pulseShapeHistMap_[*itr]->Write();
261  //TProfile* t1 = (TProfile*) (pulseShapeHistMap_[*itr]->ProfileX());
262  //t1->Write();
263  }
264 
265  file_->Write();
266  file_->Close();
267 }
268 
270  using namespace std;
271  ostringstream myStream;
272  myStream << num << flush;
273  return (myStream.str()); //returns the string form of the stringstream object
274 }
T getUntrackedParameter(std::string const &, T const &) const
std::map< int, TH1F * > ampHistMap_
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
Ecal readout channel identification [32:20] Unused (so far) [19:13] DCC id [12:6] tower [5:3] strip [...
edm::InputTag EBUncalibratedRecHitCollection_
std::vector< T >::const_iterator const_iterator
std::map< int, TH1F * > cutAmpHistMap_
EcalMGPASample sample(int i) const
Definition: EcalDataFrame.h:29
edm::InputTag EEUncalibratedRecHitCollection_
int gainId() const
get the gainId (2 bits)
EcalPulseShapeGrapher(const edm::ParameterSet &)
std::string getSliceFromFed(int)
Definition: EcalFedMap.cc:86
std::vector< int > listChannels_
int iEvent
Definition: GenABIO.cc:224
int dccId() const
get the DCC (Ecal Local DCC value not global one) id
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
std::map< int, TH1F * > firstSampleHistMap_
int hashedIndex() const
Definition: EEDetId.h:183
std::map< int, TH2F * > pulseShapeHistMap_
std::map< int, TH2F * > rawPulseShapeHistMap_
boost::transform_iterator< IterHelp, boost::counting_iterator< int > > const_iterator
void analyze(const edm::Event &, const edm::EventSetup &) override
tuple size
Write out results.
int adc() const
get the ADC sample (12 bits)