CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
trklet::PurgeDuplicate Class Reference

#include <PurgeDuplicate.h>

Inheritance diagram for trklet::PurgeDuplicate:
trklet::ProcessBase

Public Member Functions

void addInput (MemoryBase *memory, std::string input) override
 
void addOutput (MemoryBase *memory, std::string output) override
 
void execute (std::vector< Track > &outputtracks_, unsigned int iSector)
 
 PurgeDuplicate (std::string name, Settings const &settings, Globals *global)
 
 ~PurgeDuplicate () override=default
 
- Public Member Functions inherited from trklet::ProcessBase
unsigned int getISeed (const std::string &name)
 
std::string const & getName () const
 
void initLayerDisk (unsigned int pos, int &layer, int &disk)
 
void initLayerDisk (unsigned int pos, int &layer, int &disk, int &layerdisk)
 
unsigned int initLayerDisk (unsigned int pos)
 
void initLayerDisksandISeed (unsigned int &layerdisk1, unsigned int &layerdisk2, unsigned int &iSeed)
 
unsigned int nbits (unsigned int power)
 
 ProcessBase (std::string name, Settings const &settings, Globals *global)
 
virtual ~ProcessBase ()=default
 

Private Member Functions

std::pair< int, int > findLayerDisk (const Stub *) const
 
std::vector< unsigned int > findOverlapRInvBins (const Tracklet *trk) const
 
unsigned int findVarRInvBin (const Tracklet *trk) const
 
std::vector< double > getInventedCoords (unsigned int, const Stub *, const Tracklet *) const
 
std::vector< double > getInventedCoordsExtended (unsigned int, const Stub *, const Tracklet *) const
 
std::vector< const Stub * > getInventedSeedingStub (unsigned int, const Tracklet *, const std::vector< const Stub *> &) const
 
double getPhiRes (Tracklet *curTracklet, const Stub *curStub) const
 
bool isSeedingStub (int, int, int) const
 
bool isTrackInBin (const std::vector< unsigned int > &vec, unsigned int num) const
 
std::string l1tinfo (const L1TStub *, std::string) const
 

Private Attributes

std::vector< std::vector< std::pair< int, int > > > inputstubidslists_
 
std::vector< std::vector< const Stub * > > inputstublists_
 
std::vector< TrackFitMemory * > inputtrackfits_
 
std::vector< Tracklet * > inputtracklets_
 
std::vector< Track * > inputtracks_
 
std::vector< std::vector< std::pair< int, int > > > mergedstubidslists_
 
std::vector< CleanTrackMemory * > outputtracklets_
 

Additional Inherited Members

- Protected Attributes inherited from trklet::ProcessBase
Globalsglobals_
 
std::string name_
 
Settings const & settings_
 

Detailed Description

Definition at line 27 of file PurgeDuplicate.h.

Constructor & Destructor Documentation

◆ PurgeDuplicate()

PurgeDuplicate::PurgeDuplicate ( std::string  name,
Settings const &  settings,
Globals global 
)

Definition at line 25 of file PurgeDuplicate.cc.

26  : ProcessBase(name, settings, global) {}
ProcessBase(std::string name, Settings const &settings, Globals *global)
Definition: ProcessBase.cc:14

◆ ~PurgeDuplicate()

trklet::PurgeDuplicate::~PurgeDuplicate ( )
overridedefault

Member Function Documentation

◆ addInput()

void PurgeDuplicate::addInput ( MemoryBase memory,
std::string  input 
)
overridevirtual

Implements trklet::ProcessBase.

Definition at line 54 of file PurgeDuplicate.cc.

References cms::cuda::assert(), Exception, input, PixelMapPlotter::inputs, inputtrackfits_, mps_setup::memory, trklet::ProcessBase::name_, trklet::ProcessBase::settings_, createJobs::tmp, and trklet::Settings::writetrace().

54  {
55  if (settings_.writetrace()) {
56  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
57  << input;
58  }
59  unordered_set<string> inputs = {"trackin",
60  "trackin1",
61  "trackin2",
62  "trackin3",
63  "trackin4",
64  "trackin5",
65  "trackin6",
66  "trackin7",
67  "trackin8",
68  "trackin9",
69  "trackin10",
70  "trackin11",
71  "trackin12"};
72  if (inputs.find(input) != inputs.end()) {
73  auto* tmp = dynamic_cast<TrackFitMemory*>(memory);
74  assert(tmp != nullptr);
75  inputtrackfits_.push_back(tmp);
76  return;
77  }
78  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find input: " << input;
79 }
Log< level::Info, true > LogVerbatim
std::string name_
Definition: ProcessBase.h:38
Settings const & settings_
Definition: ProcessBase.h:40
bool writetrace() const
Definition: Settings.h:193
assert(be >=bs)
static std::string const input
Definition: EdmProvDump.cc:50
std::vector< TrackFitMemory * > inputtrackfits_
tmp
align.sh
Definition: createJobs.py:716

◆ addOutput()

void PurgeDuplicate::addOutput ( MemoryBase memory,
std::string  output 
)
overridevirtual

Implements trklet::ProcessBase.

Definition at line 28 of file PurgeDuplicate.cc.

References cms::cuda::assert(), Exception, mps_setup::memory, trklet::ProcessBase::name_, convertSQLitetoXML_cfg::output, PatBasicFWLiteJetAnalyzer_Selector_cfg::outputs, outputtracklets_, trklet::ProcessBase::settings_, createJobs::tmp, and trklet::Settings::writetrace().

28  {
29  if (settings_.writetrace()) {
30  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
31  << output;
32  }
33  unordered_set<string> outputs = {"trackout",
34  "trackout1",
35  "trackout2",
36  "trackout3",
37  "trackout4",
38  "trackout5",
39  "trackout6",
40  "trackout7",
41  "trackout8",
42  "trackout9",
43  "trackout10",
44  "trackout11"};
45  if (outputs.find(output) != outputs.end()) {
46  auto* tmp = dynamic_cast<CleanTrackMemory*>(memory);
47  assert(tmp != nullptr);
48  outputtracklets_.push_back(tmp);
49  return;
50  }
51  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find output: " << output;
52 }
Log< level::Info, true > LogVerbatim
std::string name_
Definition: ProcessBase.h:38
Settings const & settings_
Definition: ProcessBase.h:40
bool writetrace() const
Definition: Settings.h:193
std::vector< CleanTrackMemory * > outputtracklets_
assert(be >=bs)
tmp
align.sh
Definition: createJobs.py:716

◆ execute()

void PurgeDuplicate::execute ( std::vector< Track > &  outputtracks_,
unsigned int  iSector 
)

Definition at line 81 of file PurgeDuplicate.cc.

