CMS 3D CMS Logo

PixelTrackProducerFromSoA.cc
Go to the documentation of this file.
23 
28 
32 
33 #include "storeTracks.h"
35 
41 public:
42  using IndToEdm = std::vector<uint16_t>;
43 
44  explicit PixelTrackProducerFromSoA(const edm::ParameterSet &iConfig);
45  ~PixelTrackProducerFromSoA() override = default;
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
48 
49  // using HitModuleStart = std::array<uint32_t, gpuClustering::maxNumModules + 1>;
51 
52 private:
53  void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override;
54 
55  // Event Data tokens
60  // Event Setup tokens
63 
64  int32_t const minNumberOfHits_;
66 };
67 
69  : tBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
70  tokenTrack_(consumes<PixelTrackHeterogeneous>(iConfig.getParameter<edm::InputTag>("trackSrc"))),
71  cpuHits_(consumes<SiPixelRecHitCollectionNew>(iConfig.getParameter<edm::InputTag>("pixelRecHitLegacySrc"))),
72  hmsToken_(consumes<HMSstorage>(iConfig.getParameter<edm::InputTag>("pixelRecHitLegacySrc"))),
73  idealMagneticFieldToken_(esConsumes()),
74  ttTopoToken_(esConsumes()),
75  minNumberOfHits_(iConfig.getParameter<int>("minNumberOfHits")),
76  minQuality_(pixelTrack::qualityByName(iConfig.getParameter<std::string>("minQuality"))) {
78  throw cms::Exception("PixelTrackConfiguration")
79  << iConfig.getParameter<std::string>("minQuality") + " is not a pixelTrack::Quality";
80  }
82  throw cms::Exception("PixelTrackConfiguration")
83  << iConfig.getParameter<std::string>("minQuality") + " not supported";
84  }
85  produces<TrackingRecHitCollection>();
86  produces<reco::TrackExtraCollection>();
87  // TrackCollection refers to TrackingRechit and TrackExtra
88  // collections, need to declare its production after them to work
89  // around a rare race condition in framework scheduling
90  produces<reco::TrackCollection>();
91  produces<IndToEdm>();
92 }
93 
96  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
97  desc.add<edm::InputTag>("trackSrc", edm::InputTag("pixelTracksSoA"));
98  desc.add<edm::InputTag>("pixelRecHitLegacySrc", edm::InputTag("siPixelRecHitsPreSplittingLegacy"));
99  desc.add<int>("minNumberOfHits", 0);
100  desc.add<std::string>("minQuality", "loose");
101  descriptions.addWithDefaultLabel(desc);
102 }
103 
106  const edm::EventSetup &iSetup) const {
107  // enum class Quality : uint8_t { bad = 0, edup, dup, loose, strict, tight, highPurity };
116 
117  // std::cout << "Converting gpu helix in reco tracks" << std::endl;
118 
119  auto indToEdmP = std::make_unique<IndToEdm>();
120  auto &indToEdm = *indToEdmP;
121 
122  auto const &idealField = iSetup.getData(idealMagneticFieldToken_);
123 
125 
126  auto const &httopo = iSetup.getData(ttTopoToken_);
127 
128  const auto &bsh = iEvent.get(tBeamSpot_);
129  GlobalPoint bs(bsh.x0(), bsh.y0(), bsh.z0());
130 
131  auto const &rechits = iEvent.get(cpuHits_);
132  std::vector<TrackingRecHit const *> hitmap;
133  auto const &rcs = rechits.data();
134  auto nhits = rcs.size();
135  hitmap.resize(nhits, nullptr);
136 
137  auto const *hitsModuleStart = iEvent.get(hmsToken_).get();
138  auto fc = hitsModuleStart;
139 
140  for (auto const &h : rcs) {
141  auto const &thit = static_cast<BaseTrackerRecHit const &>(h);
142  auto detI = thit.det()->index();
143  auto const &clus = thit.firstClusterRef();
144  assert(clus.isPixel());
145  auto i = fc[detI] + clus.pixelCluster().originalId();
146  if (i >= hitmap.size())
147  hitmap.resize(i + 256, nullptr); // only in case of hit overflow in one module
148  assert(nullptr == hitmap[i]);
149  hitmap[i] = &h;
150  }
151 
152  std::vector<const TrackingRecHit *> hits;
153  hits.reserve(5);
154 
155  const auto &tsoa = *iEvent.get(tokenTrack_);
156 
157  auto const *quality = tsoa.qualityData();
158  auto const &fit = tsoa.stateAtBS;
159  auto const &hitIndices = tsoa.hitIndices;
160  auto nTracks = tsoa.nTracks();
161 
162  tracks.reserve(nTracks);
163 
164  int32_t nt = 0;
165 
166  //sort index by pt
167  std::vector<int32_t> sortIdxs(nTracks);
168  std::iota(sortIdxs.begin(), sortIdxs.end(), 0);
169  std::sort(
170  sortIdxs.begin(), sortIdxs.end(), [&](int32_t const i1, int32_t const i2) { return tsoa.pt(i1) > tsoa.pt(i2); });
171 
172  //store the index of the SoA: indToEdm[index_SoAtrack] -> index_edmTrack (if it exists)
173  indToEdm.resize(sortIdxs.size(), -1);
174  for (const auto &it : sortIdxs) {
175  auto nHits = tsoa.nHits(it);
176  assert(nHits >= 3);
177  auto q = quality[it];
178  if (q < minQuality_)
179  continue;
180  if (nHits < minNumberOfHits_)
181  continue;
182  indToEdm[it] = nt;
183  ++nt;
184 
185  hits.resize(nHits);
186  auto b = hitIndices.begin(it);
187  for (int iHit = 0; iHit < nHits; ++iHit)
188  hits[iHit] = hitmap[*(b + iHit)];
189 
190  // mind: this values are respect the beamspot!
191 
192  float chi2 = tsoa.chi2(it);
193  float phi = tsoa.phi(it);
194 
195  riemannFit::Vector5d ipar, opar;
196  riemannFit::Matrix5d icov, ocov;
197  fit.copyToDense(ipar, icov, it);
198  riemannFit::transformToPerigeePlane(ipar, icov, opar, ocov);
199 
200  LocalTrajectoryParameters lpar(opar(0), opar(1), opar(2), opar(3), opar(4), 1.);
202  for (int i = 0; i < 5; ++i)
203  for (int j = i; j < 5; ++j)
204  m(i, j) = ocov(i, j);
205 
206  float sp = std::sin(phi);
207  float cp = std::cos(phi);
208  Surface::RotationType rot(sp, -cp, 0, 0, 0, -1.f, cp, sp, 0);
209 
210  Plane impPointPlane(bs, rot);
212  impPointPlane.toGlobal(lpar.position()), impPointPlane.toGlobal(lpar.momentum()), lpar.charge(), &idealField);
213  JacobianLocalToCurvilinear jl2c(impPointPlane, lpar, idealField);
214 
215  AlgebraicSymMatrix55 mo = ROOT::Math::Similarity(jl2c.jacobian(), m);
216 
217  int ndof = 2 * hits.size() - 5;
218  chi2 = chi2 * ndof;
219  GlobalPoint vv = gp.position();
220  math::XYZPoint pos(vv.x(), vv.y(), vv.z());
221  GlobalVector pp = gp.momentum();
222  math::XYZVector mom(pp.x(), pp.y(), pp.z());
223 
224  auto track = std::make_unique<reco::Track>(chi2, ndof, pos, mom, gp.charge(), CurvilinearTrajectoryError(mo));
225 
226  // bad and edup not supported as fit not present or not reliable
227  auto tkq = recoQuality[int(q)];
228  track->setQuality(tkq);
229  // loose,tight and HP are inclusive
230  if (reco::TrackBase::highPurity == tkq) {
231  track->setQuality(reco::TrackBase::tight);
232  track->setQuality(reco::TrackBase::loose);
233  } else if (reco::TrackBase::tight == tkq) {
234  track->setQuality(reco::TrackBase::loose);
235  }
236  track->setQuality(tkq);
237  // filter???
238  tracks.emplace_back(track.release(), hits);
239  }
240  // std::cout << "processed " << nt << " good tuples " << tracks.size() << "out of " << indToEdm.size() << std::endl;
241 
242  // store tracks
243  storeTracks(iEvent, tracks, httopo);
244  iEvent.put(std::move(indToEdmP));
245 }
246 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Quality qualityByName(std::string const &name)
void storeTracks(Ev &ev, const TWH &tracksWithHits, const TrackerTopology &ttopo)
Definition: storeTracks.h:19
T z() const
Definition: PV3DBase.h:61
TrackQuality
track quality
Definition: TrackBase.h:150
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int index() const
Definition: GeomDet.h:83
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
assert(be >=bs)
Definition: Plane.h:16
const edm::EDGetTokenT< SiPixelRecHitCollectionNew > cpuHits_
const GeomDet * det() const
const edm::EDGetTokenT< HMSstorage > hmsToken_
T x() const
Definition: PV3DBase.h:59
Eigen::Matrix< double, 5, 1 > Vector5d
Definition: FitResult.h:16
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
__host__ __device__ void transformToPerigeePlane(VI5 const &ip, MI5 const &icov, VO5 &op, MO5 &ocov)
Definition: FitUtils.h:218
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttTopoToken_
LocalVector momentum() const
Momentum vector in the local frame.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double f[11][100]
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::vector< uint16_t > IndToEdm
int nt
Definition: AMPTWrapper.h:42
TrackCharge charge() const
Charge (-1, 0 or 1)
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
auto const & tracks
cannot be loose
void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const edm::EDGetTokenT< reco::BeamSpot > tBeamSpot_
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
double b
Definition: hdecay.h:118
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
PixelTrackProducerFromSoA(const edm::ParameterSet &iConfig)
~PixelTrackProducerFromSoA() override=default
Eigen::Matrix< double, 5, 5 > Matrix5d
Definition: FitResult.h:20
fixed size matrix
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenTrack_
HLT enums.
std::vector< TrackWithRecHits > TracksWithRecHits
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > idealMagneticFieldToken_
pixelTrack::Quality const minQuality_
string quality
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
LocalPoint position() const
Local x and y position coordinates.
def move(src, dest)
Definition: eostools.py:511