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 public:
39 
40  ~TotemRPUVPatternFinder() override;
41 
42  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 {
73  };
74 
76  std::map<unsigned int, RPSettings> exceptionalSettings;
77 
79 
82  double z0,
83  double threshold,
84  unsigned int planes_required,
87 };
88 
89 //----------------------------------------------------------------------------------------------------
90 
91 using namespace std;
92 using namespace edm;
93 
94 //----------------------------------------------------------------------------------------------------
95 
97  : tagRecHit(conf.getParameter<edm::InputTag>("tagRecHit")),
98  verbosity(conf.getUntrackedParameter<unsigned int>("verbosity", 0)),
99  minPlanesPerProjectionToSearch(conf.getParameter<unsigned int>("minPlanesPerProjectionToSearch")),
100  minPlanesPerProjectionToFit(conf.getParameter<unsigned int>("minPlanesPerProjectionToFit")),
101  maxHitsPerPlaneToSearch(conf.getParameter<unsigned int>("maxHitsPerPlaneToSearch")),
102  lrcgn(new FastLineRecognition(conf.getParameter<double>("clusterSize_a"),
103  conf.getParameter<double>("clusterSize_b"))),
104  threshold(conf.getParameter<double>("threshold")),
105  max_a_toFit(conf.getParameter<double>("max_a_toFit")) {
106  for (const auto &ps : conf.getParameter<vector<ParameterSet>>("exceptionalSettings")) {
107  unsigned int rpId = ps.getParameter<unsigned int>("rpId");
108 
109  RPSettings settings;
110  settings.minPlanesPerProjectionToFit_U = ps.getParameter<unsigned int>("minPlanesPerProjectionToFit_U");
111  settings.minPlanesPerProjectionToFit_V = ps.getParameter<unsigned int>("minPlanesPerProjectionToFit_V");
112  settings.threshold_U = ps.getParameter<double>("threshold_U");
113  settings.threshold_V = ps.getParameter<double>("threshold_V");
114 
115  exceptionalSettings[rpId] = settings;
116  }
117 
118  detSetVectorTotemRPRecHitToken = consumes<edm::DetSetVector<TotemRPRecHit>>(tagRecHit);
119 
120  produces<DetSetVector<TotemRPUVPattern>>();
121 }
122 
123 //----------------------------------------------------------------------------------------------------
124 
126 
127 //----------------------------------------------------------------------------------------------------
128 
130  double z0,
131  double threshold_loc,
132  unsigned int planes_required,
134  DetSet<TotemRPUVPattern> &patterns) {
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  p.setProjection(proj);
142 
143  p.setFittable(true);
144 
145  set<unsigned int> planes;
146  for (const auto &ds : p.hits())
147  planes.insert(TotemRPDetId(ds.detId()).plane());
148 
149  if (planes.size() < planes_required)
150  p.setFittable(false);
151 
152  if (fabs(p.a()) > max_a_toFit)
153  p.setFittable(false);
154 
155  patterns.push_back(p);
156  }
157 }
158 
159 //----------------------------------------------------------------------------------------------------
160 
162  if (verbosity > 5)
163  LogVerbatim("TotemRPUVPatternFinder")
164  << ">> TotemRPUVPatternFinder::produce " << event.id().run() << ":" << event.id().event();
165 
166  // geometry
168  es.get<VeryForwardRealGeometryRecord>().get(geometry);
169  if (geometryWatcher.check(es))
170  lrcgn->resetGeometry(geometry.product());
171 
172  // get input
174  event.getByToken(detSetVectorTotemRPRecHitToken, input);
175 
176  // prepare output
177  DetSetVector<TotemRPUVPattern> patternsVector;
178 
179  // split input per RP and per U/V projection
180  struct RPData {
181  DetSetVector<TotemRPRecHit> hits_U, hits_V;
182  map<uint8_t, uint16_t> planeOccupancy_U, planeOccupancy_V;
183  };
184  map<unsigned int, RPData> rpData;
185 
186  for (auto &ids : *input) {
187  TotemRPDetId detId(ids.detId());
188  unsigned int plane = detId.plane();
189  bool uDir = detId.isStripsCoordinateUDirection();
190 
191  CTPPSDetId rpId = detId.rpId();
192 
193  RPData &data = rpData[rpId];
194 
195  for (auto &h : ids) {
196  if (uDir) {
197  auto &ods = data.hits_U.find_or_insert(ids.detId());
198  ods.push_back(h);
199  data.planeOccupancy_U[plane]++;
200  } else {
201  auto &ods = data.hits_V.find_or_insert(ids.detId());
202  ods.push_back(h);
203  data.planeOccupancy_V[plane]++;
204  }
205  }
206  }
207 
208  // track recognition pot by pot
209  for (auto it : rpData) {
210  CTPPSDetId rpId(it.first);
211  RPData &data = it.second;
212 
213  // merge default and exceptional settings (if available)
214  unsigned int minPlanesPerProjectionToFit_U = minPlanesPerProjectionToFit;
215  unsigned int minPlanesPerProjectionToFit_V = minPlanesPerProjectionToFit;
216  double threshold_U = threshold;
217  double threshold_V = threshold;
218 
219  auto setIt = exceptionalSettings.find(rpId);
220  if (setIt != exceptionalSettings.end()) {
221  minPlanesPerProjectionToFit_U = setIt->second.minPlanesPerProjectionToFit_U;
222  minPlanesPerProjectionToFit_V = setIt->second.minPlanesPerProjectionToFit_V;
223  threshold_U = setIt->second.threshold_U;
224  threshold_V = setIt->second.threshold_V;
225  }
226 
227  auto &uColl = data.planeOccupancy_U;
228  auto &vColl = data.planeOccupancy_V;
229 
230  if (verbosity > 5) {
231  LogVerbatim("TotemRPUVPatternFinder")
232  << "\tRP " << rpId << "\n\t\tall planes: u = " << uColl.size() << ", v = " << vColl.size();
233  }
234 
235  // count planes with clean data (no showers, noise, ...)
236  unsigned int uPlanes = 0, vPlanes = 0;
237  for (auto pit : uColl)
238  if (pit.second <= maxHitsPerPlaneToSearch)
239  uPlanes++;
240 
241  for (auto pit : vColl)
242  if (pit.second <= maxHitsPerPlaneToSearch)
243  vPlanes++;
244 
245  if (verbosity > 5)
246  LogVerbatim("TotemRPUVPatternFinder") << "\t\tplanes with clean data: u = " << uPlanes << ", v = " << vPlanes;
247 
248  // discard RPs with too few reasonable planes
250  continue;
251 
252  // prepare data containers
253  DetSet<TotemRPUVPattern> &patterns = patternsVector.find_or_insert(rpId);
254 
255  // "typical" z0 for the RP
256  double z0 = geometry->rp(rpId)->translation().z();
257 
258  // u then v recognition
259  recognizeAndSelect(TotemRPUVPattern::projU, z0, threshold_U, minPlanesPerProjectionToFit_U, data.hits_U, patterns);
260 
261  recognizeAndSelect(TotemRPUVPattern::projV, z0, threshold_V, minPlanesPerProjectionToFit_V, data.hits_V, patterns);
262 
263  if (verbosity > 5) {
264  LogVerbatim("TotemRPUVPatternFinder") << "\t\tpatterns:";
265  for (const auto &p : patterns) {
266  unsigned int n_hits = 0;
267  for (auto &hds : p.hits())
268  n_hits += hds.size();
269 
270  LogVerbatim("TotemRPUVPatternFinder")
271  << "\t\t\tproj = " << ((p.projection() == TotemRPUVPattern::projU) ? "U" : "V") << ", a = " << p.a()
272  << ", b = " << p.b() << ", w = " << p.w() << ", fittable = " << p.fittable() << ", hits = " << n_hits;
273  }
274  }
275  }
276 
277  // save output
278  event.put(make_unique<DetSetVector<TotemRPUVPattern>>(patternsVector));
279 }
280 
281 //----------------------------------------------------------------------------------------------------
282 
285 
286  desc.add<edm::InputTag>("tagRecHit", edm::InputTag("totemRPRecHitProducer"))
287  ->setComment("input rechits collection to retrieve");
288  desc.addUntracked<unsigned int>("verbosity", 0);
289  desc.add<unsigned int>("maxHitsPerPlaneToSearch", 5)
290  ->setComment("minimum threshold of hits multiplicity to flag the pattern as dirty");
291  desc.add<unsigned int>("minPlanesPerProjectionToSearch", 3)
292  ->setComment(
293  "minimal number of reasonable (= not empty and not dirty) planes per projection and per RP, to start the "
294  "pattern search");
295  desc.add<double>("clusterSize_a", 0.02 /* rad */)
296  ->setComment("(full) cluster size (in rad) in slope-intercept space");
297  desc.add<double>("clusterSize_b", 0.3 /* mm */);
298 
299  desc.add<double>("threshold", 2.99)
300  ->setComment(
301  "minimal weight of (Hough) cluster to accept it as candidate\n"
302  " weight of cluster = sum of weights of contributing points\n"
303  " weight of point = sigma0 / sigma_of_point\n"
304  "most often: weight of point ~ 1, thus cluster weight is roughly number of contributing points");
305 
306  desc.add<unsigned int>("minPlanesPerProjectionToFit", 3)
307  ->setComment(
308  "minimal number of planes (in the recognised patterns) per projection and per RP, to tag the candidate as "
309  "fittable");
310 
311  desc.add<bool>("allowAmbiguousCombination", false)
312  ->setComment(
313  "whether to allow combination of most significant U and V pattern, in case there several of them.\n"
314  "don't set it to True, unless you have reason");
315 
316  desc.add<double>("max_a_toFit", 10.0)
317  ->setComment(
318  "maximal angle (in any projection) to mark the candidate as fittable -> controls track parallelity with "
319  "beam\n"
320  "huge value -> no constraint");
321 
322  edm::ParameterSetDescription exceptions_validator;
323  exceptions_validator.add<unsigned int>("rpId")->setComment("RP id according to CTPPSDetId");
324  exceptions_validator.add<unsigned int>("minPlanesPerProjectionToFit_U");
325  exceptions_validator.add<unsigned int>("minPlanesPerProjectionToFit_V");
326  exceptions_validator.add<double>("threshold_U");
327  exceptions_validator.add<double>("threshold_V");
328 
329  std::vector<edm::ParameterSet> exceptions_default;
330  desc.addVPSet("exceptionalSettings", exceptions_validator, exceptions_default);
331 
332  descr.add("totemRPUVPatternFinder", desc);
333 }
334 
335 //----------------------------------------------------------------------------------------------------
336 
Detector ID class for TOTEM Si strip detectors.
Definition: TotemRPDetId.h:30
T getParameter(std::string const &) const
void setComment(std::string const &value)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
TotemRPUVPatternFinder(const edm::ParameterSet &conf)
Translation translation() const
Definition: DetGeomDesc.h:65
void push_back(const T &t)
Definition: DetSet.h:67
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
FastLineRecognition * lrcgn
the line recognition algorithm
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const DetGeomDesc * rp(unsigned int id) const
returns geometry of a RP box
uint32_t plane() const
Definition: TotemRPDetId.h:46
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:48
reference find_or_insert(det_id_type id)
Definition: DetSetVector.h:234
CTPPSDetId rpId() const
Definition: CTPPSDetId.h:78
block of (exceptional) settings for 1 RP
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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 ...
void resetGeometry(const CTPPSGeometry *_g)
edm::EDGetTokenT< edm::DetSetVector< TotemRPRecHit > > detSetVectorTotemRPRecHitToken
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Class performing optimized hough transform to recognize lines.
edm::ESWatcher< VeryForwardRealGeometryRecord > geometryWatcher
void produce(edm::Event &e, const edm::EventSetup &c) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
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:79
static void fillDescriptions(edm::ConfigurationDescriptions &)
T get() const
Definition: EventSetup.h:73
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