test
CMS 3D CMS Logo

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