CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
L1TMuonBarrelKalmanRegionModule Class Reference

#include <L1TMuonBarrelKalmanRegionModule.h>

Classes

class  SeedSorter
 

Public Member Functions

 L1TMuonBarrelKalmanRegionModule (const edm::ParameterSet &, int wheel, int sector)
 
L1MuKBMTrackCollection process (L1TMuonBarrelKalmanAlgo *, const L1MuKBMTCombinedStubRefVector &stubs, int bx)
 
const int wheel () const
 
 ~L1TMuonBarrelKalmanRegionModule ()
 

Private Member Functions

L1MuKBMTrackCollection cleanHigher (const L1MuKBMTrackCollection &tracks1, const L1MuKBMTrackCollection &tracks2)
 
L1MuKBMTrackCollection cleanLower (const L1MuKBMTrackCollection &tracks1, const L1MuKBMTrackCollection &tracks2)
 
L1MuKBMTrackCollection cleanRegion (const L1MuKBMTrackCollection &, const L1MuKBMTrackCollection &, const L1MuKBMTrackCollection &)
 
L1MuKBMTrackCollection selfClean (const L1MuKBMTrackCollection &tracks)
 
L1MuKBMTrackCollection sort4 (const L1MuKBMTrackCollection &in)
 

Private Attributes

int nextSector_
 
int nextWheel_
 
int previousSector_
 
int sector_
 
int verbose_
 
int wheel_
 

Detailed Description

Definition at line 8 of file L1TMuonBarrelKalmanRegionModule.h.

Constructor & Destructor Documentation

L1TMuonBarrelKalmanRegionModule::L1TMuonBarrelKalmanRegionModule ( const edm::ParameterSet iConfig,
int  wheel,
int  sector 
)

Definition at line 3 of file L1TMuonBarrelKalmanRegionModule.cc.

References nextSector_, nextWheel_, and previousSector_.

6  : verbose_(iConfig.getParameter<int>("verbose")), sector_(sector), wheel_(wheel) {
7  if (sector == 11) {
8  nextSector_ = 0;
9  previousSector_ = 10;
10  } else if (sector == 0) {
11  nextSector_ = 1;
12  previousSector_ = 11;
13  } else {
14  nextSector_ = sector + 1;
15  previousSector_ = sector - 1;
16  }
17 
18  switch (wheel) {
19  case -2:
20  nextWheel_ = -1;
21  break;
22 
23  case -1:
24  nextWheel_ = 0;
25  break;
26 
27  case 0:
28  nextWheel_ = 999;
29  break;
30 
31  case 1:
32  nextWheel_ = 0;
33  break;
34 
35  case 2:
36  nextWheel_ = 1;
37  break;
38 
39  default:
40  nextWheel_ = 999;
41  break;
42  }
43 }
T getParameter(std::string const &) const
L1TMuonBarrelKalmanRegionModule::~L1TMuonBarrelKalmanRegionModule ( )

Definition at line 45 of file L1TMuonBarrelKalmanRegionModule.cc.

45 {}

Member Function Documentation

L1MuKBMTrackCollection L1TMuonBarrelKalmanRegionModule::cleanHigher ( const L1MuKBMTrackCollection tracks1,
const L1MuKBMTrackCollection tracks2 
)
private

Definition at line 149 of file L1TMuonBarrelKalmanRegionModule.cc.

References mps_fire::i, dqmiolumiharvest::j, keep, MillePedeFileConverter_cfg::out, and parallelization::uint.

Referenced by cleanRegion().

150  {
152 
153  for (uint i = 0; i < tracks1.size(); ++i) {
154  bool keep = true;
155 
156  for (uint j = 0; j < tracks2.size(); ++j) {
157  if (tracks1[i].overlapTrack(tracks2[j])) {
158  if (tracks1[i].rank() <= tracks2[j].rank()) {
159  keep = false;
160  }
161  }
162  }
163  if (keep)
164  out.push_back(tracks1[i]);
165  }
166 
167  return out;
168 }
const int keep
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
L1MuKBMTrackCollection L1TMuonBarrelKalmanRegionModule::cleanLower ( const L1MuKBMTrackCollection tracks1,
const L1MuKBMTrackCollection tracks2 
)
private

Definition at line 170 of file L1TMuonBarrelKalmanRegionModule.cc.

References mps_fire::i, dqmiolumiharvest::j, keep, MillePedeFileConverter_cfg::out, and parallelization::uint.

Referenced by cleanRegion().

