CMS 3D CMS Logo

PFClusterTimeSelector.cc
Go to the documentation of this file.
11 
12 #include <memory>
13 #include <vector>
14 
16 public:
18  ~PFClusterTimeSelector() override;
19 
20  void beginRun(const edm::Run& run, const edm::EventSetup& es) override;
21 
22  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
23 
24 protected:
25  struct CutInfo {
26  double depth;
27  double minE;
28  double maxE;
29  double minTime;
30  double maxTime;
31  bool endcap;
32  };
33 
34  // ----------access to event data
36  std::vector<CutInfo> cutInfo_;
37 };
38 
41 
42 using namespace std;
43 using namespace edm;
44 
46  : clusters_(consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("src"))) {
47  std::vector<edm::ParameterSet> cuts = iConfig.getParameter<std::vector<edm::ParameterSet> >("cuts");
48  for (const auto& cut : cuts) {
49  CutInfo info;
50  info.depth = cut.getParameter<double>("depth");
51  info.minE = cut.getParameter<double>("minEnergy");
52  info.maxE = cut.getParameter<double>("maxEnergy");
53  info.minTime = cut.getParameter<double>("minTime");
54  info.maxTime = cut.getParameter<double>("maxTime");
55  info.endcap = cut.getParameter<bool>("endcap");
56  cutInfo_.push_back(info);
57  }
58 
59  produces<reco::PFClusterCollection>();
60  produces<reco::PFClusterCollection>("OOT");
61 }
62 
65  iEvent.getByToken(clusters_, clusters);
66  auto out = std::make_unique<reco::PFClusterCollection>();
67  auto outOOT = std::make_unique<reco::PFClusterCollection>();
68 
69  for (const auto& cluster : *clusters) {
70  const double energy = cluster.energy();
71  const double time = cluster.time();
72  const double depth = cluster.depth();
73  const PFLayer::Layer layer = cluster.layer();
74  for (const auto& info : cutInfo_) {
75  if (energy < info.minE || energy > info.maxE)
76  continue;
77  if (depth < 0.9 * info.depth || depth > 1.1 * info.depth)
78  continue;
79  if ((info.endcap &&
81  (((!info.endcap) && (layer == PFLayer::ECAL_ENDCAP || layer == PFLayer::HCAL_ENDCAP))))
82  continue;
83 
84  if (time > info.minTime && time < info.maxTime)
85  out->push_back(cluster);
86  else
87  outOOT->push_back(cluster);
88 
89  break;
90  }
91  }
92 
93  iEvent.put(std::move(out));
94  iEvent.put(std::move(outOOT), "OOT");
95 }
96 
98 
99 // ------------ method called once each job just before starting event loop ------------
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
static const TGPicture * info(bool iBackgroundIsBlack)
std::vector< CutInfo > cutInfo_
PFClusterTimeSelector(const edm::ParameterSet &)
~PFClusterTimeSelector() override
void beginRun(const edm::Run &run, const edm::EventSetup &es) override
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
int iEvent
Definition: GenABIO.cc:224
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Layer
layer definition
Definition: PFLayer.h:29
std::vector< l1t::PFCluster > PFClusterCollection
Definition: PFCluster.h:74
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::PFClusterCollection > clusters_
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45