CMS 3D CMS Logo

Phase2ITValidateTrackingRecHit.cc
Go to the documentation of this file.
1 // Package: Phase2ITValidateTrackingRecHit
2 // Class: Phase2ITValidateTrackingRecHit
3 //
7 //
8 // Author: Shubhi Parolia, Suvankar Roy Chowdhury
9 // Date: June 2020
10 //
11 // system include files
12 #include <memory>
13 #include <map>
14 #include <vector>
15 #include <algorithm>
39 //--- for SimHit association
44 //DQM
48 // base class
51 
53 public:
56  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
57  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
58 
59 private:
60  void fillITHistos(const edm::Event& iEvent,
61  const TrackerHitAssociator& associateRecHit,
63  const std::map<unsigned int, SimTrack>& selectedSimTrackMap);
66  const double simtrackminpt_;
69  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> simHitTokens_;
70 };
71 
73  : Phase2ITValidateRecHitBase(iConfig),
74  config_(iConfig),
75  trackerHitAssociatorConfig_(iConfig, consumesCollector()),
76  simtrackminpt_(iConfig.getParameter<double>("SimTrackMinPt")),
77  tokenTracks_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
78  simTracksToken_(consumes<edm::SimTrackContainer>(iConfig.getParameter<edm::InputTag>("simTracksSrc"))) {
79  edm::LogInfo("Phase2ITValidateRecHit") << ">>> Construct Phase2ITValidateRecHit ";
80  for (const auto& itName : config_.getParameter<std::vector<std::string>>("ROUList")) {
81  simHitTokens_.push_back(consumes<std::vector<PSimHit>>(edm::InputTag("g4SimHits", itName)));
82  }
83 }
84 //
86  edm::LogInfo("Phase2ITValidateTrackingRecHit") << ">>> Destroy Phase2ITValidateTrackingRecHit ";
87 }
88 
90  std::vector<edm::Handle<edm::PSimHitContainer>> simHits;
91  for (const auto& itoken : simHitTokens_) {
92  const auto& simHitHandle = iEvent.getHandle(itoken);
93  if (!simHitHandle.isValid())
94  continue;
95  simHits.emplace_back(simHitHandle);
96  }
97  // Get the SimTracks and push them in a map of id, SimTrack
98  const auto& simTracks = iEvent.getHandle(simTracksToken_);
99 
100  std::map<unsigned int, SimTrack> selectedSimTrackMap;
101  for (const auto& simTrackIt : *simTracks) {
102  if (simTrackIt.momentum().pt() > simtrackminpt_) {
103  selectedSimTrackMap.emplace(simTrackIt.trackId(), simTrackIt);
104  }
105  }
107  fillITHistos(iEvent, associateRecHit, simHits, selectedSimTrackMap);
108 }
109 
111  const TrackerHitAssociator& associateRecHit,
113  const std::map<unsigned int, SimTrack>& selectedSimTrackMap) {
114  const auto& tracks = iEvent.getHandle(tokenTracks_);
115  if (!tracks.isValid())
116  return;
117 
118  std::map<std::string, unsigned int> nrechitLayerMap_primary;
119 
120  // loop over tracks
121  for (const auto& track : *tracks) {
122  // loop over hits
123  for (auto const& hit : track.recHits()) {
124  if (!hit->isValid())
125  continue;
126 
127  auto id = hit->geographicalId();
128  // check that we are in the pixel
129  auto subdetid = (id.subdetId());
130  if (!(subdetid == PixelSubdetector::PixelBarrel) && !(subdetid == PixelSubdetector::PixelEndcap))
131  continue;
132 
133  const GeomDetUnit* geomDetunit(tkGeom_->idToDetUnit(id));
134  if (!geomDetunit)
135  continue;
136  // determine the detector we are in
138  if (nrechitLayerMap_primary.find(key) == nrechitLayerMap_primary.end()) {
139  nrechitLayerMap_primary.emplace(key, 1);
140  } else {
141  nrechitLayerMap_primary[key] += 1;
142  }
143 
144  const SiPixelRecHit* rechit = dynamic_cast<const SiPixelRecHit*>(hit);
145  if (!rechit)
146  continue;
147 
148  const std::vector<SimHitIdpr>& matchedId = associateRecHit.associateHitId(*rechit);
149  const PSimHit* simhitClosest = nullptr;
150  float minx = 10000;
151  LocalPoint lp = rechit->localPosition();
152  for (const auto& simHitCol : simHits) {
153  for (const auto& simhitIt : *simHitCol) {
154  if (id.rawId() != simhitIt.detUnitId())
155  continue;
156  for (const auto& mId : matchedId) {
157  if (simhitIt.trackId() == mId.first) {
158  if (!simhitClosest || std::abs(simhitIt.localPosition().x() - lp.x()) < minx) {
159  minx = std::abs(simhitIt.localPosition().x() - lp.x());
160  simhitClosest = &simhitIt;
161  }
162  }
163  }
164  } //end loop over PSimhitcontainers
165  } //end loop over simHits
166 
167  if (!simhitClosest)
168  continue;
169 
170  // call the base class method to fill the plots
171  fillRechitHistos(simhitClosest, rechit, selectedSimTrackMap, nrechitLayerMap_primary);
172 
173  } // loop over tracking rechits
174  } // loop over tracks
175 
176  //fill nRecHit counter per layer
177  for (const auto& lme : nrechitLayerMap_primary) {
178  layerMEs_[lme.first].numberRecHitsprimary->Fill(nrechitLayerMap_primary[lme.first]);
179  }
180 }
181 
183  // rechitValidIT
185 
186  // call the base fillPsetDescription for the plots bookings
188 
189  //to be used in TrackerHitAssociator
190  desc.add<bool>("associatePixel", true);
191  desc.add<bool>("associateStrip", false);
192  desc.add<bool>("usePhase2Tracker", true);
193  desc.add<bool>("associateRecoTracks", false);
194  desc.add<bool>("associateHitbySimTrack", true);
195  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis", "Pixel"));
196  desc.add<std::vector<std::string>>("ROUList",
197  {
198  "TrackerHitsPixelBarrelLowTof",
199  "TrackerHitsPixelBarrelHighTof",
200  "TrackerHitsPixelEndcapLowTof",
201  "TrackerHitsPixelEndcapHighTof",
202  });
203  //
204  desc.add<edm::InputTag>("simTracksSrc", edm::InputTag("g4SimHits"));
205  desc.add<edm::InputTag>("SimVertexSource", edm::InputTag("g4SimHits"));
206  desc.add<double>("SimTrackMinPt", 2.0);
207  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"));
208  desc.add<std::string>("TopFolderName", "TrackerPhase2ITTrackingRecHitV");
209  desc.add<bool>("Verbosity", false);
210  descriptions.add("Phase2ITValidateTrackingRecHit", desc);
211 }
212 //define this as a plug-in
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::map< std::string, RecHitME > layerMEs_
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::string getITHistoId(uint32_t det_id, const TrackerTopology *tTopo)
T x() const
Definition: PV3DBase.h:59
TrackerHitAssociator::Config trackerHitAssociatorConfig_
int iEvent
Definition: GenABIO.cc:224
void fillITHistos(const edm::Event &iEvent, const TrackerHitAssociator &associateRecHit, const std::vector< edm::Handle< edm::PSimHitContainer >> &simHits, const std::map< unsigned int, SimTrack > &selectedSimTrackMap)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Phase2ITValidateTrackingRecHit(const edm::ParameterSet &)
const edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
Log< level::Info, false > LogInfo
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< reco::TrackCollection > tokenTracks_
auto const & tracks
cannot be loose
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
LocalPoint localPosition() const override
void fillRechitHistos(const PSimHit *simhitClosest, const SiPixelRecHit *rechit, const std::map< unsigned int, SimTrack > &selectedSimTrackMap, std::map< std::string, unsigned int > &nrechitLayerMap_primary)
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > simHitTokens_
std::vector< SimTrack > SimTrackContainer
Our base class.
Definition: SiPixelRecHit.h:23