CMS 3D CMS Logo

L1TkMuMantra.h
Go to the documentation of this file.
1 #ifndef L1Trigger_L1TTrackMatch_L1TKMUMANTRA_H
2 #define L1Trigger_L1TTrackMatch_L1TKMUMANTRA_H
3 
4 /*
5 ** class : GenericDataFormat
6 ** author : L.Cadamuro (UF)
7 ** date : 4/11/2019
8 ** brief : very generic structs to be used as inputs to the correlator
9 ** : to make sure that Mantra can handle muons and tracks from all the detectors
10 */
11 
12 namespace L1TkMuMantraDF {
13 
14  struct track_df {
15  double pt; // GeV
16  double eta; // rad, -inf / +inf
17  double theta; // rad, 0 -> +90-90
18  double phi; // rad, -pi / + pi
19  int nstubs; //
20  double chi2; //
21  int charge; // -1. +1
22  };
23 
24  struct muon_df {
25  double pt; // GeV
26  double eta; // rad, -inf / +inf
27  double theta; // rad, 0 -> +90-90
28  double phi; // rad, -pi / + pi
29  int charge; // -1. +1
30  };
31 } // namespace L1TkMuMantraDF
32 
33 /*
34 ** class : L1TkMuMantra
35 ** author : L.Cadamuro (UF)
36 ** date : 4/11/2019
37 ** brief : correlates muons and tracks using pre-encoded windows
38 */
39 
40 #include <iostream>
41 #include <vector>
42 #include <string>
43 #include <utility>
45 #include <cmath>
46 
47 #include "TFile.h"
50 
51 class L1TkMuMantra {
52 public:
53  L1TkMuMantra(const std::vector<double>& bounds, TFile* fIn_theta, TFile* fIn_phi, std::string name);
55 
56  // returns a vector with the same size of muons, each with an index to the matched L1 track, or -1 if no match is found
57  std::vector<int> find_match(const std::vector<L1TkMuMantraDF::track_df>& tracks,
58  const std::vector<L1TkMuMantraDF::muon_df>& muons);
59 
60  void test(double eta, double pt);
61 
62  void relax_windows(double& low, double cent, double& high); // will modify low and high
63 
64  void set_safety_factor(float sf_l, float sf_h) {
65  safety_factor_l_ = sf_l;
66  safety_factor_h_ = sf_h;
67  if (verbosity_ > 0)
68  LogTrace("L1TkMuMantra") << name_ << " safety factor LOW is " << safety_factor_l_ << std::endl;
69  if (verbosity_ > 0)
70  LogTrace("L1TkMuMantra") << name_ << " safety factor HIGH is " << safety_factor_h_ << std::endl;
71  }
72 
73  int sign(double x) {
74  if (x == 0)
75  return 1;
76  return (0 < x) - (x < 0);
77  }
78 
79  void setArbitrationType(std::string type); // MaxPt, MinDeltaPt
80 
81  // static functions, meant to be used from outside to interface with MAnTra
82  static std::vector<double> prepare_corr_bounds(std::string fname, std::string hname);
83 
84  // converters
85  static double eta_to_theta(double x) {
86  // give theta in rad
87  return (2. * atan(exp(-1. * x)));
88  }
89 
90  static double to_mpio2_pio2(double x) {
91  // put the angle in radians between -pi/2 and pi/2
92  while (x >= 0.5 * M_PI)
93  x -= M_PI;
94  while (x < -0.5 * M_PI)
95  x += M_PI;
96  return x;
97  }
98 
99 private:
100  int findBin(double val);
101 
103 
104  int nbins_; // counts the number of MuMatchWindow = bounds_.size() - 1
105  std::vector<double> bounds_; // counts the boundaries of the MuMatchWindow (in eta/theta)
106  std::vector<MuMatchWindow> wdws_theta_;
107  std::vector<MuMatchWindow> wdws_phi_;
108 
109  int min_nstubs = 4; // >= min_nstubs
110  double max_chi2 = 100; // < max_chi2
111 
112  float safety_factor_l_; // increase the lower theta/phi threshold by this fractions w.r.t. the center
113  float safety_factor_h_; // increase the upper theta/phi threshold by this fractions w.r.t. the center
114 
115  enum sortParType {
116  kMaxPt, // pick the highest pt track matched
117  kMinDeltaPt // pick the track with the smallest pt difference w.r.t the muon
118  };
119 
121 
122  int verbosity_ = 0;
123 };
124 
125 #endif // L1TKMUMANTRA_H
static double eta_to_theta(double x)
Definition: L1TkMuMantra.h:85
static double to_mpio2_pio2(double x)
Definition: L1TkMuMantra.h:90
int findBin(double val)
Definition: L1TkMuMantra.cc:82
sortParType sort_type_
Definition: L1TkMuMantra.h:120
L1TkMuMantra(const std::vector< double > &bounds, TFile *fIn_theta, TFile *fIn_phi, std::string name)
Definition: L1TkMuMantra.cc:10
double max_chi2
Definition: L1TkMuMantra.h:110
std::string name_
Definition: L1TkMuMantra.h:102
#define LogTrace(id)
void relax_windows(double &low, double cent, double &high)
static std::vector< double > prepare_corr_bounds(std::string fname, std::string hname)
std::vector< MuMatchWindow > wdws_phi_
Definition: L1TkMuMantra.h:107
int sign(double x)
Definition: L1TkMuMantra.h:73
#define M_PI
float safety_factor_h_
Definition: L1TkMuMantra.h:113
float safety_factor_l_
Definition: L1TkMuMantra.h:112
void set_safety_factor(float sf_l, float sf_h)
Definition: L1TkMuMantra.h:64
auto const & tracks
cannot be loose
string fname
main script
void test(double eta, double pt)
Definition: L1TkMuMantra.cc:99
std::vector< MuMatchWindow > wdws_theta_
Definition: L1TkMuMantra.h:106
std::vector< int > find_match(const std::vector< L1TkMuMantraDF::track_df > &tracks, const std::vector< L1TkMuMantraDF::muon_df > &muons)
void setArbitrationType(std::string type)
std::vector< double > bounds_
Definition: L1TkMuMantra.h:105