CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
DTNoiseCalibration Class Reference

#include <DTNoiseCalibration.h>

Inheritance diagram for DTNoiseCalibration:
edm::EDAnalyzer edm::EDConsumerBase

Classes

class  DTNoiseCalibration
 

Public Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c)
 
void beginJob ()
 
void beginRun (const edm::Run &run, const edm::EventSetup &setup)
 
 DTNoiseCalibration (const edm::ParameterSet &ps)
 Constructor. More...
 
void endJob ()
 
virtual ~DTNoiseCalibration ()
 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
 
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
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

std::string getChamberName (const DTChamberId &) const
 
std::string getChannelName (const DTWireId &) const
 
std::string getLayerName (const DTLayerId &) const
 
std::string getSuperLayerName (const DTSuperLayerId &) const
 

Private Attributes

std::map< DTChamberId, TH1F * > chamberOccupancyVsLumiMap_
 
std::map< DTChamberId, TH1F * > chamberOccupancyVsTimeMap_
 
std::string dbLabel_
 
int defaultTtrig_
 
edm::InputTag digiLabel_
 
edm::ESHandle< DTGeometrydtGeom_
 
TH1F * hTDCTriggerWidth_
 
unsigned int lumiMax_
 
double maximumNoiseRate_
 
int nevents_
 
bool readDB_
 
TFile * rootFile_
 
time_t runBeginTime_
 
time_t runEndTime_
 
std::map< DTLayerId, TH1F * > theHistoOccupancyMap_
 
std::map< DTWireId, TH1F * > theHistoOccupancyVsLumiMap_
 
int timeWindowOffset_
 
double triggerWidth_
 
edm::ESHandle< DTTtrigtTrigMap_
 
bool useAbsoluteRate_
 
bool useTimeWindow_
 
std::vector< DTWireIdwireIdWithHisto_
 

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

Definition at line 31 of file DTNoiseCalibration.h.

Constructor & Destructor Documentation

Constructor.

Definition at line 42 of file DTNoiseCalibration.cc.

References defaultTtrig_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), readDB_, rootFile_, DTAnalyzerDetailed_cfi::rootFileName, relativeConstraints::station, makeMuonMisalignmentScenario::wheel, and wireIdWithHisto_.

42  :
43  digiLabel_( pset.getParameter<InputTag>("digiLabel") ),
44  useTimeWindow_( pset.getParameter<bool>("useTimeWindow") ),
45  triggerWidth_( pset.getParameter<double>("triggerWidth") ),
46  timeWindowOffset_( pset.getParameter<int>("timeWindowOffset") ),
47  maximumNoiseRate_( pset.getParameter<double>("maximumNoiseRate") ),
48  useAbsoluteRate_( pset.getParameter<bool>("useAbsoluteRate") ),
49  readDB_(true), defaultTtrig_(0),
50  dbLabel_( pset.getUntrackedParameter<string>("dbLabel", "") ),
51  //fastAnalysis_( pset.getParameter<bool>("fastAnalysis", true) ),
52  wireIdWithHisto_( std::vector<DTWireId>() ),
53  lumiMax_(3000)
54  {
55 
56  // Get the debug parameter for verbose output
57  //debug = ps.getUntrackedParameter<bool>("debug");
58  /*// The analysis type
59  // The wheel & sector interested for the time-dependent analysis
60  wh = ps.getUntrackedParameter<int>("wheel", 0);
61  sect = ps.getUntrackedParameter<int>("sector", 6);*/
62 
63  if( pset.exists("defaultTtrig") ){
64  readDB_ = false;
65  defaultTtrig_ = pset.getParameter<int>("defaultTtrig");
66  }
67 
68  if( pset.exists("cellsWithHisto") ){
69  vector<string> cellsWithHisto = pset.getParameter<vector<string> >("cellsWithHisto");
70  for(vector<string>::const_iterator cell = cellsWithHisto.begin(); cell != cellsWithHisto.end(); ++cell){
71  //FIXME: Use regex to check whether format is right
72  if( (*cell) != "" && (*cell) != "None"){
73  stringstream linestr;
74  int wheel,station,sector,sl,layer,wire;
75  linestr << (*cell);
76  linestr >> wheel >> station >> sector >> sl >> layer >> wire;
77  wireIdWithHisto_.push_back(DTWireId(wheel,station,sector,sl,layer,wire));
78  }
79  }
80  }
81 
82  // The root file which will contain the histos
83  string rootFileName = pset.getUntrackedParameter<string>("rootFileName","noise.root");
84  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
85  rootFile_->cd();
86 }
edm::InputTag digiLabel_
std::vector< DTWireId > wireIdWithHisto_
DTNoiseCalibration::~DTNoiseCalibration ( )
virtual