References Reference_intrackfit_cff::barrel, newFWLiteAna::bin, trklet::Settings::extended(), spr::find(), findOverlapRInvBins(), findVarRInvBin(), dqmdumpme::first, trklet::HybridFit::Fit(), trklet::Tracklet::fit(), groupFilesInBlocks::fout, getInventedSeedingStub(), getPhiRes(), getTrack(), trklet::Tracklet::getTrack(), trklet::ProcessBase::globals_, getRunAppsInfo::grid, mps_fire::i, trklet::Track::ichisq(), inputstubidslists_, inputstublists_, inputtrackfits_, inputtracklets_, inputtracks_, createfilelist::int, trklet::Settings::inventStubs(), isTrackInBin(), dqmiolumiharvest::j, M_PI, SiStripPI::max, trklet::Settings::maxStep(), trklet::Settings::mergeComparison(), mergedstubidslists_, SiStripPI::min, trklet::Settings::minIndStubs(), trklet::ProcessBase::name_, trklet::Settings::numTracksComparedPerBin(), outputtracklets_, BeamMonitor_cff::phiBin, muonRecoAnalyzer_cfi::ptBin, trklet::Settings::removalType(), edm::second(), trklet::Tracklet::seedIndex(), trklet::Track::setSector(), trklet::Track::setStubIDprefit(), trklet::Track::setStubIDpremerge(), trklet::ProcessBase::settings_, findQualityFiles::size, jetUpdater_cfi::sort, trklet::Track::stubID(), ecalIM::uniqueIndex(), trklet::Settings::varRInvBins(), and trklet::Settings::writeMonitorData().

