CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SETPatternRecognition.cc
Go to the documentation of this file.
1 
16 #include "TMath.h"
17 
18 
19 using namespace edm;
20 using namespace std;
21 
23 : MuonSeedVPatternRecognition(parameterSet.getParameter<ParameterSet>("SETTrajBuilderParameters").getParameter<ParameterSet>("FilterParameters"))
24 {
25  const string metname = "Muon|RecoMuon|SETPatternRecognition";
26  // Parameter set for the Builder
27  ParameterSet trajectoryBuilderParameters = parameterSet.getParameter<ParameterSet>("SETTrajBuilderParameters");
28  // The inward-outward fitter (starts from seed state)
29  ParameterSet filterPSet = trajectoryBuilderParameters.getParameter<ParameterSet>("FilterParameters");
30  maxActiveChambers = filterPSet.getParameter<int>("maxActiveChambers");
31  useRPCs = filterPSet.getParameter<bool>("EnableRPCMeasurement");
32  DTRecSegmentLabel = filterPSet.getParameter<edm::InputTag>("DTRecSegmentLabel");
33  CSCRecSegmentLabel = filterPSet.getParameter<edm::InputTag>("CSCRecSegmentLabel");
34  RPCRecSegmentLabel = filterPSet.getParameter<edm::InputTag>("RPCRecSegmentLabel");
35 
36  outsideChamberErrorScale = filterPSet.getParameter<double>("OutsideChamberErrorScale");
37  minLocalSegmentAngle = filterPSet.getParameter<double>("MinLocalSegmentAngle");
38  //----
39 
40 }
41 
42 //---- it is a "cluster recognition" actually; the pattern recognition is in SETFilter
44  std::vector<MuonRecHitContainer> & segments_clusters)
45 {
46  const string metname = "Muon|RecoMuon|SETMuonSeedSeed";
47 
48  //---- Build collection of all segments
49  MuonRecHitContainer muonRecHits;
50  MuonRecHitContainer muonRecHits_DT2D_hasPhi;
51  MuonRecHitContainer muonRecHits_DT2D_hasZed;
52  MuonRecHitContainer muonRecHits_RPC;
53 
54  // ********************************************;
55  // Get the DT-Segment collection from the Event
56  // ********************************************;
57 
59  event.getByLabel(DTRecSegmentLabel, dtRecHits);
60  std::vector<DTChamberId> chambers_DT;
61  std::vector<DTChamberId>::const_iterator chIt_DT;
62  for (DTRecSegment4DCollection::const_iterator rechit = dtRecHits->begin(); rechit!=dtRecHits->end();++rechit) {
63  bool insert = true;
64  for(chIt_DT=chambers_DT.begin(); chIt_DT != chambers_DT.end(); ++chIt_DT){
65  if (
66  ((*rechit).chamberId().wheel()) == ((*chIt_DT).wheel()) &&
67  ((*rechit).chamberId().station() == (*chIt_DT).station()) &&
68  ((*rechit).chamberId().sector() == (*chIt_DT).sector())){
69  insert = false;
70  }
71  }
72  if (insert){
73  chambers_DT.push_back((*rechit).chamberId());
74  }
75  if(segmentCleaning((*rechit).geographicalId(),
76  rechit->localPosition(), rechit->localPositionError(),
77  rechit->localDirection(), rechit->localDirectionError(),
78  rechit->chi2(), rechit->degreesOfFreedom())){
79  continue;
80  }
81  if( (rechit->hasZed() && rechit->hasPhi()) ) {
82  muonRecHits.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
83  }
84  else if(rechit->hasZed()) {
85  muonRecHits_DT2D_hasZed.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
86  }
87  else if(rechit->hasPhi()) { // safeguard
88  muonRecHits_DT2D_hasPhi.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
89  }
90  else {
91  //std::cout<<"Warning in "<<metname<<": DT segment which claims to have neither phi nor Z."<<std::endl;
92  }
93  }
94  //std::cout<<"DT done"<<std::endl;
95 
96  // ********************************************;
97  // Get the CSC-Segment collection from the event
98  // ********************************************;
99 
101  event.getByLabel(CSCRecSegmentLabel, cscSegments);
102  std::vector<CSCDetId> chambers_CSC;
103  std::vector<CSCDetId>::const_iterator chIt_CSC;
104  for(CSCSegmentCollection::const_iterator rechit=cscSegments->begin(); rechit != cscSegments->end(); ++rechit) {
105  bool insert = true;
106  for(chIt_CSC=chambers_CSC.begin(); chIt_CSC != chambers_CSC.end(); ++chIt_CSC){
107  if (((*rechit).cscDetId().chamber() == (*chIt_CSC).chamber()) &&
108  ((*rechit).cscDetId().station() == (*chIt_CSC).station()) &&
109  ((*rechit).cscDetId().ring() == (*chIt_CSC).ring()) &&
110  ((*rechit).cscDetId().endcap() == (*chIt_CSC).endcap())){
111  insert = false;
112  }
113  }
114  if (insert){
115  chambers_CSC.push_back((*rechit).cscDetId().chamberId());
116  }
117  if(segmentCleaning((*rechit).geographicalId(),
118  rechit->localPosition(), rechit->localPositionError(),
119  rechit->localDirection(), rechit->localDirectionError(),
120  rechit->chi2(), rechit->degreesOfFreedom())){
121  continue;
122  }
123  muonRecHits.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
124  }
125  //std::cout<<"CSC done"<<std::endl;
126 
127  // ********************************************;
128  // Get the RPC-Hit collection from the event
129  // ********************************************;
130 
132  event.getByLabel(RPCRecSegmentLabel, rpcRecHits);
133  if(useRPCs){
134  for(RPCRecHitCollection::const_iterator rechit=rpcRecHits->begin(); rechit != rpcRecHits->end(); ++rechit) {
135  // RPCs are special
136  const LocalVector localDirection(0.,0.,1.);
137  const LocalError localDirectionError (0.,0.,0.);
138  const double chi2 = 1.;
139  const int ndf = 1;
140  if(segmentCleaning((*rechit).geographicalId(),
141  rechit->localPosition(), rechit->localPositionError(),
142  localDirection, localDirectionError,
143  chi2, ndf)){
144  continue;
145  }
146  muonRecHits_RPC.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
147  }
148  }
149  //std::cout<<"RPC done"<<std::endl;
150  //
151  if(int(chambers_DT.size() + chambers_CSC.size()) > maxActiveChambers){
152  // std::cout <<" Too many active chambers : nDT = "<<chambers_DT.size()<<
153  // " nCSC = "<<chambers_CSC.size()<<" Skip them all."<<std::endl;
154  edm::LogWarning("tooManyActiveChambers")<<" Too many active chambers : nDT = "<<chambers_DT.size()
155  <<" nCSC = "<<chambers_CSC.size()<<" Skip them all.";
156  muonRecHits.clear();
157  muonRecHits_DT2D_hasPhi.clear();
158  muonRecHits_DT2D_hasZed.clear();
159  muonRecHits_RPC.clear();
160  }
161  //---- Find "pre-clusters" from all segments; these contain potential muon candidates
162 
163  //---- From all the hits (i.e. segments; sometimes "rechits" is also used with the same meaning;
164  //---- this convention has meaning in the global reconstruction though could be misleading
165  //---- from a local reconstruction point of view; "local rechits" are used in the backward fit only)
166  //---- make clusters of hits; a trajectory could contain hits from one cluster only
167 
168  // the clustering procedure is very similar to the one used in the segment reconstruction
169 
170  bool useDT2D_hasPhi = true;
171  bool useDT2D_hasZed = true;
172  double dXclusBoxMax = 0.60; // phi - can be as large as 15 - 20 degrees for 6 GeV muons
173  double dYclusBoxMax = 0.;
174 
175  // this is the main selection criteria; the value of 0.02 rad seems wide enough to
176  // contain any hit from a passing muon and still narrow enough to remove good part of
177  // possible "junk" hits
178  // (Comment: it may be better to allow maximum difference between any two hits in a trajectory
179  // to be 0.02 or 0.04 or ...; currently the requirement below is imposed on two consecutive hits)
180 
181  dYclusBoxMax = 0.02; // theta // hardoded - remove it!
182 
183  // X and Y are distance variables - we use eta and phi here
184 
185  float dXclus = 0.0;
186  float dXclus_box = 0.0;
187  float dYclus_box = 0.0;
188 
190 
191  std::vector< MuonRecHitContainer > seeds;
192 
193  std::vector<float> running_meanX;
194  std::vector<float> running_meanY;
195 
196  std::vector<float> seed_minX;
197  std::vector<float> seed_maxX;
198  std::vector<float> seed_minY;
199  std::vector<float> seed_maxY;
200 
201  // split rechits into subvectors and return vector of vectors:
202  // Loop over rechits
203  // Create one seed per hit
204  for (MuonRecHitContainer::const_iterator it = muonRecHits.begin(); it != muonRecHits.end(); ++it ) {
205 
206  // try to avoid using 2D DT segments. We will add them later to the
207  // clusters they are most likely to belong to. Might need to add them
208  // to more than just one cluster, if we find them to be consistent with
209  // more than one. This would lead to an implicit sharing of hits between
210  // SA muon candidates.
211 
212  temp.clear();
213 
214  temp.push_back((*it));
215 
216  seeds.push_back(temp);
217 
218  // First added hit in seed defines the mean to which the next hit is compared
219  // for this seed.
220 
221  running_meanX.push_back( (*it)->globalPosition().phi() );
222  running_meanY.push_back( (*it)->globalPosition().theta() );
223 
224  // set min/max X and Y for box containing the hits in the precluster:
225  seed_minX.push_back( (*it)->globalPosition().phi() );
226  seed_maxX.push_back( (*it)->globalPosition().phi() );
227  seed_minY.push_back( (*it)->globalPosition().theta() );
228  seed_maxY.push_back( (*it)->globalPosition().theta() );
229  }
230 
231  // merge clusters that are too close
232  // measure distance between final "running mean"
233  for(unsigned int NNN = 0; NNN < seeds.size(); ++NNN) {
234 
235  for(unsigned int MMM = NNN+1; MMM < seeds.size(); ++MMM) {
236  if(running_meanX[MMM] == 999999. || running_meanX[NNN] == 999999. ) {
237  // LogDebug("CSC") << "CSCSegmentST::clusterHits: Warning: Skipping used seeds, this should happen - inform developers!\n";
238  //std::cout<<"We should never see this line now!!!"<<std::endl;
239  continue; //skip seeds that have been used
240  }
241 
242  // Some complications for using phi as a clustering variable due to wrap-around (-pi = pi)
243  // Define temporary mean, min, and max variables for the cluster which could be merged (NNN)
244  double temp_meanX = running_meanX[NNN];
245  double temp_minX = seed_minX[NNN];
246  double temp_maxX = seed_maxX[NNN];
247 
248  // check if the difference between the two phi values is greater than pi
249  // if so, need to shift temporary values by 2*pi to make a valid comparison
250  dXclus = running_meanX[NNN] - running_meanX[MMM];
251  if (dXclus > TMath::Pi()) {
252  temp_meanX = temp_meanX - 2.*TMath::Pi();
253  temp_minX = temp_minX - 2.*TMath::Pi();
254  temp_maxX = temp_maxX - 2.*TMath::Pi();
255  }
256  if (dXclus < -TMath::Pi()) {
257  temp_meanX = temp_meanX + 2.*TMath::Pi();
258  temp_minX = temp_minX + 2.*TMath::Pi();
259  temp_maxX = temp_maxX + 2.*TMath::Pi();
260  }
261 
262  // // calculate cut criteria for simple running mean distance cut:
263  // // not sure that these values are really used anywhere
264 
265  // calculate minmal distance between precluster boxes containing the hits:
266  // use the temp variables from above for phi of the NNN cluster
267  if ( temp_meanX > running_meanX[MMM] ) dXclus_box = temp_minX - seed_maxX[MMM];
268  else dXclus_box = seed_minX[MMM] - temp_maxX;
269  if ( running_meanY[NNN] > running_meanY[MMM] ) dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
270  else dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
271 
272 
273  if( dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax ) {
274  // merge clusters!
275  // merge by adding seed NNN to seed MMM and erasing seed NNN
276 
277  // calculate running mean for the merged seed:
278  // use the temp variables from above for phi of the NNN cluster
279  running_meanX[MMM] = (temp_meanX*seeds[NNN].size() + running_meanX[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
280  running_meanY[MMM] = (running_meanY[NNN]*seeds[NNN].size() + running_meanY[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
281 
282  // update min/max X and Y for box containing the hits in the merged cluster:
283  // use the temp variables from above for phi of the NNN cluster
284  if ( temp_minX <= seed_minX[MMM] ) seed_minX[MMM] = temp_minX;
285  if ( temp_maxX > seed_maxX[MMM] ) seed_maxX[MMM] = temp_maxX;
286  if ( seed_minY[NNN] <= seed_minY[MMM] ) seed_minY[MMM] = seed_minY[NNN];
287  if ( seed_maxY[NNN] > seed_maxY[MMM] ) seed_maxY[MMM] = seed_maxY[NNN];
288 
289  // now check to see if the running mean has moved outside of the allowed -pi to pi region
290  // if so, then adjust shift all values up or down by 2 * pi
291  if (running_meanX[MMM] > TMath::Pi()) {
292  running_meanX[MMM] = running_meanX[MMM] - 2.*TMath::Pi();
293  seed_minX[MMM] = seed_minX[MMM] - 2.*TMath::Pi();
294  seed_maxX[MMM] = seed_maxX[MMM] - 2.*TMath::Pi();
295  }
296  if (running_meanX[MMM] < -TMath::Pi()) {
297  running_meanX[MMM] = running_meanX[MMM] + 2.*TMath::Pi();
298  seed_minX[MMM] = seed_minX[MMM] + 2.*TMath::Pi();
299  seed_maxX[MMM] = seed_maxX[MMM] + 2.*TMath::Pi();
300  }
301 
302  // add seed NNN to MMM (lower to larger number)
303  seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
304 
305  // mark seed NNN as used (at the moment just set running mean to 999999.)
306  running_meanX[NNN] = 999999.;
307  running_meanY[NNN] = 999999.;
308  // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to
309  // next seed (NNN+1)
310  break;
311  }
312 
313  }
314  }
315  bool tooCloseClusters = false;
316  if(seeds.size()>1){
317  std::vector <double> seedTheta(seeds.size());
318  for(unsigned int iSeed = 0;iSeed<seeds.size();++iSeed){
319  seedTheta[iSeed] = seeds[iSeed][0]->globalPosition().theta();
320  if(iSeed){
321  double dTheta = fabs(seedTheta[iSeed] - seedTheta[iSeed-1]);
322  if (dTheta < 0.5){ //? should be something more clever
323  tooCloseClusters = true;
324  break;
325  }
326  }
327  }
328 
329  }
330 
331  // have formed clusters from all hits except for 2D DT segments. Now add the 2D segments to the
332  // compatible clusters. For this we compare the mean cluster postition with the
333  // 2D segment position. We should use the valid coordinate only and use the bad coordinate
334  // as a cross check.
335  for(unsigned int NNN = 0; NNN < seeds.size(); ++NNN) {
336  if(running_meanX[NNN] == 999999.) continue; //skip seeds that have been marked as used up in merging
337 
338  // We have a valid cluster - loop over all 2D segments.
339  if(useDT2D_hasZed) {
340  for (MuonRecHitContainer::const_iterator it2 = muonRecHits_DT2D_hasZed.begin(); it2 != muonRecHits_DT2D_hasZed.end(); ++it2 ) {
341  // check that global theta of 2-D segment lies within cluster box plus or minus allowed slop
342  if (((*it2)->globalPosition().theta() < seed_maxY[NNN] + dYclusBoxMax) && ((*it2)->globalPosition().theta() > seed_minY[NNN] - dYclusBoxMax)) {
343  // check that global phi of 2-D segment (assumed to be center of chamber since no phi hit info)
344  // matches with cluster box plus or minus allowed slop given that the true phi value could be
345  // anywhere within a given chamber (+/- 5 degrees ~ 0.09 radians from center)
346  if(
347  !(
348  (
349  ((*it2)->globalPosition().phi() + 0.09) < (seed_minX[NNN] - dXclusBoxMax)
350  &&
351  ((*it2)->globalPosition().phi() - 0.09) < (seed_minX[NNN] - dXclusBoxMax)
352  )
353  ||
354  (
355  ((*it2)->globalPosition().phi() + 0.09) > (seed_maxX[NNN] + dXclusBoxMax)
356  &&
357  ((*it2)->globalPosition().phi() - 0.09) > (seed_maxX[NNN] + dXclusBoxMax)
358  )
359  )
360  ) { // we have checked that the 2Dsegment is within tight theta boundaries and loose phi boundaries of the current cluster -> add it
361  seeds[NNN].push_back((*it2));
362 
363  }
364  }
365  }
366 
367  }
368 
369  // put DT hasphi loop here
370  if (useDT2D_hasPhi) {
371 
372  for (MuonRecHitContainer::const_iterator it2 = muonRecHits_DT2D_hasPhi.begin(); it2 != muonRecHits_DT2D_hasPhi.end(); ++it2 ) {
373  if (((*it2)->globalPosition().phi() < seed_maxX[NNN] + dXclusBoxMax) && ((*it2)->globalPosition().phi() > seed_minX[NNN] - dXclusBoxMax)) {
374  if(
375  !(
376  (
377  ((*it2)->globalPosition().theta() + 0.3) < (seed_minY[NNN] - dYclusBoxMax)
378  &&
379  ((*it2)->globalPosition().theta() - 0.3) < (seed_minY[NNN] - dYclusBoxMax)
380  )
381  ||
382  (
383  ((*it2)->globalPosition().theta() + 0.3) > (seed_maxY[NNN] + dYclusBoxMax)
384  &&
385  ((*it2)->globalPosition().theta() - 0.3) > (seed_maxY[NNN] + dYclusBoxMax)
386  )
387  )
388  ) { // we have checked that the 2Dsegment is within tight phi boundaries and loose theta boundaries of the current cluster -> add it
389  seeds[NNN].push_back((*it2)); // warning - neeed eta/theta switch here
390 
391  }
392  }
393  }
394  } // DT2D_hastPhi loop
395 
396  // put RPC loop here
397  int secondCh = 0;
398  DetId detId_prev;
399  if(seeds[NNN].size()>1){// actually we should check how many chambers with measurements are present
400  for(unsigned int iRH = 0 ;iRH<seeds[NNN].size() ;++iRH){
401  if( iRH && detId_prev != seeds[NNN][iRH]->hit()->geographicalId()){
402  ++secondCh;
403  break;
404  }
405  detId_prev = seeds[NNN][iRH]->hit()->geographicalId();
406  }
407  }
408 
409  if (useRPCs && !secondCh && !tooCloseClusters) {
410  for (MuonRecHitContainer::const_iterator it2 = muonRecHits_RPC.begin(); it2 != muonRecHits_RPC.end(); ++it2 ) {
411  if (((*it2)->globalPosition().phi() < seed_maxX[NNN] + dXclusBoxMax) && ((*it2)->globalPosition().phi() > seed_minX[NNN] - dXclusBoxMax)) {
412  if(
413  !(
414  (
415  ((*it2)->globalPosition().theta() + 0.3) < (seed_minY[NNN] - dYclusBoxMax)
416  &&
417  ((*it2)->globalPosition().theta() - 0.3) < (seed_minY[NNN] - dYclusBoxMax)
418  )
419  ||
420  (
421  ((*it2)->globalPosition().theta() + 0.3) > (seed_maxY[NNN] + dYclusBoxMax)
422  &&
423  ((*it2)->globalPosition().theta() - 0.3) > (seed_maxY[NNN] + dYclusBoxMax)
424  )
425  )
426  ) { // we have checked that the 2Dsegment is within tight phi boundaries and loose theta boundaries of the current cluster -> add it
427  seeds[NNN].push_back((*it2)); // warning - neeed eta/theta switch here
428 
429  }
430  }
431  }
432  } // RPC loop
433  }
434 
435  // hand over the final seeds to the output
436  // would be more elegant if we could do the above step with
437  // erasing the merged ones, rather than the
438  for(unsigned int NNN = 0; NNN < seeds.size(); ++NNN) {
439  if(running_meanX[NNN] == 999999.) continue; //skip seeds that have been marked as used up in merging
440  //std::cout<<"Next Cluster..."<<std::endl;
441  segments_clusters.push_back(seeds[NNN]);
442  }
443 }
444 
445 
447  const LocalPoint& localPosition, const LocalError& localError,
448  const LocalVector& localDirection, const LocalError& localDirectionError,
449  const double& chi2, const int& ndf){
450  // drop segments which are "bad"
451  bool dropTheSegment = true;
452  const GeomDet* geomDet = theService->trackingGeometry()->idToDet( detId );
453  // only segments whithin the boundaries of the chamber
454  bool insideCh = geomDet->surface().bounds().inside(localPosition, localError,outsideChamberErrorScale);
455 
456  // Don't use segments (nearly) parallel to the chamberi;
457  // the direction vector is normalized (R=1)
458  bool parallelSegment = fabs(localDirection.z())>minLocalSegmentAngle? false: true;
459 
460  if(insideCh && !parallelSegment){
461  dropTheSegment = false;
462  }
463  // use chi2 too? (DT, CSCs, RPCs; 2D, 4D;...)
464 
465 
466  return dropTheSegment;
467 }
const double Pi
T getParameter(std::string const &) const
virtual bool inside(const Local3DPoint &) const =0
Determine if the point is inside the bounds.
const std::string metname
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
T z() const
Definition: PV3DBase.h:63
bool segmentCleaning(const DetId &detId, const LocalPoint &localPosition, const LocalError &localError, const LocalVector &localDirection, const LocalError &localDirectionError, const double &chi2, const int &ndf)
#define end
Definition: vmac.h:38
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
MuonServiceProxy * theService
Definition: DetId.h:20
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
const Bounds & bounds() const
Definition: BoundSurface.h:89
virtual void produce(const edm::Event &event, const edm::EventSetup &eSetup, std::vector< MuonRecHitContainer > &result)
#define begin
Definition: vmac.h:31
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
SETPatternRecognition(const edm::ParameterSet &pset)
bool insert(Storage &, ItemType *, const IdTag &)
static MuonRecHitPointer specificBuild(const GeomDet *geom, const TrackingRecHit *rh)
tuple size
Write out results.