CMS 3D CMS Logo

L1TMuonBarrelKalmanLUTs.cc
Go to the documentation of this file.
3 
6  lutFile_ = new TFile(path.fullPath().c_str());
7 
8  lut_[3 * 64 + 8] = (TH1*)lutFile_->Get("gain_8_3");
9  lut_[2 * 64 + 8] = (TH1*)lutFile_->Get("gain_8_2");
10  lut_[2 * 64 + 12] = (TH1*)lutFile_->Get("gain_12_2");
11  lut_[2 * 64 + 4] = (TH1*)lutFile_->Get("gain_4_2");
12  lut_[1 * 64 + 12] = (TH1*)lutFile_->Get("gain_12_1");
13  lut_[1 * 64 + 10] = (TH1*)lutFile_->Get("gain_10_1");
14  lut_[1 * 64 + 6] = (TH1*)lutFile_->Get("gain_6_1");
15  lut_[1 * 64 + 14] = (TH1*)lutFile_->Get("gain_14_1");
16  lut_[3] = (TH1*)lutFile_->Get("gain_3_0");
17  lut_[5] = (TH1*)lutFile_->Get("gain_5_0");
18  lut_[6] = (TH1*)lutFile_->Get("gain_6_0");
19  lut_[7] = (TH1*)lutFile_->Get("gain_7_0");
20  lut_[9] = (TH1*)lutFile_->Get("gain_9_0");
21  lut_[10] = (TH1*)lutFile_->Get("gain_10_0");
22  lut_[11] = (TH1*)lutFile_->Get("gain_11_0");
23  lut_[12] = (TH1*)lutFile_->Get("gain_12_0");
24  lut_[13] = (TH1*)lutFile_->Get("gain_13_0");
25  lut_[14] = (TH1*)lutFile_->Get("gain_14_0");
26  lut_[15] = (TH1*)lutFile_->Get("gain_15_0");
27 
28  lut2HH_[3 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_3_HH");
29  lut2HH_[2 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_2_HH");
30  lut2HH_[2 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_2_HH");
31  lut2HH_[1 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_1_HH");
32  lut2HH_[1 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_1_HH");
33  lut2HH_[1 * 64 + 2] = (TH1*)lutFile_->Get("gain2_2_1_HH");
34 
35  lut2LH_[3 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_3_LH");
36  lut2LH_[2 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_2_LH");
37  lut2LH_[2 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_2_LH");
38  lut2LH_[1 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_1_LH");
39  lut2LH_[1 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_1_LH");
40  lut2LH_[1 * 64 + 2] = (TH1*)lutFile_->Get("gain2_2_1_LH");
41 
42  lut2HL_[3 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_3_HL");
43  lut2HL_[2 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_2_HL");
44  lut2HL_[2 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_2_HL");
45  lut2HL_[1 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_1_HL");
46  lut2HL_[1 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_1_HL");
47  lut2HL_[1 * 64 + 2] = (TH1*)lutFile_->Get("gain2_2_1_HL");
48 
49  lut2LL_[3 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_3_LL");
50  lut2LL_[2 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_2_LL");
51  lut2LL_[2 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_2_LL");
52  lut2LL_[1 * 64 + 8] = (TH1*)lutFile_->Get("gain2_8_1_LL");
53  lut2LL_[1 * 64 + 4] = (TH1*)lutFile_->Get("gain2_4_1_LL");
54  lut2LL_[1 * 64 + 2] = (TH1*)lutFile_->Get("gain2_2_1_LL");
55 
56  coarseEta_[3] = (TH1*)lutFile_->Get("coarseEta_3");
57  coarseEta_[5] = (TH1*)lutFile_->Get("coarseEta_5");
58  coarseEta_[6] = (TH1*)lutFile_->Get("coarseEta_6");
59  coarseEta_[7] = (TH1*)lutFile_->Get("coarseEta_7");
60  coarseEta_[9] = (TH1*)lutFile_->Get("coarseEta_9");
61 
62  coarseEta_[10] = (TH1*)lutFile_->Get("coarseEta_10");
63  coarseEta_[11] = (TH1*)lutFile_->Get("coarseEta_11");
64  coarseEta_[12] = (TH1*)lutFile_->Get("coarseEta_12");
65  coarseEta_[13] = (TH1*)lutFile_->Get("coarseEta_13");
66  coarseEta_[14] = (TH1*)lutFile_->Get("coarseEta_14");
67  coarseEta_[15] = (TH1*)lutFile_->Get("coarseEta_15");
68 }
69 
71  lutFile_->Close();
72  if (lutFile_ != nullptr)
73  delete lutFile_;
74 }
75 
76 std::vector<float> L1TMuonBarrelKalmanLUTs::trackGain(uint step, uint bitmask, uint K) {
77  std::vector<float> gain(4, 0.0);
78  const TH1* h = lut_[64 * step + bitmask];
79  gain[0] = h->GetBinContent(K + 1);
80  gain[2] = h->GetBinContent(1024 + K + 1);
81  return gain;
82 }
83 
84 std::vector<float> L1TMuonBarrelKalmanLUTs::trackGain2(uint step, uint bitmask, uint K, uint qual1, uint qual2) {
85  std::vector<float> gain(4, 0.0);
86 
87  // printf("Track gain %d %d %d\n",step,bitmask,K);
88  const TH1* h;
89  if (qual1 < 4) {
90  if (qual2 < 4)
91  h = lut2LL_[64 * step + bitmask];
92  else
93  h = lut2LH_[64 * step + bitmask];
94  } else {
95  if (qual2 < 4)
96  h = lut2HL_[64 * step + bitmask];
97  else
98  h = lut2HH_[64 * step + bitmask];
99  }
100  gain[0] = h->GetBinContent(K + 1);
101  gain[1] = h->GetBinContent(512 + K + 1);
102  gain[2] = h->GetBinContent(2 * 512 + K + 1);
103  gain[3] = h->GetBinContent(3 * 512 + K + 1);
104  return gain;
105 }
106 
107 std::pair<float, float> L1TMuonBarrelKalmanLUTs::vertexGain(uint bitmask, uint K) {
108  const TH1* h = lut_[bitmask];
109  std::pair<float, float> gain(-h->GetBinContent(K + 1), -h->GetBinContent(1024 + K + 1));
110  return gain;
111 }
112 
114  const TH1* h = coarseEta_[pattern];
115  return uint(h->GetBinContent(h->GetXaxis()->FindBin(mask)));
116 }
L1TMuonBarrelKalmanLUTs(const std::string &)
constexpr uint32_t mask
Definition: gpuClustering.h:24
std::map< uint, const TH1 * > lut2HH_
std::pair< float, float > vertexGain(uint, uint)
std::map< uint, const TH1 * > coarseEta_
std::map< uint, const TH1 * > lut2HL_
std::map< uint, const TH1 * > lut2LH_
std::map< uint, const TH1 * > lut_
std::vector< float > trackGain2(uint, uint, uint, uint, uint)
std::map< uint, const TH1 * > lut2LL_
step
Definition: StallMonitor.cc:98
std::vector< float > trackGain(uint, uint, uint)
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4