81  {
82  inputtracklets_.clear();
83  inputtracks_.clear();
84 
85  inputstubidslists_.clear();
86  inputstublists_.clear();
87  mergedstubidslists_.clear();
88 
89  if (settings_.removalType() != "merge") {
90  for (auto& inputtrackfit : inputtrackfits_) {
91  if (inputtrackfit->nTracks() == 0)
92  continue;
93  for (unsigned int j = 0; j < inputtrackfit->nTracks(); j++) {
94  Track* aTrack = inputtrackfit->getTrack(j)->getTrack();
95  aTrack->setSector(iSector);
96  inputtracks_.push_back(aTrack);
97  }
98  }
99  if (inputtracks_.empty())
100  return;
101  }
102 
103  unsigned int numTrk = inputtracks_.size();
104 
106  // Hybrid Removal //
108 #ifdef USEHYBRID
109 
110  if (settings_.removalType() == "merge") {
111  // Track seed & duplicate flag
112  std::vector<std::pair<int, bool>> trackInfo;
113  // Flag for tracks in multiple bins that get merged but are not in the correct variable bin
114  std::vector<bool> trackBinInfo;
115  // Vector to store the relative rank of the track candidate for merging, based on seed type
116  std::vector<int> seedRank;
117 
118  // Stubs on every track
119  std::vector<std::vector<const Stub*>> inputstublistsall;
120  // (layer, unique stub index within layer) of each stub on every track
121  std::vector<std::vector<std::pair<int, int>>> mergedstubidslistsall;
122  std::vector<std::vector<std::pair<int, int>>> inputstubidslistsall;
123  std::vector<Tracklet*> inputtrackletsall;
124 
125  std::vector<unsigned int> prefTracks; // Stores all the tracks that are sent to the KF from each bin
126  std::vector<int> prefTrackFit; // Stores the track seed that corresponds to the associated track in prefTracks
127 
128  for (unsigned int bin = 0; bin < settings_.varRInvBins().size() - 1; bin++) {
129  // Get vectors from TrackFit and save them
130  // inputtracklets: Tracklet objects from the FitTrack (not actually fit yet)
131  // inputstublists: L1Stubs for that track
132  // inputstubidslists: Stub stubIDs for that 3rack
133  // mergedstubidslists: the same as inputstubidslists, but will be used during duplicate removal
134  for (unsigned int i = 0; i < inputtrackfits_.size(); i++) {
135  if (inputtrackfits_[i]->nStublists() == 0)
136  continue;
137  if (inputtrackfits_[i]->nStublists() != inputtrackfits_[i]->nTracks())
138  throw "Number of stublists and tracks don't match up!";
139  for (unsigned int j = 0; j < inputtrackfits_[i]->nStublists(); j++) {
141  if (inputtracklets_.size() >= settings_.maxStep("DR"))
142  continue;
143  Tracklet* aTrack = inputtrackfits_[i]->getTrack(j);
145  std::vector<const Stub*> stublist = inputtrackfits_[i]->getStublist(j);
146  inputstublists_.push_back(stublist);
147  std::vector<std::pair<int, int>> stubidslist = inputtrackfits_[i]->getStubidslist(j);
148  inputstubidslists_.push_back(stubidslist);
149  mergedstubidslists_.push_back(stubidslist);
150 
151  // Encoding: L1L2=0, L2L3=1, L3L4=2, L5L6=3, D1D2=4, D3D4=5, L1D1=6, L2D1=7
152  // Best Guess: L1L2 > L1D1 > L2L3 > L2D1 > D1D2 > L3L4 > L5L6 > D3D4
153  // Best Rank: L1L2 > L3L4 > D3D4 > D1D2 > L2L3 > L2D1 > L5L6 > L1D1
154  // Rank-Informed Guess: L1L2 > L3L4 > L1D1 > L2L3 > L2D1 > D1D2 > L5L6 > D3D4
155  unsigned int curSeed = aTrack->seedIndex();
156  std::vector<int> ranks{1, 5, 2, 7, 4, 3, 8, 6};
157  if (settings_.extended())
158  seedRank.push_back(9);
159  else
160  seedRank.push_back(ranks[curSeed]);
161 
162  if (stublist.size() != stubidslist.size())
163  throw "Number of stubs and stubids don't match up!";
164 
165  trackInfo.emplace_back(i, false);
166  trackBinInfo.emplace_back(false);
167  } else
168  continue;
169  }
170  }
171 
172  if (inputtracklets_.empty())
173  continue;
174  const unsigned int numStublists = inputstublists_.size();
175 
176  if (settings_.inventStubs()) {
177  for (unsigned int itrk = 0; itrk < numStublists; itrk++) {
179  }
180  }
181 
182  // Initialize all-false 2D array of tracks being duplicates to other tracks
183  bool dupMap[numStublists][numStublists]; // Ends up symmetric
184  for (unsigned int itrk = 0; itrk < numStublists; itrk++) {
185  for (unsigned int jtrk = 0; jtrk < numStublists; jtrk++) {
186  dupMap[itrk][jtrk] = false;
187  }
188  }
189 
190  // Used to check if a track is in two bins, is not a duplicate in either bin, so is sent out twice
191  bool noMerge[numStublists];
192  for (unsigned int itrk = 0; itrk < numStublists; itrk++) {
193  noMerge[itrk] = false;
194  }
195 
196  // Find duplicates; Fill dupMap by looping over all pairs of "tracks"
197  // numStublists-1 since last track has no other to compare to
198  for (unsigned int itrk = 0; itrk < numStublists - 1; itrk++) {
199  for (unsigned int jtrk = itrk + 1; jtrk < numStublists; jtrk++) {
200  if (itrk >= settings_.numTracksComparedPerBin())
201  continue;
202  // Get primary track stubids = (layer, unique stub index within layer)
203  const std::vector<std::pair<int, int>>& stubsTrk1 = inputstubidslists_[itrk];
204 
205  // Get and count secondary track stubids
206  const std::vector<std::pair<int, int>>& stubsTrk2 = inputstubidslists_[jtrk];
207 
208  // Count number of layers that share stubs, and the number of UR that each track hits
209  unsigned int nShareLay = 0;
210  if (settings_.mergeComparison() == "CompareAll") {
211  bool layerArr[16];
212  for (auto& i : layerArr) {
213  i = false;
214  };
215  for (const auto& st1 : stubsTrk1) {
216  for (const auto& st2 : stubsTrk2) {
217  if (st1.first == st2.first && st1.second == st2.second) { // tracks share stub
218  // Converts layer/disk encoded in st1->first to an index in the layer array
219  int i = st1.first; // layer/disk
220  bool barrel = (i > 0 && i < 10);
221  bool endcapA = (i > 10);
222  bool endcapB = (i < 0);
223  int lay = barrel * (i - 1) + endcapA * (i - 5) - endcapB * i; // encode in range 0-15
224  if (!layerArr[lay]) {
225  nShareLay++;
226  layerArr[lay] = true;
227  }
228  }
229  }
230  }
231  } else if (settings_.mergeComparison() == "CompareBest") {
232  std::vector<const Stub*> fullStubslistsTrk1 = inputstublists_[itrk];
233  std::vector<const Stub*> fullStubslistsTrk2 = inputstublists_[jtrk];
234 
235  // Arrays to store the index of the best stub in each layer
236  int layStubidsTrk1[16];
237  int layStubidsTrk2[16];
238  for (int i = 0; i < 16; i++) {
239  layStubidsTrk1[i] = -1;
240  layStubidsTrk2[i] = -1;
241  }
242  // For each stub on the first track, find the stub with the best residual and store its index in the layStubidsTrk1 array
243  for (unsigned int stcount = 0; stcount < stubsTrk1.size(); stcount++) {
244  int i = stubsTrk1[stcount].first; // layer/disk
245  bool barrel = (i > 0 && i < 10);
246  bool endcapA = (i > 10);
247  bool endcapB = (i < 0);
248  int lay = barrel * (i - 1) + endcapA * (i - 5) - endcapB * i; // encode in range 0-15
249  double nres = getPhiRes(inputtracklets_[itrk], fullStubslistsTrk1[stcount]);
250  double ores = 0;
251  if (layStubidsTrk1[lay] != -1)
252  ores = getPhiRes(inputtracklets_[itrk], fullStubslistsTrk1[layStubidsTrk1[lay]]);
253  if (layStubidsTrk1[lay] == -1 || nres < ores) {
254  layStubidsTrk1[lay] = stcount;
255  }
256  }
257  // For each stub on the second track, find the stub with the best residual and store its index in the layStubidsTrk1 array
258  for (unsigned int stcount = 0; stcount < stubsTrk2.size(); stcount++) {
259  int i = stubsTrk2[stcount].first; // layer/disk
260  bool barrel = (i > 0 && i < 10);
261  bool endcapA = (i > 10);
262  bool endcapB = (i < 0);
263  int lay = barrel * (i - 1) + endcapA * (i - 5) - endcapB * i; // encode in range 0-15
264  double nres = getPhiRes(inputtracklets_[jtrk], fullStubslistsTrk2[stcount]);
265  double ores = 0;
266  if (layStubidsTrk2[lay] != -1)
267  ores = getPhiRes(inputtracklets_[jtrk], fullStubslistsTrk2[layStubidsTrk2[lay]]);
268  if (layStubidsTrk2[lay] == -1 || nres < ores) {
269  layStubidsTrk2[lay] = stcount;
270  }
271  }
272  // For all 16 layers (6 layers and 10 disks), count the number of layers who's best stub on both tracks are the same
273  for (int i = 0; i < 16; i++) {
274  int t1i = layStubidsTrk1[i];
275  int t2i = layStubidsTrk2[i];
276  if (t1i != -1 && t2i != -1 && stubsTrk1[t1i].first == stubsTrk2[t2i].first &&
277  stubsTrk1[t1i].second == stubsTrk2[t2i].second)
278  nShareLay++;
279  }
280  }
281  // Fill duplicate map
282  if (nShareLay >= settings_.minIndStubs()) { // For number of shared stub merge condition
283  dupMap[itrk][jtrk] = true;
284  dupMap[jtrk][itrk] = true;
285  }
286  }
287  }
288 
289  // Check to see if the track is a duplicate
290  for (unsigned int itrk = 0; itrk < numStublists; itrk++) {
291  for (unsigned int jtrk = 0; jtrk < numStublists; jtrk++) {
292  if (dupMap[itrk][jtrk]) {
293  noMerge[itrk] = true;
294  }
295  }
296  }
297 
298  // If the track isn't a duplicate, and if it's in more than one bin, and it is not in the proper varrinvbin, then mark it so it won't be sent to output
299  for (unsigned int itrk = 0; itrk < numStublists; itrk++) {
300  if (noMerge[itrk] == false) {
301  if ((findOverlapRInvBins(inputtracklets_[itrk]).size() > 1) &&
302  (findVarRInvBin(inputtracklets_[itrk]) != bin)) {
303  trackInfo[itrk].second = true;
304  }
305  }
306  }
307  // Merge duplicate tracks
308  for (unsigned int itrk = 0; itrk < numStublists - 1; itrk++) {
309  for (unsigned int jtrk = itrk + 1; jtrk < numStublists; jtrk++) {
310  // Merge a track with its first duplicate found.
311  if (dupMap[itrk][jtrk]) {
312  // Set preferred track based on seed rank
313  int preftrk;
314  int rejetrk;
315  if (seedRank[itrk] < seedRank[jtrk]) {
316  preftrk = itrk;
317  rejetrk = jtrk;
318  } else {
319  preftrk = jtrk;
320  rejetrk = itrk;
321  }
322 
323  // If the preffered track is in more than one bin, but not in the proper varrinvbin, then mark as true
324  if ((findOverlapRInvBins(inputtracklets_[preftrk]).size() > 1) &&
325  (findVarRInvBin(inputtracklets_[preftrk]) != bin)) {
326  trackBinInfo[preftrk] = true;
327  trackBinInfo[rejetrk] = true;
328  } else {
329  // Get a merged stub list
330  std::vector<const Stub*> newStubList;
331  std::vector<const Stub*> stubsTrk1 = inputstublists_[preftrk];
332  std::vector<const Stub*> stubsTrk2 = inputstublists_[rejetrk];
333  std::vector<unsigned int> stubsTrk1indices;
334  std::vector<unsigned int> stubsTrk2indices;
335  for (unsigned int stub1it = 0; stub1it < stubsTrk1.size(); stub1it++) {
336  stubsTrk1indices.push_back(stubsTrk1[stub1it]->l1tstub()->uniqueIndex());
337  }
338  for (unsigned int stub2it = 0; stub2it < stubsTrk2.size(); stub2it++) {
339  stubsTrk2indices.push_back(stubsTrk2[stub2it]->l1tstub()->uniqueIndex());
340  }
341  newStubList = stubsTrk1;
342  for (unsigned int stub2it = 0; stub2it < stubsTrk2.size(); stub2it++) {
343  if (find(stubsTrk1indices.begin(), stubsTrk1indices.end(), stubsTrk2indices[stub2it]) ==
344  stubsTrk1indices.end()) {
345  newStubList.push_back(stubsTrk2[stub2it]);
346  }
347  }
348  // Overwrite stublist of preferred track with merged list
349  inputstublists_[preftrk] = newStubList;
350 
351  std::vector<std::pair<int, int>> newStubidsList;
352  std::vector<std::pair<int, int>> stubidsTrk1 = mergedstubidslists_[preftrk];
353  std::vector<std::pair<int, int>> stubidsTrk2 = mergedstubidslists_[rejetrk];
354  newStubidsList = stubidsTrk1;
355 
356  for (unsigned int stub2it = 0; stub2it < stubsTrk2.size(); stub2it++) {
357  if (find(stubsTrk1indices.begin(), stubsTrk1indices.end(), stubsTrk2indices[stub2it]) ==
358  stubsTrk1indices.end()) {
359  newStubidsList.push_back(stubidsTrk2[stub2it]);
360  }
361  }
362  // Overwrite stubidslist of preferred track with merged list
363  mergedstubidslists_[preftrk] = newStubidsList;
364 
365  // Mark that rejected track has been merged into another track
366  trackInfo[rejetrk].second = true;
367  }
368  }
369  }
370  }
371 
372  for (unsigned int ktrk = 0; ktrk < numStublists; ktrk++) {
373  if ((trackInfo[ktrk].second != true) && (trackBinInfo[ktrk] != true)) {
374  prefTracks.push_back(ktrk);
375  prefTrackFit.push_back(trackInfo[ktrk].first);
376  inputtrackletsall.push_back(inputtracklets_[ktrk]);
377  inputstublistsall.push_back(inputstublists_[ktrk]);
378  inputstubidslistsall.push_back(inputstubidslists_[ktrk]);
379  mergedstubidslistsall.push_back(mergedstubidslists_[ktrk]);
380  }
381  }
382 
383  // Need to clear all the vectors which will be used in the next bin
384  seedRank.clear();
385  trackInfo.clear();
386  trackBinInfo.clear();
387  inputtracklets_.clear();
388  inputstublists_.clear();
389  inputstubidslists_.clear();
390  mergedstubidslists_.clear();
391  }
392 
393  // Make the final track objects, fit with KF, and send to output
394  for (unsigned int itrk = 0; itrk < prefTracks.size(); itrk++) {
395  Tracklet* tracklet = inputtrackletsall[itrk];
396  std::vector<const Stub*> trackstublist = inputstublistsall[itrk];
397 
398  // Run KF track fit
399  HybridFit hybridFitter(iSector, settings_, globals_);
400  hybridFitter.Fit(tracklet, trackstublist);
401 
402  // If the track was accepted (and thus fit), add to output
403  if (tracklet->fit()) {
404  // Add fitted Track to output (later converted to TTTrack)
405  Track* outtrack = tracklet->getTrack();
406  outtrack->setSector(iSector);
407  // Also store fitted track as more detailed Tracklet object.
408  outputtracklets_[prefTrackFit[itrk]]->addTrack(tracklet);
409 
410  // Add all tracks to standalone root file output
411  outtrack->setStubIDpremerge(inputstubidslistsall[itrk]);
412  outtrack->setStubIDprefit(mergedstubidslistsall[itrk]);
413  outputtracks_.push_back(*outtrack);
414  }
415  }
416  }
417 
418 #endif
419 
421  // Grid removal //
423  if (settings_.removalType() == "grid") {
424  // Sort tracks by ichisq/DoF so that removal will keep the lower ichisq/DoF track
425  std::sort(inputtracks_.begin(), inputtracks_.end(), [](const Track* lhs, const Track* rhs) {
426  return lhs->ichisq() / lhs->stubID().size() < rhs->ichisq() / rhs->stubID().size();
427  });
428  bool grid[35][40] = {{false}};
429 
430  for (unsigned int itrk = 0; itrk < numTrk; itrk++) {
431  if (inputtracks_[itrk]->duplicate())
432  edm::LogPrint("Tracklet") << "WARNING: Track already tagged as duplicate!!";
433 
434  double phiBin = (inputtracks_[itrk]->phi0(settings_) - 2 * M_PI / 27 * iSector) / (2 * M_PI / 9 / 50) + 9;
435  phiBin = std::max(phiBin, 0.);
436  phiBin = std::min(phiBin, 34.);
437 
438  double ptBin = 1 / inputtracks_[itrk]->pt(settings_) * 40 + 20;
439  ptBin = std::max(ptBin, 0.);
440  ptBin = std::min(ptBin, 39.);
441 
442  if (grid[(int)phiBin][(int)ptBin])
443  inputtracks_[itrk]->setDuplicate(true);
444  grid[(int)phiBin][(int)ptBin] = true;
445 
446  double phiTest = inputtracks_[itrk]->phi0(settings_) - 2 * M_PI / 27 * iSector;
447  if (phiTest < -2 * M_PI / 27)
448  edm::LogVerbatim("Tracklet") << "track phi too small!";
449  if (phiTest > 2 * 2 * M_PI / 27)
450  edm::LogVerbatim("Tracklet") << "track phi too big!";
451  }
452  } // end grid removal
453 
455  // ichi + nstub removal //
457  if (settings_.removalType() == "ichi" || settings_.removalType() == "nstub") {
458  for (unsigned int itrk = 0; itrk < numTrk - 1; itrk++) { // numTrk-1 since last track has no other to compare to
459 
460  // If primary track is a duplicate, it cannot veto any...move on
461  if (inputtracks_[itrk]->duplicate() == 1)
462  continue;
463 
464  unsigned int nStubP = 0;
465  vector<unsigned int> nStubS(numTrk);
466  vector<unsigned int> nShare(numTrk);
467  // Get and count primary stubs
468  std::map<int, int> stubsTrk1 = inputtracks_[itrk]->stubID();
469  nStubP = stubsTrk1.size();
470 
471  for (unsigned int jtrk = itrk + 1; jtrk < numTrk; jtrk++) {
472  // Skip duplicate tracks
473  if (inputtracks_[jtrk]->duplicate() == 1)
474  continue;
475 
476  // Get and count secondary stubs
477  std::map<int, int> stubsTrk2 = inputtracks_[jtrk]->stubID();
478  nStubS[jtrk] = stubsTrk2.size();
479 
480  // Count shared stubs
481  for (auto& st : stubsTrk1) {
482  if (stubsTrk2.find(st.first) != stubsTrk2.end()) {
483  if (st.second == stubsTrk2[st.first])
484  nShare[jtrk]++;
485  }
486  }
487  }
488 
489  // Tag duplicates
490  for (unsigned int jtrk = itrk + 1; jtrk < numTrk; jtrk++) {
491  // Skip duplicate tracks
492  if (inputtracks_[jtrk]->duplicate() == 1)
493  continue;
494 
495  // Chi2 duplicate removal
496  if (settings_.removalType() == "ichi") {
497  if ((nStubP - nShare[jtrk] < settings_.minIndStubs()) ||
498  (nStubS[jtrk] - nShare[jtrk] < settings_.minIndStubs())) {
499  if ((int)inputtracks_[itrk]->ichisq() / (2 * inputtracks_[itrk]->stubID().size() - 4) >
500  (int)inputtracks_[jtrk]->ichisq() / (2 * inputtracks_[itrk]->stubID().size() - 4)) {
501  inputtracks_[itrk]->setDuplicate(true);
502  } else if ((int)inputtracks_[itrk]->ichisq() / (2 * inputtracks_[itrk]->stubID().size() - 4) <=
503  (int)inputtracks_[jtrk]->ichisq() / (2 * inputtracks_[itrk]->stubID().size() - 4)) {
504  inputtracks_[jtrk]->setDuplicate(true);
505  } else {
506  edm::LogVerbatim("Tracklet") << "Error: Didn't tag either track in duplicate pair.";
507  }
508  }
509  } // end ichi removal
510 
511  // nStub duplicate removal
512  if (settings_.removalType() == "nstub") {
513  if ((nStubP - nShare[jtrk] < settings_.minIndStubs()) && (nStubP < nStubS[jtrk])) {
514  inputtracks_[itrk]->setDuplicate(true);
515  } else if ((nStubS[jtrk] - nShare[jtrk] < settings_.minIndStubs()) && (nStubS[jtrk] <= nStubP)) {
516  inputtracks_[jtrk]->setDuplicate(true);
517  } else {
518  edm::LogVerbatim("Tracklet") << "Error: Didn't tag either track in duplicate pair.";
519  }
520  } // end nstub removal
521 
522  } // end tag duplicates
523 
524  } // end loop over primary track
525 
526  } // end ichi + nstub removal
527 
528  //Add tracks to output
529  if (settings_.removalType() != "merge") {
530  for (unsigned int i = 0; i < inputtrackfits_.size(); i++) {
531  for (unsigned int j = 0; j < inputtrackfits_[i]->nTracks(); j++) {
532  if (inputtrackfits_[i]->getTrack(j)->getTrack()->duplicate() == 0) {
533  if (settings_.writeMonitorData("Seeds")) {
534  ofstream fout("seeds.txt", ofstream::app);
535  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector << " "
536  << inputtrackfits_[i]->getTrack(j)->getISeed() << endl;
537  fout.close();
538  }
539  outputtracklets_[i]->addTrack(inputtrackfits_[i]->getTrack(j));
540  }
541  //For root file:
542  outputtracks_.push_back(*inputtrackfits_[i]->getTrack(j)->getTrack());
543  }
544  }
545  }
546 }
size
Write out results.
Log< level::Info, true > LogVerbatim
unsigned int seedIndex() const
Definition: Tracklet.h:224
unsigned int maxStep(std::string module) const
Definition: Settings.h:123
std::string name_
Definition: ProcessBase.h:38
std::vector< std::vector< std::pair< int, int > > > inputstubidslists_
Settings const & settings_
Definition: ProcessBase.h:40
Globals * globals_
Definition: ProcessBase.h:41
double getPhiRes(Tracklet *curTracklet, const Stub *curStub) const
std::vector< CleanTrackMemory * > outputtracklets_
std::vector< unsigned int > findOverlapRInvBins(const Tracklet *trk) const
unsigned int minIndStubs() const
Definition: Settings.h:248
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< Track * > inputtracks_
unsigned int findVarRInvBin(const Tracklet *trk) const
Track * getTrack()
Definition: Tracklet.h:192
bool fit() const
Definition: Tracklet.h:197
U second(std::pair< T, U > const &p)
unsigned int numTracksComparedPerBin() const
Definition: Settings.h:296
std::string removalType() const
Definition: Settings.h:249
std::vector< TrackFitMemory * > inputtrackfits_
bool writeMonitorData(std::string module) const
Definition: Settings.h:116
bool inventStubs() const
Definition: Settings.h:272
#define M_PI
bool isTrackInBin(const std::vector< unsigned int > &vec, unsigned int num) const
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
void setSector(int nsec)
Definition: Track.h:34
std::vector< Tracklet * > inputtracklets_
std::vector< std::vector< std::pair< int, int > > > mergedstubidslists_
std::string mergeComparison() const
Definition: Settings.h:250
bool extended() const
Definition: Settings.h:266
std::vector< const Stub * > getInventedSeedingStub(unsigned int, const Tracklet *, const std::vector< const Stub *> &) const
int ichisq() const
Definition: Track.h:40
std::vector< std::vector< const Stub * > > inputstublists_
int uniqueIndex(int eta, int phi)
const std::vector< double > varRInvBins() const
Definition: Settings.h:298
void setStubIDprefit(std::vector< std::pair< int, int >> stubIDprefit)
Definition: Track.h:36
const std::map< int, int > & stubID() const
Definition: Track.h:42
void setStubIDpremerge(std::vector< std::pair< int, int >> stubIDpremerge)
Definition: Track.h:35

