Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Groups
Pages
CUDADataFormats
Track
interface
TrackSoAHeterogeneousT.h
Go to the documentation of this file.
1
#ifndef CUDADataFormats_Track_TrackHeterogeneousT_H
2
#define CUDADataFormats_Track_TrackHeterogeneousT_H
3
4
#include <string>
5
#include <algorithm>
6
7
#include "
CUDADataFormats/Track/interface/TrajectoryStateSoAT.h
"
8
#include "
Geometry/TrackerGeometryBuilder/interface/phase1PixelTopology.h
"
9
#include "
HeterogeneousCore/CUDAUtilities/interface/HistoContainer.h
"
10
11
#include "
CUDADataFormats/Common/interface/HeterogeneousSoA.h
"
12
13
namespace
pixelTrack {
14
enum class
Quality
: uint8_t {
bad
= 0,
edup
,
dup
,
loose
,
strict
,
tight
,
highPurity
,
notQuality
};
15
constexpr uint32_t
qualitySize
{uint8_t(
Quality::notQuality
)};
16
const
std::string
qualityName
[
qualitySize
]{
"bad"
,
"edup"
,
"dup"
,
"loose"
,
"strict"
,
"tight"
,
"highPurity"
};
17
inline
Quality
qualityByName
(
std::string
const
&
name
) {
18
auto
qp =
std::find
(
qualityName
,
qualityName
+
qualitySize
, name) -
qualityName
;
19
return
static_cast<
Quality
>
(qp);
20
}
21
}
// namespace pixelTrack
22
23
template
<
int
32_t S>
24
class
TrackSoAHeterogeneousT
{
25
public
:
26
static
constexpr int32_t
stride
() {
return
S
; }
27
28
using
Quality
=
pixelTrack::Quality
;
29
using
hindex_type
= uint32_t;
30
using
HitContainer
=
cms::cuda::OneToManyAssoc<hindex_type, S + 1, 5 * S>
;
31
32
// Always check quality is at least loose!
33
// CUDA does not support enums in __lgc ...
34
private
:
35
eigenSoA::ScalarSoA<uint8_t, S>
quality_
;
36
37
public
:
38
constexpr
Quality
quality
(int32_t
i
)
const
{
return
(
Quality
)(
quality_
(i)); }
39
constexpr
Quality
&
quality
(int32_t
i
) {
return
(
Quality
&)(
quality_
(i)); }
40
constexpr
Quality
const
*
qualityData
()
const
{
return
(
Quality
const
*)(
quality_
.
data
()); }
41
constexpr
Quality
*
qualityData
() {
return
(
Quality
*)(
quality_
.
data
()); }
42
43
// this is chi2/ndof as not necessarely all hits are used in the fit
44
eigenSoA::ScalarSoA<float, S>
chi2
;
45
46
eigenSoA::ScalarSoA<int8_t, S>
nLayers
;
47
48
constexpr
int
nHits
(
int
i
)
const
{
return
detIndices
.size(i); }
49
50
constexpr
bool
isTriplet
(
int
i
)
const
{
return
nLayers
(i) == 3; }
51
52
constexpr
int
computeNumberOfLayers
(int32_t
i
)
const
{
53
// layers are in order and we assume tracks are either forward or backward
54
auto
pdet =
detIndices
.begin(i);
55
int
nl = 1;
56
auto
ol =
phase1PixelTopology::getLayer
(*pdet);
57
for
(; pdet <
detIndices
.end(i); ++pdet) {
58
auto
il =
phase1PixelTopology::getLayer
(*pdet);
59
if
(il != ol)
60
++nl;
61
ol = il;
62
}
63
return
nl;
64
}
65
66
// State at the Beam spot
67
// phi,tip,1/pt,cotan(theta),zip
68
TrajectoryStateSoAT<S>
stateAtBS
;
69
eigenSoA::ScalarSoA<float, S>
eta
;
70
eigenSoA::ScalarSoA<float, S>
pt
;
71
constexpr
float
charge
(int32_t
i
)
const
{
return
std::copysign(1.
f
,
stateAtBS
.state(i)(2)); }
72
constexpr
float
phi
(int32_t
i
)
const
{
return
stateAtBS
.state(i)(0); }
73
constexpr
float
tip
(int32_t
i
)
const
{
return
stateAtBS
.state(i)(1); }
74
constexpr
float
zip
(int32_t
i
)
const
{
return
stateAtBS
.state(i)(4); }
75
76
// state at the detector of the outermost hit
77
// representation to be decided...
78
// not yet filled on GPU
79
// TrajectoryStateSoA<S> stateAtOuterDet;
80
81
HitContainer
hitIndices
;
82
HitContainer
detIndices
;
83
};
84
85
namespace
pixelTrack {
86
87
#ifdef GPU_SMALL_EVENTS
88
// kept for testing and debugging
89
constexpr uint32_t
maxNumber
() {
return
2 * 1024; }
90
#else
91
// tested on MC events with 55-75 pileup events
92
constexpr uint32_t
maxNumber
() {
return
32 * 1024; }
93
#endif
94
95
using
TrackSoA
=
TrackSoAHeterogeneousT<maxNumber()>
;
96
using
TrajectoryState
=
TrajectoryStateSoAT<maxNumber()>
;
97
using
HitContainer
=
TrackSoA::HitContainer
;
98
99
}
// namespace pixelTrack
100
101
#endif // CUDADataFormats_Track_TrackHeterogeneousT_H
pixelTrack::Quality::notQuality
TrackSoAHeterogeneousT::nLayers
eigenSoA::ScalarSoA< int8_t, S > nLayers
Definition:
TrackSoAHeterogeneousT.h:46
pixelTrack::Quality::strict
TrackSoAHeterogeneousT::isTriplet
constexpr bool isTriplet(int i) const
Definition:
TrackSoAHeterogeneousT.h:50
eigenSoA::ScalarSoA< uint8_t, S >
mps_fire.i
i
Definition:
mps_fire.py:428
TrackSoAHeterogeneousT::computeNumberOfLayers
constexpr int computeNumberOfLayers(int32_t i) const
Definition:
TrackSoAHeterogeneousT.h:52
pixelTrack::Quality
Quality
Definition:
TrackSoAHeterogeneousT.h:14
pixelTrack::qualityByName
Quality qualityByName(std::string const &name)
Definition:
TrackSoAHeterogeneousT.h:17
pixelTrack::qualityName
const std::string qualityName[qualitySize]
Definition:
TrackSoAHeterogeneousT.h:16
TrackSoAHeterogeneousT::qualityData
constexpr Quality * qualityData()
Definition:
TrackSoAHeterogeneousT.h:41
TrackSoAHeterogeneousT::stride
static constexpr int32_t stride()
Definition:
TrackSoAHeterogeneousT.h:26
pixelTrack::Quality::edup
TrackSoAHeterogeneousT::qualityData
constexpr Quality const * qualityData() const
Definition:
TrackSoAHeterogeneousT.h:40
pixelTrack::maxNumber
constexpr uint32_t maxNumber()
Definition:
TrackSoAHeterogeneousT.h:92
TrackSoAHeterogeneousT::zip
constexpr float zip(int32_t i) const
Definition:
TrackSoAHeterogeneousT.h:74
TrackSoAHeterogeneousT::quality_
eigenSoA::ScalarSoA< uint8_t, S > quality_
Definition:
TrackSoAHeterogeneousT.h:35
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition:
FindCaloHit.cc:19
TrajectoryStateSoAT
Definition:
TrajectoryStateSoAT.h:8
mergeVDriftHistosByStation.name
string name
Definition:
mergeVDriftHistosByStation.py:78
TrackSoAHeterogeneousT::hindex_type
uint32_t hindex_type
Definition:
TrackSoAHeterogeneousT.h:29
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
phase1PixelTopology.h
TrackSoAHeterogeneousT::charge
constexpr float charge(int32_t i) const
Definition:
TrackSoAHeterogeneousT.h:71
TrackSoAHeterogeneousT::pt
eigenSoA::ScalarSoA< float, S > pt
Definition:
TrackSoAHeterogeneousT.h:70
TrackSoAHeterogeneousT::chi2
eigenSoA::ScalarSoA< float, S > chi2
Definition:
TrackSoAHeterogeneousT.h:44
TrackSoAHeterogeneousT::tip
constexpr float tip(int32_t i) const
Definition:
TrackSoAHeterogeneousT.h:73
TrackSoAHeterogeneousT::nHits
constexpr int nHits(int i) const
Definition:
TrackSoAHeterogeneousT.h:48
pixelTrack::Quality::dup
pixelTrack::qualitySize
constexpr uint32_t qualitySize
Definition:
TrackSoAHeterogeneousT.h:15
pixelTrack::Quality::tight
TrackSoAHeterogeneousT::stateAtBS
TrajectoryStateSoAT< S > stateAtBS
Definition:
TrackSoAHeterogeneousT.h:68
phase1PixelTopology::getLayer
constexpr uint8_t getLayer(uint32_t detId)
Definition:
phase1PixelTopology.h:112
TrajectoryStateSoAT.h
HistoContainer.h
TrackSoAHeterogeneousT::HitContainer
cms::cuda::OneToManyAssoc< hindex_type, S+1, 5 *S > HitContainer
Definition:
TrackSoAHeterogeneousT.h:30
pixelTrack::Quality::highPurity
S
double S(const TLorentzVector &, const TLorentzVector &)
Definition:
Particle.cc:97
pixelTrack::Quality::loose
TrackSoAHeterogeneousT::phi
constexpr float phi(int32_t i) const
Definition:
TrackSoAHeterogeneousT.h:72
TrackSoAHeterogeneousT::detIndices
HitContainer detIndices
Definition:
TrackSoAHeterogeneousT.h:82
TrackSoAHeterogeneousT::quality
constexpr Quality quality(int32_t i) const
Definition:
TrackSoAHeterogeneousT.h:38
TrackSoAHeterogeneousT
Definition:
TrackSoAHeterogeneousT.h:24
eigenSoA::ScalarSoA::data
__host__ __device__ constexpr Scalar * data()
Definition:
eigenSoA.h:26
TrackSoAHeterogeneousT::quality
constexpr Quality & quality(int32_t i)
Definition:
TrackSoAHeterogeneousT.h:39
cms::cuda::OneToManyAssoc
Definition:
OneToManyAssoc.h:143
TrackSoAHeterogeneousT::eta
eigenSoA::ScalarSoA< float, S > eta
Definition:
TrackSoAHeterogeneousT.h:69
TrackSoAHeterogeneousT::hitIndices
HitContainer hitIndices
Definition:
TrackSoAHeterogeneousT.h:81
pixelTrack::Quality::bad
validate-o2o-wbm.f
tuple f
Definition:
validate-o2o-wbm.py:27
HeterogeneousSoA.h
Generated for CMSSW Reference Manual by
1.8.5