171  {
173 
174  for (uint i = 0; i < tracks1.size(); ++i) {
175  bool keep = true;
176 
177  for (uint j = 0; j < tracks2.size(); ++j) {
178  if (tracks1[i].overlapTrack(tracks2[j])) {
179  if (tracks1[i].rank() < tracks2[j].rank()) {
180  keep = false;
181  }
182  }
183  }
184  if (keep)
185  out.push_back(tracks1[i]);
186  }
187 
188  return out;
189 }
const int keep
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
L1MuKBMTrackCollection L1TMuonBarrelKalmanRegionModule::cleanRegion ( const L1MuKBMTrackCollection tracks2,
const L1MuKBMTrackCollection tracks3,
const L1MuKBMTrackCollection tracks4 
)
private

Definition at line 275 of file L1TMuonBarrelKalmanRegionModule.cc.

References cleanHigher(), cleanLower(), MillePedeFileConverter_cfg::out, selfClean(), and sort4().

Referenced by process().

277  {
278  L1MuKBMTrackCollection cleaned2 = selfClean(tracks2);
279  L1MuKBMTrackCollection cleaned3 = selfClean(tracks3);
280 
281  L1MuKBMTrackCollection cleaned23 = cleanHigher(cleaned2, tracks3);
282  L1MuKBMTrackCollection cleaned32 = cleanLower(cleaned3, tracks2);
283 
284  //printf("Cleaned sizes = 2=%d 3=%d 23=%d 32=%d \n",int(cleaned2.size()),int(cleaned3.size()),int(cleaned23.size()),int(cleaned32.size()));
285 
286  //merge 2,3
288 
289  if (!cleaned32.empty())
290  step1.insert(step1.end(), cleaned32.begin(), cleaned32.end());
291  if (!cleaned23.empty())
292  step1.insert(step1.end(), cleaned23.begin(), cleaned23.end());
293 
294  //take the best 2
295  L1MuKBMTrackCollection sorted23 = sort4(step1);
296  // printf("Sorted 23 =%d\n",int(sorted23.size()));
297 
298  //Now clean the tracks 4 between them
299  L1MuKBMTrackCollection cleaned4 = selfClean(tracks4);
300  // printf("Cleaned4 =%d\n",int(cleaned4.size()));
301 
302  //Now clean the 23 tracks from tracks4
303  L1MuKBMTrackCollection cleanedSorted23 = cleanHigher(sorted23, tracks4);
304 
305  // printf("Cleaned23 from 4 =%d\n",int(cleanedSorted23.size()));
306 
307  //Now clean the tracks4 from sorted 23
308  L1MuKBMTrackCollection cleanedSorted4 = cleanLower(cleaned4, sorted23);
309  // printf("Cleaned4 from 23 =%d\n",int(cleanedSorted4.size()));
310 
311  //Now merge all of those
313 
314  if (!cleanedSorted4.empty())
315  step2.insert(step2.end(), cleanedSorted4.begin(), cleanedSorted4.end());
316  if (!cleanedSorted23.empty())
317  step2.insert(step2.end(), cleanedSorted23.begin(), cleanedSorted23.end());
318 
320  return out;
321 }
L1MuKBMTrackCollection cleanLower(const L1MuKBMTrackCollection &tracks1, const L1MuKBMTrackCollection &tracks2)
L1MuKBMTrackCollection cleanHigher(const L1MuKBMTrackCollection &tracks1, const L1MuKBMTrackCollection &tracks2)
L1MuKBMTrackCollection sort4(const L1MuKBMTrackCollection &in)
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
L1MuKBMTrackCollection selfClean(const L1MuKBMTrackCollection &tracks)
L1MuKBMTrackCollection L1TMuonBarrelKalmanRegionModule::process ( L1TMuonBarrelKalmanAlgo trackMaker,
const L1MuKBMTCombinedStubRefVector stubs,
int  bx 
)

Definition at line 47 of file L1TMuonBarrelKalmanRegionModule.cc.

References l1GtPatternGenerator_cfi::bx, L1TMuonBarrelKalmanAlgo::chain(), cleanRegion(), createfilelist::int, nextSector_, nextWheel_, MillePedeFileConverter_cfg::out, previousSector_, sector_, SurveyInfoScenario_cff::seed, MCScenario_CRAFT1_22X::sorter(), verbose_, and wheel_.

Referenced by wheel().

