1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
|
// -*- C++ -*-
//
// Package: ZgMassFilter
// Class: ZgMassFilter
//
/*
Description: filter events based on the Pythia particle information
Implementation: inherits from generic EDFilter
*/
//
// Original Author: Alexey Ferapontov
// Created: Thu July 26 11:57:54 CDT 2012
//
//
#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDFilter.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/EDGetToken.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
#include "TLorentzVector.h"
#include <cmath>
#include <cstdlib>
#include <vector>
class ZgMassFilter : public edm::global::EDFilter<> {
public:
explicit ZgMassFilter(const edm::ParameterSet&);
bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
private:
const edm::EDGetTokenT<edm::HepMCProduct> token_;
const double minDileptonMass;
const double minZgMass;
};
using namespace edm;
using namespace std;
ZgMassFilter::ZgMassFilter(const edm::ParameterSet& iConfig)
: token_(consumes<edm::HepMCProduct>(
iConfig.getUntrackedParameter("moduleLabel", edm::InputTag("generator", "unsmeared")))),
minDileptonMass(iConfig.getUntrackedParameter("MinDileptonMass", 0.)),
minZgMass(iConfig.getUntrackedParameter("MinZgMass", 0.)) {}
bool ZgMassFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
bool accepted = false;
Handle<HepMCProduct> evt;
iEvent.getByToken(token_, evt);
const HepMC::GenEvent* myGenEvent = evt->GetEvent();
vector<TLorentzVector> Lepton;
Lepton.clear();
vector<TLorentzVector> Photon;
Photon.clear();
for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
++p) {
if ((*p)->status() == 1 && (abs((*p)->pdg_id()) == 11 || abs((*p)->pdg_id()) == 13 || abs((*p)->pdg_id()) == 15)) {
TLorentzVector LeptP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
Lepton.push_back(LeptP);
}
if (abs((*p)->pdg_id()) == 22 && (*p)->status() == 1) {
TLorentzVector PhotP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
Photon.push_back(PhotP);
}
}
if (Lepton.size() > 1 && (Lepton[0] + Lepton[1]).M() > minDileptonMass) {
if ((Lepton[0] + Lepton[1] + Photon[0]).M() > minZgMass) {
accepted = true;
}
}
return accepted;
}
DEFINE_FWK_MODULE(ZgMassFilter);
|