Destructor.

Definition at line 416 of file DTNoiseCalibration.cc.

References rootFile_.

416  {
417  rootFile_->Close();
418 }

Member Function Documentation

void DTNoiseCalibration::analyze ( const edm::Event e,
const edm::EventSetup c 
)

Definition at line 118 of file DTNoiseCalibration.cc.

References DTSuperLayerId::chamberId(), chamberOccupancyVsLumiMap_, chamberOccupancyVsTimeMap_, DTTimeUnits::counts, defaultTtrig_, digiLabel_, dtGeom_, Exception, spr::find(), DTTopology::firstChannel(), objects.autophobj::float, DTTtrig::get(), getChamberName(), getChannelName(), getLayerName(), hTDCTriggerWidth_, createfilelist::int, gen::k, DTAnalyzerDetailed_cfi::kFactor, DTTopology::lastChannel(), DTGeometry::layer(), LogTrace, lumiMax_, nevents_, readDB_, rootFile_, runBeginTime_, runEndTime_, edm::second(), DTLayer::specificTopology(), mps_update::status, theHistoOccupancyMap_, theHistoOccupancyVsLumiMap_, edm::EventBase::time(), timeWindowOffset_, triggerWidth_, tTrigMap_, useTimeWindow_, edm::Timestamp::value(), and wireIdWithHisto_.

