CMS 3D CMS Logo

Phase2OTValidateRecHit.cc
Go to the documentation of this file.
1 // Package: Phase2OTValidateRecHit
2 // Class: Phase2OTValidateRecHit
3 //
7 //
8 // Author: Suvankar Roy Chowdhury
9 // Date: March 2021
10 //
11 // system include files
20 
30 
44 
47 
49 public:
51  ~Phase2OTValidateRecHit() override;
52  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
53 
54  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
55 
56 private:
57  void fillOTHistos(const edm::Event& iEvent,
58  const TrackerHitAssociator& associateRecHit,
60  const std::map<unsigned int, SimTrack>& selectedSimTrackMap);
61 
63  const double simtrackminpt_;
66  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> simHitTokens_;
67 };
68 
69 //
70 // constructors
71 //
73  : Phase2OTValidateRecHitBase(iConfig),
74  trackerHitAssociatorConfig_(iConfig, consumesCollector()),
75  simtrackminpt_(iConfig.getParameter<double>("SimTrackMinPt")),
76  tokenRecHitsOT_(consumes<Phase2TrackerRecHit1DCollectionNew>(iConfig.getParameter<edm::InputTag>("rechitsSrc"))),
77  simTracksToken_(consumes<edm::SimTrackContainer>(iConfig.getParameter<edm::InputTag>("simTracksSrc"))) {
78  edm::LogInfo("Phase2OTValidateRecHit") << ">>> Construct Phase2OTValidateRecHit ";
79  for (const auto& itag : config_.getParameter<std::vector<edm::InputTag>>("PSimHitSource"))
80  simHitTokens_.push_back(consumes<edm::PSimHitContainer>(itag));
81 }
82 
83 //
84 // destructor
85 //
87  // do anything here that needs to be done at desctruction time
88  // (e.g. close files, deallocate resources etc.)
89  edm::LogInfo("Phase2OTValidateRecHit") << ">>> Destroy Phase2OTValidateRecHit ";
90 }
91 
92 //
93 // -- Analyze
94 //
96  std::vector<edm::Handle<edm::PSimHitContainer>> simHits;
97  for (const auto& itoken : simHitTokens_) {
98  const auto& simHitHandle = iEvent.getHandle(itoken);
99  if (!simHitHandle.isValid())
100  continue;
101  simHits.emplace_back(simHitHandle);
102  }
103  // Get the SimTracks and push them in a map of id, SimTrack
104  const auto& simTracks = iEvent.getHandle(simTracksToken_);
105  std::map<unsigned int, SimTrack> selectedSimTrackMap;
106  for (const auto& simTrackIt : *simTracks)
107  if (simTrackIt.momentum().pt() > simtrackminpt_) {
108  selectedSimTrackMap.insert(std::make_pair(simTrackIt.trackId(), simTrackIt));
109  }
111  fillOTHistos(iEvent, associateRecHit, simHits, selectedSimTrackMap);
112 }
113 
115  const TrackerHitAssociator& associateRecHit,
117  const std::map<unsigned int, SimTrack>& selectedSimTrackMap) {
118  // Get the RecHits
119  const auto& rechits = iEvent.getHandle(tokenRecHitsOT_);
120  if (!rechits.isValid())
121  return;
122  std::map<std::string, unsigned int> nrechitLayerMapP_primary;
123  std::map<std::string, unsigned int> nrechitLayerMapS_primary;
124  // Loop over modules
126  for (const auto& DSViter : *rechits) {
127  // Get the detector unit's id
128  unsigned int rawid(DSViter.detId());
129  DetId detId(rawid);
130  // determine the detector we are in
134  if (nrechitLayerMapP_primary.find(key) == nrechitLayerMapP_primary.end()) {
135  nrechitLayerMapP_primary.insert(std::make_pair(key, DSViter.size()));
136  } else {
137  nrechitLayerMapP_primary[key] += DSViter.size();
138  }
140  if (nrechitLayerMapS_primary.find(key) == nrechitLayerMapS_primary.end()) {
141  nrechitLayerMapS_primary.insert(std::make_pair(key, DSViter.size()));
142  } else {
143  nrechitLayerMapS_primary[key] += DSViter.size();
144  }
145  }
146  //loop over rechits for a single detId
147  for (const auto& rechit : DSViter) {
148  //GetSimHits
149  const std::vector<SimHitIdpr>& matchedId = associateRecHit.associateHitId(rechit);
150  const PSimHit* simhitClosest = nullptr;
151  LocalPoint lp = rechit.localPosition();
152  float mind = 1e4;
153  for (const auto& simHitCol : simHits) {
154  for (const auto& simhitIt : *simHitCol) {
155  if (detId.rawId() != simhitIt.detUnitId())
156  continue;
157  for (auto& mId : matchedId) {
158  if (simhitIt.trackId() == mId.first) {
159  float dx = simhitIt.localPosition().x() - lp.x();
160  float dy = simhitIt.localPosition().y() - lp.y();
161  float dist = dx * dx + dy * dy;
162  if (!simhitClosest || dist < mind) {
163  mind = dist;
164  simhitClosest = &simhitIt;
165  }
166  }
167  }
168  } //end loop over PSimhitcontainers
169  } //end loop over simHits
170  if (!simhitClosest)
171  continue;
173  simhitClosest, &rechit, selectedSimTrackMap, nrechitLayerMapP_primary, nrechitLayerMapS_primary);
174 
175  } //end loop over rechits of a detId
176  } //End loop over DetSetVector
177 
178  //fill nRecHits per event
179  //fill nRecHit counter per layer
180  for (auto& lme : nrechitLayerMapP_primary) {
181  layerMEs_[lme.first].numberRecHitsprimary_P->Fill(nrechitLayerMapP_primary[lme.first]);
182  }
183  for (auto& lme : nrechitLayerMapS_primary) {
184  layerMEs_[lme.first].numberRecHitsprimary_S->Fill(nrechitLayerMapS_primary[lme.first]);
185  }
186 }
187 
190  // call the base fillPsetDescription for the plots bookings
192 
193  //for macro-pixel sensors
195  desc.add<edm::InputTag>("SimVertexSource", edm::InputTag("g4SimHits"));
196  desc.add<bool>("associatePixel", false);
197  desc.add<std::string>("TopFolderName", "TrackerPhase2OTRecHitV");
198  desc.add<bool>("associateHitbySimTrack", true);
199  desc.add<bool>("Verbosity", false);
200  desc.add<bool>("associateStrip", true);
201  desc.add<edm::InputTag>("phase2TrackerSimLinkSrc", edm::InputTag("simSiPixelDigis", "Tracker"));
202  desc.add<bool>("associateRecoTracks", false);
203  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis", "Pixel"));
204  desc.add<bool>("usePhase2Tracker", true);
205  desc.add<edm::InputTag>("rechitsSrc", edm::InputTag("siPhase2RecHits"));
206  desc.add<edm::InputTag>("simTracksSrc", edm::InputTag("g4SimHits"));
207  desc.add<double>("SimTrackMinPt", 2.0);
208  desc.add<std::vector<edm::InputTag>>("PSimHitSource",
209  {
210  edm::InputTag("g4SimHits:TrackerHitsTIBLowTof"),
211  edm::InputTag("g4SimHits:TrackerHitsTIBHighTof"),
212  edm::InputTag("g4SimHits:TrackerHitsTIDLowTof"),
213  edm::InputTag("g4SimHits:TrackerHitsTIDHighTof"),
214  edm::InputTag("g4SimHits:TrackerHitsTOBLowTof"),
215  edm::InputTag("g4SimHits:TrackerHitsTOBHighTof"),
216  edm::InputTag("g4SimHits:TrackerHitsTECLowTof"),
217  edm::InputTag("g4SimHits:TrackerHitsTECHighTof"),
218  edm::InputTag("g4SimHits:TrackerHitsPixelBarrelLowTof"),
219  edm::InputTag("g4SimHits:TrackerHitsPixelBarrelHighTof"),
220  edm::InputTag("g4SimHits:TrackerHitsPixelEndcapLowTof"),
221  edm::InputTag("g4SimHits:TrackerHitsPixelEndcapHighTof"),
222  });
223 
224  desc.add<std::vector<std::string>>("ROUList",
225  {"TrackerHitsPixelBarrelLowTof",
226  "TrackerHitsPixelBarrelHighTof",
227  "TrackerHitsTIBLowTof",
228  "TrackerHitsTIBHighTof",
229  "TrackerHitsTIDLowTof",
230  "TrackerHitsTIDHighTof",
231  "TrackerHitsTOBLowTof",
232  "TrackerHitsTOBHighTof",
233  "TrackerHitsPixelEndcapLowTof",
234  "TrackerHitsPixelEndcapHighTof",
235  "TrackerHitsTECLowTof",
236  "TrackerHitsTECHighTof"});
237 
238  descriptions.add("Phase2OTValidateRecHit", desc);
239 }
240 
241 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void fillOTRecHitHistos(const PSimHit *simhitClosest, const Phase2TrackerRecHit1D *rechit, const std::map< unsigned int, SimTrack > &selectedSimTrackMap, std::map< std::string, unsigned int > &nrechitLayerMapP_primary, std::map< std::string, unsigned int > &nrechitLayerMapS_primary)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > simHitTokens_
TrackerHitAssociator::Config trackerHitAssociatorConfig_
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
std::string getOTHistoId(uint32_t det_id, const TrackerTopology *tTopo)
ModuleType getDetectorType(DetId) const
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Phase2OTValidateRecHit(const edm::ParameterSet &)
Log< level::Info, false > LogInfo
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
const edm::EDGetTokenT< Phase2TrackerRecHit1DCollectionNew > tokenRecHitsOT_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
void fillOTHistos(const edm::Event &iEvent, const TrackerHitAssociator &associateRecHit, const std::vector< edm::Handle< edm::PSimHitContainer >> &simHits, const std::map< unsigned int, SimTrack > &selectedSimTrackMap)
std::vector< SimTrack > SimTrackContainer
std::map< std::string, RecHitME > layerMEs_
const edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_