CMS 3D CMS Logo

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

#include <DTChamberEfficiencyTask.h>

Inheritance diagram for DTChamberEfficiencyTask:
DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks > edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::Accumulator, Args... > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup) override
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) override
 To reset the MEs. More...
 
void dqmBeginRun (const edm::Run &run, const edm::EventSetup &setup) override
 BeginRun. More...
 
 DTChamberEfficiencyTask (const edm::ParameterSet &pset)
 Constructor. More...
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) final
 
 ~DTChamberEfficiencyTask () override
 Destructor. More...
 
- Public Member Functions inherited from DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >
void accumulate (edm::Event const &event, edm::EventSetup const &setup) override
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
 DQMOneEDAnalyzer ()
 
void endRun (edm::Run const &, edm::EventSetup const &) final
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::Accumulator, Args... >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 
bool wantsGlobalLuminosityBlocks () const noexcept final
 
bool wantsGlobalRuns () const noexcept final
 
bool wantsInputProcessBlocks () const noexcept final
 
bool wantsProcessBlocks () const noexcept final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const noexcept
 
bool wantsStreamRuns () const noexcept
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
std::vector< bool > const & recordProvenanceList () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
TypeLabelList const & typeLabelList () const
 used by the fwk to register the list of products of this module More...
 
 ~ProducerBase () noexcept(false) 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
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Protected Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
- Protected Member Functions inherited from DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >
virtual void dqmEndRun (edm::Run const &, edm::EventSetup const &)
 
- Protected Member Functions inherited from edm::ProducerBase
template<Transition Tr = Transition::Event>
auto produces (std::string instanceName) noexcept
 declare what type of product will make and with which optional label More...
 
template<Transition B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<BranchType B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces ()
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces ()
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces ()
 
template<Transition Tr = Transition::Event>
auto produces () noexcept
 
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Private Member Functions

void bookHistos (DQMStore::IBooker &ibooker, DTChamberId chId)
 
const DTRecSegment4DgetBestSegment (const DTRecSegment4DCollection::range &segs) const
 
const DTRecSegment4DgetBestSegment (const DTRecSegment4D *s1, const DTRecSegment4D *s2) const
 
LocalPoint interpolate (const DTRecSegment4D &seg1, const DTRecSegment4D &seg3, const DTChamberId &MB2) const
 
bool isGoodSegment (const DTRecSegment4D &seg) const
 

Private Attributes

bool debug
 
bool detailedAnalysis
 
const DTGeometrydtGeom
 
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
 
edm::ESGetToken< DTGeometry, MuonGeometryRecordmuonGeomToken_
 
bool onlineMonitor
 
edm::ParameterSet parameters
 
edm::EDGetTokenT< DTRecSegment4DCollectionrecHits4DToken_
 
edm::Handle< DTRecSegment4DCollectionsegs
 
double theMinChi2NormSegment
 
double theMinCloseDist
 
unsigned int theMinHitsSegment
 

Additional Inherited Members

- Public Types inherited from DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
template<typename T >
using BranchAliasSetterT = ProductRegistryHelper::BranchAliasSetterT< T >
 
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex > >
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Attributes inherited from DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

DQM Analysis of 4D DT segments, it produces plots about:

Class based on the code written by S. Lacaprara : RecoLocalMuon / DTSegment / test / DTEffAnalyzer.h

Author
G. Mila - INFN Torino

Definition at line 35 of file DTChamberEfficiencyTask.h.

Constructor & Destructor Documentation

◆ DTChamberEfficiencyTask()

DTChamberEfficiencyTask::DTChamberEfficiencyTask ( const edm::ParameterSet pset)

Constructor.

Definition at line 23 of file DTChamberEfficiencyTask.cc.

References debug, detailedAnalysis, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), ProducerED_cfi::InputTag, onlineMonitor, parameters, muonDTDigis_cfi::pset, recHits4DToken_, theMinChi2NormSegment, theMinCloseDist, and theMinHitsSegment.

