CMS 3D CMS Logo

CTPPSPixelLocalTrackProducer.cc
Go to the documentation of this file.
2 
12 
16 
18 
22 
25 
26 #include <string>
27 #include <vector>
28 
30 
32 
34 
37 
39 public:
41 
43 
44  void produce(edm::Event &, const edm::EventSetup &) override;
45  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
46 
47 private:
53 
58  std::vector<uint32_t> listOfAllPlanes_;
59 
60  std::unique_ptr<RPixDetPatternFinder> patternFinder_;
61  std::unique_ptr<RPixDetTrackFinder> trackFinder_;
62 
64 };
65 
66 //------------------------------------------------------------------------------------------------//
67 
70  verbosity_ = parameterSet.getUntrackedParameter<int>("verbosity");
71  maxHitPerRomanPot_ = parameterSet.getParameter<int>("maxHitPerRomanPot");
72  maxHitPerPlane_ = parameterSet.getParameter<int>("maxHitPerPlane");
73  maxTrackPerRomanPot_ = parameterSet.getParameter<int>("maxTrackPerRomanPot");
74  maxTrackPerPattern_ = parameterSet.getParameter<int>("maxTrackPerPattern");
75  numberOfPlanesPerPot_ = parameterSet.getParameter<int>("numberOfPlanesPerPot");
76 
77  std::string patternFinderAlgorithm = parameterSet.getParameter<std::string>("patternFinderAlgorithm");
78  std::string trackFitterAlgorithm = parameterSet.getParameter<std::string>("trackFinderAlgorithm");
79 
80  // pattern algorithm selector
81  if (patternFinderAlgorithm == "RPixRoadFinder") {
82  patternFinder_ = std::unique_ptr<RPixRoadFinder>(new RPixRoadFinder(parameterSet));
83  } else {
84  throw cms::Exception("CTPPSPixelLocalTrackProducer")
85  << "Pattern finder algorithm" << patternFinderAlgorithm << " does not exist";
86  }
87 
88  listOfAllPlanes_.reserve(6);
89  for (uint32_t i = 0; i < numberOfPlanesPerPot_; ++i) {
90  listOfAllPlanes_.push_back(i);
91  }
92 
93  //tracking algorithm selector
94  if (trackFitterAlgorithm == "RPixPlaneCombinatoryTracking") {
95  trackFinder_ = std::unique_ptr<RPixPlaneCombinatoryTracking>(new RPixPlaneCombinatoryTracking(parameterSet));
96  } else {
97  throw cms::Exception("CTPPSPixelLocalTrackProducer")
98  << "Tracking fitter algorithm" << trackFitterAlgorithm << " does not exist";
99  }
100  trackFinder_->setListOfPlanes(listOfAllPlanes_);
101  trackFinder_->initialize();
102 
103  tokenCTPPSPixelRecHit_ = consumes<edm::DetSetVector<CTPPSPixelRecHit>>(edm::InputTag(inputTag_));
104 
105  produces<edm::DetSetVector<CTPPSPixelLocalTrack>>();
106 }
107 
108 //------------------------------------------------------------------------------------------------//
109 
111 
112 //------------------------------------------------------------------------------------------------//
113 
116 
117  desc.add<std::string>("label", "ctppsPixelRecHits")
118  ->setComment("label of the RecHits input for the tracking algorithm");
119  desc.add<std::string>("patternFinderAlgorithm", "RPixRoadFinder")->setComment("algorithm type for pattern finder");
120  desc.add<std::string>("trackFinderAlgorithm", "RPixPlaneCombinatoryTracking")
121  ->setComment("algorithm type for track finder");
122  desc.add<uint>("trackMinNumberOfPoints", 3)->setComment("minimum number of planes to produce a track");
123  desc.addUntracked<int>("verbosity", 0)->setComment("verbosity for track producer");
124  desc.add<double>("maximumChi2OverNDF", 5.)->setComment("maximum Chi2OverNDF for accepting the track");
125  desc.add<double>("maximumXLocalDistanceFromTrack", 0.2)
126  ->setComment("maximum x distance in mm to associate a point not used for fit to the track");
127  desc.add<double>("maximumYLocalDistanceFromTrack", 0.3)
128  ->setComment("maximum y distance in mm to associate a point not used for fit to the track");
129  desc.add<int>("maxHitPerPlane", 20)
130  ->setComment("maximum hits per plane, events with higher number will not be fitted");
131  desc.add<int>("maxHitPerRomanPot", 60)
132  ->setComment("maximum hits per roman pot, events with higher number will not be fitted");
133  desc.add<int>("maxTrackPerRomanPot", 10)
134  ->setComment("maximum tracks per roman pot, events with higher track will not be saved");
135  desc.add<int>("maxTrackPerPattern", 5)
136  ->setComment("maximum tracks per pattern, events with higher track will not be saved");
137  desc.add<int>("numberOfPlanesPerPot", 6)->setComment("number of planes per pot");
138  desc.add<double>("roadRadius", 1.0)->setComment("radius of pattern search window");
139  desc.add<int>("minRoadSize", 3)->setComment("minimum number of points in a pattern");
140  desc.add<int>("maxRoadSize", 20)->setComment("maximum number of points in a pattern");
141 
142  descriptions.add("ctppsPixelLocalTracks", desc);
143 }
144 
145 //------------------------------------------------------------------------------------------------//
146 
148  // Step A: get inputs
149 
153 
154  // get geometry
155  edm::ESHandle<CTPPSGeometry> geometryHandler;
156  iSetup.get<VeryForwardRealGeometryRecord>().get(geometryHandler);
157  const CTPPSGeometry &geometry = *geometryHandler;
158  geometryWatcher_.check(iSetup);
159 
160  std::vector<CTPPSPixelDetId> listOfPotWithHighOccupancyPlanes;
161  std::map<CTPPSPixelDetId, uint32_t> mapHitPerPot;
162 
163  for (const auto &recHitSet : recHitVector) {
164  if (verbosity_ > 2)
165  edm::LogInfo("CTPPSPixelLocalTrackProducer")
166  << "Hits found in plane = " << recHitSet.detId() << " number = " << recHitSet.size();
167  CTPPSPixelDetId tmpRomanPotId = CTPPSPixelDetId(recHitSet.detId()).rpId();
168  uint32_t hitOnPlane = recHitSet.size();
169 
170  //Get the number of hits per pot
171  if (mapHitPerPot.find(tmpRomanPotId) == mapHitPerPot.end()) {
172  mapHitPerPot[tmpRomanPotId] = hitOnPlane;
173  } else
174  mapHitPerPot[tmpRomanPotId] += hitOnPlane;
175 
176  //check is the plane occupancy is too high and save the corresponding pot
177  if (maxHitPerPlane_ >= 0 && hitOnPlane > (uint32_t)maxHitPerPlane_) {
178  if (verbosity_ > 2)
179  edm::LogInfo("CTPPSPixelLocalTrackProducer")
180  << " ---> To many hits in the plane, pot will be excluded from tracking cleared";
181  listOfPotWithHighOccupancyPlanes.push_back(tmpRomanPotId);
182  }
183  }
184 
185  //remove rechit for pot with too many hits or containing planes with too many hits
186  for (const auto &recHitSet : recHitVector) {
187  const auto tmpDetectorId = CTPPSPixelDetId(recHitSet.detId());
188  const auto tmpRomanPotId = tmpDetectorId.rpId();
189 
190  if ((maxHitPerRomanPot_ >= 0 && mapHitPerPot[tmpRomanPotId] > (uint32_t)maxHitPerRomanPot_) ||
191  find(listOfPotWithHighOccupancyPlanes.begin(), listOfPotWithHighOccupancyPlanes.end(), tmpRomanPotId) !=
192  listOfPotWithHighOccupancyPlanes.end()) {
193  edm::DetSet<CTPPSPixelRecHit> &tmpDetSet = recHitVector[tmpDetectorId];
194  tmpDetSet.clear();
195  }
196  }
197 
199 
200  // Pattern finder
201 
202  patternFinder_->clear();
203  patternFinder_->setHits(&recHitVector);
204  patternFinder_->setGeometry(&geometry);
205  patternFinder_->findPattern();
206  std::vector<RPixDetPatternFinder::Road> patternVector = patternFinder_->getPatterns();
207 
208  //loop on all the patterns
209  int numberOfTracks = 0;
210 
211  for (const auto &pattern : patternVector) {
212  CTPPSPixelDetId firstHitDetId = CTPPSPixelDetId(pattern.at(0).detId);
213  CTPPSPixelDetId romanPotId = firstHitDetId.rpId();
214 
215  std::map<CTPPSPixelDetId, std::vector<RPixDetPatternFinder::PointInPlane>>
216  hitOnPlaneMap; //hit of the pattern organized by plane
217 
218  //loop on all the hits of the pattern
219  for (const auto &hit : pattern) {
220  CTPPSPixelDetId hitDetId = CTPPSPixelDetId(hit.detId);
221  CTPPSPixelDetId tmpRomanPotId = hitDetId.rpId();
222 
223  if (tmpRomanPotId != romanPotId) { //check that the hits belong to the same tracking station
224  throw cms::Exception("CTPPSPixelLocalTrackProducer")
225  << "Hits in the pattern must belong to the same tracking station";
226  }
227 
228  if (hitOnPlaneMap.find(hitDetId) ==
229  hitOnPlaneMap.end()) { //add the plane key in case it is the first hit of that plane
230  std::vector<RPixDetPatternFinder::PointInPlane> hitOnPlane;
231  hitOnPlane.push_back(hit);
232  hitOnPlaneMap[hitDetId] = hitOnPlane;
233  } else
234  hitOnPlaneMap[hitDetId].push_back(hit); //add the hit to an existing plane key
235  }
236 
237  trackFinder_->clear();
238  trackFinder_->setRomanPotId(romanPotId);
239  trackFinder_->setHits(&hitOnPlaneMap);
240  trackFinder_->setGeometry(&geometry);
241  trackFinder_->setZ0(geometry.rpTranslation(romanPotId).z());
242  trackFinder_->findTracks(iEvent.getRun().id().run());
243  std::vector<CTPPSPixelLocalTrack> tmpTracksVector = trackFinder_->getLocalTracks();
244 
245  if (verbosity_ > 2)
246  edm::LogInfo("CTPPSPixelLocalTrackProducer") << "tmpTracksVector = " << tmpTracksVector.size();
247  if (maxTrackPerPattern_ >= 0 && tmpTracksVector.size() > (uint32_t)maxTrackPerPattern_) {
248  if (verbosity_ > 2)
249  edm::LogInfo("CTPPSPixelLocalTrackProducer") << " ---> To many tracks in the pattern, cleared";
250  continue;
251  }
252 
253  for (const auto &track : tmpTracksVector) {
254  ++numberOfTracks;
255  edm::DetSet<CTPPSPixelLocalTrack> &tmpDetSet = foundTracks.find_or_insert(romanPotId);
256  tmpDetSet.push_back(track);
257  }
258  }
259 
260  if (verbosity_ > 1)
261  edm::LogInfo("CTPPSPixelLocalTrackProducer") << "Number of tracks will be saved = " << numberOfTracks;
262 
263  for (const auto &track : foundTracks) {
264  if (verbosity_ > 1)
265  edm::LogInfo("CTPPSPixelLocalTrackProducer")
266  << "Track found in detId = " << track.detId() << " number = " << track.size();
267  if (maxTrackPerRomanPot_ >= 0 && track.size() > (uint32_t)maxTrackPerRomanPot_) {
268  if (verbosity_ > 1)
269  edm::LogInfo("CTPPSPixelLocalTrackProducer") << " ---> Too many tracks in the pot, cleared";
270  CTPPSPixelDetId tmpRomanPotId = CTPPSPixelDetId(track.detId());
271  edm::DetSet<CTPPSPixelLocalTrack> &tmpDetSet = foundTracks[tmpRomanPotId];
272  tmpDetSet.clear();
273  }
274  }
275 
276  iEvent.put(std::make_unique<edm::DetSetVector<CTPPSPixelLocalTrack>>(foundTracks));
277 
278  return;
279 }
280 
ConfigurationDescriptions.h
edm::ESWatcher::check
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
VeryForwardRealGeometryRecord
Event setup record containing the real (actual) geometry information.
Definition: VeryForwardRealGeometryRecord.h:22
edm::DetSet::push_back
void push_back(const T &t)
Definition: DetSet.h:66
edm::DetSetVector< CTPPSPixelRecHit >
Handle.h
mps_fire.i
i
Definition: mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
input
static const std::string input
Definition: EdmProvDump.cc:48
CTPPSPixelLocalTrack.h
MessageLogger.h
CTPPSGeometry
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:29
CTPPSPixelLocalTrackProducer::maxHitPerRomanPot_
int maxHitPerRomanPot_
Definition: CTPPSPixelLocalTrackProducer.cc:50
edm::ESWatcher< VeryForwardRealGeometryRecord >
ESHandle.h
edm::DetSet
Definition: DetSet.h:23
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:32
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm::DetSet::clear
void clear()
Definition: DetSet.h:71
geometry
Definition: geometry.py:1
CTPPSPixelDetId.h
edm::LogInfo
Definition: MessageLogger.h:254
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
CTPPSPixelLocalTrackProducer::trackFinder_
std::unique_ptr< RPixDetTrackFinder > trackFinder_
Definition: CTPPSPixelLocalTrackProducer.cc:61
EDProducer.h
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
year_2016_postTS2_cff.rpId
rpId
Definition: year_2016_postTS2_cff.py:23
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
CTPPSPixelLocalTrackProducer::maxTrackPerPattern_
int maxTrackPerPattern_
Definition: CTPPSPixelLocalTrackProducer.cc:52
edm::Handle
Definition: AssociativeIterator.h:50
parallelization.uint
uint
Definition: parallelization.py:124
CTPPSPixelRecHit.h
CTPPSPixelLocalTrackProducer::maxHitPerPlane_
int maxHitPerPlane_
Definition: CTPPSPixelLocalTrackProducer.cc:49
MakerMacros.h
CTPPSPixelLocalTrackProducer::~CTPPSPixelLocalTrackProducer
~CTPPSPixelLocalTrackProducer() override
Definition: CTPPSPixelLocalTrackProducer.cc:110
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
RPixDetPatternFinder.h
CTPPSPixelLocalTrackProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: CTPPSPixelLocalTrackProducer.cc:147
RPixRoadFinder
Definition: RPixRoadFinder.h:36
edm::ESHandle
Definition: DTSurvey.h:22
CTPPSPixelLocalTrackProducer::verbosity_
int verbosity_
Definition: CTPPSPixelLocalTrackProducer.cc:48
CTPPSPixelLocalTrackProducer::maxTrackPerRomanPot_
int maxTrackPerRomanPot_
Definition: CTPPSPixelLocalTrackProducer.cc:51
CTPPSPixelLocalTrackProducer::listOfAllPlanes_
std::vector< uint32_t > listOfAllPlanes_
Definition: CTPPSPixelLocalTrackProducer.cc:58
ParameterSetDescription.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
FastTrackerRecHitMaskProducer_cfi.recHits
recHits
Definition: FastTrackerRecHitMaskProducer_cfi.py:8
CTPPSPixelLocalTrackProducer::CTPPSPixelLocalTrackProducer
CTPPSPixelLocalTrackProducer(const edm::ParameterSet &parameterSet)
Definition: CTPPSPixelLocalTrackProducer.cc:68
edm::ParameterSetDescription::addUntracked
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:100
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
RPixRoadFinder.h
Event.h
iEvent
int iEvent
Definition: GenABIO.cc:224
CTPPSPixelLocalTrackProducer::patternFinder_
std::unique_ptr< RPixDetPatternFinder > patternFinder_
Definition: CTPPSPixelLocalTrackProducer.cc:60
topSingleLeptonDQM_PU_cfi.pattern
pattern
Definition: topSingleLeptonDQM_PU_cfi.py:39
edm::stream::EDProducer
Definition: EDProducer.h:38
CTPPSPixelLocalTrackProducer::numberOfPlanesPerPot_
uint32_t numberOfPlanesPerPot_
Definition: CTPPSPixelLocalTrackProducer.cc:57
edm::EventSetup
Definition: EventSetup.h:57
CTPPSPixelLocalTrackProducer::geometryWatcher_
edm::ESWatcher< VeryForwardRealGeometryRecord > geometryWatcher_
Definition: CTPPSPixelLocalTrackProducer.cc:56
DetSetVector.h
CTPPSPixelDetId
Definition: CTPPSPixelDetId.h:16
get
#define get
RPixPlaneCombinatoryTracking
Definition: RPixPlaneCombinatoryTracking.h:24
VeryForwardRealGeometryRecord.h
InputTag.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
DetId.h
Frameworkfwd.h
ESWatcher.h
Exception
Definition: hltDiff.cc:246
CTPPSPixelLocalTrackProducer::tokenCTPPSPixelRecHit_
edm::EDGetTokenT< edm::DetSetVector< CTPPSPixelRecHit > > tokenCTPPSPixelRecHit_
Definition: CTPPSPixelLocalTrackProducer.cc:55
RPixPlaneCombinatoryTracking.h
edm::parameterSet
ParameterSet const & parameterSet(Provenance const &provenance, ProcessHistory const &history)
Definition: Provenance.cc:11
CTPPSPixelLocalTrackProducer
Definition: CTPPSPixelLocalTrackProducer.cc:38
EventSetup.h
CTPPSPixelLocalTrackProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: CTPPSPixelLocalTrackProducer.cc:114
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
edm::DetSetVector::find_or_insert
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:234
DetSet.h
ParameterSet.h
edm::ParameterDescriptionNode::setComment
void setComment(std::string const &value)
Definition: ParameterDescriptionNode.cc:106
CTPPSDetId::rpId
CTPPSDetId rpId() const
Definition: CTPPSDetId.h:82
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
RPixDetTrackFinder.h
hit
Definition: SiStripHitEffFromCalibTree.cc:88
CTPPSPixelLocalTrackProducer::run
void run(const edm::DetSetVector< CTPPSPixelRecHit > &input, edm::DetSetVector< CTPPSPixelLocalTrack > &output)
CTPPSPixelLocalTrackProducer::inputTag_
edm::InputTag inputTag_
Definition: CTPPSPixelLocalTrackProducer.cc:54