CMS 3D CMS Logo

PixelTrackProducerFromSoA.cc
Go to the documentation of this file.
23 
28 
32 
33 #include "storeTracks.h"
35 
39 
44 template <typename TrackerTraits>
48 
49 public:
50  using IndToEdm = std::vector<uint32_t>;
51 
52  explicit PixelTrackProducerFromSoAT(const edm::ParameterSet &iConfig);
53  ~PixelTrackProducerFromSoAT() override = default;
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
56 
58 
59 private:
60  void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override;
61 
62  // Event Data tokens
67  // Event Setup tokens
70 
71  int32_t const minNumberOfHits_;
73 };
74 
75 template <typename TrackerTraits>
77  : tBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
78  tokenTrack_(consumes(iConfig.getParameter<edm::InputTag>("trackSrc"))),
79  cpuHits_(consumes<SiPixelRecHitCollectionNew>(iConfig.getParameter<edm::InputTag>("pixelRecHitLegacySrc"))),
80  hmsToken_(consumes<HMSstorage>(iConfig.getParameter<edm::InputTag>("pixelRecHitLegacySrc"))),
81  idealMagneticFieldToken_(esConsumes()),
82  ttTopoToken_(esConsumes()),
83  minNumberOfHits_(iConfig.getParameter<int>("minNumberOfHits")),
84  minQuality_(pixelTrack::qualityByName(iConfig.getParameter<std::string>("minQuality"))) {
86  throw cms::Exception("PixelTrackConfiguration")
87  << iConfig.getParameter<std::string>("minQuality") + " is not a pixelTrack::Quality";
88  }
90  throw cms::Exception("PixelTrackConfiguration")
91  << iConfig.getParameter<std::string>("minQuality") + " not supported";
92  }
93  produces<TrackingRecHitCollection>();
94  produces<reco::TrackExtraCollection>();
95  // TrackCollection refers to TrackingRechit and TrackExtra
96  // collections, need to declare its production after them to work
97  // around a rare race condition in framework scheduling
98  produces<reco::TrackCollection>();
99  produces<IndToEdm>();
100 }
101 
102 template <typename TrackerTraits>
105  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
106  desc.add<edm::InputTag>("trackSrc", edm::InputTag("pixelTracksSoA"));
107  desc.add<edm::InputTag>("pixelRecHitLegacySrc", edm::InputTag("siPixelRecHitsPreSplittingLegacy"));
108  desc.add<int>("minNumberOfHits", 0);
109  desc.add<std::string>("minQuality", "loose");
110  descriptions.addWithDefaultLabel(desc);
111 }
112 
113 template <typename TrackerTraits>
116  const edm::EventSetup &iSetup) const {
117  // enum class Quality : uint8_t { bad = 0, edup, dup, loose, strict, tight, highPurity };
126 
127  // std::cout << "Converting gpu helix in reco tracks" << std::endl;
128 
129  auto indToEdmP = std::make_unique<IndToEdm>();
130  auto &indToEdm = *indToEdmP;
131 
132  auto const &idealField = iSetup.getData(idealMagneticFieldToken_);
133 
135 
136  auto const &httopo = iSetup.getData(ttTopoToken_);
137 
138  const auto &bsh = iEvent.get(tBeamSpot_);
139  GlobalPoint bs(bsh.x0(), bsh.y0(), bsh.z0());
140 
141  auto const &rechits = iEvent.get(cpuHits_);
142  std::vector<TrackingRecHit const *> hitmap;
143  auto const &rcs = rechits.data();
144  auto nhits = rcs.size();
145 
146  hitmap.resize(nhits, nullptr);
147 
148  auto const *hitsModuleStart = iEvent.get(hmsToken_).get();
149  auto fc = hitsModuleStart;
150 
151  for (auto const &h : rcs) {
152  auto const &thit = static_cast<BaseTrackerRecHit const &>(h);
153  auto detI = thit.det()->index();
154  auto const &clus = thit.firstClusterRef();
155  assert(clus.isPixel());
156  auto i = fc[detI] + clus.pixelCluster().originalId();
157  if (i >= hitmap.size())
158  hitmap.resize(i + 256, nullptr); // only in case of hit overflow in one module
159 
160  assert(nullptr == hitmap[i]);
161  hitmap[i] = &h;
162  }
163 
164  std::vector<const TrackingRecHit *> hits;
165  hits.reserve(5);
166 
167  auto const &tsoa = iEvent.get(tokenTrack_);
168  auto const *quality = tsoa.view().quality();
169  auto const &hitIndices = tsoa.view().hitIndices();
170  auto nTracks = tsoa.view().nTracks();
171 
172  tracks.reserve(nTracks);
173 
174  int32_t nt = 0;
175 
176  //sort index by pt
177  std::vector<int32_t> sortIdxs(nTracks);
178  std::iota(sortIdxs.begin(), sortIdxs.end(), 0);
179  std::sort(sortIdxs.begin(), sortIdxs.end(), [&](int32_t const i1, int32_t const i2) {
180  return tsoa.view()[i1].pt() > tsoa.view()[i2].pt();
181  });
182 
183  //store the index of the SoA: indToEdm[index_SoAtrack] -> index_edmTrack (if it exists)
184  indToEdm.resize(sortIdxs.size(), -1);
185  for (const auto &it : sortIdxs) {
186  auto nHits = tracksHelpers::nHits(tsoa.view(), it);
187  assert(nHits >= 3);
188  auto q = quality[it];
189 
190  if (q < minQuality_)
191  continue;
192  if (nHits < minNumberOfHits_) //move to nLayers?
193  continue;
194  indToEdm[it] = nt;
195  ++nt;
196 
197  hits.resize(nHits);
198  auto b = hitIndices.begin(it);
199  for (int iHit = 0; iHit < nHits; ++iHit)
200  hits[iHit] = hitmap[*(b + iHit)];
201 
202  // mind: this values are respect the beamspot!
203 
204  float chi2 = tsoa.view()[it].chi2();
205  float phi = tracksHelpers::phi(tsoa.view(), it);
206 
207  riemannFit::Vector5d ipar, opar;
208  riemannFit::Matrix5d icov, ocov;
209  tracksHelpers::template copyToDense<riemannFit::Vector5d, riemannFit::Matrix5d>(tsoa.view(), ipar, icov, it);
210  riemannFit::transformToPerigeePlane(ipar, icov, opar, ocov);
211 
212  LocalTrajectoryParameters lpar(opar(0), opar(1), opar(2), opar(3), opar(4), 1.);
214  for (int i = 0; i < 5; ++i)
215  for (int j = i; j < 5; ++j)
216  m(i, j) = ocov(i, j);
217 
218  float sp = std::sin(phi);
219  float cp = std::cos(phi);
220  Surface::RotationType rot(sp, -cp, 0, 0, 0, -1.f, cp, sp, 0);
221 
222  Plane impPointPlane(bs, rot);
224  impPointPlane.toGlobal(lpar.position()), impPointPlane.toGlobal(lpar.momentum()), lpar.charge(), &idealField);
225  JacobianLocalToCurvilinear jl2c(impPointPlane, lpar, idealField);
226 
227  AlgebraicSymMatrix55 mo = ROOT::Math::Similarity(jl2c.jacobian(), m);
228 
229  int ndof = 2 * hits.size() - 5;
230  chi2 = chi2 * ndof;
231  GlobalPoint vv = gp.position();
232  math::XYZPoint pos(vv.x(), vv.y(), vv.z());
233  GlobalVector pp = gp.momentum();
234  math::XYZVector mom(pp.x(), pp.y(), pp.z());
235 
236  auto track = std::make_unique<reco::Track>(chi2, ndof, pos, mom, gp.charge(), CurvilinearTrajectoryError(mo));
237 
238  // bad and edup not supported as fit not present or not reliable
239  auto tkq = recoQuality[int(q)];
240  track->setQuality(tkq);
241  // loose,tight and HP are inclusive
242  if (reco::TrackBase::highPurity == tkq) {
243  track->setQuality(reco::TrackBase::tight);
244  track->setQuality(reco::TrackBase::loose);
245  } else if (reco::TrackBase::tight == tkq) {
246  track->setQuality(reco::TrackBase::loose);
247  }
248  track->setQuality(tkq);
249  // filter???
250  tracks.emplace_back(track.release(), hits);
251  }
252  // std::cout << "processed " << nt << " good tuples " << tracks.size() << "out of " << indToEdm.size() << std::endl;
253 
254  // store tracks
255  storeTracks(iEvent, tracks, httopo);
256  iEvent.put(std::move(indToEdmP));
257 }
258 
261 
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
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::EDGetTokenT< reco::BeamSpot > tBeamSpot_
const edm::EDGetTokenT< TrackSoAHost > tokenTrack_
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
const edm::EDGetTokenT< SiPixelRecHitCollectionNew > cpuHits_
TrackQuality
track quality
Definition: TrackBase.h:150
int index() const
Definition: GeomDet.h:83
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
~PixelTrackProducerFromSoAT() override=default
assert(be >=bs)
Definition: Plane.h:16
const GeomDet * det() const
string quality
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< MagneticField, IdealMagneticFieldRecord > idealMagneticFieldToken_
LocalVector momentum() const
Momentum vector in the local frame.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double f[11][100]
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
const edm::EDGetTokenT< HMSstorage > hmsToken_
pixelTrack::Quality const minQuality_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
void produce(edm::StreamID streamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
PixelTrackProducerFromSoAT(const edm::ParameterSet &iConfig)
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
double b
Definition: hdecay.h:118
Eigen::Matrix< double, 5, 5 > Matrix5d
Definition: FitResult.h:20
fixed size matrix
HLT enums.
std::vector< TrackWithRecHits > TracksWithRecHits
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttTopoToken_
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)