CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTTTrigCalibration Class Reference

#include <DTTTrigCalibration.h>

Inheritance diagram for DTTTrigCalibration:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 Fill the time boxes. More...
 
 DTTTrigCalibration (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob ()
 Fit the time box rising edge and write the resulting ttrig to the DB. More...
 
virtual ~DTTTrigCalibration ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void dumpTTrigMap (const DTTtrig *tTrig) const
 
std::string getOccupancyName (const DTLayerId &slId) const
 
std::string getTBoxName (const DTSuperLayerId &slId) const
 
void plotTTrig (const DTTtrig *tTrig) const
 

Private Attributes

bool checkNoisyChannels
 
bool debug
 
std::string digiLabel
 
bool doSubtractT0
 
bool findTMeanAndSigma
 
double kFactor
 
int maxDigiPerLayer
 
int maxTDCCounts
 
TFile * theFile
 
DTTimeBoxFittertheFitter
 
std::map< DTSuperLayerId, TH1F * > theHistoMap
 
std::map< DTLayerId, TH1F * > theOccupancyMap
 
DTTTrigBaseSynctheSync
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Analyzer class which fills time box plots with SL granularity for t_trig computation, fits the rising edge and write results to DB. The time boxes are written to file.

Author
G. Cerminara - INFN Torino

Definition at line 32 of file DTTTrigCalibration.h.

Constructor & Destructor Documentation

DTTTrigCalibration::DTTTrigCalibration ( const edm::ParameterSet pset)

Constructor.

Definition at line 43 of file DTTTrigCalibration.cc.

References gather_cfg::cout, debug, dtTPAnalyzer_cfg::digiLabel, reco::get(), edm::ParameterSet::getUntrackedParameter(), dtT0Analyzer_cfg::rootFileName, and interactiveExample::theFile.

43  {
44  // Get the debug parameter for verbose output
45  debug = pset.getUntrackedParameter<bool>("debug");
46 
47  // Get the label to retrieve digis from the event
48  digiLabel = pset.getUntrackedParameter<string>("digiLabel");
49 
50  // Switch on/off the DB writing
51  findTMeanAndSigma = pset.getUntrackedParameter<bool>("fitAndWrite", false);
52 
53  // The TDC time-window (ns)
54  maxTDCCounts = 5000 * pset.getUntrackedParameter<int>("tdcRescale", 1);
55  //The maximum number of digis per layer
56  maxDigiPerLayer = pset.getUntrackedParameter<int>("maxDigiPerLayer", 10);
57 
58  // The root file which will contain the histos
59  string rootFileName = pset.getUntrackedParameter<string>("rootFileName");
60  theFile = new TFile(rootFileName.c_str(), "RECREATE");
61  theFile->cd();
62  theFitter = new DTTimeBoxFitter();
63  if(debug)
65 
66  double sigmaFit = pset.getUntrackedParameter<double>("sigmaTTrigFit",10.);
67  theFitter->setFitSigma(sigmaFit);
68 
69  doSubtractT0 = pset.getUntrackedParameter<bool>("doSubtractT0","false");
70  // Get the synchronizer
71  if(doSubtractT0) {
72  theSync = DTTTrigSyncFactory::get()->create(pset.getUntrackedParameter<string>("tTrigMode"),
73  pset.getUntrackedParameter<ParameterSet>("tTrigModeConfig"));
74  } else {
75  theSync = 0;
76  }
77 
78  checkNoisyChannels = pset.getUntrackedParameter<bool>("checkNoisyChannels","false");
79 
80  // the kfactor to be uploaded in the ttrig DB
81  kFactor = pset.getUntrackedParameter<double>("kFactor",-0.7);
82 
83  if(debug)
84  cout << "[DTTTrigCalibration]Constructor called!" << endl;
85 }
T getUntrackedParameter(std::string const &, T const &) const
void setVerbosity(unsigned int lvl)
Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug.
DTTTrigBaseSync * theSync
void setFitSigma(double sigma)
tuple cout
Definition: gather_cfg.py:121
DTTimeBoxFitter * theFitter
T get(const Candidate &c)
Definition: component.h:55
DTTTrigCalibration::~DTTTrigCalibration ( )
virtual

Destructor.

Definition at line 90 of file DTTTrigCalibration.cc.

References gather_cfg::cout, debug, and interactiveExample::theFile.

90  {
91  if(debug)
92  cout << "[DTTTrigCalibration]Destructor called!" << endl;
93 
94  // // Delete all histos
95  // for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
96  // slHisto != theHistoMap.end();
97  // slHisto++) {
98  // delete (*slHisto).second;
99  // }
100 
101  theFile->Close();
102  delete theFitter;
103 }
tuple cout
Definition: gather_cfg.py:121
DTTimeBoxFitter * theFitter

Member Function Documentation

void DTTTrigCalibration::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
virtual

Fill the time boxes.

Perform the real analysis.

Implements edm::EDAnalyzer.

Definition at line 108 of file DTTTrigCalibration.cc.

References DTLayer::chamber(), DTSuperLayerId::chamberId(), gather_cfg::cout, debug, dtTPAnalyzer_cfg::digiLabel, edm::EventSetup::get(), hltrates_dqm_sourceclient-live_cfg::offset, crabStatusFromReport::statusMap, DTLayerId::superlayerId(), and interactiveExample::theFile.

108  {
109 
110  if(debug)
111  cout << "[DTTTrigCalibration] #Event: " << event.id().event() << endl;
112 
113  // Get the digis from the event
115  event.getByLabel(digiLabel, digis);
116 
118  if(checkNoisyChannels) {
119  // Get the map of noisy channels
120  eventSetup.get<DTStatusFlagRcd>().get(statusMap);
121  }
122 
123  if(doSubtractT0)
124  theSync->setES(eventSetup);
125 
126  //The chambers too noisy in this event
127  vector<DTChamberId> badChambers;
128 
129  // Iterate through all digi collections ordered by LayerId
131  for (dtLayerIt = digis->begin();
132  dtLayerIt != digis->end();
133  ++dtLayerIt){
134  // Get the iterators over the digis associated with this LayerId
135  const DTDigiCollection::Range& digiRange = (*dtLayerIt).second;
136 
137  const DTLayerId layerId = (*dtLayerIt).first;
138  const DTSuperLayerId slId = layerId.superlayerId();
139  const DTChamberId chId = slId.chamberId();
140  bool badChamber=false;
141 
142  if(debug)
143  cout<<"----------- Layer "<<layerId<<" -------------"<<endl;
144 
145  //Check if the layer is inside a noisy chamber
146  for(vector<DTChamberId>::const_iterator chamber = badChambers.begin(); chamber != badChambers.end(); ++chamber){
147  if((*chamber) == chId){
148  badChamber=true;
149  break;
150  }
151  }
152  if(badChamber) continue;
153 
154  //Check if the layer has too many digis
155  if((digiRange.second - digiRange.first) > maxDigiPerLayer){
156  if(debug)
157  cout<<"Layer "<<layerId<<"has too many digis ("<<(digiRange.second - digiRange.first)<<")"<<endl;
158  badChambers.push_back(chId);
159  continue;
160  }
161 
162  // Get the histo from the map
163  TH1F *hTBox = theHistoMap[slId];
164  if(hTBox == 0) {
165  // Book the histogram
166  theFile->cd();
167  hTBox = new TH1F(getTBoxName(slId).c_str(), "Time box (ns)", int(0.25*32.0*maxTDCCounts/25.0), 0, maxTDCCounts);
168  if(debug)
169  cout << " New Time Box: " << hTBox->GetName() << endl;
170  theHistoMap[slId] = hTBox;
171  }
172  TH1F *hO = theOccupancyMap[layerId];
173  if(hO == 0) {
174  // Book the histogram
175  theFile->cd();
176  hO = new TH1F(getOccupancyName(layerId).c_str(), "Occupancy", 100, 0, 100);
177  if(debug)
178  cout << " New Time Box: " << hO->GetName() << endl;
179  theOccupancyMap[layerId] = hO;
180  }
181 
182  // Loop over all digis in the given range
183  for (DTDigiCollection::const_iterator digi = digiRange.first;
184  digi != digiRange.second;
185  ++digi) {
186  const DTWireId wireId(layerId, (*digi).wire());
187 
188  // Check for noisy channels and skip them
189  if(checkNoisyChannels) {
190  bool isNoisy = false;
191  bool isFEMasked = false;
192  bool isTDCMasked = false;
193  bool isTrigMask = false;
194  bool isDead = false;
195  bool isNohv = false;
196  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
197  if(isNoisy) {
198  if(debug)
199  cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
200  continue;
201  }
202  }
203  theFile->cd();
204  double offset = 0;
205  if(doSubtractT0) {
206  const DTLayer* layer = 0;//fake
207  const GlobalPoint glPt;//fake
208  offset = theSync->offset(layer, wireId, glPt);
209  }
210  hTBox->Fill((*digi).time()-offset);
211  if(debug) {
212  cout << " Filling Time Box: " << hTBox->GetName() << endl;
213  cout << " offset (ns): " << offset << endl;
214  cout << " time(ns): " << (*digi).time()-offset<< endl;
215  }
216  hO->Fill((*digi).wire());
217  }
218  }
219 }
std::map< DTSuperLayerId, TH1F * > theHistoMap
DTChamberId chamberId() const
Return the corresponding ChamberId.
double offset(const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globalPos)
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
std::string getTBoxName(const DTSuperLayerId &slId) const
virtual void setES(const edm::EventSetup &setup)=0
Pass the Event Setup to the synchronization module at each event.
DTTTrigBaseSync * theSync
std::string getOccupancyName(const DTLayerId &slId) const
const T & get() const
Definition: EventSetup.h:56
std::vector< DTDigi >::const_iterator const_iterator
std::pair< const_iterator, const_iterator > Range
tuple cout
Definition: gather_cfg.py:121
std::map< DTLayerId, TH1F * > theOccupancyMap
void DTTTrigCalibration::dumpTTrigMap ( const DTTtrig tTrig) const
private

Definition at line 302 of file DTTTrigCalibration.cc.

References DTTtrig::begin(), gather_cfg::cout, and DTTtrig::end().

302  {
303  static const double convToNs = 25./32.;
304  for(DTTtrig::const_iterator ttrig = tTrig->begin();
305  ttrig != tTrig->end(); ++ttrig) {
306  cout << "Wh: " << (*ttrig).first.wheelId
307  << " St: " << (*ttrig).first.stationId
308  << " Sc: " << (*ttrig).first.sectorId
309  << " Sl: " << (*ttrig).first.slId
310  << " TTrig mean (ns): " << (*ttrig).second.tTrig * convToNs
311  << " TTrig sigma (ns): " << (*ttrig).second.tTrms * convToNs<< endl;
312  }
313 }
std::vector< std::pair< DTTtrigId, DTTtrigData > >::const_iterator const_iterator
Access methods to data.
Definition: DTTtrig.h:183
const_iterator begin() const
Definition: DTTtrig.cc:354
tuple cout
Definition: gather_cfg.py:121
const_iterator end() const
Definition: DTTtrig.cc:359
void DTTTrigCalibration::endJob ( void  )
virtual

Fit the time box rising edge and write the resulting ttrig to the DB.

Reimplemented from edm::EDAnalyzer.

Definition at line 222 of file DTTTrigCalibration.cc.

References gather_cfg::cout, debug, DTTimeUnits::ns, DTTtrig::set(), interactiveExample::theFile, and DTCalibDBUtils::writeToDB().

222  {
223  if(debug)
224  cout << "[DTTTrigCalibration]Writing histos to file!" << endl;
225 
226  // Write all time boxes to file
227  theFile->cd();
228  for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
229  slHisto != theHistoMap.end();
230  ++slHisto) {
231  (*slHisto).second->Write();
232  }
233  for(map<DTLayerId, TH1F*>::const_iterator slHisto = theOccupancyMap.begin();
234  slHisto != theOccupancyMap.end();
235  ++slHisto) {
236  (*slHisto).second->Write();
237  }
238 
239  if(findTMeanAndSigma) {
240  // Create the object to be written to DB
241  DTTtrig* tTrig = new DTTtrig();
242 
243  // Loop over the map, fit the histos and write the resulting values to the DB
244  for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
245  slHisto != theHistoMap.end();
246  ++slHisto) {
247  pair<double, double> meanAndSigma = theFitter->fitTimeBox((*slHisto).second);
248  tTrig->set((*slHisto).first,
249  meanAndSigma.first,
250  meanAndSigma.second,
251  kFactor,
253 
254  if(debug) {
255  cout << " SL: " << (*slHisto).first
256  << " mean = " << meanAndSigma.first
257  << " sigma = " << meanAndSigma.second << endl;
258  }
259  }
260 
261  // Print the ttrig map
262  dumpTTrigMap(tTrig);
263 
264  // Plot the tTrig
265  plotTTrig(tTrig);
266 
267  if(debug)
268  cout << "[DTTTrigCalibration]Writing ttrig object to DB!" << endl;
269 
270 
271  // FIXME: to be read from cfg?
272  string tTrigRecord = "DTTtrigRcd";
273 
274  // Write the object to DB
275  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
276  }
277 
278 }
std::map< DTSuperLayerId, TH1F * > theHistoMap
int set(int wheelId, int stationId, int sectorId, int slId, float tTrig, float tTrms, float kFact, DTTimeUnits::type unit)
Definition: DTTtrig.cc:248
void plotTTrig(const DTTtrig *tTrig) const
std::pair< double, double > fitTimeBox(TH1F *hTimeBox)
Fit the rising edge of the time box returning mean value and sigma (first and second respectively) ...
tuple cout
Definition: gather_cfg.py:121
std::map< DTLayerId, TH1F * > theOccupancyMap
void dumpTTrigMap(const DTTtrig *tTrig) const
DTTimeBoxFitter * theFitter
static void writeToDB(std::string record, T *payload)
string DTTTrigCalibration::getOccupancyName ( const DTLayerId slId) const
private