◆ findLayerDisk()

std::pair< int, int > PurgeDuplicate::findLayerDisk ( const Stub st) const
private

Definition at line 592 of file PurgeDuplicate.cc.

References trklet::Stub::layerdisk(), and trklet::N_LAYER.

Referenced by getInventedCoords(), getInventedCoordsExtended(), and getInventedSeedingStub().

592  {
593  std::pair<int, int> layer_disk;
594  layer_disk.first = st->layerdisk() + 1;
595  if (layer_disk.first > N_LAYER) {
596  layer_disk.first = 0;
597  }
598  layer_disk.second = st->layerdisk() - (N_LAYER - 1);
599  if (layer_disk.second < 0) {
600  layer_disk.second = 0;
601  }
602  return layer_disk;
603 }
unsigned int layerdisk() const
Definition: Stub.cc:185
constexpr int N_LAYER
Definition: Settings.h:25

◆ findOverlapRInvBins()

std::vector< unsigned int > PurgeDuplicate::findOverlapRInvBins ( const Tracklet trk) const
private

Definition at line 756 of file PurgeDuplicate.cc.

References mps_fire::i, trklet::Settings::overlapSize(), trklet::Tracklet::rinv(), trklet::ProcessBase::settings_, and trklet::Settings::varRInvBins().

