File indexing completed on 2024-09-08 23:51:46
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054 #include "FWCore/Framework/interface/Event.h"
0055 #include "FWCore/Framework/interface/Frameworkfwd.h"
0056 #include "FWCore/Framework/interface/stream/EDFilter.h"
0057 #include "FWCore/Framework/interface/MakerMacros.h"
0058 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0059 #include "FWCore/Utilities/interface/EDGetToken.h"
0060 #include "FWCore/Utilities/interface/InputTag.h"
0061 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0062
0063 #include "TH2.h"
0064
0065 #include <cmath>
0066 #include <cstdlib>
0067
0068 class HerwigMaxPtPartonFilter : public edm::stream::EDFilter<> {
0069 public:
0070 explicit HerwigMaxPtPartonFilter(const edm::ParameterSet&);
0071
0072 bool filter(edm::Event&, const edm::EventSetup&) override;
0073
0074 private:
0075 const edm::EDGetTokenT<edm::HepMCProduct> token_;
0076
0077 const double minptcut;
0078 const double maxptcut;
0079 const int processID;
0080
0081 TH2D hFSPartons_JS_PtWgting;
0082 };
0083
0084 using namespace edm;
0085 using namespace std;
0086
0087 HerwigMaxPtPartonFilter::HerwigMaxPtPartonFilter(const edm::ParameterSet& iConfig)
0088 : token_(consumes<edm::HepMCProduct>(
0089 iConfig.getUntrackedParameter("moduleLabel", edm::InputTag("generator", "unsmeared")))),
0090 minptcut(iConfig.getUntrackedParameter("MinPt", 0.)),
0091 maxptcut(iConfig.getUntrackedParameter("MaxPt", 10000.)),
0092 processID(iConfig.getUntrackedParameter("ProcessID", 0)),
0093 hFSPartons_JS_PtWgting(
0094 "hFSPartons_JS_PtWgting", "#phi-#eta Space of FS Partons (p_{T} wgt'ing)", 20, -5.205, 5.205, 32, -M_PI, M_PI) {
0095 }
0096
0097 bool HerwigMaxPtPartonFilter::filter(edm::Event& iEvent, const edm::EventSetup&) {
0098
0099 hFSPartons_JS_PtWgting.Reset();
0100
0101 bool accepted = false;
0102 bool isFSQuark = false;
0103 double maxPartonPt = 0.0;
0104
0105
0106
0107 int pos1stCluster = 0;
0108
0109
0110 long counter = 0;
0111
0112 Handle<HepMCProduct> evt;
0113 iEvent.getByToken(token_, evt);
0114
0115 const HepMC::GenEvent* myGenEvent = evt->GetEvent();
0116
0117 if (processID == 0 || processID == myGenEvent->signal_process_id()) {
0118
0119
0120 for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
0121 ++p) {
0122 if (abs((*p)->pdg_id()) == 91) {
0123 break;
0124 }
0125 pos1stCluster++;
0126 }
0127
0128
0129 for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
0130 ++p) {
0131
0132 if ((*p)->momentum().perp() > 1.0) {
0133
0134 if (abs((*p)->pdg_id()) == 1 || abs((*p)->pdg_id()) == 2 || abs((*p)->pdg_id()) == 3 ||
0135 abs((*p)->pdg_id()) == 4 || abs((*p)->pdg_id()) == 5) {
0136 if (counter < pos1stCluster && ((*p)->status() == 158 || (*p)->status() == 159)) {
0137 isFSQuark = true;
0138 }
0139 }
0140 }
0141
0142 if (isFSQuark) {
0143 hFSPartons_JS_PtWgting.Fill(
0144 (*p)->momentum().eta(), (*p)->momentum().phi(), (*p)->momentum().perp());
0145 }
0146
0147 counter++;
0148 isFSQuark = false;
0149 }
0150
0151 maxPartonPt = hFSPartons_JS_PtWgting.GetMaximum();
0152
0153
0154 if (maxPartonPt >= minptcut && maxPartonPt < maxptcut) {
0155 accepted = true;
0156
0157 }
0158 }
0159
0160 else {
0161 accepted = true;
0162 }
0163 return accepted;
0164 }
0165
0166 DEFINE_FWK_MODULE(HerwigMaxPtPartonFilter);