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 } // namespace dtsegment2dsl
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  edm::Run const &run,
60  edm::EventSetup const &setup,
61  Histograms &histograms) const {
62  // Book the histos
63  histograms.h2DHitSuperPhi = std::make_unique<HRes2DHit>("SuperPhi", booker, doall_, local_);
64  if (doall_) {
65  histograms.h2DHitEff_SuperPhi = std::make_unique<HEff2DHit>("SuperPhi", booker);
66  }
67 }
68 
69 // The real analysis
71  edm::EventSetup const &setup,
72  Histograms const &histograms) const {
73  // Get the DT Geometry
74  ESHandle<DTGeometry> dtGeom;
75  setup.get<MuonGeometryRecord>().get(dtGeom);
76 
77  // Get the SimHit collection from the event
79  event.getByToken(simHitToken_, simHits); // FIXME: second string to be removed
80 
81  // Map simHits by chamber
82  map<DTChamberId, PSimHitContainer> simHitsPerCh;
83  for (const auto &simHit : *simHits) {
84  // Create the id of the chamber (the simHits in the DT known their wireId)
85  DTChamberId chamberId = (((DTWireId(simHit.detUnitId())).layerId()).superlayerId()).chamberId();
86  // Fill the map
87  simHitsPerCh[chamberId].push_back(simHit);
88  }
89 
90  // Get the 4D rechits from the event
92  event.getByToken(segment4DToken_, segment4Ds);
93 
94  if (!segment4Ds.isValid()) {
95  if (debug_) {
96  cout << "[DTSegment2DSLPhiQuality]**Warning: no 4D Segments with label: " << segment4DLabel_
97  << " in this event, skipping!" << endl;
98  }
99  return;
100  }
101 
102  // Loop over all chambers containing a segment
104  for (chamberId = segment4Ds->id_begin(); chamberId != segment4Ds->id_end(); ++chamberId) {
105  //------------------------- simHits ---------------------------//
106  // Get simHits of each chamber
107  PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
108 
109  // Map simhits per wire
110  map<DTWireId, PSimHitContainer> simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
111  map<DTWireId, const PSimHit *> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
112  int nMuSimHit = muSimHitPerWire.size();
113  if (nMuSimHit == 0 || nMuSimHit == 1) {
114  if (debug_ && nMuSimHit == 1) {
115  cout << "[DTSegment2DSLPhiQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
116  }
117  continue; // If no or only one mu SimHit is found skip this chamber
118  }
119  if (debug_) {
120  cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
121  }
122 
123  // Find outer and inner mu SimHit to build a segment
124  pair<const PSimHit *, const PSimHit *> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
125 
126  // Find direction and position of the sim Segment in Chamber RF
127  pair<LocalVector, LocalPoint> dirAndPosSimSegm =
128  DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit, (*chamberId), &(*dtGeom));
129 
130  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
131  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
132  const DTChamber *chamber = dtGeom->chamber(*chamberId);
133  GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
134 
135  // Atan(x/z) angle and x position in SL RF
136  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
137  float posSimSeg = simSegmLocalPos.x();
138  // Position (in eta, phi coordinates) in lobal RF
139  float etaSimSeg = simSegmGlobalPos.eta();
140  float phiSimSeg = simSegmGlobalPos.phi();
141 
142  if (debug_) {
143  cout << " Simulated segment: local direction " << simSegmLocalDir << endl
144  << " local position " << simSegmLocalPos << endl
145  << " angle " << angleSimSeg << endl;
146  }
147 
148  //---------------------------- recHits --------------------------//
149  // Get the range of rechit for the corresponding chamberId
150  bool recHitFound = false;
151  DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
152  int nsegm = distance(range.first, range.second);
153  if (debug_) {
154  cout << " Chamber: " << *chamberId << " has " << nsegm << " 4D segments" << endl;
155  }
156 
157  if (nsegm != 0) {
158  // Find the best RecHit: look for the 4D RecHit with the phi angle closest
159  // to that of segment made of SimHits.
160  // RecHits must have delta alpha and delta position within 5 sigma of
161  // the residual distribution (we are looking for residuals of segments
162  // usefull to the track fit) for efficency purpose
163  const DTRecSegment2D *bestRecHit = nullptr;
164  bool bestRecHitFound = false;
165  double deltaAlpha = 99999;
166 
167  // Loop over the recHits of this chamberId
168  for (DTRecSegment4DCollection::const_iterator segment4D = range.first; segment4D != range.second; ++segment4D) {
169  // Check the dimension
170  if ((*segment4D).dimension() != 4) {
171  if (debug_) {
172  cout << "[DTSegment2DSLPhiQuality]***Error: This is not 4D "
173  "segment!!!"
174  << endl;
175  }
176  continue;
177  }
178 
179  // Get 2D superPhi segments from 4D segments
180  const DTChamberRecSegment2D *phiSegment2D = (*segment4D).phiSegment();
181  if ((*phiSegment2D).dimension() != 2) {
182  if (debug_) {
183  cout << "[DTSegment2DQuality]***Error: This is not 2D segment!!!" << endl;
184  }
185  abort();
186  }
187 
188  // Segment Local Direction and position (in Chamber RF)
189  LocalVector recSegDirection = (*phiSegment2D).localDirection();
190 
191  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
192  if (debug_) {
193  cout << " RecSegment direction: " << recSegDirection << endl
194  << " position : " << (*phiSegment2D).localPosition() << endl
195  << " alpha : " << recSegAlpha << endl;
196  }
197 
198  if (fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
199  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
200  bestRecHit = &(*phiSegment2D);
201  bestRecHitFound = true;
202  }
203  } // End of Loop over all 4D RecHits of this chambers
204 
205  if (bestRecHitFound) {
206  // Best rechit direction and position in Chamber RF
207  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
208  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
209 
210  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
211  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
212 
213  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
214  if (fabs(angleBestRHit - angleSimSeg) < 5 * sigmaResAngle_ &&
215  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5 * sigmaResPos_) {
216  recHitFound = true;
217  }
218 
219  // Fill Residual histos
220  histograms.h2DHitSuperPhi->fill(angleSimSeg,
221  angleBestRHit,
222  posSimSeg,
223  bestRecHitLocalPos.x(),
224  etaSimSeg,
225  phiSimSeg,
226  sqrt(bestRecHitLocalPosErr.xx()),
227  sqrt(bestRecHitLocalDirErr.xx()));
228  }
229  } // end of if (nsegm!= 0)
230 
231  // Fill Efficiency plot
232  if (doall_) {
233  histograms.h2DHitEff_SuperPhi->fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
234  }
235  } // End of loop over chambers
236 }
237 
238 // 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::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
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::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
identifier iterator
Definition: RangeMap.h:135
std::atomic< bool > debug
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
std::pair< const PSimHit *, const PSimHit * > findMuSimSegment(const std::map< DTWireId, const PSimHit * > &mapWireAndMuSimHit)
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)
T eta() const
Definition: PV3DBase.h:76
HLT enums.
LocalVector localDirection() const override
the local direction in SL frame
T get() const
Definition: EventSetup.h:71
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:45