Referenced by execute().

756  {
757  double rInv = trk->rinv();
758  const double overlapSize = settings_.overlapSize();
759  const std::vector<double>& varRInvBins = settings_.varRInvBins();
760  std::vector<unsigned int> chosenBins;
761  for (long unsigned int i = 0; i < varRInvBins.size() - 1; i++) {
762  if ((rInv < varRInvBins[i + 1] + overlapSize) && (rInv > varRInvBins[i] - overlapSize)) {
763  chosenBins.push_back(i);
764  }
765  }
766  return chosenBins;
767 }
double overlapSize() const
Definition: Settings.h:294
Settings const & settings_
Definition: ProcessBase.h:40
double rinv() const
Definition: Tracklet.h:120
const std::vector< double > varRInvBins() const
Definition: Settings.h:298

◆ findVarRInvBin()

unsigned int PurgeDuplicate::findVarRInvBin ( const Tracklet trk) const
private

Definition at line 736 of file PurgeDuplicate.cc.

References trigObjTnPSource_cfi::bins, HLT_2023v12_cff::distance, trklet::Tracklet::rinv(), trklet::ProcessBase::settings_, pfDeepBoostedJetPreprocessParams_cfi::upper_bound, and trklet::Settings::varRInvBins().

Referenced by execute().

