Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 10:01:40

0001 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0002 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
0003 #include "FWCore/Framework/interface/ESHandle.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/Framework/interface/stream/EDProducer.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include <memory>
0013 #include <vector>
0014 
0015 class PFClusterTimeSelector : public edm::stream::EDProducer<> {
0016 public:
0017   explicit PFClusterTimeSelector(const edm::ParameterSet&);
0018   ~PFClusterTimeSelector() override;
0019 
0020   void beginRun(const edm::Run& run, const edm::EventSetup& es) override;
0021 
0022   void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0023 
0024 protected:
0025   struct CutInfo {
0026     double depth;
0027     double minE;
0028     double maxE;
0029     double minTime;
0030     double maxTime;
0031     bool endcap;
0032   };
0033 
0034   // ----------access to event data
0035   edm::EDGetTokenT<reco::PFClusterCollection> clusters_;
0036   std::vector<CutInfo> cutInfo_;
0037 };
0038 
0039 #include "FWCore/Framework/interface/MakerMacros.h"
0040 DEFINE_FWK_MODULE(PFClusterTimeSelector);
0041 
0042 using namespace std;
0043 using namespace edm;
0044 
0045 PFClusterTimeSelector::PFClusterTimeSelector(const edm::ParameterSet& iConfig)
0046     : clusters_(consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("src"))) {
0047   std::vector<edm::ParameterSet> cuts = iConfig.getParameter<std::vector<edm::ParameterSet> >("cuts");
0048   for (const auto& cut : cuts) {
0049     CutInfo info;
0050     info.depth = cut.getParameter<double>("depth");
0051     info.minE = cut.getParameter<double>("minEnergy");
0052     info.maxE = cut.getParameter<double>("maxEnergy");
0053     info.minTime = cut.getParameter<double>("minTime");
0054     info.maxTime = cut.getParameter<double>("maxTime");
0055     info.endcap = cut.getParameter<bool>("endcap");
0056     cutInfo_.push_back(info);
0057   }
0058 
0059   produces<reco::PFClusterCollection>();
0060   produces<reco::PFClusterCollection>("OOT");
0061 }
0062 
0063 void PFClusterTimeSelector::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0064   edm::Handle<reco::PFClusterCollection> clusters;
0065   iEvent.getByToken(clusters_, clusters);
0066   auto out = std::make_unique<reco::PFClusterCollection>();
0067   auto outOOT = std::make_unique<reco::PFClusterCollection>();
0068 
0069   for (const auto& cluster : *clusters) {
0070     const double energy = cluster.energy();
0071     const double time = cluster.time();
0072     const double depth = cluster.depth();
0073     const PFLayer::Layer layer = cluster.layer();
0074     for (const auto& info : cutInfo_) {
0075       if (energy < info.minE || energy > info.maxE)
0076         continue;
0077       if (depth < 0.9 * info.depth || depth > 1.1 * info.depth)
0078         continue;
0079       if ((info.endcap &&
0080            (layer == PFLayer::ECAL_BARREL || layer == PFLayer::HCAL_BARREL1 || layer == PFLayer::HCAL_BARREL2)) ||
0081           (((!info.endcap) && (layer == PFLayer::ECAL_ENDCAP || layer == PFLayer::HCAL_ENDCAP))))
0082         continue;
0083 
0084       if (time > info.minTime && time < info.maxTime)
0085         out->push_back(cluster);
0086       else
0087         outOOT->push_back(cluster);
0088 
0089       break;
0090     }
0091   }
0092 
0093   iEvent.put(std::move(out));
0094   iEvent.put(std::move(outOOT), "OOT");
0095 }
0096 
0097 PFClusterTimeSelector::~PFClusterTimeSelector() = default;
0098 
0099 // ------------ method called once each job just before starting event loop  ------------
0100 void PFClusterTimeSelector::beginRun(const edm::Run& run, const EventSetup& es) {}