CMS 3D CMS Logo

TotemRPUVPatternFinder.cc
Go to the documentation of this file.
1 /****************************************************************************
2 *
3 * This is a part of TOTEM offline software.
4 * Authors:
5 * Jan Kašpar (jan.kaspar@gmail.com)
6 *
7 ****************************************************************************/
8 
16 
22 
25 
27 
28 //----------------------------------------------------------------------------------------------------
29 
37 {
38  public:
40 
41  virtual ~TotemRPUVPatternFinder();
42 
43  virtual void produce(edm::Event& e, const edm::EventSetup& c) override;
44 
45  private:
48 
49  unsigned int verbosity;
50 
53 
56 
59 
62 
64  double threshold;
65 
67  double max_a_toFit;
68 
70  struct RPSettings
71  {
74  };
75 
77  std::map<unsigned int, RPSettings> exceptionalSettings;
78 
80 
82  void recognizeAndSelect(TotemRPUVPattern::ProjectionType proj, double z0, double threshold,
83  unsigned int planes_required,
85 };
86 
87 //----------------------------------------------------------------------------------------------------
88 
89 using namespace std;
90 using namespace edm;
91 
92 //----------------------------------------------------------------------------------------------------
93 
95  tagRecHit(conf.getParameter<edm::InputTag>("tagRecHit")),
96  verbosity(conf.getUntrackedParameter<unsigned int>("verbosity", 0)),
97  minPlanesPerProjectionToSearch(conf.getParameter<unsigned int>("minPlanesPerProjectionToSearch")),
98  minPlanesPerProjectionToFit(conf.getParameter<unsigned int>("minPlanesPerProjectionToFit")),
99  maxHitsPerPlaneToSearch(conf.getParameter<unsigned int>("maxHitsPerPlaneToSearch")),
100  lrcgn(new FastLineRecognition(conf.getParameter<double>("clusterSize_a"), conf.getParameter<double>("clusterSize_b"))),
101  threshold(conf.getParameter<double>("threshold")),
102  max_a_toFit(conf.getParameter<double>("max_a_toFit"))
103 {
104  for (const auto &ps : conf.getParameter< vector<ParameterSet> >("exceptionalSettings"))
105  {
106  unsigned int rpId = ps.getParameter<unsigned int>("rpId");
107 
108  RPSettings settings;
109  settings.minPlanesPerProjectionToFit_U = ps.getParameter<unsigned int>("minPlanesPerProjectionToFit_U");
110  settings.minPlanesPerProjectionToFit_V = ps.getParameter<unsigned int>("minPlanesPerProjectionToFit_V");
111  settings.threshold_U = ps.getParameter<double>("threshold_U");
112  settings.threshold_V = ps.getParameter<double>("threshold_V");
113 
114  exceptionalSettings[rpId] = settings;
115  }
116 
117  detSetVectorTotemRPRecHitToken = consumes<edm::DetSetVector<TotemRPRecHit> >(tagRecHit);
118 
119  produces<DetSetVector<TotemRPUVPattern>>();
120 }
121 
122 //----------------------------------------------------------------------------------------------------
123 
125 {
126  delete lrcgn;
127 }
128 
129 //----------------------------------------------------------------------------------------------------
130 
132  double z0, double threshold_loc, unsigned int planes_required,
134 {
135  // run recognition
136  DetSet<TotemRPUVPattern> newPatterns;
137  lrcgn->getPatterns(hits, z0, threshold_loc, newPatterns);
138 
139  // set pattern properties and copy to the global pattern collection
140  for (auto &p : newPatterns)
141  {
142  p.setProjection(proj);
143 
144  p.setFittable(true);
145 
146  set<unsigned int> planes;
147  for (const auto &ds : p.getHits())
148  planes.insert(TotemRPDetId(ds.detId()).plane());
149 
150  if (planes.size() < planes_required)
151  p.setFittable(false);
152 
153  if (fabs(p.getA()) > max_a_toFit)
154  p.setFittable(false);
155 
156  patterns.push_back(p);
157  }
158 }
159 
160 //----------------------------------------------------------------------------------------------------
161 
163 {
164  if (verbosity > 5)
165  LogVerbatim("TotemRPUVPatternFinder")
166  << ">> TotemRPUVPatternFinder::produce " << event.id().run() << ":" << event.id().event();
167 
168  // geometry
170  es.get<VeryForwardRealGeometryRecord>().get(geometry);
171  if (geometryWatcher.check(es))
172  lrcgn->resetGeometry(geometry.product());
173 
174  // get input
176  event.getByToken(detSetVectorTotemRPRecHitToken, input);
177 
178  // prepare output
179  DetSetVector<TotemRPUVPattern> patternsVector;
180 
181  // split input per RP and per U/V projection
182  struct RPData
183  {
184  DetSetVector<TotemRPRecHit> hits_U, hits_V;
185  map<uint8_t, uint16_t> planeOccupancy_U, planeOccupancy_V;
186  };
187  map<unsigned int, RPData> rpData;
188 
189  for (auto &ids : *input)
190  {
191  TotemRPDetId detId(ids.detId());
192  unsigned int plane = detId.plane();
193  bool uDir = detId.isStripsCoordinateUDirection();
194 
195  CTPPSDetId rpId = detId.getRPId();
196 
197  RPData &data = rpData[rpId];
198 
199  for (auto &h : ids)
200  {
201  if (uDir)
202  {
203  auto &ods = data.hits_U.find_or_insert(ids.detId());
204  ods.push_back(h);
205  data.planeOccupancy_U[plane]++;
206  } else {
207  auto &ods = data.hits_V.find_or_insert(ids.detId());
208  ods.push_back(h);
209  data.planeOccupancy_V[plane]++;
210  }
211  }
212  }
213 
214  // track recognition pot by pot
215  for (auto it : rpData)
216  {
217  CTPPSDetId rpId(it.first);
218  RPData &data = it.second;
219 
220  // merge default and exceptional settings (if available)
221  unsigned int minPlanesPerProjectionToFit_U = minPlanesPerProjectionToFit;
222  unsigned int minPlanesPerProjectionToFit_V = minPlanesPerProjectionToFit;
223  double threshold_U = threshold;
224  double threshold_V = threshold;
225 
226  auto setIt = exceptionalSettings.find(rpId);
227  if (setIt != exceptionalSettings.end())
228  {
229  minPlanesPerProjectionToFit_U = setIt->second.minPlanesPerProjectionToFit_U;
230  minPlanesPerProjectionToFit_V = setIt->second.minPlanesPerProjectionToFit_V;
231  threshold_U = setIt->second.threshold_U;
232  threshold_V = setIt->second.threshold_V;
233  }
234 
235  auto &uColl = data.planeOccupancy_U;
236  auto &vColl = data.planeOccupancy_V;
237 
238  if (verbosity > 5)
239  {
240  LogVerbatim("TotemRPUVPatternFinder")
241  << "\tRP " << rpId
242  << "\n\t\tall planes: u = " << uColl.size() << ", v = " << vColl.size();
243  }
244 
245  // count planes with clean data (no showers, noise, ...)
246  unsigned int uPlanes = 0, vPlanes = 0;
247  for (auto pit : uColl)
248  if (pit.second <= maxHitsPerPlaneToSearch)
249  uPlanes++;
250 
251  for (auto pit : vColl)
252  if (pit.second <= maxHitsPerPlaneToSearch)
253  vPlanes++;
254 
255  if (verbosity > 5)
256  LogVerbatim("TotemRPUVPatternFinder") << "\t\tplanes with clean data: u = " << uPlanes << ", v = " << vPlanes;
257 
258  // discard RPs with too few reasonable planes
260  continue;
261 
262  // prepare data containers
263  DetSet<TotemRPUVPattern> &patterns = patternsVector.find_or_insert(rpId);
264 
265  // "typical" z0 for the RP
266  double z0 = geometry->GetRPDevice(rpId)->translation().z();
267 
268  // u then v recognition
269  recognizeAndSelect(TotemRPUVPattern::projU, z0, threshold_U, minPlanesPerProjectionToFit_U, data.hits_U, patterns);
270 
271  recognizeAndSelect(TotemRPUVPattern::projV, z0, threshold_V, minPlanesPerProjectionToFit_V, data.hits_V, patterns);
272 
273  if (verbosity > 5)
274  {
275  LogVerbatim("TotemRPUVPatternFinder") << "\t\tpatterns:";
276  for (const auto &p : patterns)
277  {
278  unsigned int n_hits = 0;
279  for (auto &hds : p.getHits())
280  n_hits += hds.size();
281 
282  LogVerbatim("TotemRPUVPatternFinder")
283  << "\t\t\tproj = " << ((p.getProjection() == TotemRPUVPattern::projU) ? "U" : "V")
284  << ", a = " << p.getA()
285  << ", b = " << p.getB()
286  << ", w = " << p.getW()
287  << ", fittable = " << p.getFittable()
288  << ", hits = " << n_hits;
289  }
290  }
291  }
292 
293  // save output
294  event.put(make_unique<DetSetVector<TotemRPUVPattern>>(patternsVector));
295 }
296 
297 //----------------------------------------------------------------------------------------------------
298 
Detector ID class for TOTEM Si strip detectors.
Definition: TotemRPDetId.h:30
T getParameter(std::string const &) const
TotemRPUVPatternFinder(const edm::ParameterSet &conf)
DetGeomDesc * GetRPDevice(unsigned int id) const
returns geometry of a RP box
void push_back(const T &t)
Definition: DetSet.h:68
FastLineRecognition * lrcgn
the line recognition algorithm
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void resetGeometry(const TotemRPGeometry *_g)
uint32_t plane() const
Definition: TotemRPDetId.h:49
void getPatterns(const edm::DetSetVector< TotemRPRecHit > &input, double _z0, double threshold, edm::DetSet< TotemRPUVPattern > &patterns)
Event setup record containing the real (actual) geometry information.
static std::string const input
Definition: EdmProvDump.cc:44
const PhiMemoryImage patterns[9]
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:254
block of (exceptional) settings for 1 RP
std::map< unsigned int, RPSettings > exceptionalSettings
exceptional settings: RP Id –> settings
double threshold
minimal weight of (Hough) cluster to accept it as candidate
void recognizeAndSelect(TotemRPUVPattern::ProjectionType proj, double z0, double threshold, unsigned int planes_required, const edm::DetSetVector< TotemRPRecHit > &hits, edm::DetSet< TotemRPUVPattern > &patterns)
executes line recognition in a projection
unsigned char minPlanesPerProjectionToFit
minimal required number of active planes per projection to mark track candidate as fittable ...
DDTranslation translation() const
Definition: DetGeomDesc.h:87
CTPPSDetId getRPId() const
Definition: CTPPSDetId.h:97
edm::EDGetTokenT< edm::DetSetVector< TotemRPRecHit > > detSetVectorTotemRPRecHitToken
Class performing optimized hough transform to recognize lines.
edm::ESWatcher< VeryForwardRealGeometryRecord > geometryWatcher
virtual void produce(edm::Event &e, const edm::EventSetup &c) override
const T & get() const
Definition: EventSetup.h:56
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
unsigned char minPlanesPerProjectionToSearch
minimal required number of active planes per projection to even start track recognition ...
Class to recognize straight line tracks, based on optimized Hough trasform.
ESHandle< TrackerGeometry > geometry
HLT enums.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
unsigned int maxHitsPerPlaneToSearch
above this limit, planes are considered noisy
T const * product() const
Definition: ESHandle.h:86
double max_a_toFit
maximal angle (in any projection) to mark candidate as fittable - controls track parallelity ...
Definition: event.py:1