736  {
737  std::vector<double> rInvBins = settings_.varRInvBins();
738 
739  //Get rinverse of track
740  double rInv = trk->rinv();
741 
742  //Check between what 2 values in rinvbins rinv is between
743  auto bins = std::upper_bound(rInvBins.begin(), rInvBins.end(), rInv);
744 
745  //return integer for bin index
746  unsigned int rIndx = std::distance(rInvBins.begin(), bins);
747  if (rIndx == std::distance(rInvBins.end(), bins))
748  return rInvBins.size() - 2;
749  else if (bins == rInvBins.begin())
750  return std::distance(rInvBins.begin(), bins);
751  else
752  return rIndx - 1;
753 }
Settings const & settings_
Definition: ProcessBase.h:40
double rinv() const
Definition: Tracklet.h:120
const std::vector< double > varRInvBins() const
Definition: Settings.h:298

◆ getInventedCoords()

std::vector< double > PurgeDuplicate::getInventedCoords ( unsigned int  iSector,
const Stub st,
const Tracklet tracklet 
) const
private

Definition at line 617 of file PurgeDuplicate.cc.

References funct::abs(), trklet::Tracklet::disk(), trklet::Settings::dphisector(), trklet::Settings::dphisectorHG(), findLayerDisk(), dqmdumpme::first, trklet::Stub::isBarrel(), trklet::Tracklet::phi0(), reco::reduceRange(), trklet::Tracklet::rinv(), trklet::Settings::rmean(), edm::second(), trklet::ProcessBase::settings_, funct::sin(), trklet::Tracklet::t(), trackerHitRTTI::vector, trklet::Tracklet::z0(), and trklet::Settings::zmean().

Referenced by getInventedSeedingStub().

619  {
620  int stubLayer = (findLayerDisk(st)).first;
621  int stubDisk = (findLayerDisk(st)).second;
622 
623  double stub_phi = -99;
624  double stub_z = -99;
625  double stub_r = -99;
626 
627  double tracklet_rinv = tracklet->rinv();
628 
629  if (st->isBarrel()) {
630  stub_r = settings_.rmean(stubLayer - 1);
631  stub_phi = tracklet->phi0() - std::asin(stub_r * tracklet_rinv / 2);
632  stub_phi = stub_phi + iSector * settings_.dphisector() - 0.5 * settings_.dphisectorHG();
633  stub_phi = reco::reduceRange(stub_phi);
634  stub_z = tracklet->z0() + 2 * tracklet->t() * 1 / tracklet_rinv * std::asin(stub_r * tracklet_rinv / 2);
635  } else {
636  stub_z = settings_.zmean(stubDisk - 1) * tracklet->disk() / abs(tracklet->disk());
637  stub_phi = tracklet->phi0() - (stub_z - tracklet->z0()) * tracklet_rinv / 2 / tracklet->t();
638  stub_phi = stub_phi + iSector * settings_.dphisector() - 0.5 * settings_.dphisectorHG();
639  stub_phi = reco::reduceRange(stub_phi);
640  stub_r = 2 / tracklet_rinv * std::sin((stub_z - tracklet->z0()) * tracklet_rinv / 2 / tracklet->t());
641  }
642 
643  std::vector invented_coords{stub_r, stub_z, stub_phi};
644  return invented_coords;
645 }
double t() const
Definition: Tracklet.h:123
int disk() const
Definition: Tracklet.cc:782
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
double dphisectorHG() const
Definition: Settings.h:314
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool isBarrel() const
Definition: Stub.h:81
Settings const & settings_
Definition: ProcessBase.h:40
double phi0() const
Definition: Tracklet.h:121
U second(std::pair< T, U > const &p)
double dphisector() const
Definition: Settings.h:323
double rmean(unsigned int iLayer) const
Definition: Settings.h:171
double rinv() const
Definition: Tracklet.h:120
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double zmean(unsigned int iDisk) const
Definition: Settings.h:174
std::pair< int, int > findLayerDisk(const Stub *) const
double z0() const
Definition: Tracklet.h:124

◆ getInventedCoordsExtended()

std::vector< double > PurgeDuplicate::getInventedCoordsExtended ( unsigned int  iSector,
const Stub st,
const Tracklet tracklet 
) const
private

Definition at line 647 of file PurgeDuplicate.cc.

References funct::abs(), HLT_2023v12_cff::beta, funct::cos(), trklet::Tracklet::d0(), trklet::Tracklet::disk(), trklet::Settings::dphisector(), trklet::Settings::dphisectorHG(), findLayerDisk(), dqmdumpme::first, trklet::Stub::isBarrel(), trklet::Stub::l1tstub(), trklet::L1TStub::phi(), trklet::Tracklet::phi0(), trklet::L1TStub::r(), reco::reduceRange(), trklet::Tracklet::rinv(), trklet::Settings::rmean(), edm::second(), fileCollector::seed, trklet::Tracklet::seedIndex(), trklet::ProcessBase::settings_, mathSSE::sqrt(), trklet::Tracklet::t(), trackerHitRTTI::vector, trklet::L1TStub::z(), trklet::Tracklet::z0(), and trklet::Settings::zmean().

Referenced by getInventedSeedingStub().

