forked from cms-sw/cmssw
-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathPixelTriplet.h
2770 lines (2496 loc) · 146 KB
/
PixelTriplet.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#ifndef PixelTriplet_cuh
#define PixelTriplet_cuh
#include "RecoTracker/LSTCore/interface/alpaka/Constants.h"
#include "RecoTracker/LSTCore/interface/alpaka/Module.h"
#include "Triplet.h"
#include "Segment.h"
#include "MiniDoublet.h"
#include "Hit.h"
#include "Quintuplet.h"
namespace SDL {
// One pixel segment, one outer tracker triplet!
struct pixelTriplets {
unsigned int* pixelSegmentIndices;
unsigned int* tripletIndices;
unsigned int* nPixelTriplets;
unsigned int* totOccupancyPixelTriplets;
float* pixelRadiusError;
float* rPhiChiSquared;
float* rPhiChiSquaredInwards;
float* rzChiSquared;
FPX* pixelRadius;
FPX* tripletRadius;
FPX* pt;
FPX* eta;
FPX* phi;
FPX* eta_pix;
FPX* phi_pix;
FPX* score;
bool* isDup;
bool* partOfPT5;
uint8_t* logicalLayers;
unsigned int* hitIndices;
uint16_t* lowerModuleIndices;
FPX* centerX;
FPX* centerY;
template <typename TBuff>
void setData(TBuff& pixelTripletsBuffer) {
pixelSegmentIndices = alpaka::getPtrNative(pixelTripletsBuffer.pixelSegmentIndices_buf);
tripletIndices = alpaka::getPtrNative(pixelTripletsBuffer.tripletIndices_buf);
nPixelTriplets = alpaka::getPtrNative(pixelTripletsBuffer.nPixelTriplets_buf);
totOccupancyPixelTriplets = alpaka::getPtrNative(pixelTripletsBuffer.totOccupancyPixelTriplets_buf);
pixelRadius = alpaka::getPtrNative(pixelTripletsBuffer.pixelRadius_buf);
tripletRadius = alpaka::getPtrNative(pixelTripletsBuffer.tripletRadius_buf);
pt = alpaka::getPtrNative(pixelTripletsBuffer.pt_buf);
eta = alpaka::getPtrNative(pixelTripletsBuffer.eta_buf);
phi = alpaka::getPtrNative(pixelTripletsBuffer.phi_buf);
eta_pix = alpaka::getPtrNative(pixelTripletsBuffer.eta_pix_buf);
phi_pix = alpaka::getPtrNative(pixelTripletsBuffer.phi_pix_buf);
score = alpaka::getPtrNative(pixelTripletsBuffer.score_buf);
isDup = alpaka::getPtrNative(pixelTripletsBuffer.isDup_buf);
partOfPT5 = alpaka::getPtrNative(pixelTripletsBuffer.partOfPT5_buf);
logicalLayers = alpaka::getPtrNative(pixelTripletsBuffer.logicalLayers_buf);
hitIndices = alpaka::getPtrNative(pixelTripletsBuffer.hitIndices_buf);
lowerModuleIndices = alpaka::getPtrNative(pixelTripletsBuffer.lowerModuleIndices_buf);
centerX = alpaka::getPtrNative(pixelTripletsBuffer.centerX_buf);
centerY = alpaka::getPtrNative(pixelTripletsBuffer.centerY_buf);
pixelRadiusError = alpaka::getPtrNative(pixelTripletsBuffer.pixelRadiusError_buf);
rPhiChiSquared = alpaka::getPtrNative(pixelTripletsBuffer.rPhiChiSquared_buf);
rPhiChiSquaredInwards = alpaka::getPtrNative(pixelTripletsBuffer.rPhiChiSquaredInwards_buf);
rzChiSquared = alpaka::getPtrNative(pixelTripletsBuffer.rzChiSquared_buf);
}
};
template <typename TDev>
struct pixelTripletsBuffer : pixelTriplets {
Buf<TDev, unsigned int> pixelSegmentIndices_buf;
Buf<TDev, unsigned int> tripletIndices_buf;
Buf<TDev, unsigned int> nPixelTriplets_buf;
Buf<TDev, unsigned int> totOccupancyPixelTriplets_buf;
Buf<TDev, FPX> pixelRadius_buf;
Buf<TDev, FPX> tripletRadius_buf;
Buf<TDev, FPX> pt_buf;
Buf<TDev, FPX> eta_buf;
Buf<TDev, FPX> phi_buf;
Buf<TDev, FPX> eta_pix_buf;
Buf<TDev, FPX> phi_pix_buf;
Buf<TDev, FPX> score_buf;
Buf<TDev, bool> isDup_buf;
Buf<TDev, bool> partOfPT5_buf;
Buf<TDev, uint8_t> logicalLayers_buf;
Buf<TDev, unsigned int> hitIndices_buf;
Buf<TDev, uint16_t> lowerModuleIndices_buf;
Buf<TDev, FPX> centerX_buf;
Buf<TDev, FPX> centerY_buf;
Buf<TDev, float> pixelRadiusError_buf;
Buf<TDev, float> rPhiChiSquared_buf;
Buf<TDev, float> rPhiChiSquaredInwards_buf;
Buf<TDev, float> rzChiSquared_buf;
template <typename TQueue, typename TDevAcc>
pixelTripletsBuffer(unsigned int maxPixelTriplets, TDevAcc const& devAccIn, TQueue& queue)
: pixelSegmentIndices_buf(allocBufWrapper<unsigned int>(devAccIn, maxPixelTriplets, queue)),
tripletIndices_buf(allocBufWrapper<unsigned int>(devAccIn, maxPixelTriplets, queue)),
nPixelTriplets_buf(allocBufWrapper<unsigned int>(devAccIn, 1, queue)),
totOccupancyPixelTriplets_buf(allocBufWrapper<unsigned int>(devAccIn, 1, queue)),
pixelRadius_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
tripletRadius_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
pt_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
eta_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
phi_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
eta_pix_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
phi_pix_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
score_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
isDup_buf(allocBufWrapper<bool>(devAccIn, maxPixelTriplets, queue)),
partOfPT5_buf(allocBufWrapper<bool>(devAccIn, maxPixelTriplets, queue)),
logicalLayers_buf(allocBufWrapper<uint8_t>(devAccIn, maxPixelTriplets * Params_pT3::kLayers, queue)),
hitIndices_buf(allocBufWrapper<unsigned int>(devAccIn, maxPixelTriplets * Params_pT3::kHits, queue)),
lowerModuleIndices_buf(allocBufWrapper<uint16_t>(devAccIn, maxPixelTriplets * Params_pT3::kLayers, queue)),
centerX_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
centerY_buf(allocBufWrapper<FPX>(devAccIn, maxPixelTriplets, queue)),
pixelRadiusError_buf(allocBufWrapper<float>(devAccIn, maxPixelTriplets, queue)),
rPhiChiSquared_buf(allocBufWrapper<float>(devAccIn, maxPixelTriplets, queue)),
rPhiChiSquaredInwards_buf(allocBufWrapper<float>(devAccIn, maxPixelTriplets, queue)),
rzChiSquared_buf(allocBufWrapper<float>(devAccIn, maxPixelTriplets, queue)) {
alpaka::memset(queue, nPixelTriplets_buf, 0u);
alpaka::memset(queue, totOccupancyPixelTriplets_buf, 0u);
alpaka::memset(queue, partOfPT5_buf, false);
alpaka::wait(queue);
}
};
ALPAKA_FN_ACC ALPAKA_FN_INLINE void addPixelTripletToMemory(struct SDL::modules& modulesInGPU,
struct SDL::miniDoublets& mdsInGPU,
struct SDL::segments& segmentsInGPU,
struct SDL::triplets& tripletsInGPU,
struct SDL::pixelTriplets& pixelTripletsInGPU,
unsigned int pixelSegmentIndex,
unsigned int tripletIndex,
float pixelRadius,
float tripletRadius,
float centerX,
float centerY,
float rPhiChiSquared,
float rPhiChiSquaredInwards,
float rzChiSquared,
unsigned int pixelTripletIndex,
float pt,
float eta,
float phi,
float eta_pix,
float phi_pix,
float score) {
pixelTripletsInGPU.pixelSegmentIndices[pixelTripletIndex] = pixelSegmentIndex;
pixelTripletsInGPU.tripletIndices[pixelTripletIndex] = tripletIndex;
pixelTripletsInGPU.pixelRadius[pixelTripletIndex] = __F2H(pixelRadius);
pixelTripletsInGPU.tripletRadius[pixelTripletIndex] = __F2H(tripletRadius);
pixelTripletsInGPU.pt[pixelTripletIndex] = __F2H(pt);
pixelTripletsInGPU.eta[pixelTripletIndex] = __F2H(eta);
pixelTripletsInGPU.phi[pixelTripletIndex] = __F2H(phi);
pixelTripletsInGPU.eta_pix[pixelTripletIndex] = __F2H(eta_pix);
pixelTripletsInGPU.phi_pix[pixelTripletIndex] = __F2H(phi_pix);
pixelTripletsInGPU.isDup[pixelTripletIndex] = false;
pixelTripletsInGPU.score[pixelTripletIndex] = __F2H(score);
pixelTripletsInGPU.centerX[pixelTripletIndex] = __F2H(centerX);
pixelTripletsInGPU.centerY[pixelTripletIndex] = __F2H(centerY);
pixelTripletsInGPU.logicalLayers[Params_pT3::kLayers * pixelTripletIndex] = 0;
pixelTripletsInGPU.logicalLayers[Params_pT3::kLayers * pixelTripletIndex + 1] = 0;
pixelTripletsInGPU.logicalLayers[Params_pT3::kLayers * pixelTripletIndex + 2] =
tripletsInGPU.logicalLayers[tripletIndex * Params_T3::kLayers];
pixelTripletsInGPU.logicalLayers[Params_pT3::kLayers * pixelTripletIndex + 3] =
tripletsInGPU.logicalLayers[tripletIndex * Params_T3::kLayers + 1];
pixelTripletsInGPU.logicalLayers[Params_pT3::kLayers * pixelTripletIndex + 4] =
tripletsInGPU.logicalLayers[tripletIndex * Params_T3::kLayers + 2];
pixelTripletsInGPU.lowerModuleIndices[Params_pT3::kLayers * pixelTripletIndex] =
segmentsInGPU.innerLowerModuleIndices[pixelSegmentIndex];
pixelTripletsInGPU.lowerModuleIndices[Params_pT3::kLayers * pixelTripletIndex + 1] =
segmentsInGPU.outerLowerModuleIndices[pixelSegmentIndex];
pixelTripletsInGPU.lowerModuleIndices[Params_pT3::kLayers * pixelTripletIndex + 2] =
tripletsInGPU.lowerModuleIndices[Params_T3::kLayers * tripletIndex];
pixelTripletsInGPU.lowerModuleIndices[Params_pT3::kLayers * pixelTripletIndex + 3] =
tripletsInGPU.lowerModuleIndices[Params_T3::kLayers * tripletIndex + 1];
pixelTripletsInGPU.lowerModuleIndices[Params_pT3::kLayers * pixelTripletIndex + 4] =
tripletsInGPU.lowerModuleIndices[Params_T3::kLayers * tripletIndex + 2];
unsigned int pixelInnerMD = segmentsInGPU.mdIndices[2 * pixelSegmentIndex];
unsigned int pixelOuterMD = segmentsInGPU.mdIndices[2 * pixelSegmentIndex + 1];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex] = mdsInGPU.anchorHitIndices[pixelInnerMD];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 1] = mdsInGPU.outerHitIndices[pixelInnerMD];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 2] = mdsInGPU.anchorHitIndices[pixelOuterMD];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 3] = mdsInGPU.outerHitIndices[pixelOuterMD];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 4] =
tripletsInGPU.hitIndices[Params_T3::kHits * tripletIndex];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 5] =
tripletsInGPU.hitIndices[Params_T3::kHits * tripletIndex + 1];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 6] =
tripletsInGPU.hitIndices[Params_T3::kHits * tripletIndex + 2];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 7] =
tripletsInGPU.hitIndices[Params_T3::kHits * tripletIndex + 3];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 8] =
tripletsInGPU.hitIndices[Params_T3::kHits * tripletIndex + 4];
pixelTripletsInGPU.hitIndices[Params_pT3::kHits * pixelTripletIndex + 9] =
tripletsInGPU.hitIndices[Params_T3::kHits * tripletIndex + 5];
pixelTripletsInGPU.rPhiChiSquared[pixelTripletIndex] = rPhiChiSquared;
pixelTripletsInGPU.rPhiChiSquaredInwards[pixelTripletIndex] = rPhiChiSquaredInwards;
pixelTripletsInGPU.rzChiSquared[pixelTripletIndex] = rzChiSquared;
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runPixelTrackletDefaultAlgopT3(TAcc const& acc,
struct SDL::modules& modulesInGPU,
struct SDL::objectRanges& rangesInGPU,
struct SDL::miniDoublets& mdsInGPU,
struct SDL::segments& segmentsInGPU,
uint16_t& pixelLowerModuleIndex,
uint16_t& outerInnerLowerModuleIndex,
uint16_t& outerOuterLowerModuleIndex,
unsigned int& innerSegmentIndex,
unsigned int& outerSegmentIndex,
float& zOut,
float& rtOut,
float& deltaPhiPos,
float& deltaPhi,
float& betaIn,
float& betaOut,
float& pt_beta,
float& zLo,
float& zHi,
float& rtLo,
float& rtHi,
float& zLoPointed,
float& zHiPointed,
float& sdlCut,
float& betaInCut,
float& betaOutCut,
float& deltaBetaCut,
float& kZ) {
zLo = -999;
zHi = -999;
rtLo = -999;
rtHi = -999;
zLoPointed = -999;
zHiPointed = -999;
kZ = -999;
betaInCut = -999;
short outerInnerLowerModuleSubdet = modulesInGPU.subdets[outerInnerLowerModuleIndex];
short outerOuterLowerModuleSubdet = modulesInGPU.subdets[outerOuterLowerModuleIndex];
unsigned int firstMDIndex = segmentsInGPU.mdIndices[Params_LS::kLayers * innerSegmentIndex];
unsigned int secondMDIndex = segmentsInGPU.mdIndices[Params_LS::kLayers * innerSegmentIndex + 1];
unsigned int thirdMDIndex = segmentsInGPU.mdIndices[Params_LS::kLayers * outerSegmentIndex];
unsigned int fourthMDIndex = segmentsInGPU.mdIndices[Params_LS::kLayers * outerSegmentIndex + 1];
if (outerInnerLowerModuleSubdet == SDL::Barrel and
(outerOuterLowerModuleSubdet == SDL::Barrel or outerOuterLowerModuleSubdet == SDL::Endcap)) {
return runTripletDefaultAlgoPPBB(acc,
modulesInGPU,
rangesInGPU,
mdsInGPU,
segmentsInGPU,
pixelLowerModuleIndex,
outerInnerLowerModuleIndex,
outerOuterLowerModuleIndex,
innerSegmentIndex,
outerSegmentIndex,
firstMDIndex,
secondMDIndex,
thirdMDIndex,
fourthMDIndex,
zOut,
rtOut,
deltaPhiPos,
deltaPhi,
betaIn,
betaOut,
pt_beta,
zLo,
zHi,
zLoPointed,
zHiPointed,
sdlCut,
betaOutCut,
deltaBetaCut);
} else if (outerInnerLowerModuleSubdet == SDL::Endcap and outerOuterLowerModuleSubdet == SDL::Endcap) {
return runTripletDefaultAlgoPPEE(acc,
modulesInGPU,
rangesInGPU,
mdsInGPU,
segmentsInGPU,
pixelLowerModuleIndex,
outerInnerLowerModuleIndex,
outerOuterLowerModuleIndex,
innerSegmentIndex,
outerSegmentIndex,
firstMDIndex,
secondMDIndex,
thirdMDIndex,
fourthMDIndex,
zOut,
rtOut,
deltaPhiPos,
deltaPhi,
betaIn,
betaOut,
pt_beta,
zLo,
rtLo,
rtHi,
sdlCut,
betaInCut,
betaOutCut,
deltaBetaCut,
kZ);
}
return false;
};
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT3RZChiSquaredCuts(struct SDL::modules& modulesInGPU,
uint16_t& lowerModuleIndex1,
uint16_t& lowerModuleIndex2,
uint16_t& lowerModuleIndex3,
float& rzChiSquared) {
const int layer1 = modulesInGPU.layers[lowerModuleIndex1] +
6 * (modulesInGPU.subdets[lowerModuleIndex1] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex1] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex1] == SDL::TwoS);
const int layer2 = modulesInGPU.layers[lowerModuleIndex2] +
6 * (modulesInGPU.subdets[lowerModuleIndex2] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex2] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex2] == SDL::TwoS);
const int layer3 = modulesInGPU.layers[lowerModuleIndex3] +
6 * (modulesInGPU.subdets[lowerModuleIndex3] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex3] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex3] == SDL::TwoS);
if (layer1 == 8 and layer2 == 9 and layer3 == 10) {
return rzChiSquared < 13.6067f;
} else if (layer1 == 8 and layer2 == 9 and layer3 == 15) {
return rzChiSquared < 5.5953f;
} else if (layer1 == 7 and layer2 == 8 and layer3 == 9) {
return rzChiSquared < 3.9263f;
}
/*
else if(layer1 == 7 and layer2 == 8 and layer3 == 14)
{
// PS+PS+2S in endcap layers 1+2+3, which is not really feasible in the current geometry,
// without skipping barrel layers 1 and 2 (not allowed by algorithm logic).
}
*/
else if (layer1 == 1 and layer2 == 2 and layer3 == 3) {
return rzChiSquared < 9.4377f;
} else if (layer1 == 1 and layer2 == 2 and layer3 == 7) {
return rzChiSquared < 9.9975f;
} else if (layer1 == 1 and layer2 == 7 and layer3 == 8) {
return rzChiSquared < 8.6369f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 7) {
return rzChiSquared < 37.945f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 12) {
return rzChiSquared < 43.0167f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 4) {
return rzChiSquared < 8.6923f;
} else if (layer1 == 2 and layer2 == 7 and layer3 == 8) {
return rzChiSquared < 11.9672f;
} else if (layer1 == 2 and layer2 == 7 and layer3 == 13) {
return rzChiSquared < 16.2133f;
}
//default - category not found!
return true;
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computeChiSquaredpT3(TAcc const& acc,
unsigned int nPoints,
float* xs,
float* ys,
float* delta1,
float* delta2,
float* slopes,
bool* isFlat,
float g,
float f,
float radius) {
//given values of (g, f, radius) and a set of points (and its uncertainties)
//compute chi squared
float c = g * g + f * f - radius * radius;
float chiSquared = 0.f;
float absArctanSlope, angleM, xPrime, yPrime, sigma2;
for (size_t i = 0; i < nPoints; i++) {
absArctanSlope = ((slopes[i] != SDL::SDL_INF) ? alpaka::math::abs(acc, alpaka::math::atan(acc, slopes[i]))
: 0.5f * float(M_PI));
if (xs[i] > 0 and ys[i] > 0) {
angleM = 0.5f * float(M_PI) - absArctanSlope;
} else if (xs[i] < 0 and ys[i] > 0) {
angleM = absArctanSlope + 0.5f * float(M_PI);
} else if (xs[i] < 0 and ys[i] < 0) {
angleM = -(absArctanSlope + 0.5f * float(M_PI));
} else if (xs[i] > 0 and ys[i] < 0) {
angleM = -(0.5f * float(M_PI) - absArctanSlope);
} else {
angleM = 0;
}
if (not isFlat[i]) {
xPrime = xs[i] * alpaka::math::cos(acc, angleM) + ys[i] * alpaka::math::sin(acc, angleM);
yPrime = ys[i] * alpaka::math::cos(acc, angleM) - xs[i] * alpaka::math::sin(acc, angleM);
} else {
xPrime = xs[i];
yPrime = ys[i];
}
sigma2 = 4 * ((xPrime * delta1[i]) * (xPrime * delta1[i]) + (yPrime * delta2[i]) * (yPrime * delta2[i]));
chiSquared += (xs[i] * xs[i] + ys[i] * ys[i] - 2 * g * xs[i] - 2 * f * ys[i] + c) *
(xs[i] * xs[i] + ys[i] * ys[i] - 2 * g * xs[i] - 2 * f * ys[i] + c) / sigma2;
}
return chiSquared;
};
//TODO: merge this one and the pT5 function later into a single function
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT3RPhiChiSquared(TAcc const& acc,
struct SDL::modules& modulesInGPU,
uint16_t* lowerModuleIndices,
float& g,
float& f,
float& radius,
float* xs,
float* ys) {
float delta1[3]{}, delta2[3]{}, slopes[3];
bool isFlat[3]{};
float chiSquared = 0;
float inv1 = widthPS / width2S;
float inv2 = pixelPSZpitch / width2S;
for (size_t i = 0; i < 3; i++) {
ModuleType moduleType = modulesInGPU.moduleType[lowerModuleIndices[i]];
short moduleSubdet = modulesInGPU.subdets[lowerModuleIndices[i]];
short moduleSide = modulesInGPU.sides[lowerModuleIndices[i]];
float drdz = modulesInGPU.drdzs[lowerModuleIndices[i]];
slopes[i] = modulesInGPU.dxdys[lowerModuleIndices[i]];
//category 1 - barrel PS flat
if (moduleSubdet == Barrel and moduleType == PS and moduleSide == Center) {
delta1[i] = inv1;
delta2[i] = inv1;
slopes[i] = -999;
isFlat[i] = true;
}
//category 2 - barrel 2S
else if (moduleSubdet == Barrel and moduleType == TwoS) {
delta1[i] = 1;
delta2[i] = 1;
slopes[i] = -999;
isFlat[i] = true;
}
//category 3 - barrel PS tilted
else if (moduleSubdet == Barrel and moduleType == PS and moduleSide != Center) {
delta1[i] = inv1;
isFlat[i] = false;
delta2[i] = (inv2 * drdz / alpaka::math::sqrt(acc, 1 + drdz * drdz));
}
//category 4 - endcap PS
else if (moduleSubdet == Endcap and moduleType == PS) {
delta1[i] = inv1;
isFlat[i] = false;
/*
despite the type of the module layer of the lower module index, all anchor
hits are on the pixel side and all non-anchor hits are on the strip side!
*/
delta2[i] = inv2;
}
//category 5 - endcap 2S
else if (moduleSubdet == Endcap and moduleType == TwoS) {
delta1[i] = 1;
delta2[i] = 500 * inv1;
isFlat[i] = false;
}
#ifdef Warnings
else {
printf("ERROR!!!!! I SHOULDN'T BE HERE!!!! subdet = %d, type = %d, side = %d\n",
moduleSubdet,
moduleType,
moduleSide);
}
#endif
}
chiSquared = computeChiSquaredpT3(acc, 3, xs, ys, delta1, delta2, slopes, isFlat, g, f, radius);
return chiSquared;
};
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT3RPhiChiSquaredInwards(
struct SDL::modules& modulesInGPU, float& g, float& f, float& r, float* xPix, float* yPix) {
float residual = (xPix[0] - g) * (xPix[0] - g) + (yPix[0] - f) * (yPix[0] - f) - r * r;
float chiSquared = residual * residual;
residual = (xPix[1] - g) * (xPix[1] - g) + (yPix[1] - f) * (yPix[1] - f) - r * r;
chiSquared += residual * residual;
chiSquared *= 0.5f;
return chiSquared;
};
//90pc threshold
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT3RPhiChiSquaredCuts(struct SDL::modules& modulesInGPU,
uint16_t& lowerModuleIndex1,
uint16_t& lowerModuleIndex2,
uint16_t& lowerModuleIndex3,
float& chiSquared) {
const int layer1 = modulesInGPU.layers[lowerModuleIndex1] +
6 * (modulesInGPU.subdets[lowerModuleIndex1] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex1] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex1] == SDL::TwoS);
const int layer2 = modulesInGPU.layers[lowerModuleIndex2] +
6 * (modulesInGPU.subdets[lowerModuleIndex2] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex2] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex2] == SDL::TwoS);
const int layer3 = modulesInGPU.layers[lowerModuleIndex3] +
6 * (modulesInGPU.subdets[lowerModuleIndex3] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex3] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex3] == SDL::TwoS);
if (layer1 == 8 and layer2 == 9 and layer3 == 10) {
return chiSquared < 7.003f;
} else if (layer1 == 8 and layer2 == 9 and layer3 == 15) {
return chiSquared < 0.5f;
} else if (layer1 == 7 and layer2 == 8 and layer3 == 9) {
return chiSquared < 8.046f;
} else if (layer1 == 7 and layer2 == 8 and layer3 == 14) {
return chiSquared < 0.575f;
} else if (layer1 == 1 and layer2 == 2 and layer3 == 7) {
return chiSquared < 5.304f;
} else if (layer1 == 1 and layer2 == 2 and layer3 == 3) {
return chiSquared < 10.6211f;
} else if (layer1 == 1 and layer2 == 7 and layer3 == 8) {
return chiSquared < 4.617f;
} else if (layer1 == 2 and layer2 == 7 and layer3 == 8) {
return chiSquared < 8.046f;
} else if (layer1 == 2 and layer2 == 7 and layer3 == 13) {
return chiSquared < 0.435f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 7) {
return chiSquared < 9.244f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 12) {
return chiSquared < 0.287f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 4) {
return chiSquared < 18.509f;
}
return true;
};
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT3RPhiChiSquaredInwardsCuts(struct SDL::modules& modulesInGPU,
uint16_t& lowerModuleIndex1,
uint16_t& lowerModuleIndex2,
uint16_t& lowerModuleIndex3,
float& chiSquared) {
const int layer1 = modulesInGPU.layers[lowerModuleIndex1] +
6 * (modulesInGPU.subdets[lowerModuleIndex1] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex1] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex1] == SDL::TwoS);
const int layer2 = modulesInGPU.layers[lowerModuleIndex2] +
6 * (modulesInGPU.subdets[lowerModuleIndex2] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex2] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex2] == SDL::TwoS);
const int layer3 = modulesInGPU.layers[lowerModuleIndex3] +
6 * (modulesInGPU.subdets[lowerModuleIndex3] == SDL::Endcap) +
5 * (modulesInGPU.subdets[lowerModuleIndex3] == SDL::Endcap and
modulesInGPU.moduleType[lowerModuleIndex3] == SDL::TwoS);
if (layer1 == 7 and layer2 == 8 and layer3 == 9) // endcap layer 1,2,3, ps
{
return chiSquared < 22016.8055f;
} else if (layer1 == 7 and layer2 == 8 and layer3 == 14) // endcap layer 1,2,3 layer3->2s
{
return chiSquared < 935179.56807f;
} else if (layer1 == 8 and layer2 == 9 and layer3 == 10) // endcap layer 2,3,4
{
return chiSquared < 29064.12959f;
} else if (layer1 == 8 and layer2 == 9 and layer3 == 15) // endcap layer 2,3,4, layer3->2s
{
return chiSquared < 935179.5681f;
} else if (layer1 == 1 and layer2 == 2 and layer3 == 3) // barrel 1,2,3
{
return chiSquared < 1370.0113195101474f;
} else if (layer1 == 1 and layer2 == 2 and layer3 == 7) // barrel 1,2 endcap 1
{
return chiSquared < 5492.110048314815f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 4) // barrel 2,3,4
{
return chiSquared < 4160.410806470067f;
} else if (layer1 == 1 and layer2 == 7 and layer3 == 8) // barrel 1, endcap 1,2
{
return chiSquared < 29064.129591225726f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 7) // barrel 2,3 endcap 1
{
return chiSquared < 12634.215376250893f;
} else if (layer1 == 2 and layer2 == 3 and layer3 == 12) // barrel 2,3, endcap 1->2s
{
return chiSquared < 353821.69361145404f;
} else if (layer1 == 2 and layer2 == 7 and layer3 == 8) // barrel2, endcap 1,2
{
return chiSquared < 33393.26076341235f;
} else if (layer1 == 2 and layer2 == 7 and layer3 == 13) //barrel 2, endcap 1, endcap2->2s
{
return chiSquared < 935179.5680742573f;
}
return true;
};
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool checkIntervalOverlappT3(const float& firstMin,
const float& firstMax,
const float& secondMin,
const float& secondMax) {
return ((firstMin <= secondMin) && (secondMin < firstMax)) || ((secondMin < firstMin) && (firstMin < secondMax));
};
/*bounds for high Pt taken from : http://uaf-10.t2.ucsd.edu/~bsathian/SDL/T5_efficiency/efficiencies/new_efficiencies/efficiencies_20210513_T5_recovering_high_Pt_efficiencies/highE_radius_matching/highE_bounds.txt */
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passRadiusCriterionBBB(TAcc const& acc,
float const& pixelRadius,
float const& pixelRadiusError,
float const& tripletRadius) {
float tripletInvRadiusErrorBound = 0.15624f;
float pixelInvRadiusErrorBound = 0.17235f;
if (pixelRadius > 2.0f * kR1GeVf) {
pixelInvRadiusErrorBound = 0.6375f;
tripletInvRadiusErrorBound = 0.6588f;
}
float tripletRadiusInvMax = (1 + tripletInvRadiusErrorBound) / tripletRadius;
float tripletRadiusInvMin = alpaka::math::max(acc, (1 - tripletInvRadiusErrorBound) / tripletRadius, 0.0f);
float pixelRadiusInvMax =
alpaka::math::max(acc, (1 + pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius - pixelRadiusError));
float pixelRadiusInvMin =
alpaka::math::min(acc, (1 - pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius + pixelRadiusError));
return checkIntervalOverlappT3(tripletRadiusInvMin, tripletRadiusInvMax, pixelRadiusInvMin, pixelRadiusInvMax);
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passRadiusCriterionBBE(TAcc const& acc,
float const& pixelRadius,
float const& pixelRadiusError,
float const& tripletRadius) {
float tripletInvRadiusErrorBound = 0.45972f;
float pixelInvRadiusErrorBound = 0.19644f;
if (pixelRadius > 2.0f * kR1GeVf) {
pixelInvRadiusErrorBound = 0.6805f;
tripletInvRadiusErrorBound = 0.8557f;
}
float tripletRadiusInvMax = (1 + tripletInvRadiusErrorBound) / tripletRadius;
float tripletRadiusInvMin = alpaka::math::max(acc, (1 - tripletInvRadiusErrorBound) / tripletRadius, 0.0f);
float pixelRadiusInvMax =
alpaka::math::max(acc, (1 + pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius - pixelRadiusError));
float pixelRadiusInvMin =
alpaka::math::min(acc, (1 - pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius + pixelRadiusError));
return checkIntervalOverlappT3(tripletRadiusInvMin, tripletRadiusInvMax, pixelRadiusInvMin, pixelRadiusInvMax);
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passRadiusCriterionBEE(TAcc const& acc,
float const& pixelRadius,
float const& pixelRadiusError,
float const& tripletRadius) {
float tripletInvRadiusErrorBound = 1.59294f;
float pixelInvRadiusErrorBound = 0.255181f;
if (pixelRadius > 2.0f * kR1GeVf) //as good as not having selections
{
pixelInvRadiusErrorBound = 2.2091f;
tripletInvRadiusErrorBound = 2.3548f;
}
float tripletRadiusInvMax = (1 + tripletInvRadiusErrorBound) / tripletRadius;
float tripletRadiusInvMin = alpaka::math::max(acc, (1 - tripletInvRadiusErrorBound) / tripletRadius, 0.0f);
float pixelRadiusInvMax =
alpaka::math::max(acc, (1 + pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius - pixelRadiusError));
float pixelRadiusInvMin =
alpaka::math::min(acc, (1 - pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius + pixelRadiusError));
pixelRadiusInvMin = alpaka::math::max(acc, pixelRadiusInvMin, 0.0f);
return checkIntervalOverlappT3(tripletRadiusInvMin, tripletRadiusInvMax, pixelRadiusInvMin, pixelRadiusInvMax);
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passRadiusCriterionEEE(TAcc const& acc,
float const& pixelRadius,
float const& pixelRadiusError,
float const& tripletRadius) {
float tripletInvRadiusErrorBound = 1.7006f;
float pixelInvRadiusErrorBound = 0.26367f;
if (pixelRadius > 2.0f * kR1GeVf) //as good as not having selections
{
pixelInvRadiusErrorBound = 2.286f;
tripletInvRadiusErrorBound = 2.436f;
}
float tripletRadiusInvMax = (1 + tripletInvRadiusErrorBound) / tripletRadius;
float tripletRadiusInvMin = alpaka::math::max(acc, (1 - tripletInvRadiusErrorBound) / tripletRadius, 0.0f);
float pixelRadiusInvMax =
alpaka::math::max(acc, (1 + pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius - pixelRadiusError));
float pixelRadiusInvMin =
alpaka::math::min(acc, (1 - pixelInvRadiusErrorBound) / pixelRadius, 1.f / (pixelRadius + pixelRadiusError));
pixelRadiusInvMin = alpaka::math::max(acc, 0.0f, pixelRadiusInvMin);
return checkIntervalOverlappT3(tripletRadiusInvMin, tripletRadiusInvMax, pixelRadiusInvMin, pixelRadiusInvMax);
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passRadiusCriterion(TAcc const& acc,
struct SDL::modules const& modulesInGPU,
float const& pixelRadius,
float const& pixelRadiusError,
float const& tripletRadius,
uint16_t const& lowerModuleIndex,
uint16_t const& middleModuleIndex,
uint16_t const& upperModuleIndex) {
if (modulesInGPU.subdets[lowerModuleIndex] == SDL::Endcap) {
return passRadiusCriterionEEE(acc, pixelRadius, pixelRadiusError, tripletRadius);
} else if (modulesInGPU.subdets[middleModuleIndex] == SDL::Endcap) {
return passRadiusCriterionBEE(acc, pixelRadius, pixelRadiusError, tripletRadius);
} else if (modulesInGPU.subdets[upperModuleIndex] == SDL::Endcap) {
return passRadiusCriterionBBE(acc, pixelRadius, pixelRadiusError, tripletRadius);
} else {
return passRadiusCriterionBBB(acc, pixelRadius, pixelRadiusError, tripletRadius);
}
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT3RZChiSquared(TAcc const& acc,
struct SDL::modules const& modulesInGPU,
const uint16_t* lowerModuleIndices,
const float* rtPix,
const float* xPix,
const float* yPix,
const float* zPix,
const float* rts,
const float* xs,
const float* ys,
const float* zs,
float pixelSegmentPt,
float pixelSegmentPx,
float pixelSegmentPy,
float pixelSegmentPz,
int pixelSegmentCharge) {
float residual = 0;
float error2 = 0;
float RMSE = 0;
float Px = pixelSegmentPx, Py = pixelSegmentPy, Pz = pixelSegmentPz;
int charge = pixelSegmentCharge;
float x1 = xPix[1] / 100;
float y1 = yPix[1] / 100;
float z1 = zPix[1] / 100;
float r1 = rtPix[1] / 100;
float a = -2.f * k2Rinv1GeVf * 100 * charge; // multiply by 100 to make the correct length units
for (size_t i = 0; i < Params_T3::kLayers; i++) {
float zsi = zs[i] / 100;
float rtsi = rts[i] / 100;
uint16_t lowerModuleIndex = lowerModuleIndices[i];
const int moduleType = modulesInGPU.moduleType[lowerModuleIndex];
const int moduleSide = modulesInGPU.sides[lowerModuleIndex];
const int moduleSubdet = modulesInGPU.subdets[lowerModuleIndex];
// calculation is detailed documented here https://indico.cern.ch/event/1185895/contributions/4982756/attachments/2526561/4345805/helix%20pT3%20summarize.pdf
float diffr, diffz;
float p = alpaka::math::sqrt(acc, Px * Px + Py * Py + Pz * Pz);
float rou = a / p;
if (moduleSubdet == SDL::Endcap) {
float s = (zsi - z1) * p / Pz;
float x = x1 + Px / a * alpaka::math::sin(acc, rou * s) - Py / a * (1 - alpaka::math::cos(acc, rou * s));
float y = y1 + Py / a * alpaka::math::sin(acc, rou * s) + Px / a * (1 - alpaka::math::cos(acc, rou * s));
diffr = alpaka::math::abs(acc, rtsi - alpaka::math::sqrt(acc, x * x + y * y)) * 100;
}
if (moduleSubdet == SDL::Barrel) {
float paraA = r1 * r1 + 2 * (Px * Px + Py * Py) / (a * a) + 2 * (y1 * Px - x1 * Py) / a - rtsi * rtsi;
float paraB = 2 * (x1 * Px + y1 * Py) / a;
float paraC = 2 * (y1 * Px - x1 * Py) / a + 2 * (Px * Px + Py * Py) / (a * a);
float A = paraB * paraB + paraC * paraC;
float B = 2 * paraA * paraB;
float C = paraA * paraA - paraC * paraC;
float sol1 = (-B + alpaka::math::sqrt(acc, B * B - 4 * A * C)) / (2 * A);
float sol2 = (-B - alpaka::math::sqrt(acc, B * B - 4 * A * C)) / (2 * A);
float solz1 = alpaka::math::asin(acc, sol1) / rou * Pz / p + z1;
float solz2 = alpaka::math::asin(acc, sol2) / rou * Pz / p + z1;
float diffz1 = alpaka::math::abs(acc, solz1 - zsi) * 100;
float diffz2 = alpaka::math::abs(acc, solz2 - zsi) * 100;
diffz = alpaka::math::min(acc, diffz1, diffz2);
}
residual = moduleSubdet == SDL::Barrel ? diffz : diffr;
//PS Modules
if (moduleType == 0) {
error2 = pixelPSZpitch * pixelPSZpitch;
} else //2S modules
{
error2 = strip2SZpitch * strip2SZpitch;
}
//special dispensation to tilted PS modules!
if (moduleType == 0 and moduleSubdet == SDL::Barrel and moduleSide != Center) {
float drdz = modulesInGPU.drdzs[lowerModuleIndex];
error2 /= (1 + drdz * drdz);
}
RMSE += (residual * residual) / error2;
}
RMSE = alpaka::math::sqrt(acc, 0.2f * RMSE); // Divided by the degree of freedom 5.
return RMSE;
};
template <typename TAcc>
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runPixelTripletDefaultAlgo(TAcc const& acc,
struct SDL::modules& modulesInGPU,
struct SDL::objectRanges& rangesInGPU,
struct SDL::miniDoublets& mdsInGPU,
struct SDL::segments& segmentsInGPU,
struct SDL::triplets& tripletsInGPU,
unsigned int& pixelSegmentIndex,
unsigned int tripletIndex,
float& pixelRadius,
float& pixelRadiusError,
float& tripletRadius,
float& centerX,
float& centerY,
float& rzChiSquared,
float& rPhiChiSquared,
float& rPhiChiSquaredInwards,
bool runChiSquaredCuts = true) {
//run pT4 compatibility between the pixel segment and inner segment, and between the pixel and outer segment of the triplet
uint16_t pixelModuleIndex = segmentsInGPU.innerLowerModuleIndices[pixelSegmentIndex];
uint16_t lowerModuleIndex = tripletsInGPU.lowerModuleIndices[Params_T3::kLayers * tripletIndex];
uint16_t middleModuleIndex = tripletsInGPU.lowerModuleIndices[Params_T3::kLayers * tripletIndex + 1];
uint16_t upperModuleIndex = tripletsInGPU.lowerModuleIndices[Params_T3::kLayers * tripletIndex + 2];
{
//placeholder
float zOut, rtOut, deltaPhiPos, deltaPhi, betaIn, betaOut, pt_beta; //temp stuff
float zLo, zHi, rtLo, rtHi, zLoPointed, zHiPointed, sdlCut, betaInCut, betaOutCut, deltaBetaCut, kZ;
// pixel segment vs inner segment of the triplet
if (not runPixelTrackletDefaultAlgopT3(acc,
modulesInGPU,
rangesInGPU,
mdsInGPU,
segmentsInGPU,
pixelModuleIndex,
lowerModuleIndex,
middleModuleIndex,
pixelSegmentIndex,
tripletsInGPU.segmentIndices[Params_LS::kLayers * tripletIndex],
zOut,
rtOut,
deltaPhiPos,
deltaPhi,
betaIn,
betaOut,
pt_beta,
zLo,
zHi,
rtLo,
rtHi,
zLoPointed,
zHiPointed,
sdlCut,
betaInCut,
betaOutCut,
deltaBetaCut,
kZ))
return false;
//pixel segment vs outer segment of triplet
if (not runPixelTrackletDefaultAlgopT3(acc,
modulesInGPU,
rangesInGPU,
mdsInGPU,
segmentsInGPU,
pixelModuleIndex,
middleModuleIndex,
upperModuleIndex,
pixelSegmentIndex,
tripletsInGPU.segmentIndices[Params_LS::kLayers * tripletIndex + 1],
zOut,
rtOut,
deltaPhiPos,
deltaPhi,
betaIn,
betaOut,
pt_beta,
zLo,
zHi,
rtLo,
rtHi,
zLoPointed,
zHiPointed,
sdlCut,
betaInCut,
betaOutCut,
deltaBetaCut,
kZ))
return false;
}
//pt matching between the pixel ptin and the triplet circle pt
unsigned int pixelSegmentArrayIndex = pixelSegmentIndex - rangesInGPU.segmentModuleIndices[pixelModuleIndex];
float pixelSegmentPt = segmentsInGPU.ptIn[pixelSegmentArrayIndex];
float pixelSegmentPtError = segmentsInGPU.ptErr[pixelSegmentArrayIndex];
float pixelSegmentPx = segmentsInGPU.px[pixelSegmentArrayIndex];
float pixelSegmentPy = segmentsInGPU.py[pixelSegmentArrayIndex];
float pixelSegmentPz = segmentsInGPU.pz[pixelSegmentArrayIndex];
int pixelSegmentCharge = segmentsInGPU.charge[pixelSegmentArrayIndex];
float pixelG = segmentsInGPU.circleCenterX[pixelSegmentArrayIndex];
float pixelF = segmentsInGPU.circleCenterY[pixelSegmentArrayIndex];
float pixelRadiusPCA = segmentsInGPU.circleRadius[pixelSegmentArrayIndex];
unsigned int pixelInnerMDIndex = segmentsInGPU.mdIndices[Params_pLS::kLayers * pixelSegmentIndex];
unsigned int pixelOuterMDIndex = segmentsInGPU.mdIndices[Params_pLS::kLayers * pixelSegmentIndex + 1];
pixelRadius = pixelSegmentPt * kR1GeVf;
pixelRadiusError = pixelSegmentPtError * kR1GeVf;
unsigned int tripletInnerSegmentIndex = tripletsInGPU.segmentIndices[2 * tripletIndex];
unsigned int tripletOuterSegmentIndex = tripletsInGPU.segmentIndices[2 * tripletIndex + 1];
unsigned int firstMDIndex = segmentsInGPU.mdIndices[2 * tripletInnerSegmentIndex];
unsigned int secondMDIndex = segmentsInGPU.mdIndices[2 * tripletInnerSegmentIndex + 1];
unsigned int thirdMDIndex = segmentsInGPU.mdIndices[2 * tripletOuterSegmentIndex + 1];
float xs[Params_T3::kLayers] = {
mdsInGPU.anchorX[firstMDIndex], mdsInGPU.anchorX[secondMDIndex], mdsInGPU.anchorX[thirdMDIndex]};
float ys[Params_T3::kLayers] = {
mdsInGPU.anchorY[firstMDIndex], mdsInGPU.anchorY[secondMDIndex], mdsInGPU.anchorY[thirdMDIndex]};
float g, f;
tripletRadius = tripletsInGPU.circleRadius[tripletIndex];
g = tripletsInGPU.circleCenterX[tripletIndex];
f = tripletsInGPU.circleCenterY[tripletIndex];
if (not passRadiusCriterion(acc,
modulesInGPU,
pixelRadius,
pixelRadiusError,
tripletRadius,
lowerModuleIndex,
middleModuleIndex,
upperModuleIndex))
return false;
uint16_t lowerModuleIndices[Params_T3::kLayers] = {lowerModuleIndex, middleModuleIndex, upperModuleIndex};
if (runChiSquaredCuts and pixelSegmentPt < 5.0f) {
float rts[Params_T3::kLayers] = {
mdsInGPU.anchorRt[firstMDIndex], mdsInGPU.anchorRt[secondMDIndex], mdsInGPU.anchorRt[thirdMDIndex]};
float zs[Params_T3::kLayers] = {
mdsInGPU.anchorZ[firstMDIndex], mdsInGPU.anchorZ[secondMDIndex], mdsInGPU.anchorZ[thirdMDIndex]};
float rtPix[Params_pLS::kLayers] = {mdsInGPU.anchorRt[pixelInnerMDIndex], mdsInGPU.anchorRt[pixelOuterMDIndex]};
float xPix[Params_pLS::kLayers] = {mdsInGPU.anchorX[pixelInnerMDIndex], mdsInGPU.anchorX[pixelOuterMDIndex]};
float yPix[Params_pLS::kLayers] = {mdsInGPU.anchorY[pixelInnerMDIndex], mdsInGPU.anchorY[pixelOuterMDIndex]};
float zPix[Params_pLS::kLayers] = {mdsInGPU.anchorZ[pixelInnerMDIndex], mdsInGPU.anchorZ[pixelOuterMDIndex]};
rzChiSquared = computePT3RZChiSquared(acc,
modulesInGPU,
lowerModuleIndices,
rtPix,
xPix,
yPix,
zPix,
rts,
xs,
ys,
zs,
pixelSegmentPt,
pixelSegmentPx,
pixelSegmentPy,
pixelSegmentPz,
pixelSegmentCharge);
if (not passPT3RZChiSquaredCuts(
modulesInGPU, lowerModuleIndex, middleModuleIndex, upperModuleIndex, rzChiSquared))
return false;
} else {
rzChiSquared = -1;
}
rPhiChiSquared =