118  {
119  ++nevents_;
120 
121  // Get the digis from the event
122  Handle<DTDigiCollection> dtdigis;
123  event.getByLabel(digiLabel_, dtdigis);
124 
125  /*TH1F *hOccupancyHisto;
126  TH2F *hEvtPerWireH;
127  string Histo2Name;*/
128 
129  //RunNumber_t runNumber = event.id().run();
130  time_t eventTime = time_t(event.time().value()>>32);
131  unsigned int lumiSection = event.luminosityBlock();
132 
133  // Loop over digis
135  for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); ++dtLayerId_It){
136  // Define time window
137  float upperLimit = 0.;
138  if(useTimeWindow_){
139  if( readDB_ ){
140  float tTrig,tTrigRMS,kFactor;
141  DTSuperLayerId slId = ((*dtLayerId_It).first).superlayerId();
142  int status = tTrigMap_->get( slId, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts );
143  if(status != 0) throw cms::Exception("DTNoiseCalibration") << "Could not find tTrig entry in DB for" << slId << endl;
144  upperLimit = tTrig - timeWindowOffset_;
145  } else {
146  upperLimit = defaultTtrig_ - timeWindowOffset_;
147  }
148  }
149 
150  double triggerWidth_s = 0.;
151  if(useTimeWindow_) triggerWidth_s = ( (upperLimit*25)/32 )/1e9;
152  else triggerWidth_s = double(triggerWidth_/1e9);
153  LogTrace("Calibration") << ((*dtLayerId_It).first).superlayerId() << " Trigger width (s): " << triggerWidth_s;
154 
155  for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
156  digiIt!=((*dtLayerId_It).second).second; ++digiIt){
157 
158  //Check the TDC trigger width
159  int tdcTime = (*digiIt).countsTDC();
160  if( !useTimeWindow_ ){
161  if( ( ((float)tdcTime*25)/32 ) > triggerWidth_ ){
162  LogError("Calibration") << "Digi has a TDC time (ns) higher than the pre-defined TDC trigger width: " << ((float)tdcTime*25)/32;
163  continue;
164  }
165  }
166 
167  hTDCTriggerWidth_->Fill(tdcTime);
168 
169  if( useTimeWindow_ && tdcTime > upperLimit) continue;
170 
171  /*LogTrace("Calibration") << "TDC time (ns): " << ((float)tdcTime*25)/32
172  <<" --- trigger width (ns): " << ((float)upperLimit*25)/32;*/
173 
174  const DTLayerId dtLId = (*dtLayerId_It).first;
175  const DTTopology& dtTopo = dtGeom_->layer(dtLId)->specificTopology();
176  const int firstWire = dtTopo.firstChannel();
177  const int lastWire = dtTopo.lastChannel();
178  //const int nWires = dtTopo.channels();
179  const int nWires = lastWire - firstWire + 1;
180 
181  // Book the occupancy histos
182  if( theHistoOccupancyMap_.find(dtLId) == theHistoOccupancyMap_.end() ){
183  string histoName = "DigiOccupancy_" + getLayerName(dtLId);
184  rootFile_->cd();
185  TH1F* hOccupancyHisto = new TH1F(histoName.c_str(), histoName.c_str(), nWires, firstWire, lastWire+1);
186  LogTrace("Calibration") << " Created occupancy Histo: " << hOccupancyHisto->GetName();
187  theHistoOccupancyMap_[dtLId] = hOccupancyHisto;
188  }
189  theHistoOccupancyMap_[dtLId]->Fill( (*digiIt).wire(), 1./triggerWidth_s );
190 
191  // Histos vs lumi section
192  const DTChamberId dtChId = dtLId.chamberId();
193  if( chamberOccupancyVsLumiMap_.find(dtChId) == chamberOccupancyVsLumiMap_.end() ){
194  string histoName = "OccupancyVsLumi_" + getChamberName(dtChId);
195  rootFile_->cd();
196  TH1F* hOccupancyVsLumiHisto = new TH1F(histoName.c_str(), histoName.c_str(), lumiMax_, 0, lumiMax_);
197  LogTrace("Calibration") << " Created occupancy histo: " << hOccupancyVsLumiHisto->GetName();
198  chamberOccupancyVsLumiMap_[dtChId] = hOccupancyVsLumiHisto;
199  }
200  chamberOccupancyVsLumiMap_[dtChId]->Fill( lumiSection, 1./triggerWidth_s );
201 
202  const DTWireId wireId(dtLId, (*digiIt).wire());
203  if( find(wireIdWithHisto_.begin(),wireIdWithHisto_.end(),wireId) != wireIdWithHisto_.end() ){
204  if( theHistoOccupancyVsLumiMap_.find(wireId) == theHistoOccupancyVsLumiMap_.end() ){
205  string histoName = "OccupancyVsLumi_" + getChannelName(wireId);
206  rootFile_->cd();
207  TH1F* hOccupancyVsLumiHisto = new TH1F(histoName.c_str(), histoName.c_str(), lumiMax_, 0, lumiMax_);
208  LogTrace("Calibration") << " Created occupancy histo: " << hOccupancyVsLumiHisto->GetName();
209  theHistoOccupancyVsLumiMap_[wireId] = hOccupancyVsLumiHisto;
210  }
211  theHistoOccupancyVsLumiMap_[wireId]->Fill( lumiSection, 1./triggerWidth_s );
212  }
213 
214  // Histos vs time
215  if( chamberOccupancyVsTimeMap_.find(dtChId) == chamberOccupancyVsTimeMap_.end() ){
216  string histoName = "OccupancyVsTime_" + getChamberName(dtChId);
217  float secPerBin = 20.0;
218  unsigned int nBins = ( (unsigned int)(runEndTime_ - runBeginTime_) )/secPerBin;
219  rootFile_->cd();
220  TH1F* hOccupancyVsTimeHisto = new TH1F(histoName.c_str(), histoName.c_str(),
221  nBins, (unsigned int)runBeginTime_,
222  (unsigned int)runEndTime_);
223  for(int k = 0; k < hOccupancyVsTimeHisto->GetNbinsX(); ++k){
224  if( k%10 == 0 ){
225  unsigned int binLowEdge = hOccupancyVsTimeHisto->GetBinLowEdge(k+1);
226  time_t timeValue = time_t(binLowEdge);
227  hOccupancyVsTimeHisto->GetXaxis()->SetBinLabel( (k+1),ctime(&timeValue) );
228  }
229  }
230  size_t lastBin = hOccupancyVsTimeHisto->GetNbinsX();
231  unsigned int binUpperEdge = hOccupancyVsTimeHisto->GetBinLowEdge(lastBin) +
232  hOccupancyVsTimeHisto->GetBinWidth(lastBin);
233  time_t timeValue = time_t(binUpperEdge);
234  hOccupancyVsTimeHisto->GetXaxis()->SetBinLabel( (lastBin),ctime(&timeValue) );
235 
236  LogTrace("Calibration") << " Created occupancy histo: " << hOccupancyVsTimeHisto->GetName();
237  chamberOccupancyVsTimeMap_[dtChId] = hOccupancyVsTimeHisto;
238  }
239  chamberOccupancyVsTimeMap_[dtChId]->Fill( (unsigned int)eventTime, 1./triggerWidth_s );
240 
241  /*// Book the digi event plot every 1000 events if the analysis is not "fast" and if is the correct sector
242  if(!fastAnalysis &&
243  dtLId.superlayerId().chamberId().wheel()==wh &&
244  dtLId.superlayerId().chamberId().sector()==sect) {
245  if(theHistoEvtPerWireMap.find(dtLId) == theHistoEvtPerWireMap.end() ||
246  (theHistoEvtPerWireMap.find(dtLId) != theHistoEvtPerWireMap.end() &&
247  skippedPlot[dtLId] != counter)){
248  skippedPlot[dtLId] = counter;
249  stringstream toAppend; toAppend << counter;
250  Histo2Name = "DigiPerWirePerEvent_" + getLayerName(dtLId) + "_" + toAppend.str();
251  theFile->cd();
252  hEvtPerWireH = new TH2F(Histo2Name.c_str(), Histo2Name.c_str(), 1000,0.5,1000.5,nWires, firstWire, lastWire+1);
253  if(hEvtPerWireH){
254  if(debug)
255  cout << " New Histo with the number of digi per evt per wire: " << hEvtPerWireH->GetName() << endl;
256  theHistoEvtPerWireMap[dtLId]=hEvtPerWireH;
257  }
258  }
259  }*/
260  }
261  }
262 
263  /*//Fill the plot of the number of digi per event per wire
264  std::map<int,int > DigiPerWirePerEvent;
265  // LOOP OVER ALL THE CHAMBERS
266  vector<DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
267  vector<DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
268  for (; ch_it != ch_end; ++ch_it) {
269  DTChamberId ch = (*ch_it)->id();
270  vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
271  vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
272  // Loop over the SLs
273  for(; sl_it != sl_end; ++sl_it) {
274  DTSuperLayerId sl = (*sl_it)->id();
275  vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin();
276  vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
277  // Loop over the Ls
278  for(; l_it != l_end; ++l_it) {
279  DTLayerId layerId = (*l_it)->id();
280 
281  // Get the number of wires
282  const DTTopology& dtTopo = dtGeom->layer(layerId)->specificTopology();
283  const int firstWire = dtTopo.firstChannel();
284  const int lastWire = dtTopo.lastChannel();
285 
286  if (theHistoEvtPerWireMap.find(layerId) != theHistoEvtPerWireMap.end() &&
287  skippedPlot[layerId] == counter) {
288 
289  for (int wire=firstWire; wire<=lastWire; wire++) {
290  DigiPerWirePerEvent[wire]= 0;
291  }
292  // loop over all the digis of the event
293  DTDigiCollection::Range layerDigi= dtdigis->get(layerId);
294  for (DTDigiCollection::const_iterator digi = layerDigi.first;
295  digi!=layerDigi.second;
296  ++digi){
297  if((cosmicRun && (*digi).countsTDC()<upperLimit) || (!cosmicRun))
298  DigiPerWirePerEvent[(*digi).wire()]+=1;
299  }
300  // fill the digi event histo
301  for (int wire=firstWire; wire<=lastWire; wire++) {
302  theFile->cd();
303  int histoEvents = nevents - (counter*1000);
304  theHistoEvtPerWireMap[layerId]->Fill(histoEvents,wire,DigiPerWirePerEvent[wire]);
305  }
306  }
307  } //Loop Ls
308  } //Loop SLs
309  } //Loop chambers
310 
311 
312  if(nevents % 1000 == 0) {
313  counter++;
314  // save the digis event plot on file
315  for(map<DTLayerId, TH2F* >::const_iterator lHisto = theHistoEvtPerWireMap.begin();
316  lHisto != theHistoEvtPerWireMap.end();
317  lHisto++) {
318  theFile->cd();
319  if((*lHisto).second)
320  (*lHisto).second->Write();
321  }
322  theHistoEvtPerWireMap.clear();
323  }*/
324 
325 }
edm::ESHandle< DTGeometry > dtGeom_
DTChamberId chamberId() const
Return the corresponding ChamberId.
edm::InputTag digiLabel_
std::vector< DTWireId > wireIdWithHisto_
std::map< DTLayerId, TH1F * > theHistoOccupancyMap_
std::string getChannelName(const DTWireId &) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::map< DTChamberId, TH1F * > chamberOccupancyVsTimeMap_
std::string getChamberName(const DTChamberId &) const
std::map< DTChamberId, TH1F * > chamberOccupancyVsLumiMap_
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:78
std::string getLayerName(const DTLayerId &) const
std::map< DTWireId, TH1F * > theHistoOccupancyVsLumiMap_
const DTLayer * layer(DTLayerId id) const
Return a layer given its id.
Definition: DTGeometry.cc:110
U second(std::pair< T, U > const &p)
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:80
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
#define LogTrace(id)
int k[5][pyjets_maxn]
int get(int wheelId, int stationId, int sectorId, int slId, float &tTrig, float &tTrms, float &kFact, DTTimeUnits::type unit) const
get content
Definition: DTTtrig.cc:85
std::vector< DTDigi >::const_iterator const_iterator
edm::ESHandle< DTTtrig > tTrigMap_
Definition: event.py:1
void DTNoiseCalibration::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 88 of file DTNoiseCalibration.cc.

