File indexing completed on 2023-03-17 11:09:43
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027
0028 #include "HLTrigger/Muon/test/MuEnrichType1Filter.h"
0029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042 MuEnrichType1Filter::MuEnrichType1Filter(const edm::ParameterSet& iConfig) : nrejected{0}, naccepted{0} {
0043 type = iConfig.getParameter<int>("type");
0044 theGenToken = consumes<edm::HepMCProduct>(edm::InputTag("VtxSmeared"));
0045 }
0046
0047 MuEnrichType1Filter::~MuEnrichType1Filter() {
0048
0049
0050 }
0051
0052
0053
0054
0055
0056
0057 bool MuEnrichType1Filter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0058 using namespace edm;
0059 using namespace std;
0060 using namespace HepMC;
0061
0062 Handle<HepMCProduct> EvtHandle;
0063 iEvent.getByToken(theGenToken, EvtHandle);
0064 const GenEvent* Evt = EvtHandle->GetEvent();
0065 if (Evt != nullptr) {
0066 edm::LogVerbatim("MuEnrichFltr") << "------------------------------";
0067 for (HepMC::GenEvent::particle_const_iterator part = Evt->particles_begin(); part != Evt->particles_end(); ++part) {
0068 if (abs((*part)->pdg_id()) == 13) {
0069 double pt = (*part)->momentum().perp();
0070 edm::LogVerbatim("MuEnrichFltr") << "Found a muon with pt" << pt;
0071 if (pt > 4 && type == 1) {
0072 nrejected++;
0073 return false;
0074 break;
0075 } else if (pt > 2 && pt < 4 && type == 2) {
0076 nrejected++;
0077 return false;
0078 break;
0079 } else if (pt > 2 && pt < 10 && type == 3) {
0080 nrejected++;
0081 return false;
0082 break;
0083 }
0084 }
0085 }
0086 }
0087 naccepted++;
0088 return true;
0089 }
0090
0091
0092 void MuEnrichType1Filter::beginJob() {}
0093
0094
0095 void MuEnrichType1Filter::endJob() {
0096 edm::LogVerbatim("MuEnrichFltr") << "Total events" << naccepted + nrejected;
0097 edm::LogVerbatim("MuEnrichFltr") << "Acccepted events" << naccepted;
0098 edm::LogVerbatim("MuEnrichFltr") << "Rejected events" << nrejected;
0099 }
0100
0101
0102 DEFINE_FWK_MODULE(MuEnrichType1Filter);