CMS 3D CMS Logo

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

#include <DTNoiseCalibration.h>

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

Public Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 
void beginJob () override
 
void beginRun (const edm::Run &run, const edm::EventSetup &setup) override
 
 DTNoiseCalibration (const edm::ParameterSet &ps)
 Constructor. More...
 
void endJob () override
 
 ~DTNoiseCalibration () 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
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
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)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
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 &)
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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

DTNoiseCalibration::DTNoiseCalibration ( const edm::ParameterSet ps)

Constructor.

Definition at line 42 of file DTNoiseCalibration.cc.

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

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

Destructor.

Definition at line 419 of file DTNoiseCalibration.cc.

References rootFile_.

419 { rootFile_->Close(); }

Member Function Documentation

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

Definition at line 116 of file DTNoiseCalibration.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 87 of file DTNoiseCalibration.cc.

References hTDCTriggerWidth_, nevents_, and triggerWidth_.

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

Definition at line 97 of file DTNoiseCalibration.cc.

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

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

Reimplemented from edm::EDAnalyzer.

Definition at line 327 of file DTNoiseCalibration.cc.

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

Definition at line 447 of file DTNoiseCalibration.cc.

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

Referenced by analyze().

447  {
448  string chamberName = "W" + std::to_string(dtChId.wheel()) + "_St" + std::to_string(dtChId.station()) + "_Sec" +
449  std::to_string(dtChId.sector());
450 
451  return chamberName;
452 }
int sector() const
Definition: DTChamberId.h:49
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
string DTNoiseCalibration::getChannelName ( const DTWireId wId) const
private

Definition at line 421 of file DTNoiseCalibration.cc.

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

Referenced by analyze().

421  {
422  string channelName = "Wh" + std::to_string(wId.wheel()) + "_St" + std::to_string(wId.station()) + "_Sec" +
423  std::to_string(wId.sector()) + "_SL" + std::to_string(wId.superlayer()) + "_L" +
424  std::to_string(wId.layer()) + "_W" + std::to_string(wId.wire());
425  return channelName;
426 }
int layer() const
Return the layer number.
Definition: DTLayerId.h:42
std::string channelName(uint32_t, BinningType _btype=kDCC)
int wire() const
Return the wire number.
Definition: DTWireId.h:42
int superlayer() const
Return the superlayer number (deprecated method name)
int sector() const
Definition: DTChamberId.h:49
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
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(), DTSuperLayerId::superlayer(), DTLayerId::superlayerId(), and DTChamberId::wheel().

Referenced by analyze().

428  {
429  const DTSuperLayerId dtSLId = lId.superlayerId();
430  const DTChamberId dtChId = dtSLId.chamberId();
431  string layerName = "W" + std::to_string(dtChId.wheel()) + "_St" + std::to_string(dtChId.station()) + "_Sec" +
432  std::to_string(dtChId.sector()) + "_SL" + std::to_string(dtSLId.superlayer()) + "_L" +
433  std::to_string(lId.layer());
434 
435  return layerName;
436 }
DTChamberId chamberId() const
Return the corresponding ChamberId.
int layer() const
Return the layer number.
Definition: DTLayerId.h:42
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:45
int superlayer() const
Return the superlayer number (deprecated method name)
int sector() const
Definition: DTChamberId.h:49
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
string DTNoiseCalibration::getSuperLayerName ( const DTSuperLayerId dtSLId) const
private

Definition at line 438 of file DTNoiseCalibration.cc.

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

438  {
439  const DTChamberId dtChId = dtSLId.chamberId();
440 
441  string superLayerName = "W" + std::to_string(dtChId.wheel()) + "_St" + std::to_string(dtChId.station()) + "_Sec" +
442  std::to_string(dtChId.sector()) + "_SL" + std::to_string(dtSLId.superlayer());
443 
444  return superLayerName;
445 }
DTChamberId chamberId() const
Return the corresponding ChamberId.
int superlayer() const
Return the superlayer number (deprecated method name)
int sector() const
Definition: DTChamberId.h:49
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39

Member Data Documentation

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

Definition at line 88 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

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

Definition at line 90 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

std::string DTNoiseCalibration::dbLabel_
private

Definition at line 65 of file DTNoiseCalibration.h.

Referenced by beginRun().

int DTNoiseCalibration::defaultTtrig_
private

Definition at line 64 of file DTNoiseCalibration.h.

Referenced by analyze(), and DTNoiseCalibration().

edm::InputTag DTNoiseCalibration::digiLabel_
private

Definition at line 52 of file DTNoiseCalibration.h.

Referenced by analyze().

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

Definition at line 76 of file DTNoiseCalibration.h.

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

TH1F* DTNoiseCalibration::hTDCTriggerWidth_
private

Definition at line 82 of file DTNoiseCalibration.h.

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

unsigned int DTNoiseCalibration::lumiMax_
private

Definition at line 68 of file DTNoiseCalibration.h.

Referenced by analyze().

double DTNoiseCalibration::maximumNoiseRate_
private

Definition at line 56 of file DTNoiseCalibration.h.

Referenced by endJob().

int DTNoiseCalibration::nevents_
private

Definition at line 70 of file DTNoiseCalibration.h.

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

bool DTNoiseCalibration::readDB_
private

Definition at line 63 of file DTNoiseCalibration.h.

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

TFile* DTNoiseCalibration::rootFile_
private

Definition at line 80 of file DTNoiseCalibration.h.

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

time_t DTNoiseCalibration::runBeginTime_
private

Definition at line 72 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginRun().

time_t DTNoiseCalibration::runEndTime_
private

Definition at line 73 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginRun().

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

Definition at line 84 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

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

Definition at line 86 of file DTNoiseCalibration.h.

Referenced by analyze(), and endJob().

int DTNoiseCalibration::timeWindowOffset_
private

Definition at line 55 of file DTNoiseCalibration.h.

Referenced by analyze().

double DTNoiseCalibration::triggerWidth_
private

Definition at line 54 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginJob().

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

Definition at line 78 of file DTNoiseCalibration.h.

Referenced by analyze(), and beginRun().

bool DTNoiseCalibration::useAbsoluteRate_
private

Definition at line 57 of file DTNoiseCalibration.h.

Referenced by endJob().

bool DTNoiseCalibration::useTimeWindow_
private

Definition at line 53 of file DTNoiseCalibration.h.

Referenced by analyze().

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

Definition at line 67 of file DTNoiseCalibration.h.

Referenced by analyze(), and DTNoiseCalibration().