CMS 3D CMS Logo

L1TMuonBarrelKalmanAlgo.cc
Go to the documentation of this file.
1 #include <cmath>
3 
4 
6  verbose_(settings.getParameter<bool>("verbose")),
7  lutService_(new L1TMuonBarrelKalmanLUTs(settings.getParameter<std::string>("lutFile"))),
8  initK_(settings.getParameter<std::vector<double> >("initialK")),
9  initK2_(settings.getParameter<std::vector<double> >("initialK2")),
10  eLoss_(settings.getParameter<std::vector<double> >("eLoss")),
11  aPhi_(settings.getParameter<std::vector<double> >("aPhi")),
12  aPhiB_(settings.getParameter<std::vector<double> >("aPhiB")),
13  aPhiBNLO_(settings.getParameter<std::vector<double> >("aPhiBNLO")),
14  bPhi_(settings.getParameter<std::vector<double> >("bPhi")),
15  bPhiB_(settings.getParameter<std::vector<double> >("bPhiB")),
16  globalChi2Cut_(settings.getParameter<unsigned int>("globalChi2Cut")),
17  chiSquare_(settings.getParameter<std::vector<double> >("chiSquare")),
18  chiSquareCutPattern_(settings.getParameter<std::vector<int> >("chiSquareCutPattern")),
19  chiSquareCutCurv_(settings.getParameter<std::vector<int> >("chiSquareCutCurvMax")),
20  chiSquareCut_(settings.getParameter<std::vector<int> >("chiSquareCut")),
21  combos4_(settings.getParameter<std::vector<int> >("combos4")),
22  combos3_(settings.getParameter<std::vector<int> >("combos3")),
23  combos2_(settings.getParameter<std::vector<int> >("combos2")),
24  combos1_(settings.getParameter<std::vector<int> >("combos1")),
25  useOfflineAlgo_(settings.getParameter<bool>("useOfflineAlgo")),
26  mScatteringPhi_(settings.getParameter<std::vector<double> >("mScatteringPhi")),
27  mScatteringPhiB_(settings.getParameter<std::vector<double> >("mScatteringPhiB")),
28  pointResolutionPhi_(settings.getParameter<double>("pointResolutionPhi")),
29  pointResolutionPhiB_(settings.getParameter<double>("pointResolutionPhiB")),
30  pointResolutionVertex_(settings.getParameter<double>("pointResolutionVertex"))
31 
32 {
33 
34 
35 
36 }
37 
38 
39 
41  for (uint i=0;i<tracks.size();++i) {
42  printf("Code=%d, track=%d\n",tracks[i].hitPattern(),mask);
43  if (tracks[i].hitPattern()==mask)
44  return std::make_pair(true,i);
45  }
46  return std::make_pair(false,0);
47 }
48 void L1TMuonBarrelKalmanAlgo::addBMTFMuon(int bx,const L1MuKBMTrack& track, std::unique_ptr<l1t::RegionalMuonCandBxCollection>& out) {
49 
50  int K = abs(track.curvatureAtVertex());
51  //calibration
52  int sign,signValid;
53 
54  if (track.curvatureAtVertex()==0) {
55  sign=0;
56  signValid=0;
57  }
58  else if (track.curvatureAtVertex()>0) {
59  sign=0;
60  signValid=1;
61  }
62  else {
63  sign=1;
64  signValid=1;
65  }
66 
67  if (K<22)
68  K=22;
69 
70  if (K>2047)
71  K=2047;
72 
73  float lsb=1.25/float(1<<13);
74  int pt = int(2*(1.0/(lsb*float(K))));
75  if (pt>511)
76  pt=511;
77 
78  int K2 = abs(track.curvatureAtMuon());
79  if (K2<22)
80  K2=22;
81 
82  if (K2>2047)
83  K2=2047;
84  int pt2 = int(1.0/(lsb*float(K2)));
85  if (pt2>254)
86  pt2=254;
87  int eta = track.hasFineEta() ? track.fineEta() : track.coarseEta();
88 
89  int phi=track.phiAtMuon();
90  phi=int((phi*M_PI/(6.0*2048.0)+15.0*M_PI/180.)/(0.625*M_PI/180.0));
91 
92  int processor=track.sector();
93  int HF = track.hasFineEta();
94 
95  int quality=rank(track)/64;
96 
97  int dxy=abs(track.dxy())>>7;
98 
99  int trackAddr;
100  std::map<int,int> addr = trackAddress(track,trackAddr);
101 
102  l1t::RegionalMuonCand muon(pt,phi,eta,sign,signValid,quality,processor,l1t::bmtf,addr);
103  muon.setHwHF(HF);
104  muon.setHwPt2(pt2);
105  muon.setHwDXY(dxy);
106 
107 
108  //nw the words!
109  uint32_t word1=pt;
110  word1=word1 | quality<<9;
111  word1=word1 | (twosCompToBits(eta))<<13;
112  word1=word1 | HF<<22;
113  word1=word1 | (twosCompToBits(phi))<<23;
114 
115  uint32_t word2=sign;
116  word2=word2 | signValid<<1;
117  word2=word2 | dxy<<2;
118  word2=word2 | trackAddr<<4;
119  word2=word2 | (twosCompToBits(track.wheel()))<<20;
120  word2=word2 | pt2<<23;
121  muon.setDataword(word2,word1);
122  out->push_back(bx,muon);
123 }
124 
125 
126 
129 
130  bool found=false;
131  uint best=0;
132  int distance=100000;
133  uint N=0;
134  for (const auto& stub :stubs) {
135  N=N+1;
136  if (stub->stNum()!=track.step())
137  continue;
138  if (abs(track.positionAngle()-correctedPhi(stub,track.sector()))<distance) {
139  distance = abs(track.positionAngle()-correctedPhi(stub,track.sector()));
140  best=N-1;
141  found=true;
142  }
143  }
144  return std::make_pair(found,best);
145 }
146 
147 
148 
150  //Promote phiB to 12 bits
151  return 8*stub->phiB();
152 
153 }
154 
156  if (stub->scNum()==sector) {
157  return stub->phi();
158  }
159  else if ((stub->scNum()==sector-1) || (stub->scNum()==11 && sector==0)) {
160  return stub->phi()-2144;
161  }
162  else if ((stub->scNum()==sector+1) || (stub->scNum()==0 && sector==11)) {
163  return stub->phi()+2144;
164  }
165  return stub->phi();
166 }
167 
168 
170  unsigned int mask = 0;
171  for (const auto& stub : track.stubs()) {
172  mask = mask+round(pow(2,stub->stNum()-1));
173  }
174  return mask;
175 }
176 
177 
178 
179 int L1TMuonBarrelKalmanAlgo::customBitmask(unsigned int bit1,unsigned int bit2,unsigned int bit3,unsigned int bit4) {
180  return bit1*1+bit2*2+bit3*4+bit4*8;
181 }
182 
183 bool L1TMuonBarrelKalmanAlgo::getBit(int bitmask,int pos) {
184  return (bitmask & ( 1 << pos )) >> pos;
185 }
186 
187 
188 
190  int K = track.curvature();
191  int phi = track.positionAngle();
192  int phiB = track.bendingAngle();
193  unsigned int step = track.step();
194 
195  int charge=1;
196  if (K!=0)
197  charge = K/abs(K);
198 
199  //energy loss term only for MU->VERTEX
200  //int offset=int(charge*eLoss_[step-1]*K*K);
201  // if (abs(offset)>4096)
202  // offset=4096*offset/abs(offset);
203  int KNew =wrapAround(int(K/(1+charge*eLoss_[step-1]*K)),8192);
204 
205  //phi propagation
206  int phiNew =wrapAround(phi+int(aPhi_[step-1]*K)-int(bPhi_[step-1]*phiB),8192);
207 
208  //phiB propagation
209  int phiBNew = wrapAround(int(aPhiB_[step-1]*K) +int(bPhiB_[step-1]*phiB),2048);
210 
211  //Only for the propagation to vertex we use the LUT for better precision and the full function
212  if (step==1)
213  phiBNew = wrapAround(int(aPhiB_[step-1]*K/(1+charge*aPhiBNLO_[step-1]*K))+int(bPhiB_[step-1]*phiB),2048);
214 
216  //Rest of the stuff is for the offline version only
217  //where we want to check what is happening in the covariaznce matrix
218 
219  //Create the transformation matrix
220  double a[9];
221  a[0] = 1.;
222  a[1] = 0.0;
223  a[2] = 0.0;
224  a[3] = aPhi_[step-1];
225  // a[3] = 0.0;
226  a[4] = 1.0;
227  a[5] = -bPhi_[step-1];
228  //a[6]=0.0;
229  a[6] = aPhiB_[step-1];
230  a[7] = 0.0;
231  a[8] = bPhiB_[step-1];
232 
233 
234  ROOT::Math::SMatrix<double,3> P(a,9);
235 
236  const std::vector<double>& covLine = track.covariance();
237  L1MuKBMTrack::CovarianceMatrix cov(covLine.begin(),covLine.end());
238  cov = ROOT::Math::Similarity(P,cov);
239 
240 
241  //Add the multiple scattering
242  double phiRMS = mScatteringPhi_[step-1]*K*K;
243  double phiBRMS = mScatteringPhiB_[step-1]*K*K;
244 
245  std::vector<double> b(6);
246  b[0] = 0;
247  b[1] = 0;
248  b[2] =phiRMS;
249  b[3] =0;
250  b[4] = 0;
251  b[5] = phiBRMS;
252 
253  reco::Candidate::CovarianceMatrix MS(b.begin(),b.end());
254 
255  cov = cov+MS;
256 
257  if (verbose_) {
258  printf("Covariance term for phiB = %f\n",cov(2,2));
259  printf("Multiple scattering term for phiB = %f\n",MS(2,2));
260  }
261 
262 
263 
264  track.setCovariance(cov);
265  track.setCoordinates(step-1,KNew,phiNew,phiBNew);
266 
267 }
268 
269 
271  updateEta(track,stub);
272  if (useOfflineAlgo_) {
273  if (mask==3 || mask ==5 || mask==9 ||mask==6|| mask==10 ||mask==12)
274  return updateOffline(track,stub);
275  else
276  return updateOffline1D(track,stub);
277 
278  }
279  else
280  return updateLUT(track,stub,mask);
281 
282 }
283 
285  int trackK = track.curvature();
286  int trackPhi = track.positionAngle();
287  int trackPhiB = track.bendingAngle();
288 
289  int phi = correctedPhi(stub,track.sector());
290  int phiB = correctedPhiB(stub);
291 
292 
293 
294  //Update eta
295  track.setCoarseEta(int((track.coarseEta()+stub->coarseEta())/2.0));
296 
297  Vector2 residual;
298  residual[0] = phi-trackPhi;
299  residual[1] = phiB-trackPhiB;
300 
301 
302 
303 
304  // if (stub->quality()<4)
305  // phiB=trackPhiB;
306 
307  Matrix23 H;
308  H(0,0)=0.0;
309  H(0,1)=1.0;
310  H(0,2)=0.0;
311  H(1,0)=0.0;
312  H(1,1)=0.0;
313  H(1,2)=1.0;
314 
315 
317  R(0,0) = pointResolutionPhi_;
318  R(0,1) = 0.0;
319  R(1,0) = 0.0;
320  R(1,1) = pointResolutionPhiB_;
321 
322  const std::vector<double>& covLine = track.covariance();
323  L1MuKBMTrack::CovarianceMatrix cov(covLine.begin(),covLine.end());
324 
325 
326  CovarianceMatrix2 S = ROOT::Math::Similarity(H,cov)+R;
327  if (!S.Invert())
328  return false;
329  Matrix32 Gain = cov*ROOT::Math::Transpose(H)*S;
330 
331  track.setKalmanGain(track.step(),abs(trackK),Gain(0,0),Gain(0,1),Gain(1,0),Gain(1,1),Gain(2,0),Gain(2,1));
332 
333  int KNew = (trackK+int(Gain(0,0)*residual(0)+Gain(0,1)*residual(1)));
334  if (abs(KNew)>8192)
335  return false;
336 
337  int phiNew = wrapAround(trackPhi+residual(0),8192);
338  int phiBNew = wrapAround(trackPhiB+int(Gain(2,0)*residual(0)+Gain(2,1)*residual(1)),2048);
339 
340  track.setResidual(stub->stNum()-1,abs(phi-phiNew)+abs(phiB-phiBNew)/8);
341 
342 
343  if (verbose_) {
344  printf(" K = %d + %f * %f + %f * %f\n",trackK,Gain(0,0),residual(0),Gain(0,1),residual(1));
345  printf(" phiB = %d + %f * %f + %f * %f\n",trackPhiB,Gain(2,0),residual(0),Gain(2,1),residual(1));
346  }
347 
348 
349  track.setCoordinates(track.step(),KNew,phiNew,phiBNew);
350  Matrix33 covNew = cov - Gain*(H*cov);
352 
353  c(0,0)=covNew(0,0);
354  c(0,1)=covNew(0,1);
355  c(0,2)=covNew(0,2);
356  c(1,0)=covNew(1,0);
357  c(1,1)=covNew(1,1);
358  c(1,2)=covNew(1,2);
359  c(2,0)=covNew(2,0);
360  c(2,1)=covNew(2,1);
361  c(2,2)=covNew(2,2);
362  if (verbose_) {
363  printf("Post Fit Covariance Matrix %f %f %f \n",cov(0,0),cov(1,1),cov(2,2));
364 
365  }
366 
367  track.setCovariance(c);
368  track.addStub(stub);
369  track.setHitPattern(hitPattern(track));
370 
371  return true;
372 }
373 
374 
376  int trackK = track.curvature();
377  int trackPhi = track.positionAngle();
378  int trackPhiB = track.bendingAngle();
379 
380 
381  int phi = correctedPhi(stub,track.sector());
382  track.setCoarseEta(int((track.coarseEta()+stub->coarseEta())/2.0));
383 
384  double residual= phi-trackPhi;
385 
386  Matrix13 H;
387  H(0,0)=0.0;
388  H(0,1)=1.0;
389  H(0,2)=0.0;
390 
391 
392  const std::vector<double>& covLine = track.covariance();
393  L1MuKBMTrack::CovarianceMatrix cov(covLine.begin(),covLine.end());
394 
395  double S = ROOT::Math::Similarity(H,cov)(0,0)+pointResolutionPhi_;
396 
397  if (S==0.0)
398  return false;
399  Matrix31 Gain = cov*ROOT::Math::Transpose(H)/S;
400 
401  track.setKalmanGain(track.step(),abs(trackK),Gain(0,0),0.0,Gain(1,0),0.0,Gain(2,0),0.0);
402 
403  int KNew = wrapAround(trackK+int(Gain(0,0)*residual),8192);
404  int phiNew = wrapAround(trackPhi+residual,8192);
405  int phiBNew = wrapAround(trackPhiB+int(Gain(2,0)*residual),2048);
406  track.setCoordinates(track.step(),KNew,phiNew,phiBNew);
407  Matrix33 covNew = cov - Gain*(H*cov);
409 
410  c(0,0)=covNew(0,0);
411  c(0,1)=covNew(0,1);
412  c(0,2)=covNew(0,2);
413  c(1,0)=covNew(1,0);
414  c(1,1)=covNew(1,1);
415  c(1,2)=covNew(1,2);
416  c(2,0)=covNew(2,0);
417  c(2,1)=covNew(2,1);
418  c(2,2)=covNew(2,2);
419  track.setCovariance(c);
420  track.addStub(stub);
421  track.setHitPattern(hitPattern(track));
422 
423  return true;
424 }
425 
426 
427 
429  int trackK = track.curvature();
430  int trackPhi = track.positionAngle();
431  int trackPhiB = track.bendingAngle();
432 
433  //Update eta
434  track.setCoarseEta(int((track.coarseEta()+stub->coarseEta())/2.0));
435 
436  int phi = correctedPhi(stub,track.sector());
437  int phiB = correctedPhiB(stub);
438  // if (stub->quality()<6)
439  // phiB=trackPhiB;
440 
441  Vector2 residual;
442  residual[0] = phi-trackPhi;
443  residual[1] = phiB-trackPhiB;
444 
445  uint absK = abs(trackK);
446  if (absK>4095)
447  absK = 4095;
448  std::vector<float> GAIN = lutService_->trackGain(track.step(),mask,absK/4);
449  track.setKalmanGain(track.step(),abs(trackK),GAIN[0],GAIN[1],1,0,GAIN[2],GAIN[3]);
450 
451  //For the three stub stuff use only gains 0 and 4
452  if (!(mask==3 || mask ==5 || mask==9 ||mask==6|| mask==10 ||mask==12)) {
453  GAIN[1]=0.0;
454  GAIN[3]=0.0;
455  }
456 
457 
458 
459 
460  int KNew = wrapAround(trackK+int(GAIN[0]*residual(0)+GAIN[1]*residual(1)),8192);
461  int phiNew = wrapAround(trackPhi+residual(0),8192);
462  int phiBNew = wrapAround(trackPhiB+int(GAIN[2]*residual(0)+GAIN[3]*residual(1)),2048);
463  track.setCoordinates(track.step(),KNew,phiNew,phiBNew);
464  track.addStub(stub);
465  track.setHitPattern(hitPattern(track));
466  return true;
467 }
468 
469 
470 
471 
473 
474  if (stub->qeta1()>=0) {
475  if (track.hasFineEta()) {
476  uint dist2=1000;
477  uint dist1 = abs(track.fineEta()-stub->eta1());
478  if (stub->qeta2()>=0)
479  dist2 = abs(track.fineEta()-stub->eta2());
480  if (dist1<dist2)
481  track.setFineEta((stub->eta1()+track.fineEta())/2);
482  else
483  track.setFineEta((stub->eta2()+track.fineEta())/2);
484  }else {
485  if (stub->qeta2()>=0)
486  track.setFineEta((stub->eta1()+stub->eta2())/2);
487  else
488  track.setFineEta(stub->eta1());
489  }
490 
491  }
492 }
493 
494 
495 
496 
497 
498 
500  if (useOfflineAlgo_)
502  else
503  vertexConstraintLUT(track);
504 
505 }
506 
507 
509  double residual = -track.dxy();
510  Matrix13 H;
511  H(0,0)=0;
512  H(0,1)=0;
513  H(0,2)=1;
514 
515  const std::vector<double>& covLine = track.covariance();
516  L1MuKBMTrack::CovarianceMatrix cov(covLine.begin(),covLine.end());
517 
518  double S = (ROOT::Math::Similarity(H,cov))(0,0)+pointResolutionVertex_;
519  S=1.0/S;
520  Matrix31 Gain = cov*(ROOT::Math::Transpose(H))*S;
521  track.setKalmanGain(track.step(),abs(track.curvature()),Gain(0,0),Gain(1,0),Gain(2,0));
522 
523  if (verbose_) {
524  printf("sigma3=%f sigma6=%f\n",cov(0,3),cov(3,3));
525  printf(" K = %d + %f * %f\n",track.curvature(),Gain(0,0),residual);
526  }
527 
528  int KNew = wrapAround(int(track.curvature()+Gain(0,0)*residual),8192);
529  int phiNew = wrapAround(int(track.positionAngle()+Gain(1,0)*residual),8192);
530  int dxyNew = wrapAround(int(track.dxy()+Gain(2,0)*residual),8192);
531  if (verbose_)
532  printf("Post fit impact parameter=%d\n",dxyNew);
533  track.setCoordinatesAtVertex(KNew,phiNew,-residual);
534  Matrix33 covNew = cov - Gain*(H*cov);
536  c(0,0)=covNew(0,0);
537  c(0,1)=covNew(0,1);
538  c(0,2)=covNew(0,2);
539  c(1,0)=covNew(1,0);
540  c(1,1)=covNew(1,1);
541  c(1,2)=covNew(1,2);
542  c(2,0)=covNew(2,0);
543  c(2,1)=covNew(2,1);
544  c(2,2)=covNew(2,2);
545  track.setCovariance(c);
546  // track.covariance = track.covariance - Gain*H*track.covariance;
547 }
548 
549 
550 
552  double residual = -track.dxy();
553  uint absK = abs(track.curvature());
554  if (absK>2047)
555  absK = 2047;
556 
557 std::pair<float,float> GAIN = lutService_->vertexGain(track.hitPattern(),absK/2);
558  track.setKalmanGain(track.step(),abs(track.curvature()),GAIN.first,GAIN.second,-1);
559 
560  int KNew = wrapAround(int(track.curvature()+GAIN.first*residual),8192);
561  int phiNew = wrapAround(int(track.positionAngle()+GAIN.second*residual),8192);
562  track.setCoordinatesAtVertex(KNew,phiNew,-residual);
563 }
564 
565 
566 
568  int K,phiINT,etaINT;
569 
570  if (track.hasFineEta())
571  etaINT=track.fineEta();
572  else
573  etaINT=track.coarseEta();
574 
575 
576  double lsb = 1.25/float(1 << 13);
577  double lsbEta = 0.010875;
578 
579 
580  if (vertex) {
581  K = track.curvatureAtVertex();
582  if (K==0)
583  track.setCharge(1);
584  else
585  track.setCharge(K/abs(K));
586 
587  phiINT = track.phiAtVertex();
588  double phi= track.sector()*M_PI/6.0+phiINT*M_PI/(6*2048.)-2*M_PI;
589  double eta = etaINT*lsbEta;
590  if (phi<-M_PI)
591  phi=phi+2*M_PI;
592  if (K==0)
593  K=1;
594 
595 
596  float FK=fabs(K);
597  FK = fabs(0.912*FK+(2.557e-5)*FK*FK-6);
598 
599  if (FK<51)
600  FK=51;
601 
602  double pt = 1.0/(lsb*(FK));
603 
604  track.setPtEtaPhi(pt,eta,phi);
605  }
606  else {
607  K=track.curvatureAtMuon();
608  if (K==0)
609  K=1;
610  if (abs(K)<46)
611  K=46*K/abs(K);
612  double pt = 1.0/(lsb*abs(K));
613  track.setPtUnconstrained(pt);
614  }
615 }
616 
617 
618 
620  L1MuKBMTrackCollection pretracks;
621  std::vector<int> combinatorics;
622  switch(seed->stNum()) {
623  case 1:
624  combinatorics=combos1_;
625  break;
626  case 2:
627  combinatorics=combos2_;
628  break;
629 
630  case 3:
631  combinatorics=combos3_;
632  break;
633 
634  case 4:
635  combinatorics=combos4_;
636  break;
637 
638  default:
639  printf("Something really bad happend\n");
640  }
641 
642  L1MuKBMTrack nullTrack(seed,correctedPhi(seed,seed->scNum()),correctedPhiB(seed));
643 
644  for( const auto& mask : combinatorics) {
645  L1MuKBMTrack track(seed,correctedPhi(seed,seed->scNum()),correctedPhiB(seed));
646  int phiB = correctedPhiB(seed);
647  if (seed->quality()<4)
648  phiB=0;
649  int charge;
650  if (phiB==0)
651  charge = 0;
652  else
653  charge=phiB/abs(phiB);
654  int initialK = int(initK_[seed->stNum()-1]*phiB/(1+initK2_[seed->stNum()-1]*charge*phiB));
655  if (initialK>8192)
656  initialK=8192;
657  if (initialK<-8192)
658  initialK=-8192;
659 
660 
661 
662  track.setCoordinates(seed->stNum(),initialK,correctedPhi(seed,seed->scNum()),phiB);
663  track.setHitPattern(hitPattern(track));
664  //Set eta coarse
665  track.setCoarseEta(seed->coarseEta());
666  //Set fine eta:
667  if (seed->qeta1()>=0) {
668  if (seed->qeta2()>=0) {
669  track.setFineEta((seed->eta1()+seed->eta2())/2);
670  }
671  else {
672  track.setFineEta(seed->eta1());
673  }
674  }
675 
676  //set covariance
677  L1MuKBMTrack::CovarianceMatrix covariance;
678 
679 
680  float DK=512*512.;
681  covariance(0,0)=DK;
682  covariance(0,1)=0;
683  covariance(0,2)=0;
684  covariance(1,0)=0;
685  covariance(1,1)=float(pointResolutionPhi_);
686  covariance(1,2)=0;
687  covariance(2,0)=0;
688  covariance(2,1)=0;
689  covariance(2,2)=float(pointResolutionPhiB_);
690  track.setCovariance(covariance);
691  //
692  if (verbose_) {
693  printf("New Kalman fit staring at step=%d, phi=%d,phiB=%d with curvature=%d\n",track.step(),track.positionAngle(),track.bendingAngle(),track.curvature());
694  printf("BITMASK:");
695  for (unsigned int i=0;i<4;++i)
696  printf("%d",getBit(mask,i));
697  printf("\n");
698  printf("------------------------------------------------------\n");
699  printf("------------------------------------------------------\n");
700  printf("------------------------------------------------------\n");
701  printf("stubs:\n");
702  for (const auto& stub: stubs)
703  printf("station=%d phi=%d phiB=%d qual=%d \n",stub->stNum(),correctedPhi(stub,seed->scNum()),correctedPhiB(stub),stub->quality());
704  printf("------------------------------------------------------\n");
705  printf("------------------------------------------------------\n");
706 
707  }
708 
709  while(track.step()>0) {
710  // muon station 1
711  if (track.step()==1) {
712  track.setCoordinatesAtMuon(track.curvature(),track.positionAngle(),track.bendingAngle());
714  if (verbose_)
715  printf ("Unconstrained PT in Muon System: pt=%f\n",track.ptUnconstrained());
716  }
717 
718  propagate(track);
719  if (verbose_)
720  printf("propagated Coordinates step:%d,phi=%d,phiB=%d,K=%d\n",track.step(),track.positionAngle(),track.bendingAngle(),track.curvature());
721 
722  if (track.step()>0)
723  if (getBit(mask,track.step()-1)) {
724  std::pair<bool,uint> bestStub = match(track,stubs);
725  if ((!bestStub.first) || (!update(track,stubs[bestStub.second],mask)))
726  break;
727  if (verbose_) {
728  printf("updated Coordinates step:%d,phi=%d,phiB=%d,K=%d\n",track.step(),track.positionAngle(),track.bendingAngle(),track.curvature());
729  }
730  }
731  if (track.step()==0) {
732  track.setCoordinatesAtVertex(track.curvature(),track.positionAngle(),track.bendingAngle());
733  if (verbose_)
734  printf(" Coordinates before vertex constraint step:%d,phi=%d,dxy=%d,K=%d\n",track.step(),track.phiAtVertex(),track.dxy(),track.curvatureAtVertex());
736  if (uint(track.approxChi2())>globalChi2Cut_)
737  break;
739  if (verbose_) {
740  printf(" Coordinates after vertex constraint step:%d,phi=%d,dxy=%d,K=%d maximum local chi2=%d\n",track.step(),track.phiAtVertex(),track.dxy(),track.curvatureAtVertex(),track.approxChi2());
741  printf("------------------------------------------------------\n");
742  printf("------------------------------------------------------\n");
743  }
745  track.setRank(rank(track));
746  if (verbose_)
747  printf ("Floating point coordinates at vertex: pt=%f, eta=%f phi=%f\n",track.pt(),track.eta(),track.phi());
748  pretracks.push_back(track);
749  }
750  }
751  }
752 
753  //Resolve eta
754  resolveEtaUnit(pretracks);
755  //Now for all the pretracks we need only one
756  L1MuKBMTrackCollection cleaned = cleanAndSort(pretracks,1);
757 
758  if (!cleaned.empty())
759  return std::make_pair(true,cleaned[0]);
760  return std::make_pair(false,nullTrack);
761 }
762 
763 
764 
765 
766 
767 
769  //here we have a simplification of the algorithm for the sake of the emulator - rsult is identical
770  // we apply cuts on the firmware as |u -u'|^2 < a+b *K^2
771  int K = track.curvatureAtMuon();
772 
773  int chi=0;
774 
775  for (const auto& stub: track.stubs()) {
776  uint delta=abs(correctedPhi(stub,track.sector())-track.phiAtMuon()+correctedPhiB(stub)-track.phiBAtMuon()-chiSquare_[stub->stNum()-1]*K);
777  chi=chi+delta;
778  }
779  chi=chi/2;
780  if (chi>511)
781  chi=511;
782  track.setApproxChi2(chi);
783 }
784 
785 
787  // int offset=0;
788  if (hitPattern(track)==customBitmask(0,0,1,1))
789  return 65;
790  // return offset+(track.stubs().size()*2+track.quality())*80-track.approxChi2();
791  return 700+(track.stubs().size())*80-track.approxChi2();
792 
793 }
794 
795 
796 
798  if (value>maximum-1)
799  return value-2*maximum;
800  if (value<-maximum)
801  return value+2*maximum;
802  return value;
803 
804 }
805 
807  int bestFineEta=0;
808  uint bestSegments=0;
809  for (const auto& track : tracks) {
810  if (track.stubs().size()>bestSegments && track.hasFineEta()) {
811  bestFineEta = track.fineEta();
812  bestSegments=track.stubs().size();
813  }
814  }
815 
816  for (auto & track :tracks) {
817  if (bestSegments!=0 && (!track.hasFineEta())) {
818  track.setFineEta(bestFineEta);
820  }
821  }
822 
823 
824 
825 }
826 
827 
828 
831 
832  if (verbose_)
833  printf(" -----Preselected Kalman Tracks-----\n");
834 
835 
836  L1MuKBMTrackCollection pretracks;
837 
838  for(const auto& track1 : tracks) {
839  if (verbose_)
840  printf("Preselected Kalman Track charge=%d pt=%f eta=%f phi=%f curvature=%d curvature STA =%d stubs=%d chi2=%d pts=%f %f\n",track1.charge(),track1.pt(),track1.eta(),track1.phi(),track1.curvatureAtVertex(),track1.curvatureAtMuon(),int(track1.stubs().size()),track1.approxChi2(),track1.pt(),track1.ptUnconstrained());
841 
842  //Clean up
843  bool veto=false;
844  for (uint i=0;i<chiSquareCutPattern_.size();++i) {
845  if (track1.hitPattern()==chiSquareCutPattern_[i] && abs(track1.curvatureAtVertex())<chiSquareCutCurv_[i] && track1.approxChi2()>chiSquareCut_[i] && track1.curvature()*track1.dxy()<0) {
846  veto=true;
847  break;
848  }
849  }
850  if (!veto)
851  pretracks.push_back(track1);
852  }
853 
854  for(const auto& track1 : pretracks) {
855  bool keep=true;
856  for(const auto& track2 : pretracks) {
857  if (track1==track2)
858  continue;
859  if (!track1.overlapTrack(track2))
860  continue;
861  if (track1.rank()<track2.rank())
862  keep=false;
863  }
864  if (keep)
865  out.push_back(track1);
866  }
867 
868  if (verbose_) {
869  printf(" -----Algo Result Kalman Tracks-----\n");
870  for (const auto& track1 :out)
871  printf("Final Kalman Track charge=%d pt=%f eta=%f phi=%f curvature=%d curvature STA =%d stubs=%d chi2=%d pts=%f %f\n",track1.charge(),track1.pt(),track1.eta(),track1.phi(),track1.curvatureAtVertex(),track1.curvatureAtMuon(),int(track1.stubs().size()),track1.approxChi2(),track1.pt(),track1.ptUnconstrained());
872  }
873 
874 
875 
877  if (!out.empty())
878  std::sort(out.begin(),out.end(),sorter);
879 
880 
881  L1MuKBMTrackCollection exported;
882  for (uint i=0;i<out.size();++i)
883  if (i<=keep)
884  exported.push_back(out[i]);
885  return exported;
886 }
887 
888 
889 int L1TMuonBarrelKalmanAlgo::encode(bool ownwheel,int sector,bool tag) {
890  if (ownwheel) {
891  if (sector==0) {
892  if (tag)
893  return 8;
894  else
895  return 9;
896  }
897  else if (sector==1) {
898  if (tag)
899  return 10;
900  else
901  return 11;
902 
903  }
904  else {
905  if (tag)
906  return 12;
907  else
908  return 13;
909  }
910 
911  }
912  else {
913  if (sector==0) {
914  if (tag)
915  return 0;
916  else
917  return 1;
918  }
919  else if (sector==1) {
920  if (tag)
921  return 2;
922  else
923  return 3;
924 
925  }
926  else {
927  if (tag)
928  return 4;
929  else
930  return 5;
931  }
932  }
933  return 15;
934 }
935 
936 
937 
938 
939 std::map<int,int> L1TMuonBarrelKalmanAlgo::trackAddress(const L1MuKBMTrack& track,int& word) {
940  std::map<int,int> out;
941  if (track.wheel()>=0)
943  else
945 
956 
957 
958  for (const auto stub: track.stubs()) {
959  bool ownwheel = stub->whNum() == track.wheel();
960  int sector=0;
961  if ((stub->scNum()==track.sector()+1) || (stub->scNum()==0 && track.sector()==11))
962  sector=+1;
963  if ((stub->scNum()==track.sector()-1) || (stub->scNum()==11 && track.sector()==0))
964  sector=-1;
965  int addr = encode(ownwheel,sector,stub->tag());
966 
967  if (stub->stNum()==4) {
968  addr=addr & 3;
970  }
971  if (stub->stNum()==3) {
973  }
974  if (stub->stNum()==2) {
976  }
977  if (stub->stNum()==1) {
979  }
980  }
981 
982  word=0;
983  word = word | out[l1t::RegionalMuonCand::kStat1]<<12;
984  word = word | out[l1t::RegionalMuonCand::kStat2]<<8;
985  word = word | out[l1t::RegionalMuonCand::kStat3]<<4;
986  word = word | out[l1t::RegionalMuonCand::kStat4];
987 
988 
989 
990  return out;
991 }
992 
993 
994 
996  if (q>=0)
997  return q;
998  else
999  return (~q)+1;
1000 
1001 
1002 }
std::vector< double > bPhiB_
int correctedPhiB(const L1MuKBMTCombinedStubRef &)
dbl * delta
Definition: mlp_gen.cc:36
int fineEta() const
void setHwPt2(int bits)
Set compressed second displaced pT as transmitted by hardware LSB = 1.0 (8 bits)
std::vector< double > aPhi_
ROOT::Math::SMatrix< double, 3, 3 > Matrix33
void setCoordinatesAtVertex(int, int, int)
std::vector< double > aPhiB_
std::pair< bool, uint > match(const L1MuKBMTrack &, const L1MuKBMTCombinedStubRefVector &)
std::vector< double > aPhiBNLO_
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepSym< double, 2 > > CovarianceMatrix2
L1TMuonBarrelKalmanAlgo(const edm::ParameterSet &settings)
int step() const
Definition: L1MuKBMTrack.cc:90
std::vector< double > mScatteringPhiB_
int bendingAngle() const
Definition: L1MuKBMTrack.cc:74
int coarseEta() const
Definition: L1MuKBMTrack.cc:78
bool update(L1MuKBMTrack &, const L1MuKBMTCombinedStubRef &, int)
void estimateChiSquare(L1MuKBMTrack &)
int phiAtVertex() const
Definition: L1MuKBMTrack.cc:58
int rank(const L1MuKBMTrack &)
void setPtUnconstrained(float)
void setHitPattern(int)
bool updateOffline1D(L1MuKBMTrack &, const L1MuKBMTCombinedStubRef &)
std::vector< edm::Ref< L1MuKBMTCombinedStubCollection > > L1MuKBMTCombinedStubRefVector
void setPtEtaPhi(double, double, double)
void vertexConstraintLUT(L1MuKBMTrack &)
void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
int phiBAtMuon() const
Definition: L1MuKBMTrack.cc:50
std::pair< bool, uint > getByCode(const L1MuKBMTrackCollection &tracks, int mask)
std::vector< double > chiSquare_
const int keep
bool updateLUT(L1MuKBMTrack &, const L1MuKBMTCombinedStubRef &, int)
int curvatureAtMuon() const
Definition: L1MuKBMTrack.cc:44
bool hasFineEta() const
std::vector< double > initK2_
void setCoordinates(int, int, int, int)
ROOT::Math::SVector< double, 2 > Vector2
int approxChi2() const
Definition: L1MuKBMTrack.cc:82
void setCovariance(const CovarianceMatrix &)
std::vector< double > mScatteringPhi_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setCoarseEta(int)
void setDataword(uint32_t msbs, uint32_t lsbs)
void setHwHF(bool bit)
Set HF (halo / fine eta) bit (EMTF: halo -> 1; BMTF: fine eta -> 1)
void vertexConstraint(L1MuKBMTrack &)
int customBitmask(unsigned int, unsigned int, unsigned int, unsigned int)
Definition: value.py:1
std::vector< double > eLoss_
void vertexConstraintOffline(L1MuKBMTrack &)
void setFineEta(int)
const L1MuKBMTCombinedStubRefVector & stubs() const
int curvature() const
Definition: L1MuKBMTrack.cc:66
int phiAtMuon() const
Definition: L1MuKBMTrack.cc:47
std::vector< double > bPhi_
std::vector< double > initK_
void updateEta(L1MuKBMTrack &, const L1MuKBMTCombinedStubRef &)
#define M_PI
int wheel() const
Definition: L1MuKBMTrack.cc:96
ROOT::Math::SMatrix< double, 2, 3 > Matrix23
void setFloatingPointValues(L1MuKBMTrack &, bool)
void setHwDXY(int bits)
Set compressed impact parameter with respect to beamspot (4 bits)
ROOT::Math::SMatrix< double, 1, 3 > Matrix13
#define N
Definition: blowfish.cc:9
int encode(bool ownwheel, int sector, bool tag)
std::vector< int > chiSquareCutPattern_
const std::vector< double > & covariance() const
std::unique_ptr< L1TMuonBarrelKalmanLUTs > lutService_
int dxy() const
Definition: L1MuKBMTrack.cc:62
std::pair< bool, L1MuKBMTrack > chain(const L1MuKBMTCombinedStubRef &, const L1MuKBMTCombinedStubRefVector &)
def uint(string)
int curvatureAtVertex() const
Definition: L1MuKBMTrack.cc:54
double b
Definition: hdecay.h:120
std::pair< OmniClusterRef, TrackingParticleRef > P
std::vector< int > chiSquareCut_
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
void addBMTFMuon(int, const L1MuKBMTrack &, std::unique_ptr< l1t::RegionalMuonCandBxCollection > &)
L1MuKBMTrackCollection cleanAndSort(const L1MuKBMTrackCollection &, uint)
int hitPattern(const L1MuKBMTrack &)
void setKalmanGain(unsigned int step, unsigned int K, float a1, float a2, float a3, float a4=0, float a5=0, float a6=0)
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:14
double a
Definition: hdecay.h:121
int correctedPhi(const L1MuKBMTCombinedStubRef &, int)
int sector() const
Definition: L1MuKBMTrack.cc:93
std::vector< int > chiSquareCutCurv_
void addStub(const L1MuKBMTCombinedStubRef &)
step
int hitPattern() const
Definition: L1MuKBMTrack.cc:86
void setResidual(uint, int)
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Candidate.h:47
ROOT::Math::SMatrix< double, 3, 1 > Matrix31
int positionAngle() const
Definition: L1MuKBMTrack.cc:70
bool updateOffline(L1MuKBMTrack &, const L1MuKBMTCombinedStubRef &)
ROOT::Math::SMatrix< double, 3, 2 > Matrix32
void resolveEtaUnit(L1MuKBMTrackCollection &)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void setApproxChi2(int)
std::map< int, int > trackAddress(const L1MuKBMTrack &, int &)