CMS 3D CMS Logo

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

#include <DTSegment2DQuality.h>

Inheritance diagram for DTSegment2DQuality:
DQMGlobalEDAnalyzer< Histograms > edm::global::EDAnalyzer< edm::RunCache< Histograms > > edm::global::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 DTSegment2DQuality (const edm::ParameterSet &pset)
 Constructor. More...
 
- Public Member Functions inherited from edm::global::EDAnalyzer< edm::RunCache< Histograms > >
 EDAnalyzer ()=default
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void bookHistograms (DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, Histograms &) const override
 Book the DQM plots. More...
 
void dqmAnalyze (edm::Event const &, edm::EventSetup const &, Histograms const &) const override
 Perform the real analysis. More...
 

Private Attributes

bool debug_
 
edm::InputTag segment2DLabel_
 
edm::EDGetTokenT< DTRecSegment2DCollectionsegment2DToken_
 
double sigmaResAngle_
 
double sigmaResPos_
 
edm::InputTag simHitLabel_
 
edm::EDGetTokenT< edm::PSimHitContainersimHitToken_
 

Additional Inherited Members

- Public Types inherited from edm::global::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::global::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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

Basic analyzer class which accesses 2D DTSegments and plot resolution comparing reconstructed and simulated quantities

Author
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 31 of file DTSegment2DQuality.h.

Constructor & Destructor Documentation

DTSegment2DQuality::DTSegment2DQuality ( const edm::ParameterSet pset)

Constructor.

Definition at line 43 of file DTSegment2DQuality.cc.

References gather_cfg::cout, DTHitQualityUtils::debug, edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

43  {
44  // get the debug parameter for verbose output
45  debug_ = pset.getUntrackedParameter<bool>("debug");
47  // the name of the simhit collection
48  simHitLabel_ = pset.getUntrackedParameter<InputTag>("simHitLabel");
49  simHitToken_ = consumes<PSimHitContainer>(pset.getUntrackedParameter<InputTag>("simHitLabel"));
50  // the name of the 2D rec hit collection
51  segment2DLabel_ = pset.getUntrackedParameter<InputTag>("segment2DLabel");
52  segment2DToken_ = consumes<DTRecSegment2DCollection>(pset.getUntrackedParameter<InputTag>("segment2DLabel"));
53 
54  // sigma resolution on position
55  sigmaResPos_ = pset.getParameter<double>("sigmaResPos");
56  // sigma resolution on angle
57  sigmaResAngle_ = pset.getParameter<double>("sigmaResAngle");
58 
59  if (debug_) {
60  cout << "[DTSegment2DQuality] Constructor called " << endl;
61 }
62 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::PSimHitContainer > simHitToken_
edm::EDGetTokenT< DTRecSegment2DCollection > segment2DToken_
edm::InputTag simHitLabel_
edm::InputTag segment2DLabel_

Member Function Documentation

void DTSegment2DQuality::bookHistograms ( DQMStore::ConcurrentBooker booker,
edm::Run const &  run,
edm::EventSetup const &  setup,
Histograms histograms 
) const
overrideprivatevirtual

Book the DQM plots.

Implements DQMGlobalEDAnalyzer< Histograms >.

Definition at line 64 of file DTSegment2DQuality.cc.

References gather_cfg::cout, Histograms::h2DHitEff_RPhi, Histograms::h2DHitEff_RZ, Histograms::h2DHitEff_RZ_W0, Histograms::h2DHitEff_RZ_W1, Histograms::h2DHitEff_RZ_W2, Histograms::h2DHitRPhi, Histograms::h2DHitRZ, Histograms::h2DHitRZ_W0, Histograms::h2DHitRZ_W1, and Histograms::h2DHitRZ_W2.

64  {
65  histograms.h2DHitRPhi = new HRes2DHit ("RPhi", booker, true, true);
66  histograms.h2DHitRZ = new HRes2DHit ("RZ", booker, true, true);
67  histograms.h2DHitRZ_W0 = new HRes2DHit ("RZ_W0", booker, true, true);
68  histograms.h2DHitRZ_W1 = new HRes2DHit ("RZ_W1", booker, true, true);
69  histograms.h2DHitRZ_W2 = new HRes2DHit ("RZ_W2", booker, true, true);
70 
71  histograms.h2DHitEff_RPhi = new HEff2DHit ("RPhi", booker);
72  histograms.h2DHitEff_RZ = new HEff2DHit ("RZ", booker);
73  histograms.h2DHitEff_RZ_W0 = new HEff2DHit ("RZ_W0", booker);
74  histograms.h2DHitEff_RZ_W1 = new HEff2DHit ("RZ_W1", booker);
75  histograms.h2DHitEff_RZ_W2 = new HEff2DHit ("RZ_W2", booker);
76  if (debug_) {
77  cout << "[DTSegment2DQuality] hitsos created " << endl;
78  }
79 }
HEff2DHit * h2DHitEff_RZ_W0
HRes2DHit * h2DHitRPhi
HRes2DHit * h2DHitRZ_W0
HEff2DHit * h2DHitEff_RZ_W2
HEff2DHit * h2DHitEff_RPhi
HRes2DHit * h2DHitRZ_W1
HEff2DHit * h2DHitEff_RZ
HEff2DHit * h2DHitEff_RZ_W1
HRes2DHit * h2DHitRZ
HRes2DHit * h2DHitRZ_W2
void DTSegment2DQuality::dqmAnalyze ( edm::Event const &  event,
edm::EventSetup const &  setup,
Histograms const &  histograms 
) const
overrideprivate

Perform the real analysis.

Definition at line 82 of file DTSegment2DQuality.cc.

References funct::abs(), gather_cfg::cout, DEFINE_FWK_MODULE, SoftLeptonByDistance_cfi::distance, PV3DBase< T, PVType, FrameType >::eta(), HRes2DHit::fill(), HEff2DHit::fill(), DTHitQualityUtils::findMuSimSegment(), DTHitQualityUtils::findMuSimSegmentDirAndPos(), DTHitQualityUtils::findSegmentAlphaAndBeta(), edm::EventSetup::get(), Histograms::h2DHitEff_RPhi, Histograms::h2DHitEff_RZ, Histograms::h2DHitEff_RZ_W0, Histograms::h2DHitEff_RZ_W1, Histograms::h2DHitEff_RZ_W2, Histograms::h2DHitRPhi, Histograms::h2DHitRZ, Histograms::h2DHitRZ_W0, Histograms::h2DHitRZ_W1, Histograms::h2DHitRZ_W2, edm::HandleBase::isValid(), DTRecSegment2D::localDirection(), DTRecSegment2D::localDirectionError(), DTRecSegment2D::localPosition(), DTRecSegment2D::localPositionError(), DTHitQualityUtils::mapMuSimHitsPerWire(), DTHitQualityUtils::mapSimHitsPerWire(), or, PV3DBase< T, PVType, FrameType >::phi(), rpcPointValidation_cfi::simHit, trackerHits::simHits, mathSSE::sqrt(), DTGeometry::superLayer(), GeomDet::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), and LocalError::xx().

