CMS 3D CMS Logo

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

#include <L1TPhase2GMTBarrelStubProcessor.h>

Public Member Functions

 L1TPhase2GMTBarrelStubProcessor ()
 
 L1TPhase2GMTBarrelStubProcessor (const edm::ParameterSet &)
 
l1t::MuonStubCollection makeStubs (const L1Phase2MuDTPhContainer *, const L1MuDTChambThContainer *)
 
 ~L1TPhase2GMTBarrelStubProcessor ()
 

Private Member Functions

l1t::MuonStub buildStub (const L1Phase2MuDTPhDigi &, const L1MuDTChambThDigi *)
 
l1t::MuonStub buildStubNoEta (const L1Phase2MuDTPhDigi &)
 
int calculateEta (uint, int, uint, uint)
 

Private Attributes

std::vector< int > coarseEta1_
 
std::vector< int > coarseEta2_
 
std::vector< int > coarseEta3_
 
std::vector< int > coarseEta4_
 
std::vector< int > eta1_
 
std::vector< int > eta2_
 
std::vector< int > eta3_
 
double etaLSB_
 
int maxBX_
 
int minBX_
 
int minPhiQuality_
 
int phiBFactor_
 
double phiLSB_
 
std::vector< int > phiOffset_
 
int verbose_
 

Detailed Description

Definition at line 13 of file L1TPhase2GMTBarrelStubProcessor.h.

Constructor & Destructor Documentation

◆ L1TPhase2GMTBarrelStubProcessor() [1/2]

L1TPhase2GMTBarrelStubProcessor::L1TPhase2GMTBarrelStubProcessor ( )

◆ L1TPhase2GMTBarrelStubProcessor() [2/2]

L1TPhase2GMTBarrelStubProcessor::L1TPhase2GMTBarrelStubProcessor ( const edm::ParameterSet iConfig)

Definition at line 9 of file L1TPhase2GMTBarrelStubProcessor.cc.