649  {
650  int stubLayer = (findLayerDisk(st)).first;
651  int stubDisk = (findLayerDisk(st)).second;
652 
653  double stub_phi = -99;
654  double stub_z = -99;
655  double stub_r = -99;
656 
657  double rho = 1 / tracklet->rinv();
658  double rho_minus_d0 = rho + tracklet->d0(); // should be -, but otherwise does not work
659 
660  // exact helix
661  if (st->isBarrel()) {
662  stub_r = settings_.rmean(stubLayer - 1);
663 
664  double sin_val = (stub_r * stub_r + rho_minus_d0 * rho_minus_d0 - rho * rho) / (2 * stub_r * rho_minus_d0);
665  stub_phi = tracklet->phi0() - std::asin(sin_val);
666  stub_phi = stub_phi + iSector * settings_.dphisector() - 0.5 * settings_.dphisectorHG();
667  stub_phi = reco::reduceRange(stub_phi);
668 
669  double beta = std::acos((rho * rho + rho_minus_d0 * rho_minus_d0 - stub_r * stub_r) / (2 * rho * rho_minus_d0));
670  stub_z = tracklet->z0() + tracklet->t() * std::abs(rho * beta);
671  } else {
672  stub_z = settings_.zmean(stubDisk - 1) * tracklet->disk() / abs(tracklet->disk());
673 
674  double beta = (stub_z - tracklet->z0()) / (tracklet->t() * std::abs(rho)); // maybe rho should be abs value
675  double r_square = -2 * rho * rho_minus_d0 * std::cos(beta) + rho * rho + rho_minus_d0 * rho_minus_d0;
676  stub_r = sqrt(r_square);
677 
678  double sin_val = (stub_r * stub_r + rho_minus_d0 * rho_minus_d0 - rho * rho) / (2 * stub_r * rho_minus_d0);
679  stub_phi = tracklet->phi0() - std::asin(sin_val);
680  stub_phi = stub_phi + iSector * settings_.dphisector() - 0.5 * settings_.dphisectorHG();
681  stub_phi = reco::reduceRange(stub_phi);
682  }
683 
684  // TMP: for displaced tracking, exclude one of the 3 seeding stubs
685  // to be discussed
686  int seed = tracklet->seedIndex();
687  if ((seed == 8 && stubLayer == 4) || (seed == 9 && stubLayer == 5) || (seed == 10 && stubLayer == 3) ||
688  (seed == 11 && abs(stubDisk) == 1)) {
689  stub_phi = st->l1tstub()->phi();
690  stub_z = st->l1tstub()->z();
691  stub_r = st->l1tstub()->r();
692  }
693 
694  std::vector invented_coords{stub_r, stub_z, stub_phi};
695  return invented_coords;
696 }
unsigned int seedIndex() const
Definition: Tracklet.h:224
double phi() const
Definition: L1TStub.h:65
double t() const
Definition: Tracklet.h:123
int disk() const
Definition: Tracklet.cc:782
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
double dphisectorHG() const
Definition: Settings.h:314
bool isBarrel() const
Definition: Stub.h:81
Settings const & settings_
Definition: ProcessBase.h:40
double phi0() const
Definition: Tracklet.h:121
double z() const
Definition: L1TStub.h:59
U second(std::pair< T, U > const &p)
double dphisector() const
Definition: Settings.h:323
double rmean(unsigned int iLayer) const
Definition: Settings.h:171
double rinv() const
Definition: Tracklet.h:120
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
L1TStub * l1tstub()
Definition: Stub.h:77
double zmean(unsigned int iDisk) const
Definition: Settings.h:174
std::pair< int, int > findLayerDisk(const Stub *) const
double d0() const
Definition: Tracklet.h:122
double z0() const
Definition: Tracklet.h:124
double r() const
Definition: L1TStub.h:60

◆ getInventedSeedingStub()

std::vector< const Stub * > PurgeDuplicate::getInventedSeedingStub ( unsigned int  iSector,
const Tracklet tracklet,
const std::vector< const Stub *> &  originalStubsList 
) const
private

Definition at line 698 of file PurgeDuplicate.cc.

References trklet::L1TStub::allStubIndex(), funct::cos(), trklet::Settings::extended(), findLayerDisk(), dqmdumpme::first, getInventedCoords(), getInventedCoordsExtended(), isSeedingStub(), trklet::Stub::l1tstub(), edm::second(), trklet::Tracklet::seedIndex(), trklet::L1TStub::setAllStubIndex(), trklet::L1TStub::setCoords(), trklet::Stub::setl1tstub(), trklet::ProcessBase::settings_, trklet::L1TStub::setUniqueIndex(), funct::sin(), and trklet::L1TStub::uniqueIndex().

Referenced by execute().

699  {
700  std::vector<const Stub*> newStubList;
701 
702  for (unsigned int stubit = 0; stubit < originalStubsList.size(); stubit++) {
703  const Stub* thisStub = originalStubsList[stubit];
704 
705  if (isSeedingStub(tracklet->seedIndex(), (findLayerDisk(thisStub)).first, (findLayerDisk(thisStub)).second)) {
706  // get a vector containing r, z, phi
707  std::vector<double> inv_r_z_phi;
708  if (!settings_.extended())
709  inv_r_z_phi = getInventedCoords(iSector, thisStub, tracklet);
710  else {
711  inv_r_z_phi = getInventedCoordsExtended(iSector, thisStub, tracklet);
712  }
713  double stub_x_invent = inv_r_z_phi[0] * std::cos(inv_r_z_phi[2]);
714  double stub_y_invent = inv_r_z_phi[0] * std::sin(inv_r_z_phi[2]);
715  double stub_z_invent = inv_r_z_phi[1];
716 
717  Stub* invent_stub_ptr = new Stub(*thisStub);
718  const L1TStub* l1stub = thisStub->l1tstub();
719  L1TStub invent_l1stub = *l1stub;
720  invent_l1stub.setCoords(stub_x_invent, stub_y_invent, stub_z_invent);
721 
722  invent_stub_ptr->setl1tstub(new L1TStub(invent_l1stub));
723  invent_stub_ptr->l1tstub()->setAllStubIndex(l1stub->allStubIndex());
724  invent_stub_ptr->l1tstub()->setUniqueIndex(l1stub->uniqueIndex());
725 
726  newStubList.push_back(invent_stub_ptr);
727 
728  } else {
729  newStubList.push_back(thisStub);
730  }
731  }
732  return newStubList;
733 }
unsigned int seedIndex() const
Definition: Tracklet.h:224
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Settings const & settings_
Definition: ProcessBase.h:40
U second(std::pair< T, U > const &p)
std::vector< double > getInventedCoords(unsigned int, const Stub *, const Tracklet *) const
std::vector< double > getInventedCoordsExtended(unsigned int, const Stub *, const Tracklet *) const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
L1TStub * l1tstub()
Definition: Stub.h:77
void setUniqueIndex(unsigned int index)
Definition: L1TStub.h:84
unsigned int allStubIndex() const
Definition: L1TStub.h:91
std::pair< int, int > findLayerDisk(const Stub *) const
void setCoords(double x, double y, double z)
Definition: L1TStub.h:85
bool extended() const
Definition: Settings.h:266
void setl1tstub(L1TStub *l1tstub)
Definition: Stub.h:79
void setAllStubIndex(unsigned int index)
Definition: L1TStub.h:83
bool isSeedingStub(int, int, int) const
unsigned int uniqueIndex() const
Definition: L1TStub.h:92

