CMS 3D CMS Logo

DTSegment2DSLPhiQuality.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi and G. Cerminara - INFN Torino
5  */
6 
7 #include <iostream>
8 #include <map>
9 
21 
23 #include "Histograms.h"
24 
25 using namespace std;
26 using namespace edm;
27 
28 namespace dtsegment2dsl {
29 struct Histograms {
30  std::unique_ptr<HRes2DHit> h2DHitSuperPhi;
31  std::unique_ptr<HEff2DHit> h2DHitEff_SuperPhi;
32 };
33 }
34 
35 using namespace dtsegment2dsl;
36 
37 // Constructor
39  // Get the debug parameter for verbose output
40  debug_ = pset.getUntrackedParameter<bool>("debug");
41  DTHitQualityUtils::debug = debug_;
42 
43  // the name of the simhit collection
44  simHitLabel_ = pset.getUntrackedParameter<InputTag>("simHitLabel");
45  simHitToken_ = consumes<PSimHitContainer>(pset.getUntrackedParameter<InputTag>("simHitLabel"));
46  // the name of the 2D rec hit collection
47  segment4DLabel_ = pset.getUntrackedParameter<InputTag>("segment4DLabel");
48  segment4DToken_ = consumes<DTRecSegment4DCollection>(pset.getUntrackedParameter<InputTag>("segment4DLabel"));
49 
50  // sigma resolution on position
51  sigmaResPos_ = pset.getParameter<double>("sigmaResPos");
52  // sigma resolution on angle
53  sigmaResAngle_ = pset.getParameter<double>("sigmaResAngle");
54  doall_ = pset.getUntrackedParameter<bool>("doall", false);
55  local_ = pset.getUntrackedParameter<bool>("local", false);
56 }
57 
59  // Book the histos
60  histograms.h2DHitSuperPhi = std::make_unique<HRes2DHit> ("SuperPhi", booker, doall_, local_);
61  if (doall_) { histograms.h2DHitEff_SuperPhi = std::make_unique<HEff2DHit> ("SuperPhi", booker);
62 }
63 }
64 
65 // The real analysis
67  // Get the DT Geometry
68  ESHandle<DTGeometry> dtGeom;
69  setup.get<MuonGeometryRecord>().get(dtGeom);
70 
71  // Get the SimHit collection from the event
73  event.getByToken(simHitToken_, simHits); // FIXME: second string to be removed
74 
75  // Map simHits by chamber
76  map<DTChamberId, PSimHitContainer > simHitsPerCh;
77  for (const auto & simHit : *simHits) {
78  // Create the id of the chamber (the simHits in the DT known their wireId)
79  DTChamberId chamberId = (((DTWireId(simHit.detUnitId())).layerId()).superlayerId()).chamberId();
80  // Fill the map
81  simHitsPerCh[chamberId].push_back(simHit);
82  }
83 
84  // Get the 4D rechits from the event
86  event.getByToken(segment4DToken_, segment4Ds);
87 
88  if (!segment4Ds.isValid()) {
89  if (debug_) {
90  cout << "[DTSegment2DSLPhiQuality]**Warning: no 4D Segments with label: " << segment4DLabel_
91  << " in this event, skipping!" << endl;
92 }
93  return;
94  }
95 
96  // Loop over all chambers containing a segment
97  DTRecSegment4DCollection::id_iterator chamberId;
98  for (chamberId = segment4Ds->id_begin();
99  chamberId != segment4Ds->id_end();
100  ++chamberId) {
101 
102  //------------------------- simHits ---------------------------//
103  // Get simHits of each chamber
104  PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
105 
106  // Map simhits per wire
107  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
108  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
109  int nMuSimHit = muSimHitPerWire.size();
110  if (nMuSimHit == 0 || nMuSimHit == 1) {
111  if (debug_ && nMuSimHit == 1) {
112  cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
113 }
114  continue; // If no or only one mu SimHit is found skip this chamber
115  }
116  if (debug_) {
117  cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
118 }
119 
120  // Find outer and inner mu SimHit to build a segment
121  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
122 
123  // Find direction and position of the sim Segment in Chamber RF
124  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
125  (*chamberId),&(*dtGeom));
126 
127  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
128  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
129  const DTChamber* chamber = dtGeom->chamber(*chamberId);
130  GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
131 
132  // Atan(x/z) angle and x position in SL RF
133  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
134  float posSimSeg = simSegmLocalPos.x();
135  // Position (in eta, phi coordinates) in lobal RF
136  float etaSimSeg = simSegmGlobalPos.eta();
137  float phiSimSeg = simSegmGlobalPos.phi();
138 
139  if (debug_) {
140  cout << " Simulated segment: local direction " << simSegmLocalDir << endl
141  << " local position " << simSegmLocalPos << endl
142  << " angle " << angleSimSeg << endl;
143 }
144 
145  //---------------------------- recHits --------------------------//
146  // Get the range of rechit for the corresponding chamberId
147  bool recHitFound = false;
148  DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
149  int nsegm = distance(range.first, range.second);
150  if (debug_) {
151  cout << " Chamber: " << *chamberId << " has " << nsegm
152  << " 4D segments" << endl;
153 }
154 
155  if (nsegm!= 0) {
156  // Find the best RecHit: look for the 4D RecHit with the phi angle closest
157  // to that of segment made of SimHits.
158  // RecHits must have delta alpha and delta position within 5 sigma of
159  // the residual distribution (we are looking for residuals of segments
160  // usefull to the track fit) for efficency purpose
161  const DTRecSegment2D* bestRecHit = nullptr;
162  bool bestRecHitFound = false;
163  double deltaAlpha = 99999;
164 
165  // Loop over the recHits of this chamberId
166  for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
167  segment4D!= range.second;
168  ++segment4D) {
169  // Check the dimension
170  if ((*segment4D).dimension() != 4) {
171  if (debug_) { cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D segment!!!" << endl;
172 }
173  continue;
174  }
175 
176  // Get 2D superPhi segments from 4D segments
177  const DTChamberRecSegment2D* phiSegment2D = (*segment4D).phiSegment();
178  if ((*phiSegment2D).dimension() != 2) {
179  if (debug_) { cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
180 }
181  abort();
182  }
183 
184  // Segment Local Direction and position (in Chamber RF)
185  LocalVector recSegDirection = (*phiSegment2D).localDirection();
186 
187  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
188  if (debug_) {
189  cout << " RecSegment direction: " << recSegDirection << endl
190  << " position : " << (*phiSegment2D).localPosition() << endl
191  << " alpha : " << recSegAlpha << endl;
192 }
193 
194  if (fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
195  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
196  bestRecHit = &(*phiSegment2D);
197  bestRecHitFound = true;
198  }
199  } // End of Loop over all 4D RecHits of this chambers
200 
201  if (bestRecHitFound) {
202  // Best rechit direction and position in Chamber RF
203  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
204  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
205 
206  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
207  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
208 
209  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
210  if (fabs(angleBestRHit - angleSimSeg) < 5*sigmaResAngle_ &&
211  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5*sigmaResPos_) {
212  recHitFound = true;
213  }
214 
215  // Fill Residual histos
216  histograms.h2DHitSuperPhi->fill(angleSimSeg,
217  angleBestRHit,
218  posSimSeg,
219  bestRecHitLocalPos.x(),
220  etaSimSeg,
221  phiSimSeg,
222  sqrt(bestRecHitLocalPosErr.xx()),
223  sqrt(bestRecHitLocalDirErr.xx())
224  );
225  }
226  } // end of if (nsegm!= 0)
227 
228  // Fill Efficiency plot
229  if (doall_) {
230  histograms.h2DHitEff_SuperPhi->fill(etaSimSeg,
231  phiSimSeg,
232  posSimSeg,
233  angleSimSeg,
234  recHitFound);
235  }
236  } // End of loop over chambers
237 }
238 
239 // declare this as a framework plugin
std::pair< double, double > findSegmentAlphaAndBeta(const LocalVector &direction)
LocalError localPositionError() const override
local position error in SL frame
T getParameter(std::string const &) const
DTSegment2DSLPhiQuality(const edm::ParameterSet &pset)
Constructor.
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, dtsegment2dsl::Histograms const &) const override
Perform the real analysis.
T getUntrackedParameter(std::string const &, T const &) const
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:117
std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::unique_ptr< HEff2DHit > h2DHitEff_SuperPhi
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
std::atomic< bool > debug
std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
T sqrt(T t)
Definition: SSEVec.h:18
bool isValid() const
Definition: HandleBase.h:74
std::unique_ptr< HRes2DHit > h2DHitSuperPhi
LocalPoint localPosition() const override
local position in SL frame
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...
T eta() const
Definition: PV3DBase.h:76
HLT enums.
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) ...
LocalVector localDirection() const override
the local direction in SL frame
T get() const
Definition: EventSetup.h:68
void bookHistograms(DQMStore::ConcurrentBooker &, edm::Run const &, edm::EventSetup const &, dtsegment2dsl::Histograms &) const override
Book the DQM plots.
std::vector< PSimHit > PSimHitContainer
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
Definition: event.py:1
Definition: Run.h:44