24  : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
25  debug = pset.getUntrackedParameter<bool>("debug", false);
26 
27  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Constructor called!";
28 
29  parameters = pset;
30 
31  // the name of the 4D rec hits collection
33  consumes<DTRecSegment4DCollection>(edm::InputTag(parameters.getUntrackedParameter<string>("recHits4DLabel")));
34 
35  // parameters to use for the segment quality check
36  theMinHitsSegment = static_cast<unsigned int>(parameters.getParameter<int>("minHitsSegment"));
37  theMinChi2NormSegment = parameters.getParameter<double>("minChi2NormSegment");
38  // parameter to use for the exstrapolated segment check
39  theMinCloseDist = parameters.getParameter<double>("minCloseDist");
40 
41  // the running modality
42  onlineMonitor = parameters.getUntrackedParameter<bool>("onlineMonitor");
43 
44  // the analysis mode
45  detailedAnalysis = parameters.getUntrackedParameter<bool>("detailedAnalysis");
46 }
Log< level::Info, true > LogVerbatim
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
T getUntrackedParameter(std::string const &, T const &) const

◆ ~DTChamberEfficiencyTask()

DTChamberEfficiencyTask::~DTChamberEfficiencyTask ( )
override

Destructor.

Definition at line 48 of file DTChamberEfficiencyTask.cc.

48  {
49  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Destructor called!";
50 }
Log< level::Info, true > LogVerbatim

Member Function Documentation

◆ analyze()

void DTChamberEfficiencyTask::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
overridevirtual

Reimplemented from DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >.

Definition at line 145 of file DTChamberEfficiencyTask.cc.

References DTGeometry::chamber(), DTGeometry::chambers(), detailedAnalysis, dtGeom, getBestSegment(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), mergeVDriftHistosByStation::histos, histosPerCh, interpolate(), isGoodSegment(), DTRecSegment4D::localPosition(), mag(), recHits4DToken_, nano_mu_digi_cff::sector, DTChamberId::sector(), segs, DTChamberId::station(), relativeConstraints::station, theMinCloseDist, DTChamberId::wheel(), makeMuonMisalignmentScenario::wheel, x, PV3DBase< T, PVType, FrameType >::x(), y, and PV3DBase< T, PVType, FrameType >::y().