◆ getPhiRes()

double PurgeDuplicate::getPhiRes ( Tracklet curTracklet,
const Stub curStub 
) const
private

Definition at line 548 of file PurgeDuplicate.cc.

References funct::abs(), isSeedingStub(), trklet::Stub::l1tstub(), trklet::Stub::layerdisk(), trklet::N_LAYER, trklet::L1TStub::phi(), trklet::Tracklet::proj(), and trklet::Tracklet::seedIndex().

Referenced by execute().

548  {
549  double phiproj;
550  double stubphi;
551  double phires;
552  // Get phi position of stub
553  stubphi = curStub->l1tstub()->phi();
554  // Get layer that the stub is in (Layer 1->6, Disk 1->5)
555  int Layer = curStub->layerdisk() + 1;
556  if (Layer > N_LAYER) {
557  Layer = 0;
558  }
559  int Disk = curStub->layerdisk() - (N_LAYER - 1);
560  if (Disk < 0) {
561  Disk = 0;
562  }
563  // Get phi projection of tracklet
564  int seedindex = curTracklet->seedIndex();
565  // If this stub is a seed stub, set projection=phi, so that res=0
566  if (isSeedingStub(seedindex, Layer, Disk)) {
567  phiproj = stubphi;
568  // Otherwise, get projection of tracklet
569  } else {
570  phiproj = curTracklet->proj(curStub->layerdisk()).phiproj();
571  }
572  // Calculate residual
573  phires = std::abs(stubphi - phiproj);
574  return phires;
575 }
unsigned int seedIndex() const
Definition: Tracklet.h:224
double phi() const
Definition: L1TStub.h:65
Projection & proj(int layerdisk)
Definition: Tracklet.h:87
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
Definition: BoundDisk.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
L1TStub * l1tstub()
Definition: Stub.h:77
unsigned int layerdisk() const
Definition: Stub.cc:185
bool isSeedingStub(int, int, int) const
constexpr int N_LAYER
Definition: Settings.h:25

◆ isSeedingStub()

bool PurgeDuplicate::isSeedingStub ( int  seedindex,
int  Layer,
int  Disk 
) const
private

Definition at line 577 of file PurgeDuplicate.cc.

References funct::abs().

Referenced by getInventedSeedingStub(), and getPhiRes().

577  {
578  if ((seedindex == 0 && (Layer == 1 || Layer == 2)) || (seedindex == 1 && (Layer == 2 || Layer == 3)) ||
579  (seedindex == 2 && (Layer == 3 || Layer == 4)) || (seedindex == 3 && (Layer == 5 || Layer == 6)) ||
580  (seedindex == 4 && (abs(Disk) == 1 || abs(Disk) == 2)) ||
581  (seedindex == 5 && (abs(Disk) == 3 || abs(Disk) == 4)) || (seedindex == 6 && (Layer == 1 || abs(Disk) == 1)) ||
582  (seedindex == 7 && (Layer == 2 || abs(Disk) == 1)) ||
583  (seedindex == 8 && (Layer == 2 || Layer == 3 || Layer == 4)) ||
584  (seedindex == 9 && (Layer == 4 || Layer == 5 || Layer == 6)) ||
585  (seedindex == 10 && (Layer == 2 || Layer == 3 || abs(Disk) == 1)) ||
586  (seedindex == 11 && (Layer == 2 || abs(Disk) == 1 || abs(Disk) == 2)))
587  return true;
588 
589  return false;
590 }
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
Definition: BoundDisk.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ isTrackInBin()

bool PurgeDuplicate::isTrackInBin ( const std::vector< unsigned int > &  vec,
unsigned int  num 
) const
private

Definition at line 770 of file PurgeDuplicate.cc.

References spr::find(), newFWLiteAna::found, EgammaValidation_cff::num, and mps_fire::result.

Referenced by execute().

770  {
771  auto result = std::find(vec.begin(), vec.end(), num);
772  bool found = (result != vec.end());
773  return found;
774 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19

◆ l1tinfo()

std::string PurgeDuplicate::l1tinfo ( const L1TStub l1stub,
std::string  str = "" 
) const
private

Definition at line 605 of file PurgeDuplicate.cc.

References trklet::L1TStub::bend(), trklet::L1TStub::iphi(), trklet::L1TStub::layerdisk(), trklet::L1TStub::phi(), trklet::L1TStub::r(), str, AlCaHLTBitMon_QueryRunRegistry::string, and trklet::L1TStub::z().

605  {
606  std::string thestr = Form("\t %s stub info: r/z/phi:\t%f\t%f\t%f\t%d\t%f\t%d",
607  str.c_str(),
608  l1stub->r(),
609  l1stub->z(),
610  l1stub->phi(),
611  l1stub->iphi(),
612  l1stub->bend(),
613  l1stub->layerdisk());
614  return thestr;
615 }
double phi() const
Definition: L1TStub.h:65
double bend() const
Definition: L1TStub.h:63
double z() const
Definition: L1TStub.h:59
int layerdisk() const
Definition: L1TStub.h:119
double r() const
Definition: L1TStub.h:60
unsigned int iphi() const
Definition: L1TStub.h:67
#define str(s)

Member Data Documentation

◆ inputstubidslists_

std::vector<std::vector<std::pair<int, int> > > trklet::PurgeDuplicate::inputstubidslists_
private

Definition at line 60 of file PurgeDuplicate.h.

Referenced by execute().

◆ inputstublists_

std::vector<std::vector<const Stub*> > trklet::PurgeDuplicate::inputstublists_
private

Definition at line 59 of file PurgeDuplicate.h.

Referenced by execute().

◆ inputtrackfits_

std::vector<TrackFitMemory*> trklet::PurgeDuplicate::inputtrackfits_
private

Definition at line 62 of file PurgeDuplicate.h.

Referenced by addInput(), and execute().

◆ inputtracklets_

std::vector<Tracklet*> trklet::PurgeDuplicate::inputtracklets_
private

Definition at line 63 of file PurgeDuplicate.h.

Referenced by execute().

◆ inputtracks_

std::vector<Track*> trklet::PurgeDuplicate::inputtracks_
private

Definition at line 58 of file PurgeDuplicate.h.

Referenced by execute().

◆ mergedstubidslists_

std::vector<std::vector<std::pair<int, int> > > trklet::PurgeDuplicate::mergedstubidslists_
private

Definition at line 61 of file PurgeDuplicate.h.

Referenced by execute().

◆ outputtracklets_

std::vector<CleanTrackMemory*> trklet::PurgeDuplicate::outputtracklets_
private

Definition at line 64 of file PurgeDuplicate.h.

Referenced by addOutput(), and execute().