CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZgMassFilter.cc
Go to the documentation of this file.
3 #include <iostream>
4 #include "TLorentzVector.h"
5 
6 using namespace edm;
7 using namespace std;
8 
10 token_(consumes<edm::HepMCProduct>(iConfig.getUntrackedParameter("moduleLabel",edm::InputTag("generator","unsmeared")))),
11 minDileptonMass(iConfig.getUntrackedParameter("MinDileptonMass", 0.)),
12 minZgMass(iConfig.getUntrackedParameter("MinZgMass", 0.))
13 {
14 }
15 
17 {
18 }
19 
20 // ------------ method called to skim the data ------------
22 {
23  using namespace edm;
24  bool accepted = false;
26  iEvent.getByToken(token_, evt);
27  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
28 
29  vector<TLorentzVector> Lepton; Lepton.clear();
30  vector<TLorentzVector> Photon; Photon.clear();
31 
32  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
33  p != myGenEvent->particles_end(); ++p ) {
34  if ((*p)->status() == 1 && (abs((*p)->pdg_id()) == 11 || abs((*p)->pdg_id()) == 13 || abs((*p)->pdg_id()) == 15)) {
35  TLorentzVector LeptP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
36  Lepton.push_back(LeptP);
37  }
38  if ( abs((*p)->pdg_id()) == 22 && (*p)->status() == 1) {
39  TLorentzVector PhotP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
40  Photon.push_back(PhotP);
41  }
42  }
43 
44  if (Lepton.size() > 1 && (Lepton[0]+Lepton[1]).M() > minDileptonMass ) {
45  if ((Lepton[0]+Lepton[1]+Photon[0]).M() > minZgMass) {
46  accepted = true;
47  }
48  }
49 
50  return accepted;
51 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
double minDileptonMass
Definition: ZgMassFilter.h:56
ZgMassFilter(const edm::ParameterSet &)
Definition: ZgMassFilter.cc:9
edm::EDGetTokenT< edm::HepMCProduct > token_
Definition: ZgMassFilter.h:55
int iEvent
Definition: GenABIO.cc:230
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double minZgMass
Definition: ZgMassFilter.h:57
virtual bool filter(edm::Event &, const edm::EventSetup &)
Definition: ZgMassFilter.cc:21