CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 //
17 
19 
20 // system include files
21 
22 using namespace reco;
23 using namespace std;
24 using namespace edm;
25 
27  : usePixel(iConfig.getParameter<bool>("usePixel")),
28  useStrip(iConfig.getParameter<bool>("useStrip")),
29  MeVperADCPixel(iConfig.getParameter<double>("MeVperADCPixel")),
30  MeVperADCStrip(iConfig.getParameter<double>("MeVperADCStrip")),
31  minTrackHits(iConfig.getParameter<unsigned>("minTrackHits")),
32  minTrackPt(iConfig.getParameter<double>("minTrackPt")),
33  minTrackPtPrescale(iConfig.getParameter<double>("minTrackPtPrescale")),
34  maxTrackEta(iConfig.getParameter<double>("maxTrackEta")),
35  m_calibrationPath(iConfig.getParameter<string>("calibrationPath")),
36  useCalibration(iConfig.getParameter<bool>("useCalibration")),
37  shapetest(iConfig.getParameter<bool>("shapeTest")),
38  lowPtTracksPrescalePass(iConfig.getParameter<uint32_t>("lowPtTracksPrescalePass")),
39  lowPtTracksPrescaleFail(iConfig.getParameter<uint32_t>("lowPtTracksPrescaleFail")),
40  lowPtTracksEstimator(iConfig.getParameter<edm::ParameterSet>("lowPtTracksEstimatorParameters")),
41  lowPtTracksDeDxThreshold(iConfig.getParameter<double>("lowPtTracksDeDxThreshold")) {
42  produces<reco::DeDxHitInfoCollection>();
43  produces<reco::DeDxHitInfoAss>();
44  produces<edm::ValueMap<int>>("prescale");
45 
46  m_tracksTag = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"));
47  m_tkGeomToken = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>();
48 
49  if (!usePixel && !useStrip)
50  edm::LogError("DeDxHitsProducer") << "No Pixel Hits NOR Strip Hits will be saved. Running this module is useless";
51 }
52 
54 
55 // ------------ method called once each job just before starting event loop ------------
57  tkGeom = iSetup.getHandle(m_tkGeomToken);
58  if (useCalibration && calibGains.empty()) {
59  m_off = tkGeom->offsetDU(GeomDetEnumerators::PixelBarrel); //index start at the first pixel
60 
62  }
63 }
64 
66  edm::Handle<reco::TrackCollection> trackCollectionHandle;
67  iEvent.getByToken(m_tracksTag, trackCollectionHandle);
68  const TrackCollection& trackCollection(*trackCollectionHandle.product());
69 
70  // creates the output collection
71  auto resultdedxHitColl = std::make_unique<reco::DeDxHitInfoCollection>();
72 
73  std::vector<int> indices;
74  std::vector<int> prescales;
75  uint64_t state[2] = {iEvent.id().event(), iEvent.id().luminosityBlock()};
76  for (unsigned int j = 0; j < trackCollection.size(); j++) {
78 
79  //track selection
80  bool passPt = (track.pt() >= minTrackPt), passLowDeDx = false, passHighDeDx = false, pass = passPt;
81  if (!pass && (track.pt() >= minTrackPtPrescale)) {
82  if (lowPtTracksPrescalePass > 0) {
83  passHighDeDx = ((xorshift128p(state) % lowPtTracksPrescalePass) == 0);
84  }
85  if (lowPtTracksPrescaleFail > 0) {
86  passLowDeDx = ((xorshift128p(state) % lowPtTracksPrescaleFail) == 0);
87  }
88  pass = passHighDeDx || passLowDeDx;
89  }
90  if (!pass || std::abs(track.eta()) > maxTrackEta || track.numberOfValidHits() < minTrackHits) {
91  indices.push_back(-1);
92  continue;
93  }
94 
95  reco::DeDxHitInfo hitDeDxInfo;
96  auto const& trajParams = track.extra()->trajParams();
97  auto hb = track.recHitsBegin();
98  for (unsigned int h = 0; h < track.recHitsSize(); h++) {
99  auto recHit = *(hb + h);
100  if (!trackerHitRTTI::isFromDet(*recHit))
101  continue;
102 
103  auto trackDirection = trajParams[h].direction();
104  float cosine = trackDirection.z() / trackDirection.mag();
105 
106  processHit(recHit, track.p(), cosine, hitDeDxInfo, trajParams[h].position());
107  }
108 
109  if (!passPt) {
110  std::vector<DeDxHit> hits;
111  hits.reserve(hitDeDxInfo.size());
112  for (unsigned int i = 0, n = hitDeDxInfo.size(); i < n; ++i) {
113  if (hitDeDxInfo.detId(i).subdetId() <= 2) {
114  hits.push_back(DeDxHit(hitDeDxInfo.charge(i) / hitDeDxInfo.pathlength(i) * MeVperADCPixel, 0, 0, 0));
115  } else {
116  if (shapetest && !DeDxTools::shapeSelection(*hitDeDxInfo.stripCluster(i)))
117  continue;
118  hits.push_back(DeDxHit(hitDeDxInfo.charge(i) / hitDeDxInfo.pathlength(i) * MeVperADCStrip, 0, 0, 0));
119  }
120  }
121  std::sort(hits.begin(), hits.end(), std::less<DeDxHit>());
123  if (passLowDeDx) {
124  prescales.push_back(lowPtTracksPrescaleFail);
125  } else {
126  indices.push_back(-1);
127  continue;
128  }
129  } else {
130  if (passHighDeDx) {
131  prescales.push_back(lowPtTracksPrescalePass);
132  } else {
133  indices.push_back(-1);
134  continue;
135  }
136  }
137  } else {
138  prescales.push_back(1);
139  }
140  indices.push_back(resultdedxHitColl->size());
141  resultdedxHitColl->push_back(hitDeDxInfo);
142  }
144 
145  edm::OrphanHandle<reco::DeDxHitInfoCollection> dedxHitCollHandle = iEvent.put(std::move(resultdedxHitColl));
146 
147  //create map passing the handle to the matched collection
148  auto dedxMatch = std::make_unique<reco::DeDxHitInfoAss>(dedxHitCollHandle);
149  reco::DeDxHitInfoAss::Filler filler(*dedxMatch);
150  filler.insert(trackCollectionHandle, indices.begin(), indices.end());
151  filler.fill();
152  iEvent.put(std::move(dedxMatch));
153 
154  auto dedxPrescale = std::make_unique<edm::ValueMap<int>>();
155  edm::ValueMap<int>::Filler pfiller(*dedxPrescale);
156  pfiller.insert(dedxHitCollHandle, prescales.begin(), prescales.end());
157  pfiller.fill();
158  iEvent.put(std::move(dedxPrescale), "prescale");
159 }
160 
162  const float trackMomentum,
163  const float cosine,
164  reco::DeDxHitInfo& hitDeDxInfo,
165  const LocalPoint& hitLocalPos) {
166  auto const& thit = static_cast<BaseTrackerRecHit const&>(*recHit);
167  if (!thit.isValid())
168  return;
169 
170  float cosineAbs = std::max(0.00000001f, std::abs(cosine)); //make sure cosine is not 0
171 
172  auto const& clus = thit.firstClusterRef();
173  if (!clus.isValid())
174  return;
175 
176  if (clus.isPixel()) {
177  if (!usePixel)
178  return;
179 
180  const auto* detUnit = recHit->detUnit();
181  if (detUnit == nullptr)
182  detUnit = tkGeom->idToDet(thit.geographicalId());
183  float pathLen = detUnit->surface().bounds().thickness() / cosineAbs;
184  float chargeAbs = clus.pixelCluster().charge();
185  hitDeDxInfo.addHit(chargeAbs, pathLen, thit.geographicalId(), hitLocalPos, clus.pixelCluster());
186  } else if (clus.isStrip() && !thit.isMatched()) {
187  if (!useStrip)
188  return;
189 
190  const auto* detUnit = recHit->detUnit();
191  if (detUnit == nullptr)
192  detUnit = tkGeom->idToDet(thit.geographicalId());
193  int NSaturating = 0;
194  float pathLen = detUnit->surface().bounds().thickness() / cosineAbs;
195  float chargeAbs = DeDxTools::getCharge(&(clus.stripCluster()), NSaturating, *detUnit, calibGains, m_off);
196  hitDeDxInfo.addHit(chargeAbs, pathLen, thit.geographicalId(), hitLocalPos, clus.stripCluster());
197  } else if (clus.isStrip() && thit.isMatched()) {
198  if (!useStrip)
199  return;
200  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(recHit);
201  if (!matchedHit)
202  return;
203 
204  const auto* detUnitM = matchedHit->monoHit().detUnit();
205  if (detUnitM == nullptr)
206  detUnitM = tkGeom->idToDet(matchedHit->monoHit().geographicalId());
207  int NSaturating = 0;
208  float pathLen = detUnitM->surface().bounds().thickness() / cosineAbs;
209  float chargeAbs =
210  DeDxTools::getCharge(&(matchedHit->monoHit().stripCluster()), NSaturating, *detUnitM, calibGains, m_off);
211  hitDeDxInfo.addHit(chargeAbs, pathLen, thit.geographicalId(), hitLocalPos, matchedHit->monoHit().stripCluster());
212 
213  const auto* detUnitS = matchedHit->stereoHit().detUnit();
214  if (detUnitS == nullptr)
215  detUnitS = tkGeom->idToDet(matchedHit->stereoHit().geographicalId());
216  NSaturating = 0;
217  pathLen = detUnitS->surface().bounds().thickness() / cosineAbs;
218  chargeAbs =
219  DeDxTools::getCharge(&(matchedHit->stereoHit().stripCluster()), NSaturating, *detUnitS, calibGains, m_off);
220  hitDeDxInfo.addHit(chargeAbs, pathLen, thit.geographicalId(), hitLocalPos, matchedHit->stereoHit().stripCluster());
221  }
222 }
223 
224 //define this as a plug-in
double p() const
momentum vector magnitude
Definition: TrackBase.h:631
EventNumber_t event() const
Definition: EventID.h:40
edm::EDGetTokenT< reco::TrackCollection > m_tracksTag
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const GeomDetUnit * detUnit() const override
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:139
bool isFromDet(TrackingRecHit const &hit)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
float charge(size_t i) const
Definition: DeDxHitInfo.h:42
void makeCalibrationMap(const std::string &m_calibrationPath, const TrackerGeometry &tkGeom, std::vector< std::vector< float >> &calibGains, const unsigned int &m_off)
GenericTruncatedAverageDeDxEstimator lowPtTracksEstimator
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< std::vector< float > > calibGains
void produce(edm::Event &, const edm::EventSetup &) override
Log< level::Error, false > LogError
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:39
const unsigned int lowPtTracksPrescalePass
tuple MeVperADCPixel
int iEvent
Definition: GenABIO.cc:224
const float lowPtTracksDeDxThreshold
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
DeDxHitInfoProducer(const edm::ParameterSet &)
tuple minTrackPt
Definition: align_cfg.py:22
bool shapeSelection(const SiStripCluster &ampls)
Definition: DeDxTools.cc:12
void processHit(const TrackingRecHit *recHit, const float trackMomentum, const float cosine, reco::DeDxHitInfo &hitDeDxInfo, const LocalPoint &hitLocalPos)
const std::string m_calibrationPath
double pt() const
track transverse momentum
Definition: TrackBase.h:637
def move
Definition: eostools.py:511
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void addHit(const float charge, const float pathlength, const DetId &detId, const LocalPoint &pos, const SiStripCluster &stripCluster)
Definition: DeDxHitInfo.h:89
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
const unsigned int minTrackHits
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:88
const unsigned int lowPtTracksPrescaleFail
int getCharge(const SiStripCluster *cluster, int &nSatStrip, const GeomDetUnit &detUnit, const std::vector< std::vector< float >> &calibGains, const unsigned int &m_off)
uint64_t xorshift128p(uint64_t state[2])
tuple trackCollection
SiStripRecHit2D stereoHit() const
T const * product() const
Definition: Handle.h:70
unsigned long long uint64_t
Definition: Time.h:13
SiStripCluster const & stripCluster() const
size_t size() const
Definition: DeDxHitInfo.h:41
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
SiStripRecHit2D monoHit() const
tuple MeVperADCStrip
edm::EventID id() const
Definition: EventBase.h:59
virtual const GeomDetUnit * detUnit() const
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_tkGeomToken
float pathlength(size_t i) const
Definition: DeDxHitInfo.h:43
std::pair< float, float > dedx(const reco::DeDxHitCollection &Hits) override
DetId geographicalId() const
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
void beginRun(edm::Run const &run, const edm::EventSetup &) override
list indices
Definition: dqmdumpme.py:50
edm::ESHandle< TrackerGeometry > tkGeom
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
const SiStripCluster * stripCluster(size_t i) const
Definition: DeDxHitInfo.h:66
DetId detId(size_t i) const
Definition: DeDxHitInfo.h:44
const float minTrackPtPrescale
Definition: Run.h:45