CMS 3D CMS Logo

TtSemiLepSignalSel.cc
Go to the documentation of this file.
2 #include "TVector3.h"
3 
5 
6 TtSemiLepSignalSel::TtSemiLepSignalSel(const std::vector<pat::Jet>& topJets,
7  const math::XYZTLorentzVector& lepton,
8  const edm::View<pat::MET>& MET) { //function
9 
10  unsigned int nJetsMax = topJets.size();
11 
12  var_MET = MET.begin()->et();
13  var_sumEt = 0.;
14 
15  math::XYZTLorentzVector Jetsum(0., 0., 0., 0.);
16 
17  for (unsigned int i = 0; i < nJetsMax; i++) {
18  math::XYZTLorentzVector aJet = topJets[i].p4();
19  Jetsum += aJet;
20  var_sumEt += topJets[i].et();
21  }
22  massalljets = Jetsum.M();
23 
24  var_lepeta = lepton.Eta();
25 
26  math::XYZTLorentzVector Met = MET.begin()->p4();
27  const math::XYZTLorentzVector& Lep = lepton;
28  double Etjet[4];
29  double Jetjet[6];
30  double dijetmass;
31  var_mindijetmass = 99999.;
32  var_maxdijetmass = -1.;
33  int counter = 0;
34  for (int i = 0; i < 4; i++) {
35  math::XYZTLorentzVector aJet = topJets[i].p4();
36  Etjet[i] = aJet.Et();
37  for (int j = i + 1; j < 4; j++) {
38  math::XYZTLorentzVector asecJet = topJets[j].p4();
39  dijetmass = (aJet + asecJet).M();
40  if (dijetmass < var_mindijetmass)
41  var_mindijetmass = dijetmass;
42  if (dijetmass > var_maxdijetmass)
43  var_maxdijetmass = dijetmass;
44  counter++;
45  }
46  }
47 
48  var_Et1 = Etjet[0];
49 
51 
52  counter = 0;
53  for (int i = 0; i < 4; i++) {
54  math::XYZTLorentzVector aJet = topJets[i].p4();
55  for (int j = i + 1; j < 4; j++) {
56  math::XYZTLorentzVector asecJet = topJets[j].p4();
57  Jetjet[counter] = fabs(aJet.Eta() - asecJet.Eta());
58  counter++;
59  }
60  }
61 
62  var_detajet2jet3 = Jetjet[3];
63  var_detajet3jet4 = Jetjet[5];
64 
65  double Lepjet[4];
66  var_mindRjetlepton = 99999.;
67  for (int i = 0; i < 4; i++) {
68  math::XYZTLorentzVector aJet = topJets[i].p4();
69  Lepjet[i] = DeltaR(Lep, aJet);
70  if (Lepjet[i] < var_mindRjetlepton)
71  var_mindRjetlepton = Lepjet[i];
72  }
73 }
74 
76  double dPhi = fabs(v1.Phi() - v2.Phi());
77  if (dPhi > TMath::Pi())
78  dPhi = 2 * TMath::Pi() - dPhi;
79  return dPhi;
80 }
81 
83  double dPhi = DeltaPhi(v1, v2);
84  double dR = TMath::Sqrt((v1.Eta() - v2.Eta()) * (v1.Eta() - v2.Eta()) + dPhi * dPhi);
85  return dR;
86 }
87 
counter
Definition: counter.py:1
mps_fire.i
i
Definition: mps_fire.py:355
TtSemiLepSignalSel::var_maxdijetmass
double var_maxdijetmass
Definition: TtSemiLepSignalSel.h:49
TtSemiLepSignalSel::var_lepeta
double var_lepeta
Definition: TtSemiLepSignalSel.h:40
TtSemiLepSignalSel::massalljets
double massalljets
Definition: TtSemiLepSignalSel.h:53
TtSemiLepSignalSel::TtSemiLepSignalSel
TtSemiLepSignalSel()
Definition: TtSemiLepSignalSel.cc:4
TtSemiLepSignalSel::var_Et1
double var_Et1
Definition: TtSemiLepSignalSel.h:39
TtSemiLepSignalSel::var_MET
double var_MET
Definition: TtSemiLepSignalSel.h:41
TtSemiLepSignalSel::var_detajet3jet4
double var_detajet3jet4
Definition: TtSemiLepSignalSel.h:46
TtSemiLepSignalSel::~TtSemiLepSignalSel
~TtSemiLepSignalSel()
Definition: TtSemiLepSignalSel.cc:88
HLT_2018_cff.dPhi
dPhi
Definition: HLT_2018_cff.py:12290
edm::View
Definition: CaloClusterFwd.h:14
TtSemiLepSignalSel::var_mindijetmass
double var_mindijetmass
Definition: TtSemiLepSignalSel.h:48
TtSemiLepSignalSel.h
counter
static std::atomic< unsigned int > counter
Definition: SharedResourceNames.cc:15
TtFullLepDaughter::Lep
static const std::string Lep
Definition: TtFullLeptonicEvent.h:10
TtSemiLepSignalSel::var_dphiMETlepton
double var_dphiMETlepton
Definition: TtSemiLepSignalSel.h:43
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
MET
Pi
const double Pi
Definition: CosmicMuonParameters.h:18
TtSemiLepSignalSel::DeltaR
double DeltaR(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
Definition: TtSemiLepSignalSel.cc:82
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
TtSemiLepSignalSel::var_detajet2jet3
double var_detajet2jet3
Definition: TtSemiLepSignalSel.h:45
TtSemiLepSignalSel::var_sumEt
double var_sumEt
Definition: TtSemiLepSignalSel.h:38
TtSemiLepSignalSel::DeltaPhi
double DeltaPhi(const math::XYZTLorentzVector &, const math::XYZTLorentzVector &)
Definition: TtSemiLepSignalSel.cc:75
TtSemiLepSignalSel::var_mindRjetlepton
double var_mindRjetlepton
Definition: TtSemiLepSignalSel.h:51