145  {
146  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask")
147  << "[DTChamberEfficiencyTask] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event();
148 
149  // Get the 4D rechit collection from the event
150  event.getByToken(recHits4DToken_, segs);
151 
152  int bottom = 0, top = 0;
153 
154  // Loop over all the chambers
155  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
156  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
157  for (; ch_it != ch_end; ++ch_it) {
158  DTChamberId ch = (*ch_it)->id();
159  int wheel = ch.wheel();
160  int sector = ch.sector();
161  int station = ch.station();
162 
163  DTChamberId MidId(wheel, station, sector);
164 
165  // get efficiency for MB1 using MB2 and MB3
166  if (station == 1) {
167  bottom = 2;
168  top = 3;
169  }
170 
171  // get efficiency for MB2 using MB1 and MB3
172  if (station == 2) {
173  bottom = 1;
174  top = 3;
175  }
176 
177  // get efficiency for MB2 using MB2 and MB4
178  if (station == 3) {
179  bottom = 2;
180  top = 4;
181  }
182 
183  // get efficiency for MB4 using MB2 and MB3
184  if (station == 4) {
185  bottom = 2;
186  top = 3;
187  }
188 
189  // Select events with (good) segments in Bot and Top
190  DTChamberId BotId(wheel, bottom, sector);
191  DTChamberId TopId(wheel, top, sector);
192 
193  // Get segments in the bottom chambers (if any)
194  DTRecSegment4DCollection::range segsBot = segs->get(BotId);
195  int nSegsBot = segsBot.second - segsBot.first;
196  // check if any segments is there
197  if (nSegsBot == 0)
198  continue;
199 
200  vector<MonitorElement*> histos = histosPerCh[MidId];
201 
202  // Get segments in the top chambers (if any)
203  DTRecSegment4DCollection::range segsTop = segs->get(TopId);
204  int nSegsTop = segsTop.second - segsTop.first;
205 
206  // Select one segment for the bottom chamber
207  const DTRecSegment4D& bestBotSeg = getBestSegment(segsBot);
208 
209  // Select one segment for the top chamber
210  DTRecSegment4D* pBestTopSeg = nullptr;
211  if (nSegsTop > 0)
212  pBestTopSeg = const_cast<DTRecSegment4D*>(&getBestSegment(segsTop));
213  //if top chamber is MB4 sector 10, consider also sector 14
214  if (TopId.station() == 4 && TopId.sector() == 10) {
215  DTChamberId TopId14(wheel, top, 14);
216  DTRecSegment4DCollection::range segsTop14 = segs->get(TopId14);
217  int nSegsTop14 = segsTop14.second - segsTop14.first;
218  nSegsTop += nSegsTop14;
219  if (nSegsTop14) {
220  DTRecSegment4D* pBestTopSeg14 = const_cast<DTRecSegment4D*>(&getBestSegment(segsTop14));
221  // get best between sector 10 and 14
222  pBestTopSeg = const_cast<DTRecSegment4D*>(getBestSegment(pBestTopSeg, pBestTopSeg14));
223  }
224  }
225  if (!pBestTopSeg)
226  continue;
227  const DTRecSegment4D& bestTopSeg = *pBestTopSeg;
228 
229  DTRecSegment4DCollection::range segsMid = segs->get(MidId);
230  int nSegsMid = segsMid.second - segsMid.first;
231 
232  if (detailedAnalysis) {
233  // very trivial efficiency, just count segments
234  histos[3]->Fill(0);
235  if (nSegsMid > 0)
236  histos[3]->Fill(1);
237  }
238 
239  // get position at Mid by interpolating the position (not direction) of best
240  // segment in Bot and Top to Mid surface
241  LocalPoint posAtMid = interpolate(bestBotSeg, bestTopSeg, MidId);
242 
243  // is best segment good enough?
244  if (isGoodSegment(bestBotSeg) && isGoodSegment(bestTopSeg)) {
245  if (detailedAnalysis)
246  histos[4]->Fill(posAtMid.x(), posAtMid.y());
247  //check if interpolation fall inside middle chamber
248  if ((dtGeom->chamber(MidId))->surface().bounds().inside(posAtMid)) {
249  histos[0]->Fill(posAtMid.x(), posAtMid.y());
250  if (nSegsMid > 0) {
251  if (detailedAnalysis) {
252  histos[3]->Fill(2);
253  histos[5]->Fill(posAtMid.x(), posAtMid.y());
254  }
255 
256  const DTRecSegment4D& bestMidSeg = getBestSegment(segsMid);
257  // check if middle segments is good enough
258  if (isGoodSegment(bestMidSeg)) {
259  if (detailedAnalysis)
260  histos[6]->Fill(posAtMid.x(), posAtMid.y());
261  LocalPoint midSegPos = bestMidSeg.localPosition();
262 
263  // check if middle segments is also close enough
264  double dist;
265  if (bestMidSeg.hasPhi()) {
266  if (bestTopSeg.hasZed() && bestBotSeg.hasZed() && bestMidSeg.hasZed()) {
267  dist = (midSegPos - posAtMid).mag();
268  } else {
269  dist = fabs((midSegPos - posAtMid).x());
270  }
271  } else {
272  dist = fabs((midSegPos - posAtMid).y());
273  }
274  if (dist < theMinCloseDist) {
275  histos[1]->Fill(posAtMid.x(), posAtMid.y());
276  }
277  if (detailedAnalysis)
278  histos[2]->Fill(dist);
279  }
280  }
281  }
282  }
283  } // loop over stations
284 }
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:45
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
bool hasPhi() const
Does it have the Phi projection?
LocalPoint interpolate(const DTRecSegment4D &seg1, const DTRecSegment4D &seg3, const DTChamberId &MB2) const
LocalPoint localPosition() const override
Local position in Chamber frame.
edm::Handle< DTRecSegment4DCollection > segs
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const DTRecSegment4D & getBestSegment(const DTRecSegment4DCollection::range &segs) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
bool isGoodSegment(const DTRecSegment4D &seg) const
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
int sector() const
Definition: DTChamberId.h:52
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
bool hasZed() const
Does it have the Z projection?
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90