Definition at line 292 of file DTTTrigCalibration.cc.

References DTLayerId::layer(), DTChamberId::sector(), DTChamberId::station(), DTSuperLayerId::superlayer(), and DTChamberId::wheel().

292  {
293  string histoName;
294  stringstream theStream;
295  theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
296  << "_SL" << slId.superlayer() << "_L"<< slId.layer() <<"_Occupancy";
297  theStream >> histoName;
298  return histoName;
299 }
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
int superlayer() const
Return the superlayer number (deprecated method name)
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
string DTTTrigCalibration::getTBoxName ( const DTSuperLayerId slId) const
private

Definition at line 283 of file DTTTrigCalibration.cc.

References DTChamberId::sector(), DTChamberId::station(), DTSuperLayerId::superlayer(), and DTChamberId::wheel().

283  {
284  string histoName;
285  stringstream theStream;
286  theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
287  << "_SL" << slId.superlayer() << "_hTimeBox";
288  theStream >> histoName;
289  return histoName;
290 }
int superlayer() const
Return the superlayer number (deprecated method name)
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTTTrigCalibration::plotTTrig ( const DTTtrig tTrig) const
private

Definition at line 316 of file DTTTrigCalibration.cc.

References DTTtrig::begin(), and DTTtrig::end().

316  {
317 
318  TH1F* tTrig_YB1_Se10 = new TH1F("tTrig_YB1_Se10","tTrig YB1_Se10",15,1,16);
319  TH1F* tTrig_YB2_Se10 = new TH1F("tTrig_YB2_Se10","tTrig YB2_Se10",15,1,16);
320  TH1F* tTrig_YB2_Se11 = new TH1F("tTrig_YB2_Se11","tTrig YB2_Se11",12,1,13);
321 
322  static const double convToNs = 25./32.;
323  for(DTTtrig::const_iterator ttrig = tTrig->begin();
324  ttrig != tTrig->end(); ++ttrig) {
325 
326  // avoid to have wired numbers in the plot
327  float tTrigValue=0;
328  float tTrmsValue=0;
329  if ((*ttrig).second.tTrig * convToNs > 0 &&
330  (*ttrig).second.tTrig * convToNs < 32000 ) {
331  tTrigValue = (*ttrig).second.tTrig * convToNs;
332  tTrmsValue = (*ttrig).second.tTrms * convToNs;
333  }
334 
335  int binx;
336  string binLabel;
337  stringstream binLabelStream;
338  if ((*ttrig).first.sectorId != 14) {
339  binx = ((*ttrig).first.stationId-1)*3 + (*ttrig).first.slId;
340  binLabelStream << "MB"<<(*ttrig).first.stationId<<"_SL"<<(*ttrig).first.slId;
341  }
342  else {
343  binx = 12 + (*ttrig).first.slId;
344  binLabelStream << "MB14_SL"<<(*ttrig).first.slId;
345  }
346  binLabelStream >> binLabel;
347 
348  if ((*ttrig).first.wheelId == 2) {
349  if ((*ttrig).first.sectorId == 10 || (*ttrig).first.sectorId == 14) {
350  tTrig_YB2_Se10->Fill( binx,tTrigValue);
351  tTrig_YB2_Se10->SetBinError( binx, tTrmsValue);
352  tTrig_YB2_Se10->GetXaxis()->SetBinLabel(binx,binLabel.c_str());
353  tTrig_YB2_Se10->GetYaxis()->SetTitle("ns");
354  }
355  else {
356  tTrig_YB2_Se11->Fill( binx,tTrigValue);
357  tTrig_YB2_Se11->SetBinError( binx,tTrmsValue);
358  tTrig_YB2_Se11->GetXaxis()->SetBinLabel(binx,binLabel.c_str());
359  tTrig_YB2_Se11->GetYaxis()->SetTitle("ns");
360  }
361  }
362  else {
363  tTrig_YB1_Se10->Fill( binx,tTrigValue);
364  tTrig_YB1_Se10->SetBinError( binx,tTrmsValue);
365  tTrig_YB1_Se10->GetXaxis()->SetBinLabel(binx,binLabel.c_str());
366  tTrig_YB1_Se10->GetYaxis()->SetTitle("ns");
367  }
368  }
369 
370  tTrig_YB1_Se10->Write();
371  tTrig_YB2_Se10->Write();
372  tTrig_YB2_Se11->Write();
373 
374 }
std::vector< std::pair< DTTtrigId, DTTtrigData > >::const_iterator const_iterator
Access methods to data.
Definition: DTTtrig.h:183
const_iterator begin() const
Definition: DTTtrig.cc:354
const_iterator end() const
Definition: DTTtrig.cc:359

