CMS 3D CMS Logo

DTT0Calibration.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * $Date: 2012/05/11 17:17:17 $
5  * $Revision: 1.6 $
6  * \author S. Bolognesi - INFN Torino
7  * 06/08/2008 Mofified by Antonio.Vilela.Pereira@cern.ch
8  */
9 
12 
15 
18 
21 
24 
25 #include "TKey.h"
26 #include "TF1.h"
27 
28 #include <cassert>
29 
30 using namespace std;
31 using namespace edm;
32 
33 // Constructor
35  : debug(pset.getUntrackedParameter<bool>("debug")),
36  digiToken(consumes<DTDigiCollection>(pset.getUntrackedParameter<string>("digiLabel"))),
37  theFile(pset.getUntrackedParameter<string>("rootFileName", "DTT0PerLayer.root").c_str(), "RECREATE"),
38  nevents(0),
39  eventsForLayerT0(pset.getParameter<unsigned int>("eventsForLayerT0")),
40  eventsForWireT0(pset.getParameter<unsigned int>("eventsForWireT0")),
41  tpPeakWidth(pset.getParameter<double>("tpPeakWidth")),
42  tpPeakWidthPerLayer(pset.getParameter<double>("tpPeakWidthPerLayer")),
43  rejectDigiFromPeak(pset.getParameter<unsigned int>("rejectDigiFromPeak")),
44  hLayerPeaks("hLayerPeaks", "", 3000, 0, 3000),
45  spectrum(20),
46  dtGeomToken_(esConsumes()) {
47  // Get the debug parameter for verbose output
48  if (debug)
49  cout << "[DTT0Calibration]Constructor called!" << endl;
50 
52  pset.getUntrackedParameter<string>("calibWheel", "All"); //FIXME amke a vector of integer instead of a string
53  if (theCalibWheel != "All") {
54  stringstream linestr;
55  int selWheel;
56  linestr << theCalibWheel;
57  linestr >> selWheel;
58  cout << "[DTT0CalibrationPerLayer] chosen wheel " << selWheel << endl;
59  }
60 
61  // Sector/s to calibrate
63  pset.getUntrackedParameter<string>("calibSector", "All"); //FIXME amke a vector of integer instead of a string
64  if (theCalibSector != "All") {
65  stringstream linestr;
66  int selSector;
67  linestr << theCalibSector;
68  linestr >> selSector;
69  cout << "[DTT0CalibrationPerLayer] chosen sector " << selSector << endl;
70  }
71 
72  vector<string> defaultCell;
73  const auto& cellsWithHistos = pset.getUntrackedParameter<vector<string> >("cellsWithHisto", defaultCell);
74  for (const auto& cell : cellsWithHistos) {
75  stringstream linestr;
76  int wheel, sector, station, sl, layer, wire;
77  linestr << cell;
78  linestr >> wheel >> sector >> station >> sl >> layer >> wire;
80  }
81 }
82 
83 // Destructor
85  if (debug)
86  cout << "[DTT0Calibration]Destructor called!" << endl;
87 
88  theFile.Close();
89 }
90 
93  nevents++;
94 
95  // Get the digis from the event
97  event.getByToken(digiToken, digis);
98 
99  // Get the DT Geometry
100  if (nevents == 1)
101  dtGeom = eventSetup.getHandle(dtGeomToken_);
102 
103  // Iterate through all digi collections ordered by LayerId
104  for (const auto& digis_per_layer : *digis) {
105  //std::cout << __LINE__ << std::endl;
106  // Get the iterators over the digis associated with this LayerId
107  const DTDigiCollection::Range& digiRange = digis_per_layer.second;
108 
109  // Get the layerId
110  const DTLayerId layerId = digis_per_layer.first;
111  //const DTChamberId chamberId = layerId.superlayerId().chamberId();
112 
113  if ((theCalibWheel != "All") && (layerId.superlayerId().chamberId().wheel() != selWheel))
114  continue;
115  if ((theCalibSector != "All") && (layerId.superlayerId().chamberId().sector() != selSector))
116  continue;
117 
118  // Loop over all digis in the given layer
119  for (DTDigiCollection::const_iterator digi = digiRange.first; digi != digiRange.second; ++digi) {
120  const double t0 = digi->countsTDC();
121  const DTWireId wireIdtmp(layerId, (*digi).wire());
122 
123  // Use first bunch of events to fill t0 per layer
124  if (nevents <= eventsForLayerT0) {
125  // If it doesn't exist, book it
126  if (not theHistoLayerMap.count(layerId)) {
127  theHistoLayerMap[layerId] = TH1I(getHistoName(layerId).c_str(),
128  "T0 from pulses by layer (TDC counts, 1 TDC count = 0.781 ns)",
129  3000,
130  0,
131  3000);
132  if (debug)
133  cout << " New T0 per Layer Histo: " << theHistoLayerMap[layerId].GetName() << endl;
134  }
135  theHistoLayerMap[layerId].Fill(t0);
136  }
137 
138  // Use all the remaining events to compute t0 per wire
139  if (nevents > eventsForLayerT0) {
140  // Get the wireId
141  const DTWireId wireId(layerId, (*digi).wire());
142  if (debug) {
143  cout << " Wire: " << wireId << endl << " time (TDC counts): " << (*digi).countsTDC() << endl;
144  }
145 
146  //Fill the histos per wire for the chosen cells
147  if (std::find(layerIdWithWireHistos.begin(), layerIdWithWireHistos.end(), layerId) !=
149  std::find(wireIdWithHistos.begin(), wireIdWithHistos.end(), wireId) != wireIdWithHistos.end()) {
150  //If it doesn't exist, book it
151  if (theHistoWireMap.count(wireId) == 0) {
152  theHistoWireMap[wireId] = TH1I(getHistoName(wireId).c_str(),
153  "T0 from pulses by wire (TDC counts, 1 TDC count = 0.781 ns)",
154  7000,
155  0,
156  7000);
157  if (debug)
158  cout << " New T0 per wire Histo: " << theHistoWireMap[wireId].GetName() << endl;
159  }
160  theHistoWireMap[wireId].Fill(t0);
161  }
162 
163  //Select per layer
164  if (fabs(theTPPeakMap[layerId] - t0) > rejectDigiFromPeak) {
165  if (debug)
166  cout << "digi skipped because t0 too far from peak " << theTPPeakMap[layerId] << endl;
167  continue;
168  }
169 
170  //Use second bunch of events to compute a t0 reference per wire
172  if (!nDigiPerWire_ref[wireId]) {
173  mK_ref[wireId] = 0;
174  }
175  nDigiPerWire_ref[wireId] = nDigiPerWire_ref[wireId] + 1;
176  mK_ref[wireId] = mK_ref[wireId] + (t0 - mK_ref[wireId]) / nDigiPerWire_ref[wireId];
177  }
178  //Use last all the remaining events to compute the mean and sigma t0 per wire
179  else if (nevents > (eventsForLayerT0 + eventsForWireT0)) {
180  if (abs(t0 - mK_ref[wireId]) > tpPeakWidth)
181  continue;
182  if (!nDigiPerWire[wireId]) {
183  theAbsoluteT0PerWire[wireId] = 0;
184  qK[wireId] = 0;
185  mK[wireId] = 0;
186  }
187  nDigiPerWire[wireId] = nDigiPerWire[wireId] + 1;
188  theAbsoluteT0PerWire[wireId] = theAbsoluteT0PerWire[wireId] + t0;
189  qK[wireId] =
190  qK[wireId] + ((nDigiPerWire[wireId] - 1) * (t0 - mK[wireId]) * (t0 - mK[wireId]) / nDigiPerWire[wireId]);
191  mK[wireId] = mK[wireId] + (t0 - mK[wireId]) / nDigiPerWire[wireId];
192  }
193  } //end if(nevents>1000)
194  } //end loop on digi
195  } //end loop on layer
196 
197  //Use the t0 per layer histos to have an indication about the t0 position
198  if (nevents == eventsForLayerT0) {
199  for (const auto& lHisto : theHistoLayerMap) {
200  const auto& layerId = lHisto.first;
201  const auto& hist = lHisto.second;
202  if (debug)
203  cout << "Reading histogram " << hist.GetName() << " with mean " << hist.GetMean() << " and RMS "
204  << hist.GetRMS() << endl;
205 
206  //Find peaks
207  int npeaks = spectrum.Search(&hist, (tpPeakWidthPerLayer / 2.), "", 0.3);
208 
209  double* peaks = spectrum.GetPositionX();
210  //Put in a std::vector<float>
211  vector<double> peakMeans(peaks, peaks + npeaks);
212  //Sort the peaks in ascending order
213  sort(peakMeans.begin(), peakMeans.end());
214 
215  if (peakMeans.empty()) {
216  theTPPeakMap[layerId] = hist.GetMaximumBin();
217  std::cout << "No peaks found by peakfinder in layer " << layerId << ". Taking maximum bin at "
218  << theTPPeakMap[layerId] << ". Please check!" << std::endl;
219  layerIdWithWireHistos.push_back(layerId);
220  } else if (fabs(hist.GetXaxis()->FindBin(peakMeans.front()) - hist.GetXaxis()->FindBin(peakMeans.back())) <
222  theTPPeakMap[layerId] = peakMeans[peakMeans.size() / 2];
223  } else {
224  bool peak_set = false;
225  for (const auto& peak : peakMeans) {
226  // Skip if at low edge
227  if (peak - tpPeakWidthPerLayer <= 0)
228  continue;
229  // Get integral of peak
230  double sum = 0;
231  for (int ibin = peak - tpPeakWidthPerLayer; ibin < peak + tpPeakWidthPerLayer; ibin++) {
232  sum += hist.GetBinContent(ibin);
233  }
234  // Skip if peak too small
235  if (sum < hist.GetMaximum() / 2)
236  continue;
237 
238  // Passed all cuts
239  theTPPeakMap[layerId] = peak;
240  peak_set = true;
241  break;
242  }
243  if (peak_set) {
244  std::cout << "Peaks to far away from each other in layer " << layerId
245  << ". Maybe cross talk? Taking first good peak at " << theTPPeakMap[layerId] << ". Please check!"
246  << std::endl;
247  layerIdWithWireHistos.push_back(layerId);
248  } else {
249  theTPPeakMap[layerId] = hist.GetMaximumBin();
250  std::cout << "Peaks to far away from each other in layer " << layerId
251  << " and no good peak found. Taking maximum bin at " << theTPPeakMap[layerId] << ". Please check!"
252  << std::endl;
253  layerIdWithWireHistos.push_back(layerId);
254  }
255  }
256  if (peakMeans.size() > 5) {
257  std::cout << "Found more than 5 peaks in layer " << layerId << ". Please check!" << std::endl;
258  if (std::find(layerIdWithWireHistos.begin(), layerIdWithWireHistos.end(), layerId) ==
259  layerIdWithWireHistos.end())
260  layerIdWithWireHistos.push_back(layerId);
261  }
262  // Check for noise
263  int nspikes = 0;
264  for (int ibin = 0; ibin < hist.GetNbinsX(); ibin++) {
265  if (hist.GetBinContent(ibin + 1) > hist.GetMaximum() * 0.001)
266  nspikes++;
267  }
268  if (nspikes > 50) {
269  std::cout << "Found a lot of (>50) small spikes in layer " << layerId
270  << ". Please check if all wires are functioning as expected!" << std::endl;
271  if (std::find(layerIdWithWireHistos.begin(), layerIdWithWireHistos.end(), layerId) ==
272  layerIdWithWireHistos.end())
273  layerIdWithWireHistos.push_back(layerId);
274  }
275  hLayerPeaks.Fill(theTPPeakMap[layerId]);
276  }
277  }
278 }
279 
281  std::cout << "Analyzed " << nevents << " events" << std::endl;
282 
283  DTT0* t0sWRTChamber = new DTT0();
284 
285  if (debug)
286  cout << "[DTT0CalibrationPerLayer]Writing histos to file!" << endl;
287 
288  theFile.cd();
289  //hT0SectorHisto->Write();
290  hLayerPeaks.Write();
291  for (const auto& wHisto : theHistoWireMap) {
292  wHisto.second.Write();
293  }
294  for (const auto& lHisto : theHistoLayerMap) {
295  lHisto.second.Write();
296  }
297 
298  if (debug)
299  cout << "[DTT0Calibration] Compute and store t0 and sigma per wire" << endl;
300 
301  // Calculate uncertainties per wire (counting experiment)
302  for (auto& wiret0 : theAbsoluteT0PerWire) {
303  auto& wireId = wiret0.first;
304  if (nDigiPerWire[wireId] > 1)
305  theSigmaT0PerWire[wireId] = qK[wireId] / (nDigiPerWire[wireId] - 1);
306  else
307  theSigmaT0PerWire[wireId] = 999.; // Only one measurement: uncertainty -> infinity
308  // syst uncert
309  //theSigmaT0PerWire[wireId] += pow(0.5, 2));
310  // Every time the same measurement. Use Laplace estimator as estimation how propable it is to measure another value due to limited size of sample
311  if (theSigmaT0PerWire[wireId] == 0) {
312  theSigmaT0PerWire[wireId] += pow(1. / (nDigiPerWire[wireId] + 1), 2);
313  }
314  }
315 
316  // function to calculate unweighted means
317  auto unweighted_mean_function = [](const std::list<double>& values, const std::list<double>& sigmas) {
318  double mean = 0;
319  for (auto& value : values) {
320  mean += value;
321  }
322  mean /= values.size();
323 
324  double uncertainty = 0;
325  for (auto& value : values) {
326  uncertainty += pow(value - mean, 2);
327  }
328  uncertainty /= values.size();
329  uncertainty = sqrt(uncertainty);
330  return std::make_pair(mean, uncertainty);
331  };
332 
333  // correct for odd-even effect in each super layer
334  std::map<DTSuperLayerId, std::pair<double, double> > mean_sigma_even;
335  std::map<DTSuperLayerId, std::pair<double, double> > mean_sigma_odd;
336  for (const auto& superlayer : dtGeom->superLayers()) {
337  const auto superlayer_id = superlayer->id();
338  std::list<double> values_even;
339  std::list<double> sigmas_even;
340  std::list<double> values_odd;
341  std::list<double> sigmas_odd;
342 
343  for (const auto& wiret0 : theAbsoluteT0PerWire) {
344  const auto& wireId = wiret0.first;
345  if (wireId.layerId().superlayerId() == superlayer_id) {
346  const auto& t0 = wiret0.second / nDigiPerWire[wireId];
347  if (wireId.layerId().layer() % 2) {
348  values_odd.push_back(t0);
349  sigmas_odd.push_back(sqrt(theSigmaT0PerWire[wireId]));
350  } else {
351  values_even.push_back(t0);
352  sigmas_even.push_back(sqrt(theSigmaT0PerWire[wireId]));
353  }
354  }
355  }
356  // get mean and uncertainty
357  mean_sigma_even.emplace(superlayer_id, unweighted_mean_function(values_even, sigmas_even));
358  mean_sigma_odd.emplace(superlayer_id, unweighted_mean_function(values_odd, sigmas_odd));
359  }
360 
361  // filter outliers
362  for (const auto& superlayer : dtGeom->superLayers()) {
363  const auto superlayer_id = superlayer->id();
364  std::list<double> values_even;
365  std::list<double> sigmas_even;
366  std::list<double> values_odd;
367  std::list<double> sigmas_odd;
368 
369  for (const auto& wiret0 : theAbsoluteT0PerWire) {
370  const auto& wireId = wiret0.first;
371  if (wireId.layerId().superlayerId() == superlayer_id) {
372  const auto& t0 = wiret0.second / nDigiPerWire[wireId];
373  if (wireId.layerId().layer() % 2 and
374  abs(t0 - mean_sigma_odd[superlayer_id].first) < 2 * mean_sigma_odd[superlayer_id].second) {
375  values_odd.push_back(t0);
376  sigmas_odd.push_back(sqrt(theSigmaT0PerWire[wireId]));
377  } else {
378  if (abs(t0 - mean_sigma_even[superlayer_id].first) < 2 * mean_sigma_even[superlayer_id].second) {
379  values_even.push_back(t0);
380  sigmas_even.push_back(sqrt(theSigmaT0PerWire[wireId]));
381  }
382  }
383  }
384  }
385  // get mean and uncertainty
386  mean_sigma_even[superlayer_id] = unweighted_mean_function(values_even, sigmas_even);
387  mean_sigma_odd[superlayer_id] = unweighted_mean_function(values_odd, sigmas_odd);
388  }
389 
390  // apply correction
391  for (auto& wiret0 : theAbsoluteT0PerWire) {
392  const auto& wire_id = wiret0.first;
393  const auto& superlayer_id = wiret0.first.layerId().superlayerId();
394  const auto& layer = wiret0.first.layerId().layer();
395  auto& t0 = wiret0.second;
396  t0 /= nDigiPerWire[wire_id];
397  if (not layer % 2)
398  continue;
399  // t0 is reference. Changing it changes the map
400  t0 += mean_sigma_even[superlayer_id].first - mean_sigma_odd[superlayer_id].first;
401  theSigmaT0PerWire[wire_id] +=
402  pow(mean_sigma_odd[superlayer_id].second, 2) + pow(mean_sigma_even[superlayer_id].second, 2);
403  }
404 
405  // get chamber mean
406  std::map<DTChamberId, std::list<double> > values_per_chamber;
407  std::map<DTChamberId, std::list<double> > sigmas_per_chamber;
408  for (const auto& wire_t0 : theAbsoluteT0PerWire) {
409  const auto& wire_id = wire_t0.first;
410  const auto& chamber_id = wire_id.chamberId();
411  const auto& t0 = wire_t0.second;
412  values_per_chamber[chamber_id].push_back(t0);
413  sigmas_per_chamber[chamber_id].push_back(sqrt(theSigmaT0PerWire[wire_id]));
414  }
415 
416  std::map<DTChamberId, std::pair<double, double> > mean_per_chamber;
417  for (const auto& chamber_mean : values_per_chamber) {
418  const auto& chamber_id = chamber_mean.first;
419  const auto& means = chamber_mean.second;
420  const auto& sigmas = sigmas_per_chamber[chamber_id];
421  mean_per_chamber.emplace(chamber_id, unweighted_mean_function(means, sigmas));
422  }
423 
424  // calculate relative values
425  for (const auto& wire_t0 : theAbsoluteT0PerWire) {
426  const auto& wire_id = wire_t0.first;
427  const auto& chamber_id = wire_id.chamberId();
428  const auto& t0 = wire_t0.second;
429  theRelativeT0PerWire.emplace(wire_id, t0 - mean_per_chamber[chamber_id].first);
430  cout << "[DTT0Calibration] Wire " << wire_id << " has t0 " << theRelativeT0PerWire[wire_id]
431  << " (relative, after even-odd layer corrections) "
432  << " sigma " << sqrt(theSigmaT0PerWire[wire_id]) << endl;
433  }
434 
435  for (const auto& wire_t0 : theRelativeT0PerWire) {
436  const auto& wire_id = wire_t0.first;
437  const auto& t0 = wire_t0.second;
438  t0sWRTChamber->set(wire_id, t0, sqrt(theSigmaT0PerWire[wire_id]), DTTimeUnits::counts);
439  }
440 
442  if (debug)
443  cout << "[DTT0Calibration]Writing values in DB!" << endl;
444  // FIXME: to be read from cfg?
445  string t0Record = "DTT0Rcd";
446  // Write the t0 map to DB
447  DTCalibDBUtils::writeToDB(t0Record, t0sWRTChamber);
448  delete t0sWRTChamber;
449 }
450 
451 string DTT0Calibration::getHistoName(const DTWireId& wId) const {
452  string histoName;
453  stringstream theStream;
454  theStream << "Ch_" << wId.wheel() << "_" << wId.station() << "_" << wId.sector() << "_SL" << wId.superlayer() << "_L"
455  << wId.layer() << "_W" << wId.wire() << "_hT0Histo";
456  theStream >> histoName;
457  return histoName;
458 }
459 
460 string DTT0Calibration::getHistoName(const DTLayerId& lId) const {
461  string histoName;
462  stringstream theStream;
463  theStream << "Ch_" << lId.wheel() << "_" << lId.station() << "_" << lId.sector() << "_SL" << lId.superlayer() << "_L"
464  << lId.layer() << "_hT0Histo";
465  theStream >> histoName;
466  return histoName;
467 }
std::vector< DTLayerId > layerIdWithWireHistos
int set(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, float t0mean, float t0rms, DTTimeUnits::type unit)
Definition: DTT0.cc:97
int station() const
Return the station number.
Definition: DTChamberId.h:45
double tpPeakWidthPerLayer
std::map< DTLayerId, TH1I > theHistoLayerMap
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::map< DTWireId, double > theAbsoluteT0PerWire
static void writeToDB(std::string record, const T &payload)
void endJob() override
Compute the mean and the RMS of the t0 from the maps and write them to the DB with channel granularit...
unsigned int nevents
int wire() const
Return the wire number.
Definition: DTWireId.h:45
std::string getHistoName(const DTWireId &wId) const
std::vector< DTWireId > wireIdWithHistos
std::string theCalibWheel
std::map< DTWireId, double > mK_ref
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::ESHandle< DTGeometry > dtGeom
unsigned int eventsForLayerT0
U second(std::pair< T, U > const &p)
Definition: DTT0.h:48
DTChamberId chamberId() const
Return the corresponding ChamberId.
std::map< DTLayerId, double > theTPPeakMap
std::map< DTWireId, double > theSigmaT0PerWire
unsigned int rejectDigiFromPeak
T sqrt(T t)
Definition: SSEVec.h:23
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
~DTT0Calibration() override
Destructor.
std::string theCalibSector
Definition: value.py:1
edm::EDGetTokenT< DTDigiCollection > digiToken
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
Fill the maps with t0 (by channel)
std::map< DTWireId, int > nDigiPerWire
std::map< DTWireId, double > mK
int superlayer() const
Return the superlayer number (deprecated method name)
#define debug
Definition: HDRShower.cc:19
std::map< DTWireId, int > nDigiPerWire_ref
std::pair< const_iterator, const_iterator > Range
std::vector< DigiType >::const_iterator const_iterator
int layer() const
Return the layer number.
Definition: DTLayerId.h:45
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
HLT enums.
int sector() const
Definition: DTChamberId.h:52
std::map< DTWireId, double > theRelativeT0PerWire
std::map< DTWireId, double > qK
std::map< DTWireId, TH1I > theHistoWireMap
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:48
DTT0Calibration(const edm::ParameterSet &pset)
Constructor.
TSpectrum spectrum
unsigned int eventsForWireT0
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
Definition: event.py:1
const std::vector< const DTSuperLayer * > & superLayers() const
Return a vector of all SuperLayer.
Definition: DTGeometry.cc:86