File indexing completed on 2024-04-06 12:13:39
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/ConfigurationDescriptions.h"
0025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0026 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0027 #include "FWCore/Utilities/interface/EDGetToken.h"
0028 #include "FWCore/Utilities/interface/InputTag.h"
0029 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0030
0031 #include <cmath>
0032 #include <cstdlib>
0033 #include <string>
0034
0035 class DJpsiFilter : public edm::global::EDFilter<> {
0036 public:
0037 explicit DJpsiFilter(const edm::ParameterSet&);
0038
0039 static void fillDescriptions(edm::ConfigurationDescriptions&);
0040
0041 private:
0042 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0043
0044 const edm::EDGetTokenT<edm::HepMCProduct> token_;
0045 const double minPt;
0046 const double maxY;
0047 const double maxPt;
0048 const double minY;
0049 const int status;
0050 const int particleID;
0051 };
0052
0053 DJpsiFilter::DJpsiFilter(const edm::ParameterSet& iConfig)
0054 : token_(consumes<edm::HepMCProduct>(
0055 edm::InputTag(iConfig.getUntrackedParameter("moduleLabel", std::string("generator")), "unsmeared"))),
0056 minPt(iConfig.getUntrackedParameter("MinPt", 0.)),
0057 maxY(iConfig.getUntrackedParameter("MaxY", 10.)),
0058 maxPt(iConfig.getUntrackedParameter("MaxPt", 1000.)),
0059 minY(iConfig.getUntrackedParameter("MinY", 0.)),
0060 status(iConfig.getUntrackedParameter("Status", 0)),
0061 particleID(iConfig.getUntrackedParameter("ParticleID", 0)) {}
0062
0063 bool DJpsiFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0064 bool accepted = false;
0065 int n2jpsi = 0;
0066 edm::Handle<edm::HepMCProduct> evt;
0067 iEvent.getByToken(token_, evt);
0068 const HepMC::GenEvent* myGenEvent = evt->GetEvent();
0069
0070 for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
0071 ++p) {
0072 if ((*p)->status() != status)
0073 continue;
0074 double energy = (*p)->momentum().e();
0075 double pz = (*p)->momentum().pz();
0076 double momentumY = 0.5 * std::log((energy + pz) / (energy - pz));
0077 if ((*p)->momentum().perp() > minPt && std::fabs(momentumY) < maxY && (*p)->momentum().perp() < maxPt &&
0078 std::fabs(momentumY) > minY) {
0079 if (std::abs((*p)->pdg_id()) == particleID)
0080 n2jpsi++;
0081 }
0082 if (n2jpsi >= 2) {
0083 accepted = true;
0084 break;
0085 }
0086 }
0087 return accepted;
0088 }
0089
0090
0091 void DJpsiFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0092
0093
0094 edm::ParameterSetDescription desc;
0095 desc.setUnknown();
0096 descriptions.addDefault(desc);
0097 }
0098
0099 DEFINE_FWK_MODULE(DJpsiFilter);