10  : minPhiQuality_(iConfig.getParameter<int>("minPhiQuality")),
11  minBX_(iConfig.getParameter<int>("minBX")),
12  maxBX_(iConfig.getParameter<int>("maxBX")),
13  eta1_(iConfig.getParameter<std::vector<int> >("eta_1")),
14  eta2_(iConfig.getParameter<std::vector<int> >("eta_2")),
15  eta3_(iConfig.getParameter<std::vector<int> >("eta_3")),
16  coarseEta1_(iConfig.getParameter<std::vector<int> >("coarseEta_1")),
17  coarseEta2_(iConfig.getParameter<std::vector<int> >("coarseEta_2")),
18  coarseEta3_(iConfig.getParameter<std::vector<int> >("coarseEta_3")),
19  coarseEta4_(iConfig.getParameter<std::vector<int> >("coarseEta_4")),
20  phiOffset_(iConfig.getParameter<std::vector<int> >("phiOffset")),
21  phiBFactor_(iConfig.getParameter<int>("phiBDivider")),
22  verbose_(iConfig.getParameter<int>("verbose")),
23  phiLSB_(iConfig.getParameter<double>("phiLSB")),
24  etaLSB_(iConfig.getParameter<double>("etaLSB")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307

◆ ~L1TPhase2GMTBarrelStubProcessor()

L1TPhase2GMTBarrelStubProcessor::~L1TPhase2GMTBarrelStubProcessor ( )

Definition at line 26 of file L1TPhase2GMTBarrelStubProcessor.cc.

26 {}

Member Function Documentation

◆ buildStub()

l1t::MuonStub L1TPhase2GMTBarrelStubProcessor::buildStub ( const L1Phase2MuDTPhDigi phiS,
const L1MuDTChambThDigi etaS 
)
private

Definition at line 28 of file L1TPhase2GMTBarrelStubProcessor.cc.

References buildStubNoEta(), calculateEta(), HLT_2023v12_cff::eta1, HLT_2023v12_cff::eta2, etaLSB_, mps_fire::i, l1t::MuonStub::offline_coord1(), l1t::MuonStub::offline_coord2(), L1MuDTChambThDigi::position(), L1MuDTChambThDigi::quality(), L1MuDTChambThDigi::scNum(), l1t::MuonStub::setEta(), l1t::MuonStub::setOfflineQuantities(), L1MuDTChambThDigi::stNum(), parallelization::uint, and L1MuDTChambThDigi::whNum().

Referenced by makeStubs().

29  {
30  l1t::MuonStub stub = buildStubNoEta(phiS);
31 
32  //Now full eta
33  int qeta1 = -16384;
34  int qeta2 = -16384;
35  int eta1 = -16384;
36  int eta2 = -16384;
37 
38  bool hasEta = false;
39  for (uint i = 0; i < 7; ++i) {
40  if (etaS->position(i) == 0)
41  continue;
42  if (!hasEta) {
43  hasEta = true;
44  eta1 = calculateEta(i, etaS->whNum(), etaS->scNum(), etaS->stNum());
45 
46  if (etaS->quality(i) == 1)
47  qeta1 = 2;
48  else
49  qeta1 = 1;
50  } else {
51  eta2 = calculateEta(i, etaS->whNum(), etaS->scNum(), etaS->stNum());
52  if (etaS->quality(i) == 1)
53  qeta2 = 2;
54  else
55  qeta2 = 1;
56  }
57  }
58 
59  if (qeta2 > 0) { //both stubs->average
60  stub.setEta(eta1, eta2, 3);
62 
63  } else if (qeta1 > 0) { //Good single stub
64  stub.setEta(eta1, 0, 1);
65  stub.setOfflineQuantities(stub.offline_coord1(), stub.offline_coord2(), eta1 * etaLSB_, 0.0);
66  }
67 
68  return stub;
69 }
void setOfflineQuantities(double coord1, double coord2, double eta1, double eta2)
Definition: MuonStub.h:101
double offline_coord1() const
Definition: MuonStub.h:96
int position(const int i) const
l1t::MuonStub buildStubNoEta(const L1Phase2MuDTPhDigi &)
double offline_coord2() const
Definition: MuonStub.h:97
void setEta(int eta1, int eta2, int etaQ)
Definition: MuonStub.h:107
int quality(const int i) const

◆ buildStubNoEta()

l1t::MuonStub L1TPhase2GMTBarrelStubProcessor::buildStubNoEta ( const L1Phase2MuDTPhDigi phiS)
private

Definition at line 71 of file L1TPhase2GMTBarrelStubProcessor.cc.

References nano_mu_digi_cff::bx, L1Phase2MuDTPhDigi::bxNum(), coarseEta1_, coarseEta2_, coarseEta3_, coarseEta4_, PVValHelper::eta, etaLSB_, L1Phase2MuDTPhDigi::index(), createfilelist::int, M_PI, phi, L1Phase2MuDTPhDigi::phi(), L1Phase2MuDTPhDigi::phiBend(), phiBFactor_, phiLSB_, phiOffset_, quality, L1Phase2MuDTPhDigi::scNum(), nano_mu_digi_cff::sector, l1t::MuonStub::setOfflineQuantities(), Validation_hcalonly_cfi::sign, relativeConstraints::station, L1Phase2MuDTPhDigi::stNum(), makeGlobalPositionRcd_cfg::tag, parallelization::uint, makeMuonMisalignmentScenario::wheel, and L1Phase2MuDTPhDigi::whNum().

Referenced by buildStub(), and makeStubs().

71  {
72  int wheel = phiS.whNum();
73  int abswheel = fabs(phiS.whNum());
74  int sign = wheel > 0 ? 1 : -1;
75  int sector = phiS.scNum();
76  int station = phiS.stNum();
77  double globalPhi = (sector * 30) + phiS.phi() * 30. / 65535.;
78  if (globalPhi < -180)
79  globalPhi += 360;
80  if (globalPhi > 180)
81  globalPhi -= 360;
82  globalPhi = globalPhi * M_PI / 180.;
83  int phi = int(globalPhi / phiLSB_) + phiOffset_[station - 1];
84  int phiB = phiS.phiBend() / phiBFactor_;
85  uint tag = phiS.index();
86  int bx = phiS.bxNum() - 20;
87  int quality = 3;
88  uint tfLayer = phiS.stNum() - 1;
89  int eta = -16384;
90  if (station == 1) {
91  eta = coarseEta1_[abswheel];
92  } else if (station == 2) {
93  eta = coarseEta2_[abswheel];
94  } else if (station == 3) {
95  eta = coarseEta3_[abswheel];
96  } else if (station == 4) {
97  eta = coarseEta4_[abswheel];
98  }
99 
100  //override!!!
101  // eta=abswheel;
102 
103  //Now full eta
104 
105  eta = eta * sign;
106  l1t::MuonStub stub(wheel, sector, station, tfLayer, phi, phiB, tag, bx, quality, eta, 0, 0, 1);
107  stub.setOfflineQuantities(globalPhi, float(phiB), eta * etaLSB_, 0.0);
108  return stub;
109 }
string quality
#define M_PI

◆ calculateEta()

int L1TPhase2GMTBarrelStubProcessor::calculateEta ( uint  i,
int  wheel,
uint  sector,
uint  station 
)
private

Definition at line 167 of file L1TPhase2GMTBarrelStubProcessor.cc.

References PVValHelper::eta, eta1_, eta2_, eta3_, mps_fire::i, nano_mu_digi_cff::sector, relativeConstraints::station, and makeMuonMisalignmentScenario::wheel.

Referenced by buildStub().

167  {
168  int eta = 0;
169  if (wheel > 0) {
170  eta = 7 * wheel + 3 - i;
171  } else if (wheel < 0) {
172  eta = 7 * wheel + i - 3;
173  } else {
174  if (sector == 0 || sector == 3 || sector == 4 || sector == 7 || sector == 8 || sector == 11)
175  eta = i - 3;
176  else
177  eta = 3 - i;
178  }
179 
180  if (station == 1)
181  eta = eta1_[eta + 17];
182  else if (station == 2)
183  eta = eta2_[eta + 17];
184  else
185  eta = eta3_[eta + 17];
186 
187  return eta;
188 }

◆ makeStubs()

l1t::MuonStubCollection L1TPhase2GMTBarrelStubProcessor::makeStubs ( const L1Phase2MuDTPhContainer phiContainer,
const L1MuDTChambThContainer etaContainer 
)

Definition at line 111 of file L1TPhase2GMTBarrelStubProcessor.cc.

References buildStub(), buildStubNoEta(), nano_mu_digi_cff::bx, L1MuDTChambThContainer::chThetaSegm(), L1Phase2MuDTPhContainer::getContainer(), maxBX_, minBX_, minPhiQuality_, MillePedeFileConverter_cfg::out, nano_mu_digi_cff::sector, relativeConstraints::station, verbose_, and makeMuonMisalignmentScenario::wheel.

Referenced by Phase2L1TGMTStubProducer::produce().

112  {
114  for (int bx = minBX_; bx <= maxBX_; bx++) {
115  for (int wheel = -2; wheel <= 2; wheel++) {
116  for (int sector = 0; sector < 12; sector++) {
117  for (int station = 1; station < 5; station++) {
118  bool hasEta = false;
119  const L1MuDTChambThDigi* tseta = etaContainer->chThetaSegm(wheel, station, sector, bx);
120  if (tseta != nullptr) {
121  hasEta = true;
122  }
123 
124  for (const auto& phiDigi : *phiContainer->getContainer()) {
125  if ((phiDigi.bxNum() - 20) != bx || phiDigi.whNum() != wheel || phiDigi.scNum() != sector ||
126  phiDigi.stNum() != station)
127  continue;
128  if (phiDigi.quality() < minPhiQuality_)
129  continue;
130  if (hasEta) {
131  out.push_back(buildStub(phiDigi, tseta));
132  } else {
133  out.push_back(buildStubNoEta(phiDigi));
134  }
135  }
136  }
137  }
138  }
139  }
140 
141  if (verbose_) {
142  printf("Barrel Stubs\n");
143  for (const auto& stub : out)
144  printf(
145  "Barrel Stub bx=%d TF=%d etaRegion=%d phiRegion=%d depthRegion=%d coord1=%f,%d coord2=%f,%d eta1=%f,%d "
146  "eta2=%f,%d quality=%d etaQuality=%d\n",
147  stub.bxNum(),
148  stub.tfLayer(),
149  stub.etaRegion(),
150  stub.phiRegion(),
151  stub.depthRegion(),
152  stub.offline_coord1(),
153  stub.coord1(),
154  stub.offline_coord2(),
155  stub.coord2(),
156  stub.offline_eta1(),
157  stub.eta1(),
158  stub.offline_eta2(),
159  stub.eta2(),
160  stub.quality(),
161  stub.etaQuality());
162  }
163 
164  return out;
165 }
L1MuDTChambThDigi const * chThetaSegm(int wheel, int stat, int sect, int bx) const
std::vector< MuonStub > MuonStubCollection
Definition: MuonStub.h:40
l1t::MuonStub buildStub(const L1Phase2MuDTPhDigi &, const L1MuDTChambThDigi *)
l1t::MuonStub buildStubNoEta(const L1Phase2MuDTPhDigi &)
Segment_Container const * getContainer() const

Member Data Documentation

◆ coarseEta1_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::coarseEta1_
private

Definition at line 35 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ coarseEta2_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::coarseEta2_
private

Definition at line 36 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ coarseEta3_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::coarseEta3_
private

Definition at line 37 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ coarseEta4_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::coarseEta4_
private

Definition at line 38 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ eta1_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::eta1_
private

Definition at line 32 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by calculateEta().

◆ eta2_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::eta2_
private

Definition at line 33 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by calculateEta().

◆ eta3_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::eta3_
private

Definition at line 34 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by calculateEta().

◆ etaLSB_

double L1TPhase2GMTBarrelStubProcessor::etaLSB_
private

Definition at line 43 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStub(), and buildStubNoEta().

◆ maxBX_

int L1TPhase2GMTBarrelStubProcessor::maxBX_
private

Definition at line 30 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by makeStubs().

◆ minBX_

int L1TPhase2GMTBarrelStubProcessor::minBX_
private

Definition at line 29 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by makeStubs().

◆ minPhiQuality_

int L1TPhase2GMTBarrelStubProcessor::minPhiQuality_
private

Definition at line 27 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by makeStubs().

◆ phiBFactor_

int L1TPhase2GMTBarrelStubProcessor::phiBFactor_
private

Definition at line 40 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ phiLSB_

double L1TPhase2GMTBarrelStubProcessor::phiLSB_
private

Definition at line 42 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ phiOffset_

std::vector<int> L1TPhase2GMTBarrelStubProcessor::phiOffset_
private

Definition at line 39 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by buildStubNoEta().

◆ verbose_

int L1TPhase2GMTBarrelStubProcessor::verbose_
private

Definition at line 41 of file L1TPhase2GMTBarrelStubProcessor.h.

Referenced by makeStubs().