CMS 3D CMS Logo

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