CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
JetVertexMain Class Reference

#include <JetVertexMain.h>

Public Member Functions

 JetVertexMain (const edm::ParameterSet &parameters)
 
std::pair< double, bool > Main (const reco::CaloJet &jet, edm::Handle< reco::TrackCollection > tracks, double SIGNAL_V_Z, double SIGNAL_V_Z_Error)
 
 ~JetVertexMain ()
 

Private Member Functions

double DeltaR (double eta1, double eta2, double phi1, double phi2)
 
double Track_Pt (double px, double py)
 

Private Attributes

int Algo
 
double cone_size
 
double cutDeltaZ
 
double cutSigmaZ
 
std::string cutType
 
bool discriminator
 
double threshold
 

Detailed Description

Definition at line 13 of file JetVertexMain.h.

Constructor & Destructor Documentation

JetVertexMain::JetVertexMain ( const edm::ParameterSet parameters)

Definition at line 12 of file JetVertexMain.cc.

References edm::ParameterSet::getParameter(), and dtT0WireCalibration_cfg::threshold.

12  {
13 
14  cutSigmaZ = parameters.getParameter<double>("JV_sigmaZ");
15  cutDeltaZ = parameters.getParameter<double>("JV_deltaZ");
16  threshold = parameters.getParameter<double>("JV_alpha_threshold");
17  cone_size = parameters.getParameter<double>("JV_cone_size");
18  Algo = parameters.getParameter<int>("JV_type_Algo");
19  cutType = parameters.getParameter<std::string>("JV_cutType");
20 
21 }
T getParameter(std::string const &) const
double cutSigmaZ
Definition: JetVertexMain.h:30
double cone_size
Definition: JetVertexMain.h:33
std::string cutType
Definition: JetVertexMain.h:35
double threshold
Definition: JetVertexMain.h:32
Definition: fakeMenu.h:4
double cutDeltaZ
Definition: JetVertexMain.h:31
JetVertexMain::~JetVertexMain ( )
inline

Definition at line 19 of file JetVertexMain.h.

19 {};

Member Function Documentation

double JetVertexMain::DeltaR ( double  eta1,
double  eta2,
double  phi1,
double  phi2 
)
private

Definition at line 89 of file JetVertexMain.cc.

References M_PI, and mathSSE::sqrt().

89  {
90 
91  double dphi = fabs(phi1-phi2);
92  if(dphi > M_PI) dphi = 2*M_PI - dphi;
93  double deta = fabs(eta1-eta2);
94  return sqrt(dphi*dphi + deta*deta);
95 
96 }
T sqrt(T t)
Definition: SSEVec.h:28
#define M_PI
Definition: BFit3D.cc:3
std::pair< double, bool > JetVertexMain::Main ( const reco::CaloJet jet,
edm::Handle< reco::TrackCollection tracks,
double  SIGNAL_V_Z,
double  SIGNAL_V_Z_Error 
)

Definition at line 24 of file JetVertexMain.cc.

References gather_cfg::cout, reco::LeafCandidate::et(), reco::LeafCandidate::eta(), reco::LeafCandidate::phi(), mathSSE::sqrt(), and dtT0WireCalibration_cfg::threshold.

Referenced by cms::JetVertexAssociation::produce().

25  {
26 
27  std::pair<double, bool> parameter;
28 
29  double jet_et = jet.et();
30  double jet_phi = jet.phi();
31  double jet_eta = jet.eta();
32 
33  // cout<<"JET: "<<jet_et<<endl;
34  double Pt_jets_X = 0. ;
35  double Pt_jets_Y = 0. ;
36  double Pt_jets_X_tot = 0. ;
37  double Pt_jets_Y_tot = 0. ;
38 
39  TrackCollection::const_iterator track = tracks->begin ();
40 
41  if (tracks->size() > 0 ) {
42  for (; track != tracks->end (); track++) {
43  double Vertex_Z = track->vz();
44  double Vertex_Z_Error = track->dzError();
45  double track_eta = track->eta();
46  double track_phi = track->phi();
47 
48  if (DeltaR(track_eta,jet_eta, track_phi, jet_phi) < cone_size) {
49 
50  double DeltaZ = Vertex_Z-signal_vert_Z;
51  double DeltaZ_Error = sqrt((Vertex_Z_Error*Vertex_Z_Error)+(signal_vert_z_error*signal_vert_z_error));
52  Pt_jets_X_tot += track->px();
53  Pt_jets_Y_tot += track->py();
54  if (cutType == "sig") discriminator = (fabs(DeltaZ)/DeltaZ_Error) <= cutSigmaZ;
55  else discriminator = fabs(DeltaZ) < cutDeltaZ;
56 
57  if (discriminator){
58 
59  Pt_jets_X += track->px();
60  Pt_jets_Y += track->py();
61 
62  }
63 
64  }
65  }
66 }
67  double Var = -1;
68 
69  if (Algo == 1) Var = Track_Pt(Pt_jets_X, Pt_jets_Y)/jet_et;
70  else if (Algo == 2) {
71  if (Track_Pt(Pt_jets_X_tot, Pt_jets_Y_tot)!=0) Var = Track_Pt(Pt_jets_X, Pt_jets_Y)/Track_Pt(Pt_jets_X_tot, Pt_jets_Y_tot);
72  else std::cout << "[Jets] JetVertexAssociation: Warning! problems for Algo = 2: possible division by zero .." << std::endl;
73  }
74  else {
75 
76  Var = Track_Pt(Pt_jets_X, Pt_jets_Y)/jet_et;
77  std::cout << "[Jets] JetVertexAssociation: Warning! Algo = " << Algo << " not found; using Algo = 1" << std::endl;
78  }
79 
80  // cout<<"Var = "<<Var<<endl;
81 
82  if (Var >= threshold) parameter = std::pair<double, bool>(Var, true);
83  else parameter = std::pair<double, bool>(Var, false);
84 
85  return parameter;
86 
87 }
double cutSigmaZ
Definition: JetVertexMain.h:30
virtual double et() const
transverse energy
double DeltaR(double eta1, double eta2, double phi1, double phi2)
double cone_size
Definition: JetVertexMain.h:33
virtual double eta() const
momentum pseudorapidity
double Track_Pt(double px, double py)
std::string cutType
Definition: JetVertexMain.h:35
double threshold
Definition: JetVertexMain.h:32
T sqrt(T t)
Definition: SSEVec.h:28
tuple cout
Definition: gather_cfg.py:41
virtual double phi() const
momentum azimuthal angle
Definition: fakeMenu.h:4
double cutDeltaZ
Definition: JetVertexMain.h:31
double JetVertexMain::Track_Pt ( double  px,
double  py 
)
private

Definition at line 98 of file JetVertexMain.cc.

References mathSSE::sqrt().

98  {
99 
100  return sqrt(px*px+py*py);
101 
102 }
T sqrt(T t)
Definition: SSEVec.h:28

Member Data Documentation

int JetVertexMain::Algo
private

Definition at line 34 of file JetVertexMain.h.

double JetVertexMain::cone_size
private

Definition at line 33 of file JetVertexMain.h.

double JetVertexMain::cutDeltaZ
private

Definition at line 31 of file JetVertexMain.h.

double JetVertexMain::cutSigmaZ
private

Definition at line 30 of file JetVertexMain.h.

std::string JetVertexMain::cutType
private

Definition at line 35 of file JetVertexMain.h.

bool JetVertexMain::discriminator
private

Definition at line 36 of file JetVertexMain.h.

double JetVertexMain::threshold
private

Definition at line 32 of file JetVertexMain.h.