CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZgammaMassFilter.cc
Go to the documentation of this file.
3 #include <iostream>
4 #include "TLorentzVector.h"
5 
6 // order std::vector of TLorentzVector elements
7 class orderByPt {
8 public:
9  bool operator()(TLorentzVector const & a, TLorentzVector const & b) {
10  if (a.Pt() == b.Pt() ) {
11  return a.Pt() < b.Pt();
12  } else {
13  return a.Pt() > b.Pt() ;
14  }
15  }
16 };
17 
18 
19 using namespace edm;
20 using namespace std;
21 
23 {
24  token_ =consumes<edm::HepMCProduct>(iConfig.getParameter<InputTag>("HepMCProduct"));
25  minPhotonPt =iConfig.getParameter<double>("minPhotonPt");
26  minLeptonPt =iConfig.getParameter<double>("minLeptonPt");
27  minPhotonEta =iConfig.getParameter<double>("minPhotonEta");
28  minLeptonEta =iConfig.getParameter<double>("minLeptonEta");
29  maxPhotonEta =iConfig.getParameter<double>("maxPhotonEta");
30  maxLeptonEta =iConfig.getParameter<double>("maxLeptonEta");
31  minDileptonMass =iConfig.getParameter<double>("minDileptonMass");
32  minZgMass =iConfig.getParameter<double>("minZgMass");
33 }
34 
36 {
37 }
38 
39 // ------------ method called to skim the data ------------
41 {
42  using namespace edm;
43 
44  bool accepted = false;
46  iEvent.getByToken(token_, evt);
47  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
48 
49  vector<TLorentzVector> Lepton; Lepton.clear();
50  vector<TLorentzVector> Photon; Photon.clear();
51  vector<float> Charge; Charge.clear();
52 
53  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
54  p != myGenEvent->particles_end(); ++p )
55  {
56  if ((*p)->status() == 1 && (abs((*p)->pdg_id()) == 11 || abs((*p)->pdg_id()) == 13 || abs((*p)->pdg_id()) == 15)) {
57  TLorentzVector LeptP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
58  if (LeptP.Pt() > minLeptonPt) {
59  Lepton.push_back(LeptP);
60  } // if pt
61  }// if lepton
62 
63  if ( abs((*p)->pdg_id()) == 22 && (*p)->status() == 1) {
64  TLorentzVector PhotP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
65  if (PhotP.Pt() > minPhotonPt) {
66  Photon.push_back(PhotP);
67  }// if pt
68  }// if photon
69 
70  }// loop over particles
71 
72 
73  // std::cout << "\n" << "Photon size: " << Photon.size() << std::endl;
74  // for (unsigned int u=0; u<Photon.size(); u++){
75  // std::cout << "BEF photon PT: " << Photon[u].Pt() << std::endl;
76  // }
77  // std::cout << "\n" << "Lepton size: " << Lepton.size() << std::endl;
78  // for (unsigned int u=0; u<Lepton.size(); u++){
79  // std::cout << "BEF lepton PT: " << Lepton[u].Pt() << std::endl;
80  // }
81 
82  // order Lepton and Photon according to Pt
83  std::stable_sort(Photon.begin(), Photon.end(), orderByPt());
84  std::stable_sort(Lepton.begin(), Lepton.end(), orderByPt());
85 
86  // std::cout << "\n" << std::endl;
87  // std::cout << "\n" << "Photon size: " << Photon.size() << std::endl;
88  // for (unsigned int u=0; u<Photon.size(); u++){
89  // std::cout << "AFT photon PT: " << Photon[u].Pt() << std::endl;
90  // }
91  // std::cout << "\n" << "Lepton size: " << Lepton.size() << std::endl;
92  // for (unsigned int u=0; u<Lepton.size(); u++){
93  // std::cout << "AFT lepton PT: " << Lepton[u].Pt() << std::endl;
94  // }
95  // std::cout << "\n" << std::endl;
96 
97  if (
98  Photon.size() > 0 && Lepton.size() > 1 &&
99  Photon[0].Pt() > minPhotonPt &&
100  Lepton[0].Pt() > minLeptonPt &&
101  Lepton[1].Pt() > minLeptonPt &&
102  Photon[0].Eta() > minPhotonEta &&
103  Lepton[0].Eta() > minLeptonEta &&
104  Lepton[1].Eta() > minLeptonEta &&
105  Photon[0].Eta() < maxPhotonEta &&
106  Lepton[0].Eta() < maxLeptonEta &&
107  Lepton[1].Eta() < maxLeptonEta &&
108  (Lepton[0]+Lepton[1]).M() > minDileptonMass &&
109  (Lepton[0]+Lepton[1]+Photon[0]).M() > minZgMass
110  )
111  { // satisfy molteplicity, kinematics, and ll llg minimum mass
112  accepted = true;
113  }
114 
115  // std::cout << "++ returning: " << accepted << "\n" << std::endl;
116 
117  return accepted;
118 }
T getParameter(std::string const &) const
bool operator()(TLorentzVector const &a, TLorentzVector const &b)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
int iEvent
Definition: GenABIO.cc:230
virtual bool filter(edm::Event &, const edm::EventSetup &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
ZgammaMassFilter(const edm::ParameterSet &)