CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Phase2OTValidateTrackingRecHit.cc
Go to the documentation of this file.
1 // Package: Phase2OTValidateTrackingRecHit
2 // Class: Phase2OTValidateTrackingRecHit
3 //
7 //
8 // Author: Suvankar Roy Chowdhury
9 // Date: March 2021
10 //
11 // system include files
20 
30 
40 
47 
50 
52 public:
55  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
56 
57  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
58 
59 private:
60  void fillOTHistos(const edm::Event& iEvent,
61  const TrackerHitAssociator& associateRecHit,
63  const std::map<unsigned int, SimTrack>& selectedSimTrackMap);
64 
67  const double simtrackminpt_;
70  std::vector<edm::EDGetTokenT<edm::PSimHitContainer>> simHitTokens_;
71 };
72 
73 //
74 // constructors
75 //
77  : Phase2OTValidateRecHitBase(iConfig),
78  config_(iConfig),
79  trackerHitAssociatorConfig_(iConfig, consumesCollector()),
80  simtrackminpt_(iConfig.getParameter<double>("SimTrackMinPt")),
81  tokenTracks_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
82  simTracksToken_(consumes<edm::SimTrackContainer>(iConfig.getParameter<edm::InputTag>("simTracksSrc"))) {
83  edm::LogInfo("Phase2OTValidateTrackingRecHit") << ">>> Construct Phase2OTValidateTrackingRecHit ";
84  for (const auto& itag : config_.getParameter<std::vector<edm::InputTag>>("PSimHitSource"))
85  simHitTokens_.push_back(consumes<edm::PSimHitContainer>(itag));
86 }
87 
88 //
89 // destructor
90 //
92  // do anything here that needs to be done at desctruction time
93  // (e.g. close files, deallocate resources etc.)
94  edm::LogInfo("Phase2OTValidateTrackingRecHit") << ">>> Destroy Phase2OTValidateTrackingRecHit ";
95 }
96 
97 //
98 // -- Analyze
99 //
101  std::vector<edm::Handle<edm::PSimHitContainer>> simHits;
102  for (const auto& itoken : simHitTokens_) {
103  const auto& simHitHandle = iEvent.getHandle(itoken);
104  if (!simHitHandle.isValid())
105  continue;
106  simHits.emplace_back(simHitHandle);
107  }
108  // Get the SimTracks and push them in a map of id, SimTrack
109  const auto& simTracks = iEvent.getHandle(simTracksToken_);
110  std::map<unsigned int, SimTrack> selectedSimTrackMap;
111  for (edm::SimTrackContainer::const_iterator simTrackIt(simTracks->begin()); simTrackIt != simTracks->end();
112  ++simTrackIt) {
113  if (simTrackIt->momentum().pt() > simtrackminpt_) {
114  selectedSimTrackMap.insert(std::make_pair(simTrackIt->trackId(), *simTrackIt));
115  }
116  }
117  TrackerHitAssociator associateRecHit(iEvent, trackerHitAssociatorConfig_);
118  fillOTHistos(iEvent, associateRecHit, simHits, selectedSimTrackMap);
119 }
120 
122  const TrackerHitAssociator& associateRecHit,
124  const std::map<unsigned int, SimTrack>& selectedSimTrackMap) {
125  const auto& tracks = iEvent.getHandle(tokenTracks_);
126  if (!tracks.isValid())
127  return;
128  std::map<std::string, unsigned int> nrechitLayerMapP_primary;
129  std::map<std::string, unsigned int> nrechitLayerMapS_primary;
130  // Loop over tracks
131  for (const auto& track : *tracks) {
132  for (const auto& hit : track.recHits()) {
133  // Get the detector unit's id
134  if (!hit->isValid())
135  continue;
136  auto detId = hit->geographicalId();
137  // check that we are in the pixel
138  auto subdetid = (detId.subdetId());
139  if (subdetid == PixelSubdetector::PixelBarrel || subdetid == PixelSubdetector::PixelEndcap)
140  continue;
141 
142  // determine the detector we are in
146  if (nrechitLayerMapP_primary.find(key) == nrechitLayerMapP_primary.end()) {
147  nrechitLayerMapP_primary.emplace(key, 1);
148  } else {
149  nrechitLayerMapP_primary[key] += 1;
150  }
152  if (nrechitLayerMapS_primary.find(key) == nrechitLayerMapS_primary.end()) {
153  nrechitLayerMapS_primary.emplace(key, 1);
154  } else {
155  nrechitLayerMapS_primary[key] += 1;
156  }
157  }
158  //GetSimHits
159  const Phase2TrackerRecHit1D* rechit = dynamic_cast<const Phase2TrackerRecHit1D*>(hit);
160  if (!rechit) {
161  edm::LogError("Phase2OTValidateTrackingRecHit")
162  << "Cannot cast tracking rechit to Phase2TrackerRecHit1D!" << std::endl;
163  continue;
164  }
165  const std::vector<SimHitIdpr>& matchedId = associateRecHit.associateHitId(*rechit);
166  const PSimHit* simhitClosest = nullptr;
167  LocalPoint lp = rechit->localPosition();
168  float mind = 1e4;
169  for (const auto& simHitCol : simHits) {
170  for (const auto& simhitIt : *simHitCol) {
171  if (detId.rawId() != simhitIt.detUnitId())
172  continue;
173  for (auto& mId : matchedId) {
174  if (simhitIt.trackId() == mId.first) {
175  float dx = simhitIt.localPosition().x() - lp.x();
176  float dy = simhitIt.localPosition().y() - lp.y();
177  float dist = dx * dx + dy * dy;
178  if (!simhitClosest || dist < mind) {
179  mind = dist;
180  simhitClosest = &simhitIt;
181  }
182  }
183  }
184  } //end loop over PSimhitcontainers
185  } //end loop over simHits
186  if (!simhitClosest)
187  continue;
189  simhitClosest, rechit, selectedSimTrackMap, nrechitLayerMapP_primary, nrechitLayerMapS_primary);
190 
191  } //end loop over rechits of a track
192  } //End loop over tracks
193 
194  //fill nRecHits per event
195  //fill nRecHit counter per layer
196  for (auto& lme : nrechitLayerMapP_primary) {
197  layerMEs_[lme.first].numberRecHitsprimary_P->Fill(nrechitLayerMapP_primary[lme.first]);
198  }
199  for (auto& lme : nrechitLayerMapS_primary) {
200  layerMEs_[lme.first].numberRecHitsprimary_S->Fill(nrechitLayerMapS_primary[lme.first]);
201  }
202 }
203 
206  // call the base fillPsetDescription for the plots bookings
208 
209  //for macro-pixel sensors
211  desc.add<edm::InputTag>("SimVertexSource", edm::InputTag("g4SimHits"));
212  desc.add<bool>("associatePixel", false);
213  desc.add<bool>("associateHitbySimTrack", true);
214  desc.add<bool>("Verbosity", false);
215  desc.add<bool>("associateStrip", true);
216  desc.add<edm::InputTag>("phase2TrackerSimLinkSrc", edm::InputTag("simSiPixelDigis", "Tracker"));
217  desc.add<bool>("associateRecoTracks", false);
218  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis", "Pixel"));
219  desc.add<bool>("usePhase2Tracker", true);
220  desc.add<edm::InputTag>("rechitsSrc", edm::InputTag("siPhase2RecHits"));
221  desc.add<edm::InputTag>("simTracksSrc", edm::InputTag("g4SimHits"));
222  desc.add<double>("SimTrackMinPt", 2.0);
223  desc.add<std::vector<edm::InputTag>>("PSimHitSource",
224  {
225  edm::InputTag("g4SimHits:TrackerHitsTIBLowTof"),
226  edm::InputTag("g4SimHits:TrackerHitsTIBHighTof"),
227  edm::InputTag("g4SimHits:TrackerHitsTIDLowTof"),
228  edm::InputTag("g4SimHits:TrackerHitsTIDHighTof"),
229  edm::InputTag("g4SimHits:TrackerHitsTOBLowTof"),
230  edm::InputTag("g4SimHits:TrackerHitsTOBHighTof"),
231  edm::InputTag("g4SimHits:TrackerHitsTECLowTof"),
232  edm::InputTag("g4SimHits:TrackerHitsTECHighTof"),
233  edm::InputTag("g4SimHits:TrackerHitsPixelBarrelLowTof"),
234  edm::InputTag("g4SimHits:TrackerHitsPixelBarrelHighTof"),
235  edm::InputTag("g4SimHits:TrackerHitsPixelEndcapLowTof"),
236  edm::InputTag("g4SimHits:TrackerHitsPixelEndcapHighTof"),
237  });
238 
239  desc.add<std::vector<std::string>>("ROUList",
240  {"TrackerHitsPixelBarrelLowTof",
241  "TrackerHitsPixelBarrelHighTof",
242  "TrackerHitsTIBLowTof",
243  "TrackerHitsTIBHighTof",
244  "TrackerHitsTIDLowTof",
245  "TrackerHitsTIDHighTof",
246  "TrackerHitsTOBLowTof",
247  "TrackerHitsTOBHighTof",
248  "TrackerHitsPixelEndcapLowTof",
249  "TrackerHitsPixelEndcapHighTof",
250  "TrackerHitsTECLowTof",
251  "TrackerHitsTECHighTof"});
252  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"));
253  desc.add<std::string>("TopFolderName", "TrackerPhase2OTTrackingRecHitV");
254  descriptions.add("Phase2OTValidateTrackingRecHit", desc);
255 }
256 
257 //define this as a plug-in
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 fillPSetDescription(edm::ParameterSetDescription &desc)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
T y() const
Definition: PV3DBase.h:60
auto const & tracks
cannot be loose
std::vector< edm::EDGetTokenT< edm::PSimHitContainer > > simHitTokens_
Log< level::Error, false > LogError
void fillOTHistos(const edm::Event &iEvent, const TrackerHitAssociator &associateRecHit, const std::vector< edm::Handle< edm::PSimHitContainer >> &simHits, const std::map< unsigned int, SimTrack > &selectedSimTrackMap)
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::string getOTHistoId(uint32_t det_id, const TrackerTopology *tTopo)
tuple key
prepare the HTCondor submission files and eventually submit them
TrackerHitAssociator::Config trackerHitAssociatorConfig_
const edm::EDGetTokenT< reco::TrackCollection > tokenTracks_
ModuleType getDetectorType(DetId) const
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Info, false > LogInfo
Phase2OTValidateTrackingRecHit(const edm::ParameterSet &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
tuple simHits
Definition: trackerHits.py:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
LocalPoint localPosition() const override
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
T x() const
Definition: PV3DBase.h:59
std::vector< SimTrack > SimTrackContainer
std::map< std::string, RecHitME > layerMEs_