CMS 3D CMS Logo

Phase2ITValidateRecHit.cc
Go to the documentation of this file.
1 // Package: Phase2ITValidateRecHit
2 // Class: Phase2ITValidateRecHit
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>
36 //--- for SimHit association
41 //DQM
45 // base class
48 
50 public:
52  ~Phase2ITValidateRecHit() override;
53  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56 
57 private:
58  void fillITHistos(const edm::Event& iEvent,
59  const TrackerHitAssociator& associateRecHit,
61  const std::map<unsigned int, SimTrack>& selectedSimTrackMap);
62 
65  const double simtrackminpt_;
68  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> simHitTokens_;
69 };
70 
72  : Phase2ITValidateRecHitBase(iConfig),
73  config_(iConfig),
74  trackerHitAssociatorConfig_(iConfig, consumesCollector()),
75  simtrackminpt_(iConfig.getParameter<double>("SimTrackMinPt")),
76  tokenRecHitsIT_(consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("rechitsSrc"))),
77  simTracksToken_(consumes<edm::SimTrackContainer>(iConfig.getParameter<edm::InputTag>("simTracksSrc"))) {
78  edm::LogInfo("Phase2ITValidateRecHit") << ">>> Construct Phase2ITValidateRecHit ";
79  for (const auto& itName : config_.getParameter<std::vector<std::string>>("ROUList")) {
80  simHitTokens_.push_back(consumes<std::vector<PSimHit>>(edm::InputTag("g4SimHits", itName)));
81  }
82 }
83 //
85  edm::LogInfo("Phase2ITValidateRecHit") << ">>> Destroy Phase2ITValidateRecHit ";
86 }
87 
89  std::vector<edm::Handle<edm::PSimHitContainer>> simHits;
90  for (const auto& itoken : simHitTokens_) {
91  const auto& simHitHandle = iEvent.getHandle(itoken);
92  if (!simHitHandle.isValid())
93  continue;
94  simHits.emplace_back(simHitHandle);
95  }
96  // Get the SimTracks and push them in a map of id, SimTrack
97  const auto& simTracks = iEvent.getHandle(simTracksToken_);
98 
99  std::map<unsigned int, SimTrack> selectedSimTrackMap;
100  for (const auto& simTrackIt : *simTracks) {
101  if (simTrackIt.momentum().pt() > simtrackminpt_) {
102  selectedSimTrackMap.emplace(simTrackIt.trackId(), simTrackIt);
103  }
104  }
106  fillITHistos(iEvent, associateRecHit, simHits, selectedSimTrackMap);
107 }
108 
110  const TrackerHitAssociator& associateRecHit,
112  const std::map<unsigned int, SimTrack>& selectedSimTrackMap) {
113  // Get the RecHits
114  const auto& rechits = iEvent.getHandle(tokenRecHitsIT_);
115  if (!rechits.isValid())
116  return;
117  std::map<std::string, unsigned int> nrechitLayerMap_primary;
118  // Loop over modules
119  for (const auto& DSViter : *rechits) {
120  // Get the detector unit's id
121  unsigned int rawid(DSViter.detId());
122  DetId detId(rawid);
123  // determine the detector we are in
125  if (nrechitLayerMap_primary.find(key) == nrechitLayerMap_primary.end()) {
126  nrechitLayerMap_primary.emplace(key, DSViter.size());
127  } else {
128  nrechitLayerMap_primary[key] += DSViter.size();
129  }
130  //loop over rechits for a single detId
131  for (const auto& rechit : DSViter) {
132  //GetSimHits
133  const std::vector<SimHitIdpr>& matchedId = associateRecHit.associateHitId(rechit);
134  const PSimHit* simhitClosest = nullptr;
135  float minx = 10000;
136  LocalPoint lp = rechit.localPosition();
137  for (const auto& simHitCol : simHits) {
138  for (const auto& simhitIt : *simHitCol) {
139  if (detId.rawId() != simhitIt.detUnitId())
140  continue;
141  for (const auto& mId : matchedId) {
142  if (simhitIt.trackId() == mId.first) {
143  if (!simhitClosest || abs(simhitIt.localPosition().x() - lp.x()) < minx) {
144  minx = std::abs(simhitIt.localPosition().x() - lp.x());
145  simhitClosest = &simhitIt;
146  }
147  }
148  }
149  } //end loop over PSimhitcontainers
150  } //end loop over simHits
151  if (!simhitClosest)
152  continue;
153 
154  // call the base class method to fill the plots
155  fillRechitHistos(simhitClosest, &rechit, selectedSimTrackMap, nrechitLayerMap_primary);
156 
157  } //end loop over rechits of a detId
158  } //End loop over DetSetVector
159 
160  //fill nRecHit counter per layer
161  for (const auto& lme : nrechitLayerMap_primary) {
162  layerMEs_[lme.first].numberRecHitsprimary->Fill(nrechitLayerMap_primary[lme.first]);
163  }
164 }
165 
167  // rechitValidIT
169 
170  // call the base fillPsetDescription for the plots bookings
172 
173  //to be used in TrackerHitAssociator
174  desc.add<bool>("associatePixel", true);
175  desc.add<bool>("associateStrip", false);
176  desc.add<bool>("usePhase2Tracker", true);
177  desc.add<bool>("associateRecoTracks", false);
178  desc.add<bool>("associateHitbySimTrack", true);
179  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis", "Pixel"));
180  desc.add<std::vector<std::string>>("ROUList",
181  {
182  "TrackerHitsPixelBarrelLowTof",
183  "TrackerHitsPixelBarrelHighTof",
184  "TrackerHitsPixelEndcapLowTof",
185  "TrackerHitsPixelEndcapHighTof",
186  });
187  //
188  desc.add<edm::InputTag>("simTracksSrc", edm::InputTag("g4SimHits"));
189  desc.add<edm::InputTag>("SimVertexSource", edm::InputTag("g4SimHits"));
190  desc.add<double>("SimTrackMinPt", 2.0);
191  desc.add<edm::InputTag>("rechitsSrc", edm::InputTag("siPixelRecHits"));
192  desc.add<std::string>("TopFolderName", "TrackerPhase2ITRecHitV");
193  desc.add<bool>("Verbosity", false);
194  descriptions.add("Phase2ITValidateRecHit", desc);
195 }
196 //define this as a plug-in
Phase2ITValidateRecHit(const edm::ParameterSet &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > simHitTokens_
std::map< std::string, RecHitME > layerMEs_
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
TrackerHitAssociator::Config trackerHitAssociatorConfig_
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
std::string getITHistoId(uint32_t det_id, const TrackerTopology *tTopo)
T x() const
Definition: PV3DBase.h:59
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillPSetDescription(edm::ParameterSetDescription &desc)
const edm::EDGetTokenT< SiPixelRecHitCollection > tokenRecHitsIT_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Log< level::Info, false > LogInfo
Definition: DetId.h:17
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void fillITHistos(const edm::Event &iEvent, const TrackerHitAssociator &associateRecHit, const std::vector< edm::Handle< edm::PSimHitContainer >> &simHits, const std::map< unsigned int, SimTrack > &selectedSimTrackMap)
HLT enums.
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< SimTrack > SimTrackContainer