CMS 3D CMS Logo

DeDxHitInfoProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: DeDxHitInfoProducer
4 // Class: DeDxHitInfoProducer
5 //
13 //
14 // Original Author: loic Quertenmont (querten)
15 // Created: Mon Nov 21 14:09:02 CEST 2014
16 // Modifications: Tamas Almos Vami (2022)
17 
18 // system include files
19 #include <memory>
20 
42 
43 using namespace reco;
44 using namespace std;
45 using namespace edm;
46 
48 public:
49  explicit DeDxHitInfoProducer(const edm::ParameterSet&);
50  ~DeDxHitInfoProducer() override;
51 
52 private:
53  void beginRun(edm::Run const& run, const edm::EventSetup&) override;
54  void produce(edm::Event&, const edm::EventSetup&) override;
55 
56  void makeCalibrationMap(const TrackerGeometry& tkGeom_);
57  void processRec(reco::DeDxHitInfo&, const SiStripRecHit2D&, const LocalPoint&, const LocalVector&, const float&);
58  void processHit(const TrackingRecHit* recHit,
59  const float trackMomentum,
60  const LocalVector& trackDirection,
61  reco::DeDxHitInfo& hitDeDxInfo,
62  const LocalPoint& hitLocalPos);
63 
64  // ----------member data ---------------------------
65  const bool usePixel_;
66  const bool useStrip_;
67  const float theMeVperADCPixel_;
68  const float theMeVperADCStrip_;
69 
70  const unsigned int minTrackHits_;
71  const float minTrackPt_;
72  const float minTrackPtPrescale_;
73  const float maxTrackEta_;
74 
76  const bool useCalibration_;
77  const bool doShapeTest_;
78  const bool usePixelShape_;
79 
80  const unsigned int lowPtTracksPrescalePass_, lowPtTracksPrescaleFail_;
84 
92 
93  std::vector<std::vector<float>> calibGains_;
94  unsigned int offsetDU_;
95 
97  uint64_t x = state[0];
98  uint64_t const y = state[1];
99  state[0] = y;
100  x ^= x << 23; // a
101  state[1] = x ^ y ^ (x >> 17) ^ (y >> 26); // b, c
102  return state[1] + y;
103  }
104 };
105 
107  : usePixel_(iConfig.getParameter<bool>("usePixel")),
108  useStrip_(iConfig.getParameter<bool>("useStrip")),
109  theMeVperADCPixel_(iConfig.getParameter<double>("MeVperADCPixel")),
110  theMeVperADCStrip_(iConfig.getParameter<double>("MeVperADCStrip")),
111  minTrackHits_(iConfig.getParameter<unsigned>("minTrackHits")),
112  minTrackPt_(iConfig.getParameter<double>("minTrackPt")),
113  minTrackPtPrescale_(iConfig.getParameter<double>("minTrackPtPrescale")),
114  maxTrackEta_(iConfig.getParameter<double>("maxTrackEta")),
115  calibrationPath_(iConfig.getParameter<string>("calibrationPath")),
116  useCalibration_(iConfig.getParameter<bool>("useCalibration")),
117  doShapeTest_(iConfig.getParameter<bool>("shapeTest")),
118  usePixelShape_(not iConfig.getParameter<edm::InputTag>("clusterShapeCache").label().empty()),
119  lowPtTracksPrescalePass_(iConfig.getParameter<uint32_t>("lowPtTracksPrescalePass")),
120  lowPtTracksPrescaleFail_(iConfig.getParameter<uint32_t>("lowPtTracksPrescaleFail")),
121  lowPtTracksEstimator_(iConfig.getParameter<edm::ParameterSet>("lowPtTracksEstimatorParameters")),
122  lowPtTracksDeDxThreshold_(iConfig.getParameter<double>("lowPtTracksDeDxThreshold")),
123  usePixelForPrescales_(iConfig.getParameter<bool>("usePixelForPrescales")),
124  tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"))),
125  pixShapeCacheToken_(consumes<SiPixelClusterShapeCache>(iConfig.getParameter<edm::InputTag>("clusterShapeCache"))),
127  clShapeToken_(
128  esConsumes<ClusterShapeHitFilter, CkfComponentsRecord>(edm::ESInputTag("", "ClusterShapeHitFilter"))) {
129  produces<reco::DeDxHitInfoCollection>();
130  produces<reco::DeDxHitInfoAss>();
131  produces<edm::ValueMap<int>>("prescale");
132 
133  if (!usePixel_ && !useStrip_)
134  edm::LogError("DeDxHitsProducer") << "No Pixel Hits NOR Strip Hits will be saved. Running this module is useless";
135 }
136 
138 
139 // ------------ method called once each job just before starting event loop ------------
140 void DeDxHitInfoProducer::beginRun(edm::Run const& iRun, const edm::EventSetup& iSetup) {
141  tkGeom_ = iSetup.getHandle(tkGeomToken_);
142  if (useCalibration_ && calibGains_.empty()) {
143  offsetDU_ = tkGeom_->offsetDU(GeomDetEnumerators::PixelBarrel); //index start at the first pixel
144 
146  }
147 }
148 
150  edm::Handle<reco::TrackCollection> trackCollectionHandle;
151  iEvent.getByToken(tracksToken_, trackCollectionHandle);
152  const TrackCollection& trackCollection(*trackCollectionHandle.product());
153 
154  clShape_ = iSetup.getHandle(clShapeToken_);
155  if (usePixelShape_)
157 
158  // creates the output collection
159  auto resultdedxHitColl = std::make_unique<reco::DeDxHitInfoCollection>();
160 
161  std::vector<int> indices;
162  std::vector<int> prescales;
163  uint64_t state[2] = {iEvent.id().event(), iEvent.id().luminosityBlock()};
164  for (unsigned int j = 0; j < trackCollection.size(); j++) {
166 
167  //track selection
168  bool passPt = (track.pt() >= minTrackPt_), passLowDeDx = false, passHighDeDx = false, pass = passPt;
169  if (!pass && (track.pt() >= minTrackPtPrescale_)) {
170  if (lowPtTracksPrescalePass_ > 0) {
171  passHighDeDx = ((xorshift128p(state) % lowPtTracksPrescalePass_) == 0);
172  }
173  if (lowPtTracksPrescaleFail_ > 0) {
174  passLowDeDx = ((xorshift128p(state) % lowPtTracksPrescaleFail_) == 0);
175  }
176  pass = passHighDeDx || passLowDeDx;
177  }
178  if (!pass || std::abs(track.eta()) > maxTrackEta_ || track.numberOfValidHits() < minTrackHits_) {
179  indices.push_back(-1);
180  continue;
181  }
182 
183  reco::DeDxHitInfo hitDeDxInfo;
184  auto const& trajParams = track.extra()->trajParams();
185  auto hb = track.recHitsBegin();
186  for (unsigned int h = 0; h < track.recHitsSize(); h++) {
187  auto recHit = *(hb + h);
189  continue;
190 
191  processHit(recHit, track.p(), trajParams[h].direction(), hitDeDxInfo, trajParams[h].position());
192  }
193 
194  if (!passPt) {
195  std::vector<DeDxHit> hits;
196  hits.reserve(hitDeDxInfo.size());
197  for (unsigned int i = 0, n = hitDeDxInfo.size(); i < n; ++i) {
198  if (hitDeDxInfo.detId(i).subdetId() <= 2 && usePixelForPrescales_) {
199  hits.push_back(DeDxHit(hitDeDxInfo.charge(i) / hitDeDxInfo.pathlength(i) * theMeVperADCPixel_, 0, 0, 0));
200  } else if (hitDeDxInfo.detId(i).subdetId() > 2) {
201  if (doShapeTest_ && !deDxTools::shapeSelection(*hitDeDxInfo.stripCluster(i)))
202  continue;
203  hits.push_back(DeDxHit(hitDeDxInfo.charge(i) / hitDeDxInfo.pathlength(i) * theMeVperADCStrip_, 0, 0, 0));
204  }
205  }
206 
207  // In case we have a pixel only track, but usePixelForPrescales_ is false
208  if (hits.empty()) {
209  indices.push_back(-1);
210  continue;
211  }
212  std::sort(hits.begin(), hits.end(), std::less<DeDxHit>());
214  if (passLowDeDx) {
216  } else {
217  indices.push_back(-1);
218  continue;
219  }
220  } else {
221  if (passHighDeDx) {
223  } else {
224  indices.push_back(-1);
225  continue;
226  }
227  }
228  } else {
229  prescales.push_back(1);
230  }
231  indices.push_back(resultdedxHitColl->size());
232  resultdedxHitColl->push_back(hitDeDxInfo);
233  }
235 
236  edm::OrphanHandle<reco::DeDxHitInfoCollection> dedxHitCollHandle = iEvent.put(std::move(resultdedxHitColl));
237 
238  //create map passing the handle to the matched collection
239  auto dedxMatch = std::make_unique<reco::DeDxHitInfoAss>(dedxHitCollHandle);
241  filler.insert(trackCollectionHandle, indices.begin(), indices.end());
242  filler.fill();
243  iEvent.put(std::move(dedxMatch));
244 
245  auto dedxPrescale = std::make_unique<edm::ValueMap<int>>();
246  edm::ValueMap<int>::Filler pfiller(*dedxPrescale);
247  pfiller.insert(dedxHitCollHandle, prescales.begin(), prescales.end());
248  pfiller.fill();
249  iEvent.put(std::move(dedxPrescale), "prescale");
250 }
251 
253  const SiStripRecHit2D& recHit,
254  const LocalPoint& lpos,
255  const LocalVector& ldir,
256  const float& cos) {
257  uint8_t type(0);
258  int meas;
259  float pred;
260  const auto& usable = clShape_->getSizes(recHit, {}, ldir, meas, pred);
261  if (usable && meas <= int(std::abs(pred)) + 4)
263  if (clShape_->isCompatible(recHit, ldir))
265  if (usable)
267 
268  int NSaturating(0);
269  const auto& detId = recHit.geographicalId();
270  const auto* detUnit = recHit.detUnit();
271  if (detUnit == nullptr)
272  detUnit = tkGeom_->idToDet(detId);
273  const auto pathLen = detUnit->surface().bounds().thickness() / cos;
274  float chargeAbs = deDxTools::getCharge(&(recHit.stripCluster()), NSaturating, *detUnit, calibGains_, offsetDU_);
275  hitDeDxInfo.addHit(chargeAbs, pathLen, detId, lpos, type, recHit.stripCluster());
276 }
277 
279  const float trackMomentum,
280  const LocalVector& trackDirection,
281  reco::DeDxHitInfo& hitDeDxInfo,
282  const LocalPoint& hitLocalPos) {
283  auto const& thit = static_cast<BaseTrackerRecHit const&>(*recHit);
284  if (!thit.isValid())
285  return;
286 
287  //make sure cosine is not 0
288  float cosine = trackDirection.z() / trackDirection.mag();
289  float cosineAbs = std::max(0.00000001f, std::abs(cosine));
290 
291  auto const& clus = thit.firstClusterRef();
292  if (!clus.isValid())
293  return;
294 
295  const auto* detUnit = recHit->detUnit();
296  if (detUnit == nullptr) {
297  detUnit = tkGeom_->idToDet(thit.geographicalId());
298  }
299  float pathLen = detUnit->surface().bounds().thickness() / cosineAbs;
300 
301  if (clus.isPixel()) {
302  if (!usePixel_)
303  return;
304 
305  uint8_t type(0);
306  const auto& pixelDet = *dynamic_cast<const PixelGeomDetUnit*>(detUnit);
307  const auto& pixelRecHit = *dynamic_cast<const SiPixelRecHit*>(recHit);
309  ClusterShape().determineShape(pixelDet, clus.pixelCluster(), data);
310  if (data.isComplete)
312  if (usePixelShape_ && clShape_->isCompatible(pixelRecHit, trackDirection, *pixShapeCache_))
314  if (data.isComplete && data.isStraight && data.hasBigPixelsOnlyInside)
316 
317  float chargeAbs = clus.pixelCluster().charge();
318  hitDeDxInfo.addHit(chargeAbs, pathLen, thit.geographicalId(), hitLocalPos, type, clus.pixelCluster());
319  } else if (clus.isStrip() && !thit.isMatched()) {
320  if (!useStrip_)
321  return;
322 
323  processRec(hitDeDxInfo, {thit.geographicalId(), clus}, hitLocalPos, trackDirection, cosineAbs);
324  } else if (clus.isStrip() && thit.isMatched()) {
325  if (!useStrip_)
326  return;
327  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(recHit);
328  if (!matchedHit)
329  return;
330 
331  processRec(hitDeDxInfo, matchedHit->monoHit(), hitLocalPos, trackDirection, cosineAbs);
332  processRec(hitDeDxInfo, matchedHit->stereoHit(), hitLocalPos, trackDirection, cosineAbs);
333  }
334 }
335 
336 //define this as a plug-in
float pathlength(size_t i) const
Definition: DeDxHitInfo.h:47
const unsigned int minTrackHits_
const std::string calibrationPath_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::EDGetTokenT< SiPixelClusterShapeCache > pixShapeCacheToken_
SiStripRecHit2D stereoHit() const
bool isFromDet(TrackingRecHit const &hit)
T z() const
Definition: PV3DBase.h:61
void processRec(reco::DeDxHitInfo &, const SiStripRecHit2D &, const LocalPoint &, const LocalVector &, const float &)
unsigned int offsetDU(SubDetector sid) const
T const * product() const
Definition: Handle.h:70
const unsigned int lowPtTracksPrescaleFail_
void processHit(const TrackingRecHit *recHit, const float trackMomentum, const LocalVector &trackDirection, reco::DeDxHitInfo &hitDeDxInfo, const LocalPoint &hitLocalPos)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::ESHandle< ClusterShapeHitFilter > clShape_
static constexpr int Compatible
Definition: DeDxHitInfo.h:40
void produce(edm::Event &, const edm::EventSetup &) override
int getCharge(const SiStripCluster *cluster, int &nSatStrip, const GeomDetUnit &detUnit, const std::vector< std::vector< float >> &calibGains, const unsigned int &m_off)
Log< level::Error, false > LogError
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=nullptr) const
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
virtual float thickness() const =0
char const * label
int iEvent
Definition: GenABIO.cc:224
const unsigned int lowPtTracksPrescalePass_
DeDxHitInfoProducer(const edm::ParameterSet &)
void addHit(const float charge, const float pathlength, const DetId &detId, const LocalPoint &pos, const uint8_t &type, const SiStripCluster &stripCluster)
Definition: DeDxHitInfo.h:94
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T mag() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Transition
Definition: Transition.h:12
const SiStripCluster * stripCluster(size_t i) const
Definition: DeDxHitInfo.h:71
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
edm::ESHandle< TrackerGeometry > tkGeom_
static constexpr int Calibration
Definition: DeDxHitInfo.h:40
DetId detId(size_t i) const
Definition: DeDxHitInfo.h:48
trackCollection
Definition: JetHT_cfg.py:51
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const TrackerGeomDet * idToDet(DetId) const override
void makeCalibrationMap(const std::string &m_calibrationPath, const TrackerGeometry &tkGeom, std::vector< std::vector< float >> &calibGains, const unsigned int &m_off)
uint64_t xorshift128p(uint64_t state[2])
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
bool getSizes(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, int &part, ClusterData::ArrayType &meas, std::pair< float, float > &predr, PixelData const *pd=nullptr) const
unsigned long long uint64_t
Definition: Time.h:13
bool shapeSelection(const SiStripCluster &ampls)
Definition: DeDxTools.cc:13
size_t size() const
Definition: DeDxHitInfo.h:45
const edm::ESGetToken< ClusterShapeHitFilter, CkfComponentsRecord > clShapeToken_
~DeDxHitInfoProducer() override
static constexpr int Complete
Definition: DeDxHitInfo.h:40
fixed size matrix
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
GenericTruncatedAverageDeDxEstimator lowPtTracksEstimator_
edm::Handle< SiPixelClusterShapeCache > pixShapeCache_
float x
SiStripRecHit2D monoHit() const
std::pair< float, float > dedx(const reco::DeDxHitCollection &Hits) override
float charge(size_t i) const
Definition: DeDxHitInfo.h:46
std::vector< std::vector< float > > calibGains_
void beginRun(edm::Run const &run, const edm::EventSetup &) override
const float lowPtTracksDeDxThreshold_
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
Our base class.
Definition: SiPixelRecHit.h:23
const Bounds & bounds() const
Definition: Surface.h:87