References hTDCTriggerWidth_, nevents_, and triggerWidth_.

88  {
89  LogVerbatim("Calibration") << "[DTNoiseCalibration]: Begin job";
90 
91  nevents_ = 0;
92 
93  TH1::SetDefaultSumw2(true);
94  int numBin = (triggerWidth_*32/25)/50;
95  hTDCTriggerWidth_ = new TH1F("TDC_Time_Distribution", "TDC_Time_Distribution", numBin, 0, triggerWidth_*32/25);
96 }
void DTNoiseCalibration::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)

Definition at line 98 of file DTNoiseCalibration.cc.

References edm::RunBase::beginTime(), dbLabel_, dtGeom_, edm::RunBase::endTime(), edm::EventSetup::get(), readDB_, runBeginTime_, runEndTime_, tTrigMap_, and edm::Timestamp::value().

98  {
99 
100  // Get the DT Geometry
101  setup.get<MuonGeometryRecord>().get(dtGeom_);
102 
103  // tTrig
104  if( readDB_ ) setup.get<DTTtrigRcd>().get(dbLabel_,tTrigMap_);
105 
106  runBeginTime_ = time_t(run.beginTime().value()>>32);
107  runEndTime_ = time_t(run.endTime().value()>>32);
108  /*
109  nevents = 0;
110  counter = 0;
111 
112  // TDC time distribution
113  int numBin = (triggerWidth_*(32/25))/50;
114  hTDCTriggerWidth = new TH1F("TDC_Time_Distribution", "TDC_Time_Distribution", numBin, 0, triggerWidth_*(32/25));*/
115 
116 }
Timestamp const & endTime() const
Definition: RunBase.h:42
edm::ESHandle< DTGeometry > dtGeom_
const T & get() const
Definition: EventSetup.h:56
Timestamp const & beginTime() const
Definition: RunBase.h:41
TimeValue_t value() const
Definition: Timestamp.h:56
edm::ESHandle< DTTtrig > tTrigMap_
void DTNoiseCalibration::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 327 of file DTNoiseCalibration.cc.