◆ beginLuminosityBlock()

void DTChamberEfficiencyTask::beginLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  context 
)
override

To reset the MEs.

Definition at line 52 of file DTChamberEfficiencyTask.cc.

References edm::ParameterSet::getUntrackedParameter(), timingPdfMaker::histo, histosPerCh, mps_fire::i, edm::LuminosityBlockBase::id(), edm::LuminosityBlockID::luminosityBlock(), genParticles_cff::map, onlineMonitor, parameters, and findQualityFiles::size.

52  {
53  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask")
54  << "[DTChamberEfficiencyTask]: Begin of LS transition";
55 
56  if (lumiSeg.id().luminosityBlock() % parameters.getUntrackedParameter<int>("ResetCycle", 3) == 0 && onlineMonitor) {
57  for (map<DTChamberId, vector<MonitorElement*> >::const_iterator histo = histosPerCh.begin();
58  histo != histosPerCh.end();
59  histo++) {
60  int size = (*histo).second.size();
61  for (int i = 0; i < size; i++) {
62  (*histo).second[i]->Reset();
63  }
64  }
65  }
66 }
size
Write out results.
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh

◆ bookHistograms()

void DTChamberEfficiencyTask::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  context 
)
overrideprotectedvirtual

Implements DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >.

Definition at line 73 of file DTChamberEfficiencyTask.cc.

References bookHistos(), DTGeometry::chambers(), dtGeom, and dqm::implementation::NavigatorBase::setCurrentFolder().

75  {
76  ibooker.setCurrentFolder("DT/DTChamberEfficiencyTask");
77 
78  // Loop over all the chambers
79  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
80  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
81  for (; ch_it != ch_end; ++ch_it) {
82  // histo booking
83  bookHistos(ibooker, (*ch_it)->id());
84  }
85 }
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
void bookHistos(DQMStore::IBooker &ibooker, DTChamberId chId)

◆ bookHistos()

void DTChamberEfficiencyTask::bookHistos ( DQMStore::IBooker ibooker,
DTChamberId  chId 
)
private

Definition at line 88 of file DTChamberEfficiencyTask.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2D(), detailedAnalysis, mergeVDriftHistosByStation::histos, histosPerCh, nano_mu_digi_cff::sector, DTChamberId::sector(), dqm::implementation::NavigatorBase::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by bookHistograms().

