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

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
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

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
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

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.

Date:
2007/07/11 12:21:00
Revision:
1.1
Author
G. Cerminara - INFN Torino

Definition at line 34 of file DTTTrigCalibration.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 45 of file DTTTrigCalibration.cc.

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

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

Destructor.

Definition at line 92 of file DTTTrigCalibration.cc.

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

92  {
93  if(debug)
94  cout << "[DTTTrigCalibration]Destructor called!" << endl;
95 
96  // // Delete all histos
97  // for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
98  // slHisto != theHistoMap.end();
99  // slHisto++) {
100  // delete (*slHisto).second;
101  // }
102 
103  theFile->Close();
104  delete theFitter;
105 }
tuple cout
Definition: gather_cfg.py:41
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 110 of file DTTTrigCalibration.cc.

References DTLayer::chamber(), DTSuperLayerId::chamberId(), gather_cfg::cout, debug, dtT0WireCalibration_cfg::digiLabel, edm::EventSetup::get(), evf::evtn::offset(), crabStatusFromReport::statusMap, DTLayerId::superlayerId(), and interactiveExample::theFile.

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

Definition at line 304 of file DTTTrigCalibration.cc.

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

304  {
305  static const double convToNs = 25./32.;
306  for(DTTtrig::const_iterator ttrig = tTrig->begin();
307  ttrig != tTrig->end(); ttrig++) {
308  cout << "Wh: " << (*ttrig).first.wheelId
309  << " St: " << (*ttrig).first.stationId
310  << " Sc: " << (*ttrig).first.sectorId
311  << " Sl: " << (*ttrig).first.slId
312  << " TTrig mean (ns): " << (*ttrig).second.tTrig * convToNs
313  << " TTrig sigma (ns): " << (*ttrig).second.tTrms * convToNs<< endl;
314  }
315 }
std::vector< std::pair< DTTtrigId, DTTtrigData > >::const_iterator const_iterator
Access methods to data.
Definition: DTTtrig.h:176
const_iterator begin() const
Definition: DTTtrig.cc:377
tuple cout
Definition: gather_cfg.py:41
const_iterator end() const
Definition: DTTtrig.cc:382
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 224 of file DTTTrigCalibration.cc.

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

224  {
225  if(debug)
226  cout << "[DTTTrigCalibration]Writing histos to file!" << endl;
227 
228  // Write all time boxes to file
229  theFile->cd();
230  for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
231  slHisto != theHistoMap.end();
232  slHisto++) {
233  (*slHisto).second->Write();
234  }
235  for(map<DTLayerId, TH1F*>::const_iterator slHisto = theOccupancyMap.begin();
236  slHisto != theOccupancyMap.end();
237  slHisto++) {
238  (*slHisto).second->Write();
239  }
240 
241  if(findTMeanAndSigma) {
242  // Create the object to be written to DB
243  DTTtrig* tTrig = new DTTtrig();
244 
245  // Loop over the map, fit the histos and write the resulting values to the DB
246  for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
247  slHisto != theHistoMap.end();
248  slHisto++) {
249  pair<double, double> meanAndSigma = theFitter->fitTimeBox((*slHisto).second);
250  tTrig->set((*slHisto).first,
251  meanAndSigma.first,
252  meanAndSigma.second,
253  kFactor,
255 
256  if(debug) {
257  cout << " SL: " << (*slHisto).first
258  << " mean = " << meanAndSigma.first
259  << " sigma = " << meanAndSigma.second << endl;
260  }
261  }
262 
263  // Print the ttrig map
264  dumpTTrigMap(tTrig);
265 
266  // Plot the tTrig
267  plotTTrig(tTrig);
268 
269  if(debug)
270  cout << "[DTTTrigCalibration]Writing ttrig object to DB!" << endl;
271 
272 
273  // FIXME: to be read from cfg?
274  string tTrigRecord = "DTTtrigRcd";
275 
276  // Write the object to DB
277  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
278  }
279 
280 }
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:262
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:41
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 294 of file DTTTrigCalibration.cc.

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

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

Definition at line 285 of file DTTTrigCalibration.cc.

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

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

Definition at line 318 of file DTTTrigCalibration.cc.

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

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

Member Data Documentation

bool DTTTrigCalibration::checkNoisyChannels
private

Definition at line 86 of file DTTTrigCalibration.h.

bool DTTTrigCalibration::debug
private

Definition at line 66 of file DTTTrigCalibration.h.

std::string DTTTrigCalibration::digiLabel
private

Definition at line 69 of file DTTTrigCalibration.h.

bool DTTTrigCalibration::doSubtractT0
private

Definition at line 84 of file DTTTrigCalibration.h.

bool DTTTrigCalibration::findTMeanAndSigma
private

Definition at line 88 of file DTTTrigCalibration.h.

double DTTTrigCalibration::kFactor
private

Definition at line 90 of file DTTTrigCalibration.h.

int DTTTrigCalibration::maxDigiPerLayer
private

Definition at line 74 of file DTTTrigCalibration.h.

int DTTTrigCalibration::maxTDCCounts
private

Definition at line 72 of file DTTTrigCalibration.h.

TFile* DTTTrigCalibration::theFile
private

Definition at line 77 of file DTTTrigCalibration.h.

DTTimeBoxFitter* DTTTrigCalibration::theFitter
private

Definition at line 93 of file DTTTrigCalibration.h.

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

Definition at line 80 of file DTTTrigCalibration.h.

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

Definition at line 81 of file DTTTrigCalibration.h.

DTTTrigBaseSync* DTTTrigCalibration::theSync
private

Definition at line 95 of file DTTTrigCalibration.h.