Member Data Documentation

bool DTTTrigCalibration::checkNoisyChannels
private

Definition at line 84 of file DTTTrigCalibration.h.

bool DTTTrigCalibration::debug
private
std::string DTTTrigCalibration::digiLabel
private

Definition at line 67 of file DTTTrigCalibration.h.

bool DTTTrigCalibration::doSubtractT0
private

Definition at line 82 of file DTTTrigCalibration.h.

bool DTTTrigCalibration::findTMeanAndSigma
private

Definition at line 86 of file DTTTrigCalibration.h.

double DTTTrigCalibration::kFactor
private

Definition at line 88 of file DTTTrigCalibration.h.

int DTTTrigCalibration::maxDigiPerLayer
private

Definition at line 72 of file DTTTrigCalibration.h.

int DTTTrigCalibration::maxTDCCounts
private

Definition at line 70 of file DTTTrigCalibration.h.

TFile* DTTTrigCalibration::theFile
private

Definition at line 75 of file DTTTrigCalibration.h.

DTTimeBoxFitter* DTTTrigCalibration::theFitter
private

Definition at line 91 of file DTTTrigCalibration.h.

std::map<DTSuperLayerId, TH1F*> DTTTrigCalibration::theHistoMap
private

Definition at line 78 of file DTTTrigCalibration.h.

std::map<DTLayerId, TH1F*> DTTTrigCalibration::theOccupancyMap
private

Definition at line 79 of file DTTTrigCalibration.h.

DTTTrigBaseSync* DTTTrigCalibration::theSync
private

Definition at line 93 of file DTTTrigCalibration.h.