88  {
89  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << " Booking histos for CH : " << chId;
90 
91  // Compose the chamber name
92  stringstream wheel;
93  wheel << chId.wheel();
94  stringstream station;
95  station << chId.station();
96  stringstream sector;
97  sector << chId.sector();
98 
99  string HistoName = "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
100 
101  ibooker.setCurrentFolder("DT/01-DTChamberEfficiency/Task/Wheel" + wheel.str() + "/Sector" + sector.str() +
102  "/Station" + station.str());
103 
104  // Create the monitor elements
105  vector<MonitorElement*> histos;
106 
107  //efficiency selection cuts
108  // a- number of segments of the top chamber > 0 && number of segments of the bottom chamber > 0
109  // b- number of segments of the middle chamber > 0
110  // c- check of the top and bottom segment quality
111  // d- check if interpolation falls inside the middle chamber
112  // e- check of the middle segment quality
113  // f- check if the distance between the reconstructed and the exstrapolated segments is ok
114 
115  // histo for efficiency with cuts a-/c-/d-
116  histos.push_back(ibooker.book2D(
117  "hEffGoodSegVsPosDen" + HistoName, "Eff vs local position (good) ", 25, -250., 250., 25, -250., 250.));
118  // histo for efficiency with cuts a-/b-/c-/d-/e-/f-
119  histos.push_back(ibooker.book2D("hEffGoodCloseSegVsPosNum" + HistoName,
120  "Eff vs local position (good and close segs) ",
121  25,
122  -250.,
123  250.,
124  25,
125  -250.,
126  250.));
127  if (detailedAnalysis) {
128  histos.push_back(
129  ibooker.book1D("hDistSegFromExtrap" + HistoName, "Distance segments from extrap position ", 200, 0., 200.));
130  // histo for efficiency from segment counting
131  histos.push_back(ibooker.book1D("hNaiveEffSeg" + HistoName, "Naive eff ", 10, 0., 10.));
132  // histo for efficiency with cuts a-/c-
133  histos.push_back(ibooker.book2D(
134  "hEffSegVsPosDen" + HistoName, "Eff vs local position (all) ", 25, -250., 250., 25, -250., 250.));
135  // histo for efficiency with cuts a-/b-/c-/d-
136  histos.push_back(
137  ibooker.book2D("hEffSegVsPosNum" + HistoName, "Eff vs local position ", 25, -250., 250., 25, -250., 250.));
138  // histo for efficiency with cuts a-/b-/c-/d-/e-
139  histos.push_back(ibooker.book2D(
140  "hEffGoodSegVsPosNum" + HistoName, "Eff vs local position (good segs) ", 25, -250., 250., 25, -250., 250.));
141  }
142  histosPerCh[chId] = histos;
143 }
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:45
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
std::string HistoName
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
int sector() const
Definition: DTChamberId.h:52
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98

◆ dqmBeginRun()

void DTChamberEfficiencyTask::dqmBeginRun ( const edm::Run run,
const edm::EventSetup setup 
)
overridevirtual

BeginRun.

Reimplemented from DQMOneEDAnalyzer< edm::one::WatchLuminosityBlocks >.

Definition at line 68 of file DTChamberEfficiencyTask.cc.

References dtGeom, muonGeomToken_, and singleTopDQM_cfi::setup.

68  {
69  // Get the DT Geometry
70  dtGeom = &setup.getData(muonGeomToken_);
71 }
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_

◆ endLuminosityBlock()

void DTChamberEfficiencyTask::endLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  context 
)
inlinefinal

Definition at line 48 of file DTChamberEfficiencyTask.h.

48 {}

◆ getBestSegment() [1/2]

const DTRecSegment4D & DTChamberEfficiencyTask::getBestSegment ( const DTRecSegment4DCollection::range segs) const
private

Definition at line 287 of file DTChamberEfficiencyTask.cc.

References nHits, and segs.

Referenced by analyze().

287  {
289  unsigned int nHitBest = 0;
290  double chi2Best = 99999.;
291  for (DTRecSegment4DCollection::const_iterator seg = segs.first; seg != segs.second; ++seg) {
292  unsigned int nHits = ((*seg).hasPhi() ? (*seg).phiSegment()->recHits().size() : 0);
293  nHits += ((*seg).hasZed() ? (*seg).zSegment()->recHits().size() : 0);
294 
295  if (nHits == nHitBest) {
296  if ((*seg).chi2() / (*seg).degreesOfFreedom() < chi2Best) {
297  chi2Best = (*seg).chi2() / (*seg).degreesOfFreedom();
298  bestIter = seg;
299  }
300  } else if (nHits > nHitBest) {
301  nHitBest = nHits;
302  bestIter = seg;
303  }
304  }
305  return *bestIter;
306 }
edm::Handle< DTRecSegment4DCollection > segs
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits

◆ getBestSegment() [2/2]

const DTRecSegment4D * DTChamberEfficiencyTask::getBestSegment ( const DTRecSegment4D s1,
const DTRecSegment4D s2 
) const
private