49  {
52  L1MuKBMTrackCollection pretracks2;
53  L1MuKBMTrackCollection pretracks3;
54  L1MuKBMTrackCollection pretracks4;
55  for (const auto& stub : stubsAll) {
56  if (stub->bxNum() != bx)
57  continue;
58 
59  if ((stub->scNum() == nextSector_ && stub->phi() >= -112) ||
60  (stub->scNum() == previousSector_ && stub->phi() <= 111))
61  continue;
62 
63  if (stub->whNum() == wheel_ && stub->scNum() == sector_) {
64  seeds.push_back(stub);
65  stubs.push_back(stub);
66  } else if (stub->whNum() == wheel_ && (stub->scNum() == nextSector_ || stub->scNum() == previousSector_)) {
67  stubs.push_back(stub);
68  } else if (stub->whNum() == nextWheel_ &&
69  (stub->scNum() == nextSector_ || stub->scNum() == previousSector_ || stub->scNum() == sector_)) {
70  stubs.push_back(stub);
71  }
72  }
73 
74  //Sort the seeds by tag so that the emulator is aligned like the firmware
75 
76  SeedSorter sorter;
77  if (seeds.size() > 1) {
78  std::sort(seeds.begin(), seeds.end(), sorter);
79  }
80 
81  for (const auto seed : seeds) {
82  std::pair<bool, L1MuKBMTrack> trackInfo = trackMaker->chain(seed, stubs);
83  if (trackInfo.first) {
84  if (seed->stNum() == 2)
85  pretracks2.push_back(trackInfo.second);
86  if (seed->stNum() == 3)
87  pretracks3.push_back(trackInfo.second);
88  if (seed->stNum() == 4)
89  pretracks4.push_back(trackInfo.second);
90  }
91  }
92  //for (const auto& track1 :pretracks2)
93  // printf("SEED=2 Kalman Track charge=%d pt=%f eta=%f phi=%f curvature=%d curvature STA =%d stubs=%d chi2=%d pts=%f %f pattern=%d\n",track1.charge(),track1.pt(),track1.eta(),track1.phi(),track1.curvatureAtVertex(),track1.curvatureAtMuon(),int(track1.stubs().size()),track1.approxChi2(),track1.pt(),track1.ptUnconstrained(),track1.hitPattern());
94 
95  //for (const auto& track1 :pretracks3)
96  // printf("SEED=3 Kalman Track charge=%d pt=%f eta=%f phi=%f curvature=%d curvature STA =%d stubs=%d chi2=%d pts=%f %f pattern=%d\n",track1.charge(),track1.pt(),track1.eta(),track1.phi(),track1.curvatureAtVertex(),track1.curvatureAtMuon(),int(track1.stubs().size()),track1.approxChi2(),track1.pt(),track1.ptUnconstrained(),track1.hitPattern());
97 
98  //for (const auto& track1 :pretracks4)
99  // printf("SEED=4 Kalman Track charge=%d pt=%f eta=%f phi=%f curvature=%d curvature STA =%d stubs=%d chi2=%d pts=%f %f pattern=%d\n",track1.charge(),track1.pt(),track1.eta(),track1.phi(),track1.curvatureAtVertex(),track1.curvatureAtMuon(),int(track1.stubs().size()),track1.approxChi2(),track1.pt(),track1.ptUnconstrained(),track1.hitPattern());
100 
101  // trackMaker->resolveEtaUnit(pretracks);
102  L1MuKBMTrackCollection out = cleanRegion(pretracks2, pretracks3, pretracks4);
103  if (verbose_) {
104  printf(" -----Sector Processor Kalman Tracks-----\n");
105  for (const auto& track1 : out)
106  printf("Kalman Track charge=%d pt=%f eta=%f phi=%f curvature=%d curvature STA =%d stubs=%d chi2=%d pts=%f %f\n",
107  track1.charge(),
108  track1.pt(),
109  track1.eta(),
110  track1.phi(),
111  track1.curvatureAtVertex(),
112  track1.curvatureAtMuon(),
113  int(track1.stubs().size()),
114  track1.approxChi2(),
115  track1.pt(),
116  track1.ptUnconstrained());
117  }
118 
119  return out;
120 }
std::vector< edm::Ref< L1MuKBMTCombinedStubCollection > > L1MuKBMTCombinedStubRefVector
L1MuKBMTrackCollection cleanRegion(const L1MuKBMTrackCollection &, const L1MuKBMTrackCollection &, const L1MuKBMTrackCollection &)
std::pair< bool, L1MuKBMTrack > chain(const L1MuKBMTCombinedStubRef &, const L1MuKBMTCombinedStubRefVector &)
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
L1MuKBMTrackCollection L1TMuonBarrelKalmanRegionModule::selfClean ( const L1MuKBMTrackCollection tracks)
private

Definition at line 122 of file L1TMuonBarrelKalmanRegionModule.cc.

References mps_fire::i, dqmiolumiharvest::j, keep, MillePedeFileConverter_cfg::out, and parallelization::uint.

Referenced by cleanRegion().