References stringResolutionProvider_cfi::bin, chamberOccupancyVsLumiMap_, chamberOccupancyVsTimeMap_, dtGeom_, DTTopology::firstChannel(), hTDCTriggerWidth_, DTTopology::lastChannel(), DTGeometry::layer(), LogTrace, maximumNoiseRate_, nevents_, PostProcessor_cff::normalization, record, rootFile_, DTStatusFlag::setCellNoise(), DTLayer::specificTopology(), theHistoOccupancyMap_, theHistoOccupancyVsLumiMap_, and useAbsoluteRate_.

Referenced by o2olib.O2ORunMgr::executeJob().

327  {
328 
329  //LogVerbatim("Calibration") << "[DTNoiseCalibration] endjob called!";
330  LogVerbatim("Calibration") << "[DTNoiseCalibration] Total number of events analyzed: " << nevents_;
331 
332  // Save the TDC digi plot
333  rootFile_->cd();
334  hTDCTriggerWidth_->Write();
335 
336  double normalization = 1./double(nevents_);
337 
338  for(map<DTWireId, TH1F*>::const_iterator wHisto = theHistoOccupancyVsLumiMap_.begin();
339  wHisto != theHistoOccupancyVsLumiMap_.end(); ++wHisto){
340  (*wHisto).second->Scale(normalization);
341  (*wHisto).second->Write();
342  }
343 
344  for(map<DTChamberId, TH1F*>::const_iterator chHisto = chamberOccupancyVsLumiMap_.begin();
345  chHisto != chamberOccupancyVsLumiMap_.end(); ++chHisto){
346  (*chHisto).second->Scale(normalization);
347  (*chHisto).second->Write();
348  }
349 
350  for(map<DTChamberId, TH1F*>::const_iterator chHisto = chamberOccupancyVsTimeMap_.begin();
351  chHisto != chamberOccupancyVsTimeMap_.end(); ++chHisto){
352  (*chHisto).second->Scale(normalization);
353  (*chHisto).second->Write();
354  }
355 
356  // Save on file the occupancy histos and write the list of noisy cells
357  DTStatusFlag *statusMap = new DTStatusFlag();
358  for(map<DTLayerId, TH1F*>::const_iterator lHisto = theHistoOccupancyMap_.begin();
359  lHisto != theHistoOccupancyMap_.end();
360  ++lHisto){
361  /*double triggerWidth_s = 0.;
362  if( useTimeWindow_ ){
363  double triggerWidth_ns = 0.;
364  if( readDB_ ){
365  float tTrig, tTrigRMS, kFactor;
366  DTSuperLayerId slId = ((*lHisto).first).superlayerId();
367  int status = tTrigMap_->get( slId, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts );
368  if(status != 0) throw cms::Exception("DTNoiseCalibration") << "Could not find tTrig entry in DB for" << slId << endl;
369  triggerWidth_ns = tTrig - timeWindowOffset_;
370  } else{
371  triggerWidth_ns = defaultTtrig_ - timeWindowOffset_;
372  }
373  triggerWidth_ns = (triggerWidth_ns*25)/32;
374  triggerWidth_s = triggerWidth_ns/1e9;
375  } else{
376  triggerWidth_s = double(triggerWidth_/1e9);
377  }
378  LogTrace("Calibration") << (*lHisto).second->GetName() << " trigger width (s): " << triggerWidth_s;*/
379 
380  //double normalization = 1./(nevents_*triggerWidth_s);
381  if((*lHisto).second){
382  (*lHisto).second->Scale(normalization);
383  rootFile_->cd();
384  (*lHisto).second->Write();
385  const DTTopology& dtTopo = dtGeom_->layer((*lHisto).first)->specificTopology();
386  const int firstWire = dtTopo.firstChannel();
387  const int lastWire = dtTopo.lastChannel();
388  //const int nWires = dtTopo.channels();
389  const int nWires = lastWire - firstWire + 1;
390  // Find average in layer
391  double averageRate = 0.;
392  for(int bin = 1; bin <= (*lHisto).second->GetNbinsX(); ++bin)
393  averageRate += (*lHisto).second->GetBinContent(bin);
394 
395  if(nWires) averageRate /= nWires;
396  LogTrace("Calibration") << " Average rate = " << averageRate;
397 
398  for(int i_wire = firstWire; i_wire <= lastWire; ++i_wire){
399  // From definition of "noisy cell"
400  int bin = i_wire - firstWire + 1;
401  double channelRate = (*lHisto).second->GetBinContent(bin);
402  double rateOffset = (useAbsoluteRate_) ? 0. : averageRate;
403  if( (channelRate - rateOffset) > maximumNoiseRate_ ){
404  DTWireId wireID((*lHisto).first, i_wire);
405  statusMap->setCellNoise(wireID,1);
406  LogVerbatim("Calibration") << ">>> Channel noisy: " << wireID;
407  }
408  }
409  }
410  }
411  LogVerbatim("Calibration") << "Writing noise map object to DB";
412  string record = "DTStatusFlagRcd";
413  DTCalibDBUtils::writeToDB<DTStatusFlag>(record, statusMap);
414 }
JetCorrectorParameters::Record record
Definition: classes.h:7
edm::ESHandle< DTGeometry > dtGeom_
std::map< DTLayerId, TH1F * > theHistoOccupancyMap_
std::map< DTChamberId, TH1F * > chamberOccupancyVsTimeMap_
std::map< DTChamberId, TH1F * > chamberOccupancyVsLumiMap_
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:78
std::map< DTWireId, TH1F * > theHistoOccupancyVsLumiMap_
const DTLayer * layer(DTLayerId id) const
Return a layer given its id.
Definition: DTGeometry.cc:110
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:80
const DTTopology & specificTopology() const
Definition: DTLayer.cc:42
int setCellNoise(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, bool flag)
#define LogTrace(id)
bin
set the eta bin as selection string.
string DTNoiseCalibration::getChamberName ( const DTChamberId dtChId) const
private

