CMS 3D CMS Logo

CTPPSPixelLocalTrackProducer.cc
Go to the documentation of this file.
2 
12 
17 
19 
23 
26 
27 #include <memory>
28 
29 #include <string>
30 #include <vector>
31 
33 
35 
37 
40 
43 
44 namespace {
45  constexpr int rocMask = 0xE000;
46  constexpr int rocOffset = 13;
47  constexpr int rocSizeInPixels = 4160;
48 } // namespace
49 
51 public:
53 
55 
56  void produce(edm::Event &, const edm::EventSetup &) override;
57  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
58 
59 private:
65 
70 
72 
74  std::vector<uint32_t> listOfAllPlanes_;
75 
76  std::unique_ptr<RPixDetPatternFinder> patternFinder_;
77  std::unique_ptr<RPixDetTrackFinder> trackFinder_;
78 
79  // void run(const edm::DetSetVector<CTPPSPixelRecHit> &input, edm::DetSetVector<CTPPSPixelLocalTrack> &output);
80 };
81 
82 //------------------------------------------------------------------------------------------------//
83 
86  verbosity_ = parameterSet.getUntrackedParameter<int>("verbosity");
87  maxHitPerRomanPot_ = parameterSet.getParameter<int>("maxHitPerRomanPot");
88  maxHitPerPlane_ = parameterSet.getParameter<int>("maxHitPerPlane");
89  maxTrackPerRomanPot_ = parameterSet.getParameter<int>("maxTrackPerRomanPot");
90  maxTrackPerPattern_ = parameterSet.getParameter<int>("maxTrackPerPattern");
91  numberOfPlanesPerPot_ = parameterSet.getParameter<int>("numberOfPlanesPerPot");
92 
94  std::string trackFitterAlgorithm = parameterSet.getParameter<std::string>("trackFinderAlgorithm");
95 
96  // pattern algorithm selector
97  if (patternFinderAlgorithm == "RPixRoadFinder") {
98  patternFinder_ = std::make_unique<RPixRoadFinder>(parameterSet);
99  } else {
100  throw cms::Exception("CTPPSPixelLocalTrackProducer")
101  << "Pattern finder algorithm" << patternFinderAlgorithm << " does not exist";
102  }
103 
104  listOfAllPlanes_.reserve(6);
105  for (uint32_t i = 0; i < numberOfPlanesPerPot_; ++i) {
106  listOfAllPlanes_.push_back(i);
107  }
108 
109  //tracking algorithm selector
110  if (trackFitterAlgorithm == "RPixPlaneCombinatoryTracking") {
111  trackFinder_ = std::make_unique<RPixPlaneCombinatoryTracking>(parameterSet);
112  } else {
113  throw cms::Exception("CTPPSPixelLocalTrackProducer")
114  << "Tracking fitter algorithm" << trackFitterAlgorithm << " does not exist";
115  }
116  trackFinder_->setListOfPlanes(listOfAllPlanes_);
117  trackFinder_->initialize();
118 
119  tokenCTPPSPixelRecHit_ = consumes<edm::DetSetVector<CTPPSPixelRecHit>>(inputTag_);
120  tokenCTPPSGeometry_ = esConsumes<CTPPSGeometry, VeryForwardRealGeometryRecord>();
121  tokenCTPPSPixelAnalysisMask_ = esConsumes<CTPPSPixelAnalysisMask, CTPPSPixelAnalysisMaskRcd>();
122 
123  produces<edm::DetSetVector<CTPPSPixelLocalTrack>>();
124 }
125 
126 //------------------------------------------------------------------------------------------------//
127 
129 
130 //------------------------------------------------------------------------------------------------//
131 
134 
135  desc.add<edm::InputTag>("tag", edm::InputTag("ctppsPixelRecHits"))
136  ->setComment("inputTag of the RecHits input for the tracking algorithm");
137  desc.add<std::string>("patternFinderAlgorithm", "RPixRoadFinder")->setComment("algorithm type for pattern finder");
138  desc.add<std::string>("trackFinderAlgorithm", "RPixPlaneCombinatoryTracking")
139  ->setComment("algorithm type for track finder");
140  desc.add<uint>("trackMinNumberOfPoints", 3)->setComment("minimum number of planes to produce a track");
141  desc.addUntracked<int>("verbosity", 0)->setComment("verbosity for track producer");
142  desc.add<double>("maximumChi2OverNDF", 5.)->setComment("maximum Chi2OverNDF for accepting the track");
143  desc.add<double>("maximumXLocalDistanceFromTrack", 0.2)
144  ->setComment("maximum x distance in mm to associate a point not used for fit to the track");
145  desc.add<double>("maximumYLocalDistanceFromTrack", 0.3)
146  ->setComment("maximum y distance in mm to associate a point not used for fit to the track");
147  desc.add<int>("maxHitPerPlane", 20)
148  ->setComment("maximum hits per plane, events with higher number will not be fitted");
149  desc.add<int>("maxHitPerRomanPot", 60)
150  ->setComment("maximum hits per roman pot, events with higher number will not be fitted");
151  desc.add<int>("maxTrackPerRomanPot", 10)
152  ->setComment("maximum tracks per roman pot, events with higher track will not be saved");
153  desc.add<int>("maxTrackPerPattern", 5)
154  ->setComment("maximum tracks per pattern, events with higher track will not be saved");
155  desc.add<int>("numberOfPlanesPerPot", 6)->setComment("number of planes per pot");
156  desc.add<double>("roadRadius", 1.0)->setComment("radius of pattern search window");
157  desc.add<int>("minRoadSize", 3)->setComment("minimum number of points in a pattern");
158  desc.add<int>("maxRoadSize", 20)->setComment("maximum number of points in a pattern");
159  //parameters for bad pot reconstruction patch 45-220-fr 2022
160  desc.add<double>("roadRadiusBadPot", 0.5)->setComment("radius of pattern search window for bad Pot");
161  // desc.add<bool>("isBadPot", true)->setComment("flag to enable road search for bad pot");
162 
163  descriptions.add("ctppsPixelLocalTracks", desc);
164 }
165 
166 //------------------------------------------------------------------------------------------------//
167 
169  // Step A: get inputs
170 
174 
175  // get geometry
177  const CTPPSGeometry &geometry = *geometryHandler;
178  geometryWatcher_.check(iSetup);
179 
180  // get mask
181  bool isBadPot_45_220 = false;
182  if (!recHits->empty()) {
183  const auto &mask = iSetup.getData(tokenCTPPSPixelAnalysisMask_);
184 
185  // Read Mask checking if 45-220-far is masked as bad and needs special treatment
186  std::map<uint32_t, CTPPSPixelROCAnalysisMask> const &maschera = mask.analysisMask;
187 
188  bool mask_45_220[6][6] = {{false}};
189  for (auto const &det : maschera) {
190  CTPPSPixelDetId detId(det.first);
191  unsigned int rocNum = (det.first & rocMask) >> rocOffset;
192  if (rocNum > 5 || detId.plane() > 5)
193  throw cms::Exception("InvalidRocOrPlaneNumber") << "roc number from mask: " << rocNum;
194 
195  if (detId.arm() == 0 && detId.station() == 2 && detId.rp() == 3) { // pot 45-220-far
196  if (det.second.maskedPixels.size() == rocSizeInPixels) { // roc fully masked
197  mask_45_220[detId.plane()][rocNum] = true;
198  }
199  }
200  }
201 
202  // search for specific pattern that requires special reconstruction (isBadPot)
203  isBadPot_45_220 = mask_45_220[1][4] && mask_45_220[1][5] && mask_45_220[2][4] && mask_45_220[2][5] &&
204  mask_45_220[3][4] && mask_45_220[3][5] && mask_45_220[4][4] && mask_45_220[4][5];
205  }
206  std::vector<CTPPSPixelDetId> listOfPotWithHighOccupancyPlanes;
207  std::map<CTPPSPixelDetId, uint32_t> mapHitPerPot;
208 
209  for (const auto &recHitSet : recHitVector) {
210  if (verbosity_ > 2)
211  edm::LogInfo("CTPPSPixelLocalTrackProducer")
212  << "Hits found in plane = " << recHitSet.detId() << " number = " << recHitSet.size();
213  CTPPSPixelDetId tmpRomanPotId = CTPPSPixelDetId(recHitSet.detId()).rpId();
214  uint32_t hitOnPlane = recHitSet.size();
215 
216  //Get the number of hits per pot
217  if (mapHitPerPot.find(tmpRomanPotId) == mapHitPerPot.end()) {
218  mapHitPerPot[tmpRomanPotId] = hitOnPlane;
219  } else
220  mapHitPerPot[tmpRomanPotId] += hitOnPlane;
221 
222  //check is the plane occupancy is too high and save the corresponding pot
223  if (maxHitPerPlane_ >= 0 && hitOnPlane > (uint32_t)maxHitPerPlane_) {
224  if (verbosity_ > 2)
225  edm::LogInfo("CTPPSPixelLocalTrackProducer")
226  << " ---> To many hits in the plane, pot will be excluded from tracking cleared";
227  listOfPotWithHighOccupancyPlanes.push_back(tmpRomanPotId);
228  }
229  }
230 
231  //remove rechit for pot with too many hits or containing planes with too many hits
232  for (const auto &recHitSet : recHitVector) {
233  const auto tmpDetectorId = CTPPSPixelDetId(recHitSet.detId());
234  const auto tmpRomanPotId = tmpDetectorId.rpId();
235 
236  if ((maxHitPerRomanPot_ >= 0 && mapHitPerPot[tmpRomanPotId] > (uint32_t)maxHitPerRomanPot_) ||
237  find(listOfPotWithHighOccupancyPlanes.begin(), listOfPotWithHighOccupancyPlanes.end(), tmpRomanPotId) !=
238  listOfPotWithHighOccupancyPlanes.end()) {
239  edm::DetSet<CTPPSPixelRecHit> &tmpDetSet = recHitVector[tmpDetectorId];
240  tmpDetSet.clear();
241  }
242  }
243 
245 
246  // Pattern finder
247 
248  patternFinder_->clear();
249  patternFinder_->setHits(&recHitVector);
250  patternFinder_->setGeometry(&geometry);
251  patternFinder_->findPattern(isBadPot_45_220);
252  std::vector<RPixDetPatternFinder::Road> patternVector = patternFinder_->getPatterns();
253 
254  //loop on all the patterns
255  int numberOfTracks = 0;
256 
257  for (const auto &pattern : patternVector) {
258  CTPPSPixelDetId firstHitDetId = CTPPSPixelDetId(pattern.at(0).detId);
259  CTPPSPixelDetId romanPotId = firstHitDetId.rpId();
260 
261  std::map<CTPPSPixelDetId, std::vector<RPixDetPatternFinder::PointInPlane>>
262  hitOnPlaneMap; //hit of the pattern organized by plane
263 
264  //loop on all the hits of the pattern
265  for (const auto &hit : pattern) {
266  CTPPSPixelDetId hitDetId = CTPPSPixelDetId(hit.detId);
267  CTPPSPixelDetId tmpRomanPotId = hitDetId.rpId();
268 
269  if (tmpRomanPotId != romanPotId) { //check that the hits belong to the same tracking station
270  throw cms::Exception("CTPPSPixelLocalTrackProducer")
271  << "Hits in the pattern must belong to the same tracking station";
272  }
273 
274  if (hitOnPlaneMap.find(hitDetId) ==
275  hitOnPlaneMap.end()) { //add the plane key in case it is the first hit of that plane
276  std::vector<RPixDetPatternFinder::PointInPlane> hitOnPlane;
277  hitOnPlane.push_back(hit);
278  hitOnPlaneMap[hitDetId] = hitOnPlane;
279  } else
280  hitOnPlaneMap[hitDetId].push_back(hit); //add the hit to an existing plane key
281  }
282 
283  trackFinder_->clear();
284  trackFinder_->setRomanPotId(romanPotId);
285  trackFinder_->setHits(&hitOnPlaneMap);
286  trackFinder_->setGeometry(&geometry);
287  trackFinder_->setZ0(geometry.rpTranslation(romanPotId).z());
288  trackFinder_->findTracks(iEvent.getRun().id().run());
289  std::vector<CTPPSPixelLocalTrack> tmpTracksVector = trackFinder_->getLocalTracks();
290 
291  if (verbosity_ > 2)
292  edm::LogInfo("CTPPSPixelLocalTrackProducer") << "tmpTracksVector = " << tmpTracksVector.size();
293  if (maxTrackPerPattern_ >= 0 && tmpTracksVector.size() > (uint32_t)maxTrackPerPattern_) {
294  if (verbosity_ > 2)
295  edm::LogInfo("CTPPSPixelLocalTrackProducer") << " ---> To many tracks in the pattern, cleared";
296  continue;
297  }
298 
299  for (const auto &track : tmpTracksVector) {
300  ++numberOfTracks;
301  edm::DetSet<CTPPSPixelLocalTrack> &tmpDetSet = foundTracks.find_or_insert(romanPotId);
302  tmpDetSet.push_back(track);
303  }
304  }
305 
306  if (verbosity_ > 1)
307  edm::LogInfo("CTPPSPixelLocalTrackProducer") << "Number of tracks will be saved = " << numberOfTracks;
308 
309  for (const auto &track : foundTracks) {
310  if (verbosity_ > 1)
311  edm::LogInfo("CTPPSPixelLocalTrackProducer")
312  << "Track found in detId = " << track.detId() << " number = " << track.size();
313  if (maxTrackPerRomanPot_ >= 0 && track.size() > (uint32_t)maxTrackPerRomanPot_) {
314  if (verbosity_ > 1)
315  edm::LogInfo("CTPPSPixelLocalTrackProducer") << " ---> Too many tracks in the pot, cleared";
316  CTPPSPixelDetId tmpRomanPotId = CTPPSPixelDetId(track.detId());
317  edm::DetSet<CTPPSPixelLocalTrack> &tmpDetSet = foundTracks[tmpRomanPotId];
318  tmpDetSet.clear();
319  }
320  }
321 
322  iEvent.put(std::make_unique<edm::DetSetVector<CTPPSPixelLocalTrack>>(foundTracks));
323 
324  return;
325 }
326 
edm::ESWatcher< VeryForwardRealGeometryRecord > geometryWatcher_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< CTPPSGeometry, VeryForwardRealGeometryRecord > tokenCTPPSGeometry_
void push_back(const T &t)
Definition: DetSet.h:66
uint32_t arm() const
Definition: CTPPSDetId.h:51
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &, const edm::EventSetup &) override
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
ParameterSet const & parameterSet(StableProvenance const &provenance, ProcessHistory const &history)
Definition: Provenance.cc:11
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:234
T getUntrackedParameter(std::string const &, T const &) const
int iEvent
Definition: GenABIO.cc:224
std::unique_ptr< RPixDetTrackFinder > trackFinder_
CTPPSDetId rpId() const
Definition: CTPPSDetId.h:78
uint32_t plane() const
bool getData(T &iHolder) const
Definition: EventSetup.h:122
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:30
Log< level::Info, false > LogInfo
uint32_t rp() const
Definition: CTPPSDetId.h:65
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelRecHit > > tokenCTPPSPixelRecHit_
uint32_t station() const
Definition: CTPPSDetId.h:58
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
std::unique_ptr< RPixDetPatternFinder > patternFinder_
void clear()
Definition: DetSet.h:71
CTPPSPixelLocalTrackProducer(const edm::ParameterSet &parameterSet)
edm::ESGetToken< CTPPSPixelAnalysisMask, CTPPSPixelAnalysisMaskRcd > tokenCTPPSPixelAnalysisMask_