CMS 3D CMS Logo

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) override
 Fill the time boxes. More...
 
 DTTTrigCalibration (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob () override
 Fit the time box rising edge and write the resulting ttrig to the DB. More...
 
 ~DTTTrigCalibration () override
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

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

Destructor.

Definition at line 91 of file DTTTrigCalibration.cc.

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

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

Member Function Documentation

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

Fill the time boxes.

Perform the real analysis.

Definition at line 109 of file DTTTrigCalibration.cc.

References DTStatusFlag::cellStatus(), relativeConstraints::chamber, DTSuperLayerId::chamberId(), gather_cfg::cout, debug, CastorSimpleReconstructor_cfi::digiLabel, edm::EventSetup::get(), PFRecoTauDiscriminationByIsolation_cfi::offset, DTLayerId::superlayerId(), and interactiveExample::theFile.

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

Definition at line 303 of file DTTTrigCalibration.cc.

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

303  {
304  static const double convToNs = 25./32.;
305  for(DTTtrig::const_iterator ttrig = tTrig->begin();
306  ttrig != tTrig->end(); ++ttrig) {
307  cout << "Wh: " << (*ttrig).first.wheelId
308  << " St: " << (*ttrig).first.stationId
309  << " Sc: " << (*ttrig).first.sectorId
310  << " Sl: " << (*ttrig).first.slId
311  << " TTrig mean (ns): " << (*ttrig).second.tTrig * convToNs
312  << " TTrig sigma (ns): " << (*ttrig).second.tTrms * convToNs<< endl;
313  }
314 }
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
void DTTTrigCalibration::endJob ( void  )
overridevirtual

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

Reimplemented from edm::EDAnalyzer.

Definition at line 223 of file DTTTrigCalibration.cc.

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

Referenced by o2olib.O2ORunMgr::executeJob().

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

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

293  {
294  string histoName;
295  stringstream theStream;
296  theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
297  << "_SL" << slId.superlayer() << "_L"<< slId.layer() <<"_Occupancy";
298  theStream >> histoName;
299  return histoName;
300 }
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 284 of file DTTTrigCalibration.cc.

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

284  {
285  string histoName;
286  stringstream theStream;
287  theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
288  << "_SL" << slId.superlayer() << "_hTimeBox";
289  theStream >> histoName;
290  return histoName;
291 }
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 317 of file DTTTrigCalibration.cc.

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

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