Definition at line 465 of file DTNoiseCalibration.cc.

References DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by analyze().

465  {
466 
467  stringstream wheel; wheel << dtChId.wheel();
468  stringstream station; station << dtChId.station();
469  stringstream sector; sector << dtChId.sector();
470 
471  string chamberName =
472  "W" + wheel.str()
473  + "_St" + station.str()
474  + "_Sec" + sector.str();
475 
476  return chamberName;
477 }
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 DTNoiseCalibration::getChannelName ( const DTWireId wId) const
private

Definition at line 420 of file DTNoiseCalibration.cc.

References ecaldqm::binning::channelName(), DTLayerId::layer(), DTChamberId::sector(), DTChamberId::station(), DTSuperLayerId::superlayer(), DTChamberId::wheel(), and DTWireId::wire().

Referenced by analyze().

420  {
421  stringstream channelName;
422  channelName << "Wh" << wId.wheel() << "_St" << wId.station() << "_Sec" << wId.sector()
423  << "_SL" << wId.superlayer() << "_L" << wId.layer() << "_W"<< wId.wire();
424 
425  return channelName.str();
426 }
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
std::string channelName(uint32_t, BinningType _btype=kDCC)
int wire() const
Return the wire number.
Definition: DTWireId.h:56
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 DTNoiseCalibration::getLayerName ( const DTLayerId lId) const
private

