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