122  {
124 
125  for (uint i = 0; i < tracks.size(); ++i) {
126  bool keep = true;
127 
128  for (uint j = 0; j < tracks.size(); ++j) {
129  if (i == j)
130  continue;
131 
132  if (tracks[i].overlapTrack(tracks[j])) {
133  if (tracks[i].rank() < tracks[j].rank()) {
134  keep = false;
135  } else if (tracks[i].rank() == tracks[j].rank()) { //if same rank prefer seed that is high
136  if (!tracks[j].stubs()[0]->tag())
137  keep = false;
138  }
139  }
140  }
141 
142  if (keep)
143  out.push_back(tracks[i]);
144  }
145 
146  return out;
147 }
const int keep
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
L1MuKBMTrackCollection L1TMuonBarrelKalmanRegionModule::sort4 ( const L1MuKBMTrackCollection in)
private

Definition at line 191 of file L1TMuonBarrelKalmanRegionModule.cc.

References MillePedeFileConverter_cfg::out, DiDispStaMuonMonitor_cfi::pt, and reco::LeafCandidate::pt().

Referenced by cleanRegion().

191  {
193  //partial sort like in firmwarE (bitonic)
194  if (in.size() <= 2)
195  return in;
196  else if (in.size() == 3) {
197  //Step 1
198  L1MuKBMTrack s2_1;
199  L1MuKBMTrack s2_3;
200  if (in[2].pt() >= in[0].pt()) {
201  s2_1 = in[2];
202  s2_3 = in[0];
203  } else {
204  s2_1 = in[0];
205  s2_3 = in[2];
206  }
207 
208  L1MuKBMTrack s2_2 = in[1];
209  //Step 2;
210  L1MuKBMTrack s3_1 = s2_1;
211  L1MuKBMTrack s3_2;
212  L1MuKBMTrack s3_3;
213 
214  if (s2_3.pt() >= s2_2.pt()) {
215  s3_2 = s2_3;
216  s3_3 = s2_2;
217  } else {
218  s3_2 = s2_2;
219  s3_3 = s2_3;
220  }
221 
222  out.push_back(s3_1);
223  out.push_back(s3_2);
224 
225  } else {
226  //Step 1
227  L1MuKBMTrack s2_1;
228  L1MuKBMTrack s2_2;
229  L1MuKBMTrack s2_3;
230  L1MuKBMTrack s2_4;
231 
232  if (in[2].pt() >= in[0].pt()) {
233  s2_1 = in[2];
234  s2_3 = in[0];
235  } else {
236  s2_1 = in[0];
237  s2_3 = in[2];
238  }
239  if (in[3].pt() >= in[1].pt()) {
240  s2_2 = in[3];
241  s2_4 = in[1];
242  } else {
243  s2_2 = in[1];
244  s2_4 = in[3];
245  }
246  //Step 2
247  L1MuKBMTrack s3_1;
248  L1MuKBMTrack s3_2;
249  L1MuKBMTrack s3_3;
250  L1MuKBMTrack s3_4;
251 
252  if (s2_4.pt() >= s2_1.pt()) {
253  s3_1 = s2_4;
254  s3_4 = s2_1;
255  } else {
256  s3_1 = s2_1;
257  s3_4 = s2_4;
258  }
259 
260  if (s2_3.pt() >= s2_2.pt()) {
261  s3_2 = s2_3;
262  s3_3 = s2_2;
263  } else {
264  s3_2 = s2_2;
265  s3_3 = s2_3;
266  }
267 
268  out.push_back(s3_1);
269  out.push_back(s3_2);
270  }
271 
272  return out;
273 }
double pt() const final
transverse momentum
std::vector< L1MuKBMTrack > L1MuKBMTrackCollection
Definition: L1MuKBMTrack.h:15
const int L1TMuonBarrelKalmanRegionModule::wheel ( ) const
inline

Member Data Documentation

int L1TMuonBarrelKalmanRegionModule::nextSector_
private

Definition at line 21 of file L1TMuonBarrelKalmanRegionModule.h.

Referenced by L1TMuonBarrelKalmanRegionModule(), and process().

int L1TMuonBarrelKalmanRegionModule::nextWheel_
private

Definition at line 23 of file L1TMuonBarrelKalmanRegionModule.h.

Referenced by L1TMuonBarrelKalmanRegionModule(), and process().

int L1TMuonBarrelKalmanRegionModule::previousSector_
private

Definition at line 22 of file L1TMuonBarrelKalmanRegionModule.h.

Referenced by L1TMuonBarrelKalmanRegionModule(), and process().

int L1TMuonBarrelKalmanRegionModule::sector_
private

Definition at line 19 of file L1TMuonBarrelKalmanRegionModule.h.

Referenced by process().

int L1TMuonBarrelKalmanRegionModule::verbose_
private

Definition at line 18 of file L1TMuonBarrelKalmanRegionModule.h.

Referenced by process().

int L1TMuonBarrelKalmanRegionModule::wheel_
private

Definition at line 20 of file L1TMuonBarrelKalmanRegionModule.h.

Referenced by process(), and wheel().