Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:13:40

0001 // -*- C++ -*-
0002 //
0003 // Package:    LHEVpTFilter
0004 // Class:      LHEVpTFilter
0005 //
0006 /*
0007 
0008  Description: Filter to select events with V pT in a given range.
0009  (Based on LHEGenericFilter)
0010 
0011 */
0012 //
0013 
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "FWCore/Framework/interface/Frameworkfwd.h"
0016 #include "FWCore/Framework/interface/global/EDFilter.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 #include "FWCore/Utilities/interface/EDGetToken.h"
0021 #include "FWCore/Utilities/interface/InputTag.h"
0022 
0023 #include "Math/Vector4D.h"
0024 #include "Math/Vector4Dfwd.h"
0025 
0026 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
0027 
0028 #include <cmath>
0029 #include <cstdlib>
0030 #include <vector>
0031 
0032 //
0033 // class declaration
0034 //
0035 
0036 class LHEVpTFilter : public edm::global::EDFilter<> {
0037 public:
0038   explicit LHEVpTFilter(const edm::ParameterSet&);
0039 
0040 private:
0041   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0042 
0043   // ----------member data ---------------------------
0044 
0045   const edm::EDGetTokenT<LHEEventProduct> src_;
0046   const double vptMin_;  // number of particles required to pass filter
0047   const double vptMax_;  // number of particles required to pass filter
0048 };
0049 
0050 LHEVpTFilter::LHEVpTFilter(const edm::ParameterSet& iConfig)
0051     : src_(consumes<LHEEventProduct>(iConfig.getParameter<edm::InputTag>("src"))),
0052       vptMin_(iConfig.getParameter<double>("VpTMin")),
0053       vptMax_(iConfig.getParameter<double>("VpTMax")) {}
0054 
0055 // ------------ method called to skim the data  ------------
0056 bool LHEVpTFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0057   edm::Handle<LHEEventProduct> EvtHandle;
0058   iEvent.getByToken(src_, EvtHandle);
0059 
0060   std::vector<lhef::HEPEUP::FiveVector> const& lheParticles = EvtHandle->hepeup().PUP;
0061 
0062   std::vector<ROOT::Math::PxPyPzEVector> lepCands;
0063   for (unsigned int i = 0; i < lheParticles.size(); ++i) {
0064     if (EvtHandle->hepeup().ISTUP[i] != 1) {  // keep only outgoing particles
0065       continue;
0066     }
0067     unsigned absPdgId = std::abs(EvtHandle->hepeup().IDUP[i]);
0068     if (absPdgId >= 11 && absPdgId <= 16) {
0069       lepCands.push_back(
0070           ROOT::Math::PxPyPzEVector(lheParticles[i][0], lheParticles[i][1], lheParticles[i][2], lheParticles[i][3]));
0071     }
0072   }
0073   double vpt_ = -1;
0074   if (lepCands.size() == 2) {
0075     vpt_ = (lepCands[0] + lepCands[1]).pt();
0076   }
0077   if (vpt_ < vptMax_ && vpt_ >= vptMin_) {
0078     return true;
0079   } else {
0080     return false;
0081   }
0082 }
0083 
0084 //define this as a plug-in
0085 DEFINE_FWK_MODULE(LHEVpTFilter);