Definition at line 428 of file DTNoiseCalibration.cc.

References DTSuperLayerId::chamberId(), DTLayerId::layer(), DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, DTSuperLayerId::superlayer(), DTLayerId::superlayerId(), DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by analyze().

428  {
429 
430  const DTSuperLayerId dtSLId = lId.superlayerId();
431  const DTChamberId dtChId = dtSLId.chamberId();
432  stringstream Layer; Layer << lId.layer();
433  stringstream superLayer; superLayer << dtSLId.superlayer();
434  stringstream wheel; wheel << dtChId.wheel();
435  stringstream station; station << dtChId.station();
436  stringstream sector; sector << dtChId.sector();
437 
438  string layerName =
439  "W" + wheel.str()
440  + "_St" + station.str()
441  + "_Sec" + sector.str()
442  + "_SL" + superLayer.str()
443  + "_L" + Layer.str();
444 
445  return layerName;
446 }
DTChamberId chamberId() const
Return the corresponding ChamberId.
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
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 DTNoiseCalibration::getSuperLayerName ( const DTSuperLayerId dtSLId) const
private

Definition at line 448 of file DTNoiseCalibration.cc.

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

448  {
449 
450  const DTChamberId dtChId = dtSLId.chamberId();
451  stringstream superLayer; superLayer << dtSLId.superlayer();
452  stringstream wheel; wheel << dtChId.wheel();
453  stringstream station; station << dtChId.station();
454  stringstream sector; sector << dtChId.sector();
455 
456  string superLayerName =
457  "W" + wheel.str()
458  + "_St" + station.str()
459  + "_Sec" + sector.str()
460  + "_SL" + superLayer.str();
461 
462  return superLayerName;
463 }
DTChamberId chamberId() const
Return the corresponding ChamberId.
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

