File indexing completed on 2024-04-06 12:13:42
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include "DataFormats/Common/interface/Handle.h"
0020 #include "FWCore/Framework/interface/Event.h"
0021 #include "FWCore/Framework/interface/Frameworkfwd.h"
0022 #include "FWCore/Framework/interface/global/EDFilter.h"
0023 #include "FWCore/Framework/interface/MakerMacros.h"
0024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0025 #include "FWCore/Utilities/interface/EDGetToken.h"
0026 #include "FWCore/Utilities/interface/InputTag.h"
0027 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0028
0029 #include "TLorentzVector.h"
0030
0031 #include <cmath>
0032 #include <cstdlib>
0033 #include <vector>
0034
0035 class ZgMassFilter : public edm::global::EDFilter<> {
0036 public:
0037 explicit ZgMassFilter(const edm::ParameterSet&);
0038
0039 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0040
0041 private:
0042 const edm::EDGetTokenT<edm::HepMCProduct> token_;
0043 const double minDileptonMass;
0044 const double minZgMass;
0045 };
0046
0047 using namespace edm;
0048 using namespace std;
0049
0050 ZgMassFilter::ZgMassFilter(const edm::ParameterSet& iConfig)
0051 : token_(consumes<edm::HepMCProduct>(
0052 iConfig.getUntrackedParameter("moduleLabel", edm::InputTag("generator", "unsmeared")))),
0053 minDileptonMass(iConfig.getUntrackedParameter("MinDileptonMass", 0.)),
0054 minZgMass(iConfig.getUntrackedParameter("MinZgMass", 0.)) {}
0055
0056 bool ZgMassFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0057 bool accepted = false;
0058 Handle<HepMCProduct> evt;
0059 iEvent.getByToken(token_, evt);
0060 const HepMC::GenEvent* myGenEvent = evt->GetEvent();
0061
0062 vector<TLorentzVector> Lepton;
0063 Lepton.clear();
0064 vector<TLorentzVector> Photon;
0065 Photon.clear();
0066
0067 for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
0068 ++p) {
0069 if ((*p)->status() == 1 && (abs((*p)->pdg_id()) == 11 || abs((*p)->pdg_id()) == 13 || abs((*p)->pdg_id()) == 15)) {
0070 TLorentzVector LeptP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
0071 Lepton.push_back(LeptP);
0072 }
0073 if (abs((*p)->pdg_id()) == 22 && (*p)->status() == 1) {
0074 TLorentzVector PhotP((*p)->momentum().px(), (*p)->momentum().py(), (*p)->momentum().pz(), (*p)->momentum().e());
0075 Photon.push_back(PhotP);
0076 }
0077 }
0078
0079 if (Lepton.size() > 1 && (Lepton[0] + Lepton[1]).M() > minDileptonMass) {
0080 if ((Lepton[0] + Lepton[1] + Photon[0]).M() > minZgMass) {
0081 accepted = true;
0082 }
0083 }
0084
0085 return accepted;
0086 }
0087
0088 DEFINE_FWK_MODULE(ZgMassFilter);