CMS 3D CMS Logo

JetVertexMain.cc
Go to the documentation of this file.
1 
7 #include <cmath>
8 #include <string>
9 using namespace reco;
10 using namespace edm;
11 
13  cutSigmaZ = parameters.getParameter<double>("JV_sigmaZ");
14  cutDeltaZ = parameters.getParameter<double>("JV_deltaZ");
15  threshold = parameters.getParameter<double>("JV_alpha_threshold");
16  cone_size = parameters.getParameter<double>("JV_cone_size");
17  Algo = parameters.getParameter<int>("JV_type_Algo");
18  cutType = parameters.getParameter<std::string>("JV_cutType");
19 }
20 
21 std::pair<double, bool> JetVertexMain::Main(const reco::CaloJet& jet,
23  double signal_vert_Z,
24  double signal_vert_z_error) const {
25  std::pair<double, bool> parameter;
26 
27  double jet_et = jet.et();
28  double jet_phi = jet.phi();
29  double jet_eta = jet.eta();
30 
31  // cout<<"JET: "<<jet_et<<endl;
32  double Pt_jets_X = 0.;
33  double Pt_jets_Y = 0.;
34  double Pt_jets_X_tot = 0.;
35  double Pt_jets_Y_tot = 0.;
36 
37  TrackCollection::const_iterator track = tracks->begin();
38 
39  if (!tracks->empty()) {
40  for (; track != tracks->end(); track++) {
41  double Vertex_Z = track->vz();
42  double Vertex_Z_Error = track->dzError();
43  double track_eta = track->eta();
44  double track_phi = track->phi();
45 
46  if (DeltaR(track_eta, jet_eta, track_phi, jet_phi) < cone_size) {
47  double DeltaZ = Vertex_Z - signal_vert_Z;
48  double DeltaZ_Error = sqrt((Vertex_Z_Error * Vertex_Z_Error) + (signal_vert_z_error * signal_vert_z_error));
49  Pt_jets_X_tot += track->px();
50  Pt_jets_Y_tot += track->py();
51  bool discriminator;
52  if (cutType == "sig")
53  discriminator = (fabs(DeltaZ) / DeltaZ_Error) <= cutSigmaZ;
54  else
55  discriminator = fabs(DeltaZ) < cutDeltaZ;
56 
57  if (discriminator) {
58  Pt_jets_X += track->px();
59  Pt_jets_Y += track->py();
60  }
61  }
62  }
63  }
64  double Var = -1;
65 
66  if (Algo == 1)
67  Var = Track_Pt(Pt_jets_X, Pt_jets_Y) / jet_et;
68  else if (Algo == 2) {
69  if (Track_Pt(Pt_jets_X_tot, Pt_jets_Y_tot) != 0)
70  Var = Track_Pt(Pt_jets_X, Pt_jets_Y) / Track_Pt(Pt_jets_X_tot, Pt_jets_Y_tot);
71  else
72  std::cout << "[Jets] JetVertexAssociation: Warning! problems for Algo = 2: possible division by zero .."
73  << std::endl;
74  } else {
75  Var = Track_Pt(Pt_jets_X, Pt_jets_Y) / jet_et;
76  std::cout << "[Jets] JetVertexAssociation: Warning! Algo = " << Algo << " not found; using Algo = 1" << std::endl;
77  }
78 
79  // cout<<"Var = "<<Var<<endl;
80 
81  if (Var >= threshold)
82  parameter = std::pair<double, bool>(Var, true);
83  else
84  parameter = std::pair<double, bool>(Var, false);
85 
86  return parameter;
87 }
88 
89 double JetVertexMain::DeltaR(double eta1, double eta2, double phi1, double phi2) const {
90  double dphi = fabs(phi1 - phi2);
91  if (dphi > M_PI)
92  dphi = 2 * M_PI - dphi;
93  double deta = fabs(eta1 - eta2);
94  return sqrt(dphi * dphi + deta * deta);
95 }
96 
97 double JetVertexMain::Track_Pt(double px, double py) const { return sqrt(px * px + py * py); }
Jets made from CaloTowers.
Definition: CaloJet.h:27
def Var(expr, valtype, doc=None, precision=-1, lazyEval=False)
Definition: common_cff.py:17
double Track_Pt(double px, double py) const
T sqrt(T t)
Definition: SSEVec.h:19
double DeltaR(double eta1, double eta2, double phi1, double phi2) const
JetVertexMain(const edm::ParameterSet &parameters)
#define M_PI
std::pair< double, bool > Main(const reco::CaloJet &jet, edm::Handle< reco::TrackCollection > tracks, double SIGNAL_V_Z, double SIGNAL_V_Z_Error) const
fixed size matrix
HLT enums.
Definition: fakeMenu.h:6