Member Data Documentation

std::map<DTChamberId, TH1F*> DTNoiseCalibration::chamberOccupancyVsLumiMap_
private

Definition at line 89 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

std::map<DTChamberId, TH1F*> DTNoiseCalibration::chamberOccupancyVsTimeMap_
private

Definition at line 91 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

std::string DTNoiseCalibration::dbLabel_
private

Definition at line 66 of file DTNoiseCalibration.h.

Referenced by beginRun().

int DTNoiseCalibration::defaultTtrig_
private

Definition at line 65 of file DTNoiseCalibration.h.

Referenced by analyze(), and DTNoiseCalibration().

edm::InputTag DTNoiseCalibration::digiLabel_
private

Definition at line 53 of file DTNoiseCalibration.h.

Referenced by analyze().

edm::ESHandle<DTGeometry> DTNoiseCalibration::dtGeom_
private

Definition at line 77 of file DTNoiseCalibration.h.

Referenced by analyze(), beginRun(), and endJob().

TH1F* DTNoiseCalibration::hTDCTriggerWidth_
private

Definition at line 83 of file DTNoiseCalibration.h.

Referenced by analyze(), beginJob(), and endJob().

unsigned int DTNoiseCalibration::lumiMax_
private

Definition at line 69 of file DTNoiseCalibration.h.

Referenced by analyze().

double DTNoiseCalibration::maximumNoiseRate_
private

Definition at line 57 of file DTNoiseCalibration.h.

Referenced by endJob().

int DTNoiseCalibration::nevents_
private

Definition at line 71 of file DTNoiseCalibration.h.

Referenced by analyze(), beginJob(), and endJob().

bool DTNoiseCalibration::readDB_
private

Definition at line 64 of file DTNoiseCalibration.h.

Referenced by analyze(), beginRun(), and DTNoiseCalibration().

TFile* DTNoiseCalibration::rootFile_
private

Definition at line 81 of file DTNoiseCalibration.h.

Referenced by analyze(), DTNoiseCalibration(), endJob(), and ~DTNoiseCalibration().

time_t DTNoiseCalibration::runBeginTime_
private

Definition at line 73 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginRun().

time_t DTNoiseCalibration::runEndTime_
private

Definition at line 74 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginRun().

std::map<DTLayerId, TH1F*> DTNoiseCalibration::theHistoOccupancyMap_
private

Definition at line 85 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

std::map<DTWireId, TH1F*> DTNoiseCalibration::theHistoOccupancyVsLumiMap_
private

Definition at line 87 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

int DTNoiseCalibration::timeWindowOffset_
private

Definition at line 56 of file DTNoiseCalibration.h.

Referenced by analyze().

double DTNoiseCalibration::triggerWidth_
private

Definition at line 55 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginJob().

edm::ESHandle<DTTtrig> DTNoiseCalibration::tTrigMap_
private

Definition at line 79 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginRun().

bool DTNoiseCalibration::useAbsoluteRate_
private

Definition at line 58 of file DTNoiseCalibration.h.

Referenced by endJob().

bool DTNoiseCalibration::useTimeWindow_
private

Definition at line 54 of file DTNoiseCalibration.h.

Referenced by analyze().

std::vector<DTWireId> DTNoiseCalibration::wireIdWithHisto_
private

Definition at line 68 of file DTNoiseCalibration.h.

Referenced by analyze(), and DTNoiseCalibration().