Definition at line 308 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::chi2(), DTRecSegment4D::degreesOfFreedom(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), DTRecSegment4D::phiSegment(), DTRecSegment2D::recHits(), and DTRecSegment4D::zSegment().

309  {
310  if (!s1)
311  return s2;
312  if (!s2)
313  return s1;
314  unsigned int nHits1 = (s1->hasPhi() ? s1->phiSegment()->recHits().size() : 0);
315  nHits1 += (s1->hasZed() ? s1->zSegment()->recHits().size() : 0);
316 
317  unsigned int nHits2 = (s2->hasPhi() ? s2->phiSegment()->recHits().size() : 0);
318  nHits2 += (s2->hasZed() ? s2->zSegment()->recHits().size() : 0);
319 
320  if (nHits1 == nHits2) {
321  if (s1->chi2() / s1->degreesOfFreedom() < s2->chi2() / s2->degreesOfFreedom())
322  return s1;
323  else
324  return s2;
325  } else if (nHits1 > nHits2)
326  return s1;
327  return s2;
328 }
bool hasPhi() const
Does it have the Phi projection?
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
int degreesOfFreedom() const override
Degrees of freedom of the segment fit.
bool hasZed() const
Does it have the Z projection?
double chi2() const override
Chi2 of the segment fit.

◆ interpolate()

LocalPoint DTChamberEfficiencyTask::interpolate ( const DTRecSegment4D seg1,
const DTRecSegment4D seg3,
const DTChamberId MB2 
) const
private

Definition at line 330 of file DTChamberEfficiencyTask.cc.