82  {
83  // Get the DT Geometry
84  ESHandle<DTGeometry> dtGeom;
85  setup.get<MuonGeometryRecord>().get(dtGeom);
86 
87  // Get the SimHit collection from the event
89  event.getByToken(simHitToken_, simHits); // FIXME: second string to be removed
90 
91  // Map simHits by sl
92  map<DTSuperLayerId, PSimHitContainer > simHitsPerSl;
93  for (const auto & simHit : *simHits) {
94  // Create the id of the sl (the simHits in the DT known their wireId)
95  DTSuperLayerId slId = ((DTWireId(simHit.detUnitId())).layerId()).superlayerId();
96  // Fill the map
97  simHitsPerSl[slId].push_back(simHit);
98  }
99 
100  // Get the 2D rechits from the event
102  event.getByToken(segment2DToken_, segment2Ds);
103 
104  if (not segment2Ds.isValid()) {
105  if (debug_) { cout << "[DTSegment2DQuality]**Warning: no 2DSegments with label: " << segment2DLabel_
106  << " in this event, skipping !" << endl;
107 }
108  return;
109  }
110 
111  // Loop over all superlayers containing a segment
112  DTRecSegment2DCollection::id_iterator slId;
113  for (slId = segment2Ds->id_begin();
114  slId != segment2Ds->id_end();
115  ++slId) {
116 
117  //------------------------- simHits ---------------------------//
118  // Get simHits of each superlayer
119  PSimHitContainer simHits = simHitsPerSl[(*slId)];
120 
121  // Map simhits per wire
122  map<DTWireId, PSimHitContainer> simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
123  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
124  int nMuSimHit = muSimHitPerWire.size();
125  if (nMuSimHit == 0 or nMuSimHit == 1) {
126  if (debug_ and nMuSimHit == 1) {
127  cout << "[DTSegment2DQuality] Only " << nMuSimHit << " mu SimHit in this SL, skipping !" << endl;
128 }
129  continue; // If no or only one mu SimHit is found skip this SL
130  }
131  if (debug_) {
132  cout << "=== SL " << (*slId) << " has " << nMuSimHit << " SimHits" << endl;
133 }
134 
135  // Find outer and inner mu SimHit to build a segment
136  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
137  // Check that outermost and innermost SimHit are not the same
138  if (inAndOutSimHit.first == inAndOutSimHit.second ) {
139  cout << "[DTHitQualityUtils]***Warning: outermost and innermost SimHit are the same !" << endl;
140  continue;
141  }
142 
143  // Find direction and position of the sim Segment in SL RF
144  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
145  (*slId),&(*dtGeom));
146 
147  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
148  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
149  if (debug_) {
150  cout << " Simulated segment: local direction " << simSegmLocalDir << endl
151  << " local position " << simSegmLocalPos << endl;
152 }
153  const DTSuperLayer* superLayer = dtGeom->superLayer(*slId);
154  GlobalPoint simSegmGlobalPos = superLayer->toGlobal(simSegmLocalPos);
155 
156  // Atan(x/z) angle and x position in SL RF
157  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
158  float posSimSeg = simSegmLocalPos.x();
159  // Position (in eta, phi coordinates) in the global RF
160  float etaSimSeg = simSegmGlobalPos.eta();
161  float phiSimSeg = simSegmGlobalPos.phi();
162 
163  //---------------------------- recHits --------------------------//
164  // Get the range of rechit for the corresponding slId
165  bool recHitFound = false;
166  DTRecSegment2DCollection::range range = segment2Ds->get(*slId);
167  int nsegm = distance(range.first, range.second);
168  if (debug_) {
169  cout << " Sl: " << *slId << " has " << nsegm
170  << " 2D segments" << endl;
171 }
172 
173  if (nsegm != 0) {
174  // Find the best RecHit: look for the 2D RecHit with the angle closest
175  // to that of segment made of SimHits.
176  // RecHits must have delta alpha and delta position within 5 sigma of
177  // the residual distribution (we are looking for residuals of segments
178  // usefull to the track fit) for efficency purpose
179  const DTRecSegment2D* bestRecHit = nullptr;
180  bool bestRecHitFound = false;
181  double deltaAlpha = 99999;
182 
183  // Loop over the recHits of this slId
184  for (DTRecSegment2DCollection::const_iterator segment2D = range.first;
185  segment2D != range.second;
186  ++segment2D) {
187  // Check the dimension
188  if ((*segment2D).dimension() != 2) {
189  if (debug_) { cout << "[DTSegment2DQuality]***Error: This is not 2D segment !!!" << endl;
190 }
191  abort();
192  }
193  // Segment Local Direction and position (in SL RF)
194  LocalVector recSegDirection = (*segment2D).localDirection();
195  LocalPoint recSegPosition = (*segment2D).localPosition();
196 
197  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
198  if (debug_) {
199  cout << " RecSegment direction: " << recSegDirection << endl
200  << " position : " << recSegPosition << endl
201  << " alpha : " << recSegAlpha << endl;
202 }
203 
204  if (fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
205  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
206  bestRecHit = &(*segment2D);
207  bestRecHitFound = true;
208  }
209  } // End of Loop over all 2D RecHits
210 
211  if (bestRecHitFound) {
212  // Best rechit direction and position in SL RF
213  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
214  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
215 
216  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
217  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
218 
219  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
220 
221  if (fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle_ and
222  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos_) {
223  recHitFound = true;
224  }
225 
226  // Fill Residual histos
227  HRes2DHit *hRes = nullptr;
228  if ((*slId).superlayer() == 1 or (*slId).superlayer() == 3) { // RPhi SL
229  hRes = histograms.h2DHitRPhi;
230  } else if ((*slId).superlayer() == 2) { // RZ SL
231  histograms.h2DHitRZ->fill(angleSimSeg,
232  angleBestRHit,
233  posSimSeg,
234  bestRecHitLocalPos.x(),
235  etaSimSeg,
236  phiSimSeg,
237  sqrt(bestRecHitLocalPosErr.xx()),
238  sqrt(bestRecHitLocalDirErr.xx())
239  );
240  if (abs((*slId).wheel()) == 0) {
241  hRes = histograms.h2DHitRZ_W0;
242  } else if (abs((*slId).wheel()) == 1) {
243  hRes = histograms.h2DHitRZ_W1;
244  } else if (abs((*slId).wheel()) == 2) {
245  hRes = histograms.h2DHitRZ_W2;
246 }
247  }
248  hRes->fill(angleSimSeg,
249  angleBestRHit,
250  posSimSeg,
251  bestRecHitLocalPos.x(),
252  etaSimSeg,
253  phiSimSeg,
254  sqrt(bestRecHitLocalPosErr.xx()),
255  sqrt(bestRecHitLocalDirErr.xx())
256  );
257  }
258  } // end of if (nsegm != 0)
259 
260  // Fill Efficiency plot
261  HEff2DHit *hEff = nullptr;
262  if ((*slId).superlayer() == 1 or (*slId).superlayer() == 3) { // RPhi SL
263  hEff = histograms.h2DHitEff_RPhi;
264  } else if ((*slId).superlayer() == 2) { // RZ SL
265  histograms.h2DHitEff_RZ->fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
266  if (abs((*slId).wheel()) == 0) {
267  hEff = histograms.h2DHitEff_RZ_W0;
268  } else if (abs((*slId).wheel()) == 1) {
269  hEff = histograms.h2DHitEff_RZ_W1;
270  } else if (abs((*slId).wheel()) == 2) {
271  hEff = histograms.h2DHitEff_RZ_W2;
272  }
273  }
274  hEff->fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
275  } // End of loop over superlayers
276 }
LocalError localPositionError() const override
local position error in SL frame
void fill(float angleSimSegment, float angleRecSegment, float posSimSegment, float posRecSegment, float etaSimSegment, float phiSimSegment, float sigmaPos, float sigmaAngle)
Definition: Histograms.h:245
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
edm::EDGetTokenT< edm::PSimHitContainer > simHitToken_
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
void fill(float etaSimSegm, float phiSimSegm, float posSimSegm, float angleSimSegm, bool fillRecHit)
Definition: Histograms.h:309
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
static std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
static std::pair< const PSimHit *, const PSimHit * > findMuSimSegment(const std::map< DTWireId, const PSimHit * > &mapWireAndMuSimHit)
Find Innermost and outermost SimHit from Mu in a SL (they identify a simulated segment) ...
T sqrt(T t)
Definition: SSEVec.h:18
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< DTRecSegment2DCollection > segment2DToken_
bool isValid() const
Definition: HandleBase.h:74
static std::pair< double, double > findSegmentAlphaAndBeta(const LocalVector &direction)
static std::pair< LocalVector, LocalPoint > findMuSimSegmentDirAndPos(const std::pair< const PSimHit *, const PSimHit * > &inAndOutSimHit, const DetId detId, const DTGeometry *muonGeom)
Find direction and position of a segment (in local RF) from outer and inner mu SimHit in the RF of ob...
LocalPoint localPosition() const override
local position in SL frame
T eta() const
Definition: PV3DBase.h:76
LocalVector localDirection() const override
the local direction in SL frame
std::vector< PSimHit > PSimHitContainer
static std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
LocalError localDirectionError() const override
the local direction error (xx,xy,yy) in SL frame: only xx is not 0.
T x() const
Definition: PV3DBase.h:62
edm::InputTag segment2DLabel_
const DTSuperLayer * superLayer(const DTSuperLayerId &id) const
Return a DTSuperLayer given its id.
Definition: DTGeometry.cc:104

Member Data Documentation

bool DTSegment2DQuality::debug_
private

Definition at line 57 of file DTSegment2DQuality.h.

edm::InputTag DTSegment2DQuality::segment2DLabel_
private

Definition at line 46 of file DTSegment2DQuality.h.

edm::EDGetTokenT<DTRecSegment2DCollection> DTSegment2DQuality::segment2DToken_
private

Definition at line 48 of file DTSegment2DQuality.h.

double DTSegment2DQuality::sigmaResAngle_
private

Definition at line 54 of file DTSegment2DQuality.h.

double DTSegment2DQuality::sigmaResPos_
private

Definition at line 51 of file DTSegment2DQuality.h.

edm::InputTag DTSegment2DQuality::simHitLabel_
private

Definition at line 45 of file DTSegment2DQuality.h.

edm::EDGetTokenT<edm::PSimHitContainer> DTSegment2DQuality::simHitToken_
private

Definition at line 47 of file DTSegment2DQuality.h.