Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    LHEGenericFilter
0004 // Class:      LHEGenericFilter
0005 //
0006 /*
0007 
0008  Description: Filter to select events with an arbitrary number of given particle(s).
0009 
0010  Implementation: derived from MCSingleParticleFilter
0011 
0012 */
0013 //
0014 // Original Author:  Roberto Covarelli
0015 //         Created:  Wed Feb 29 04:22:16 CST 2012
0016 //
0017 //
0018 
0019 #include "DataFormats/Common/interface/Handle.h"
0020 #include "FWCore/Framework/interface/Frameworkfwd.h"
0021 #include "FWCore/Framework/interface/global/EDFilter.h"
0022 #include "FWCore/Framework/interface/Event.h"
0023 #include "FWCore/Framework/interface/MakerMacros.h"
0024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0026 #include "FWCore/Utilities/interface/EDGetToken.h"
0027 #include "FWCore/Utilities/interface/InputTag.h"
0028 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
0029 
0030 #include <cmath>
0031 #include <cstdlib>
0032 #include <string>
0033 #include <vector>
0034 
0035 //
0036 // class declaration
0037 //
0038 
0039 class LHEGenericFilter : public edm::global::EDFilter<> {
0040 public:
0041   explicit LHEGenericFilter(const edm::ParameterSet&);
0042 
0043 private:
0044   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0045 
0046   // ----------member data ---------------------------
0047 
0048   const edm::EDGetTokenT<LHEEventProduct> src_;
0049   const int numRequired_;              // number of particles required to pass filter
0050   const std::vector<int> particleID_;  // vector of particle IDs to look for
0051   enum Logic { LT, GT, EQ, NE };
0052   Logic whichlogic;
0053 };
0054 
0055 using namespace std;
0056 
0057 LHEGenericFilter::LHEGenericFilter(const edm::ParameterSet& iConfig)
0058     : src_(consumes<LHEEventProduct>(iConfig.getParameter<edm::InputTag>("src"))),
0059       numRequired_(iConfig.getParameter<int>("NumRequired")),
0060       particleID_(iConfig.getParameter<std::vector<int> >("ParticleID")) {
0061   // LT  meaning <
0062   // GT          >
0063   // EQ          =
0064   // NE          !=
0065   std::string acceptLogic = iConfig.getParameter<std::string>("AcceptLogic");
0066   if (acceptLogic == "LT")
0067     whichlogic = LT;
0068   else if (acceptLogic == "GT")
0069     whichlogic = GT;
0070   else if (acceptLogic == "EQ")
0071     whichlogic = EQ;
0072   else if (acceptLogic == "NE")
0073     whichlogic = NE;
0074   else {
0075     edm::LogError("cat_A") << "wrong input for AcceptLogic string";
0076     // at least initialize it to something reproducible
0077     whichlogic = LT;
0078   }
0079 }
0080 
0081 // ------------ method called to skim the data  ------------
0082 bool LHEGenericFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0083   edm::Handle<LHEEventProduct> EvtHandle;
0084   iEvent.getByToken(src_, EvtHandle);
0085 
0086   int nFound = 0;
0087 
0088   for (int i = 0; i < EvtHandle->hepeup().NUP; ++i) {
0089     if (EvtHandle->hepeup().ISTUP[i] != 1) {  // keep only outgoing particles
0090       continue;
0091     }
0092     for (unsigned int j = 0; j < particleID_.size(); ++j) {
0093       if (particleID_[j] == 0 || abs(particleID_[j]) == abs(EvtHandle->hepeup().IDUP[i])) {
0094         nFound++;
0095         break;  // only match a given particle once!
0096       }
0097     }  // loop over targets
0098 
0099   }  // loop over particles
0100 
0101   // event accept/reject logic
0102   if ((whichlogic == LT && nFound < numRequired_) || (whichlogic == GT && nFound > numRequired_) ||
0103       (whichlogic == EQ && nFound == numRequired_) || (whichlogic == NE && nFound != numRequired_)) {
0104     return true;
0105   } else {
0106     return false;
0107   }
0108 }
0109 
0110 //define this as a plug-in
0111 DEFINE_FWK_MODULE(LHEGenericFilter);