References DTGeometry::chamber(), DTRecSegment4D::chamberId(), DeadROC_duringRun::dir, dtGeom, DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), globals_cff::id2, DTRecSegment4D::localPosition(), toLocal(), unit(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze().

332  {
333  // Get GlobalPoition of Seg in MB1
334  GlobalPoint gpos1 = (dtGeom->chamber(seg1.chamberId()))->toGlobal(seg1.localPosition());
335 
336  // Get GlobalPoition of Seg in MB3
337  GlobalPoint gpos3 = (dtGeom->chamber(seg3.chamberId()))->toGlobal(seg3.localPosition());
338 
339  // interpolate
340  // get all in MB2 frame
341  LocalPoint pos1 = (dtGeom->chamber(id2))->toLocal(gpos1);
342  LocalPoint pos3 = (dtGeom->chamber(id2))->toLocal(gpos3);
343 
344  // case 1: 1 and 3 has both projection. No problem
345 
346  // case 2: one projection is missing for one of the segments. Keep the other's segment position
347  if (!seg1.hasZed())
348  pos1 = LocalPoint(pos1.x(), pos3.y(), pos1.z());
349  if (!seg3.hasZed())
350  pos3 = LocalPoint(pos3.x(), pos1.y(), pos3.z());
351 
352  if (!seg1.hasPhi())
353  pos1 = LocalPoint(pos3.x(), pos1.y(), pos1.z());
354  if (!seg3.hasPhi())
355  pos3 = LocalPoint(pos1.x(), pos3.y(), pos3.z());
356 
357  // direction
358  LocalVector dir = (pos3 - pos1).unit(); // z points inward!
359  LocalPoint pos2 = pos1 + dir * pos1.z() / (-dir.z());
360 
361  return pos2;
362 }
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
bool hasPhi() const
Does it have the Phi projection?
T z() const
Definition: PV3DBase.h:61
virtual DTChamberId chamberId() const
The (specific) DetId of the chamber on which the segment resides.
LocalPoint localPosition() const override
Local position in Chamber frame.
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
Basic3DVector unit() const
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
bool hasZed() const
Does it have the Z projection?
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90

◆ isGoodSegment()

bool DTChamberEfficiencyTask::isGoodSegment ( const DTRecSegment4D seg) const
private

Definition at line 364 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::chamberId(), DTRecSegment4D::chi2(), DTRecSegment4D::degreesOfFreedom(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), nHits, DTRecSegment4D::phiSegment(), DTRecSegment2D::recHits(), DTChamberId::station(), theMinChi2NormSegment, theMinHitsSegment, and DTRecSegment4D::zSegment().

Referenced by analyze().

364  {
365  if (seg.chamberId().station() != 4 && !seg.hasZed())
366  return false;
367  unsigned int nHits = (seg.hasPhi() ? seg.phiSegment()->recHits().size() : 0);
368  nHits += (seg.hasZed() ? seg.zSegment()->recHits().size() : 0);
369  return (nHits >= theMinHitsSegment && seg.chi2() / seg.degreesOfFreedom() < theMinChi2NormSegment);
370 }
int station() const
Return the station number.
Definition: DTChamberId.h:45
bool hasPhi() const
Does it have the Phi projection?
virtual DTChamberId chamberId() const
The (specific) DetId of the chamber on which the segment resides.
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
int degreesOfFreedom() const override
Degrees of freedom of the segment fit.
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
bool hasZed() const
Does it have the Z projection?
double chi2() const override
Chi2 of the segment fit.

Member Data Documentation

◆ debug

bool DTChamberEfficiencyTask::debug
private

◆ detailedAnalysis

bool DTChamberEfficiencyTask::detailedAnalysis
private

Definition at line 70 of file DTChamberEfficiencyTask.h.

Referenced by analyze(), bookHistos(), and DTChamberEfficiencyTask().

◆ dtGeom

const DTGeometry* DTChamberEfficiencyTask::dtGeom
private

Definition at line 85 of file DTChamberEfficiencyTask.h.

Referenced by analyze(), bookHistograms(), dqmBeginRun(), and interpolate().

◆ histosPerCh

std::map<DTChamberId, std::vector<MonitorElement*> > DTChamberEfficiencyTask::histosPerCh
private

Definition at line 77 of file DTChamberEfficiencyTask.h.

Referenced by analyze(), beginLuminosityBlock(), and bookHistos().

◆ muonGeomToken_

edm::ESGetToken<DTGeometry, MuonGeometryRecord> DTChamberEfficiencyTask::muonGeomToken_
private

Definition at line 84 of file DTChamberEfficiencyTask.h.

Referenced by dqmBeginRun().

◆ onlineMonitor

bool DTChamberEfficiencyTask::onlineMonitor
private

Definition at line 68 of file DTChamberEfficiencyTask.h.

Referenced by beginLuminosityBlock(), and DTChamberEfficiencyTask().

◆ parameters

edm::ParameterSet DTChamberEfficiencyTask::parameters
private

Definition at line 75 of file DTChamberEfficiencyTask.h.

Referenced by beginLuminosityBlock(), and DTChamberEfficiencyTask().

◆ recHits4DToken_

edm::EDGetTokenT<DTRecSegment4DCollection> DTChamberEfficiencyTask::recHits4DToken_
private

Definition at line 73 of file DTChamberEfficiencyTask.h.

Referenced by analyze(), and DTChamberEfficiencyTask().

◆ segs

edm::Handle<DTRecSegment4DCollection> DTChamberEfficiencyTask::segs
private

Definition at line 86 of file DTChamberEfficiencyTask.h.

Referenced by analyze(), and getBestSegment().

◆ theMinChi2NormSegment

double DTChamberEfficiencyTask::theMinChi2NormSegment
private

Definition at line 80 of file DTChamberEfficiencyTask.h.

Referenced by DTChamberEfficiencyTask(), and isGoodSegment().

◆ theMinCloseDist

double DTChamberEfficiencyTask::theMinCloseDist
private

Definition at line 81 of file DTChamberEfficiencyTask.h.

Referenced by analyze(), and DTChamberEfficiencyTask().

◆ theMinHitsSegment

unsigned int DTChamberEfficiencyTask::theMinHitsSegment
private

Definition at line 79 of file DTChamberEfficiencyTask.h.

Referenced by DTChamberEfficiencyTask(), and isGoodSegment().