CMS 3D CMS Logo

CSCSegAlgoShowering.cc
Go to the documentation of this file.
1 
10 
14 
17 
18 #include <algorithm>
19 #include <cmath>
20 #include <iostream>
21 #include <string>
22 
23 /* Constructor
24  *
25  */
27  // debug = ps.getUntrackedParameter<bool>("CSCSegmentDebug");
28  dRPhiFineMax = ps.getParameter<double>("dRPhiFineMax");
29  dPhiFineMax = ps.getParameter<double>("dPhiFineMax");
30  tanThetaMax = ps.getParameter<double>("tanThetaMax");
31  tanPhiMax = ps.getParameter<double>("tanPhiMax");
32  maxRatioResidual = ps.getParameter<double>("maxRatioResidualPrune");
33  // maxDR = ps.getParameter<double>("maxDR");
34  maxDTheta = ps.getParameter<double>("maxDTheta");
35  maxDPhi = ps.getParameter<double>("maxDPhi");
36 }
37 
38 /* Destructor:
39  *
40  */
42 
43 /* showerSeg
44  *
45  */
47  theChamber = aChamber;
48  // Initialize parameters
49  std::vector<float> x, y, gz;
50  std::vector<int> n;
51 
52  for (int i = 0; i < 6; ++i) {
53  x.push_back(0.);
54  y.push_back(0.);
55  gz.push_back(0.);
56  n.push_back(0);
57  }
58 
59  // Loop over hits to find center-of-mass position in each layer
60  for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); ++it) {
61  const CSCRecHit2D& hit = (**it);
62  const CSCDetId id = hit.cscDetId();
63  int l_id = id.layer();
64  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
65  GlobalPoint gp = layer->toGlobal(hit.localPosition());
67 
68  ++n[l_id - 1];
69  x[l_id - 1] += lp.x();
70  y[l_id - 1] += lp.y();
71  gz[l_id - 1] += gp.z();
72  }
73 
74  // Determine center of mass for each layer and average center of mass for chamber
75  float avgChamberX = 0.;
76  float avgChamberY = 0.;
77  int n_lay = 0;
78 
79  for (unsigned i = 0; i < 6; ++i) {
80  if (n[i] < 1)
81  continue;
82 
83  x[i] = x[i] / n[i];
84  y[i] = y[i] / n[i];
85  avgChamberX += x[i];
86  avgChamberY += y[i];
87  n_lay++;
88  }
89 
90  if (n_lay > 0) {
91  avgChamberX = avgChamberX / n_lay;
92  avgChamberY = avgChamberY / n_lay;
93  }
94 
95  // Create a FakeSegment origin that points back to the IP
96  // Keep all this in global coordinates until last minute to avoid screwing up +/- signs !
97 
98  LocalPoint lpCOM(avgChamberX, avgChamberY, 0.);
99  GlobalPoint gpCOM = theChamber->toGlobal(lpCOM);
100  GlobalVector gvCOM(gpCOM.x(), gpCOM.y(), gpCOM.z());
101 
102  float Gdxdz = gpCOM.x() / gpCOM.z();
103  float Gdydz = gpCOM.y() / gpCOM.z();
104 
105  // Figure out the intersection of this vector with each layer of the chamber
106  // by projecting the vector
107  std::vector<LocalPoint> layerPoints;
108 
109  for (size_t i = 0; i != 6; ++i) {
110  // Get the layer z coordinates in global frame
111  const CSCLayer* layer = theChamber->layer(i + 1);
112  LocalPoint temp(0., 0., 0.);
113  GlobalPoint gp = layer->toGlobal(temp);
114  float layer_Z = gp.z();
115 
116  // Then compute interesection of vector with that plane
117  float layer_X = Gdxdz * layer_Z;
118  float layer_Y = Gdydz * layer_Z;
119  GlobalPoint Gintersect(layer_X, layer_Y, layer_Z);
120  LocalPoint Lintersect = theChamber->toLocal(Gintersect);
121 
122  float layerX = Lintersect.x();
123  float layerY = Lintersect.y();
124  float layerZ = Lintersect.z();
125  LocalPoint layerPoint(layerX, layerY, layerZ);
126  layerPoints.push_back(layerPoint);
127  }
128 
129  std::vector<float> r_closest;
130  std::vector<int> id;
131  for (size_t i = 0; i != 6; ++i) {
132  id.push_back(-1);
133  r_closest.push_back(9999.);
134  }
135 
136  int idx = 0;
137 
138  // Loop over all hits and find hit closest to com for that layer.
139  for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); ++it) {
140  const CSCRecHit2D& hit = (**it);
141  int layId = hit.cscDetId().layer();
142 
143  const CSCLayer* layer = theChamber->layer(layId);
144  GlobalPoint gp = layer->toGlobal(hit.localPosition());
146 
147  float d_x = lp.x() - layerPoints[layId - 1].x();
148  float d_y = lp.y() - layerPoints[layId - 1].y();
149 
150  LocalPoint diff(d_x, d_y, 0.);
151 
152  if (fabs(diff.mag()) < r_closest[layId - 1]) {
153  r_closest[layId - 1] = fabs(diff.mag());
154  id[layId - 1] = idx;
155  }
156  ++idx;
157  }
158 
159  // Now fill vector of rechits closest to center of mass:
160  protoSegment.clear();
161  idx = 0;
162 
163  // Loop over all hits and find hit closest to com for that layer.
164  for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); ++it) {
165  const CSCRecHit2D& hit = (**it);
166  int layId = hit.cscDetId().layer();
167 
168  if (idx == id[layId - 1])
169  protoSegment.push_back(*it);
170 
171  ++idx;
172  }
173 
174  // Reorder hits in protosegment
175  if (gz[0] > 0.) {
176  if (gz[0] > gz[5]) {
177  reverse(protoSegment.begin(), protoSegment.end());
178  }
179  } else if (gz[0] < 0.) {
180  if (gz[0] < gz[5]) {
181  reverse(protoSegment.begin(), protoSegment.end());
182  }
183  }
184 
185  // Fit the protosegment
187 
188  // If there is one very bad hit on segment, remove it and refit
189  if (protoSegment.size() > 4)
191 
192  // If any hit on a layer is closer to segment than original, replace it and refit
193  for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); it++) {
194  const CSCRecHit2D* h = *it;
195  int layer = h->cscDetId().layer();
196  if (isHitNearSegment(h))
198  }
199 
200  // Check again for a bad hit, and remove and refit if necessary
201  if (sfit_->nhits() > 5)
203 
204  // Does the fitted line point to the IP?
205  // If it doesn't, the algorithm has probably failed i.e. that's life!
206 
207  GlobalVector protoGlobalDir = theChamber->toGlobal(sfit_->localdir());
208  double protoTheta = protoGlobalDir.theta();
209  double protoPhi = protoGlobalDir.phi();
210  double simTheta = gvCOM.theta();
211  double simPhi = gvCOM.phi();
212 
213  float dTheta = fabs(protoTheta - simTheta);
214  float dPhi = fabs(protoPhi - simPhi);
215  // float dR = sqrt(dEta*dEta + dPhi*dPhi);
216 
217  // Flag the segment with chi2=-1 of the segment isn't pointing toward origin
218  // i.e. flag that the algorithm has probably just failed (I presume we expect
219  // a segment to point to the IP if the algorithm is successful!)
220 
221  double theFlag = -1.;
222  if (dTheta > maxDTheta || dPhi > maxDPhi) {
223  } else {
224  theFlag = sfit_->chi2(); // If it points to IP, just pass fit chi2 as usual
225  }
226 
227  // Create an actual CSCSegment - retrieve all info from the fit
229  delete sfit_;
230  sfit_ = nullptr;
231 
232  return temp;
233 }
234 
235 /* isHitNearSegment
236  *
237  * Compare rechit with expected position from proto_segment
238  */
240  const CSCLayer* layer = theChamber->layer(hit->cscDetId().layer());
241 
242  // hit phi position in global coordinates
243  GlobalPoint Hgp = layer->toGlobal(hit->localPosition());
244  double Hphi = Hgp.phi();
245  if (Hphi < 0.)
246  Hphi += 2. * M_PI;
247  LocalPoint Hlp = theChamber->toLocal(Hgp);
248  double z = Hlp.z();
249 
250  double LocalX = sfit_->xfit(z);
251  double LocalY = sfit_->yfit(z);
252  LocalPoint Slp(LocalX, LocalY, z);
253  GlobalPoint Sgp = theChamber->toGlobal(Slp);
254  double Sphi = Sgp.phi();
255  if (Sphi < 0.)
256  Sphi += 2. * M_PI;
257  double R = sqrt(Sgp.x() * Sgp.x() + Sgp.y() * Sgp.y());
258 
259  double deltaPhi = Sphi - Hphi;
260  if (deltaPhi > 2. * M_PI)
261  deltaPhi -= 2. * M_PI;
262  if (deltaPhi < -2. * M_PI)
263  deltaPhi += 2. * M_PI;
264  if (deltaPhi < 0.)
265  deltaPhi = -deltaPhi;
266 
267  double RdeltaPhi = R * deltaPhi;
268 
269  if (RdeltaPhi < dRPhiFineMax && deltaPhi < dPhiFineMax)
270  return true;
271 
272  return false;
273 }
274 
275 /* Method addHit
276  *
277  * Test if can add hit to proto segment. If so, try to add it.
278  *
279  */
281  // Return true if hit was added successfully
282  // Return false if there is already a hit on the same layer
283 
284  bool ok = true;
285 
286  // Test that we are not trying to add the same hit again
287  for (ChamberHitContainer::const_iterator it = protoSegment.begin(); it != protoSegment.end(); it++)
288  if (aHit == (*it))
289  return false;
290 
291  protoSegment.push_back(aHit);
292 
293  return ok;
294 }
295 
296 /* Method compareProtoSegment
297  *
298  * For hit coming from the same layer of an existing hit within the proto segment
299  * test if achieve better chi^2 by using this hit than the other
300  *
301  */
303  // Try adding the hit to existing segment, and removing one from the same layer
304  ChamberHitContainer::iterator it;
305  for (it = protoSegment.begin(); it != protoSegment.end();) {
306  if ((*it)->cscDetId().layer() == layer) {
307  it = protoSegment.erase(it);
308  } else {
309  ++it;
310  }
311  }
312  bool ok = addHit(h, layer);
313  if (ok) {
314  CSCSegFit* newfit = new CSCSegFit(theChamber, protoSegment);
315  newfit->fit();
316  if (newfit->chi2() > sfit_->chi2()) {
317  // new fit is worse: revert to old fit
318  delete newfit;
319  } else {
320  // new fit is better
321  delete sfit_;
322  sfit_ = newfit;
323  }
324  }
325 }
326 
327 // Look for a hit with a large deviation from fit
328 
330  //@@ THIS FUNCTION HAS 3 RETURNS PATHS!
331 
332  // Only prune if have at least 5 hits
333  if (protoSegment.size() < 5)
334  return;
335 
336  // Now Study residuals
337 
338  float maxResidual = 0.;
339  float sumResidual = 0.;
340  int nHits = 0;
341 
342  ChamberHitContainer::iterator ih;
343  ChamberHitContainer::iterator ibad = protoSegment.end();
344 
345  for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
346  const CSCRecHit2D& hit = (**ih);
347  const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
348  GlobalPoint gp = layer->toGlobal(hit.localPosition());
350 
351  double u = lp.x();
352  double v = lp.y();
353  double z = lp.z();
354 
355  // double du = segfit->xdev( u, z );
356  // double dv = segfit->ydev( v, z );
357 
358  float residual = sfit_->Rdev(u, v, z); // == sqrt(du*du + dv*dv)
359 
360  sumResidual += residual;
361  ++nHits;
362  if (residual > maxResidual) {
363  maxResidual = residual;
364  ibad = ih;
365  }
366  }
367 
368  float corrAvgResidual = (sumResidual - maxResidual) / (nHits - 1);
369 
370  // Keep all hits
371  if (maxResidual / corrAvgResidual < maxRatioResidual)
372  return;
373 
374  // Drop worst hit
375  if (ibad != protoSegment.end())
376  protoSegment.erase(ibad);
377 
378  // Make a new fit
380 
381  return;
382 }
383 
385  // Create fit for the hits in the protosegment & run it
386  delete sfit_;
388  sfit_->fit();
389 }
void fit(void)
Definition: CSCSegFit.cc:13
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const CSCLayer * layer(CSCDetId id) const
Return the layer corresponding to the given id.
Definition: CSCChamber.cc:30
double chi2(void) const
Definition: CSCSegFit.h:82
void compareProtoSegment(const CSCRecHit2D *h, int layer)
virtual ~CSCSegAlgoShowering()
Destructor.
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
CSCSegment showerSeg(const CSCChamber *aChamber, const ChamberHitContainer &rechits)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
LocalPoint intercept() const
Definition: CSCSegFit.h:84
float yfit(float z) const
Definition: CSCSegFit.cc:432
bool addHit(const CSCRecHit2D *hit, int layer)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T sqrt(T t)
Definition: SSEVec.h:19
CSCSegAlgoShowering(const edm::ParameterSet &ps)
Constructor.
std::vector< const CSCRecHit2D * > ChamberHitContainer
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
#define M_PI
AlgebraicSymMatrix covarianceMatrix(void)
Definition: CSCSegFit.cc:352
float xfit(float z) const
Definition: CSCSegFit.cc:426
LocalVector localdir() const
Definition: CSCSegFit.h:85
CSCSetOfHits hits(void) const
Definition: CSCSegFit.h:79
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
size_t nhits(void) const
Definition: CSCSegFit.h:81
bool isHitNearSegment(const CSCRecHit2D *h) const
Utility functions.
float Rdev(float x, float y, float z) const
Definition: CSCSegFit.cc:438
const CSCChamber * theChamber
ChamberHitContainer protoSegment
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72