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  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
24 
25 protected:
26  struct CutInfo {
27  double depth;
28  double minE;
29  double maxE;
30  double minTime;
31  double maxTime;
32  bool endcap;
33  };
34 
35  // ----------access to event data
37  std::vector<CutInfo> cutInfo_;
38 };
39 
42 
45  desc.add<edm::InputTag>("src", {"particleFlowClusterECALWithTimeUncorrected"});
46  {
47  std::vector<edm::ParameterSet> vpset;
48  vpset.reserve(10);
49  {
51  pset.addParameter<double>("depth", 1.0);
52  pset.addParameter<double>("minEnergy", 0.0);
53  pset.addParameter<double>("maxEnergy", 1.0);
54  pset.addParameter<bool>("endcap", false);
55  pset.addParameter<double>("minTime", -12.);
56  pset.addParameter<double>("maxTime", 12.);
57  vpset.emplace_back(pset);
58  }
59  {
61  pset.addParameter<double>("depth", 1.0);
62  pset.addParameter<double>("minEnergy", 0.0);
63  pset.addParameter<double>("maxEnergy", 1.0);
64  pset.addParameter<bool>("endcap", true);
65  pset.addParameter<double>("minTime", -31.5);
66  pset.addParameter<double>("maxTime", 31.5);
67  vpset.emplace_back(pset);
68  }
69  {
71  pset.addParameter<double>("depth", 1.0);
72  pset.addParameter<double>("minEnergy", 1.0);
73  pset.addParameter<double>("maxEnergy", 2.0);
74  pset.addParameter<bool>("endcap", false);
75  pset.addParameter<double>("minTime", -6.);
76  pset.addParameter<double>("maxTime", 6.);
77  vpset.emplace_back(pset);
78  }
79  {
81  pset.addParameter<double>("depth", 1.0);
82  pset.addParameter<double>("minEnergy", 1.0);
83  pset.addParameter<double>("maxEnergy", 2.0);
84  pset.addParameter<bool>("endcap", true);
85  pset.addParameter<double>("minTime", -20.5);
86  pset.addParameter<double>("maxTime", 20.5);
87  vpset.emplace_back(pset);
88  }
89  {
91  pset.addParameter<double>("depth", 1.0);
92  pset.addParameter<double>("minEnergy", 2.0);
93  pset.addParameter<double>("maxEnergy", 5.0);
94  pset.addParameter<bool>("endcap", false);
95  pset.addParameter<double>("minTime", -4.);
96  pset.addParameter<double>("maxTime", 4.);
97  vpset.emplace_back(pset);
98  }
99  {
101  pset.addParameter<double>("depth", 1.0);
102  pset.addParameter<double>("minEnergy", 2.0);
103  pset.addParameter<double>("maxEnergy", 5.0);
104  pset.addParameter<bool>("endcap", true);
105  pset.addParameter<double>("minTime", -12.);
106  pset.addParameter<double>("maxTime", 12.);
107  vpset.emplace_back(pset);
108  }
109  {
111  pset.addParameter<double>("depth", 1.0);
112  pset.addParameter<double>("minEnergy", 5.0);
113  pset.addParameter<double>("maxEnergy", 20.0);
114  pset.addParameter<bool>("endcap", false);
115  pset.addParameter<double>("minTime", -4.);
116  pset.addParameter<double>("maxTime", 4.);
117  vpset.emplace_back(pset);
118  }
119  {
121  pset.addParameter<double>("depth", 1.0);
122  pset.addParameter<double>("minEnergy", 5.0);
123  pset.addParameter<double>("maxEnergy", 20.0);
124  pset.addParameter<bool>("endcap", true);
125  pset.addParameter<double>("minTime", -5.);
126  pset.addParameter<double>("maxTime", 5.);
127  vpset.emplace_back(pset);
128  }
129  {
131  pset.addParameter<double>("depth", 1.0);
132  pset.addParameter<double>("minEnergy", 20.0);
133  pset.addParameter<double>("maxEnergy", 1e24);
134  pset.addParameter<bool>("endcap", false);
135  pset.addParameter<double>("minTime", -4.);
136  pset.addParameter<double>("maxTime", 4.);
137  vpset.emplace_back(pset);
138  }
139  {
141  pset.addParameter<double>("depth", 1.0);
142  pset.addParameter<double>("minEnergy", 20.0);
143  pset.addParameter<double>("maxEnergy", 1e24);
144  pset.addParameter<bool>("endcap", true);
145  pset.addParameter<double>("minTime", -5.);
146  pset.addParameter<double>("maxTime", 5.);
147  vpset.emplace_back(pset);
148  }
150  psd.add<double>("depth", 1.0);
151  psd.add<double>("minEnergy", 0.0);
152  psd.add<double>("maxEnergy", 1e24);
153  psd.add<bool>("endcap", false);
154  psd.add<double>("minTime", -50.);
155  psd.add<double>("maxTime", 50.);
156  desc.addVPSet("cuts", psd, vpset);
157  }
158  descriptions.add("particleFlowClusterECALTimeSelected", desc);
159 }
160 
161 using namespace std;
162 using namespace edm;
163 
165  : clusters_(consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("src"))) {
166  std::vector<edm::ParameterSet> cuts = iConfig.getParameter<std::vector<edm::ParameterSet> >("cuts");
167  for (const auto& cut : cuts) {
168  CutInfo info;
169  info.depth = cut.getParameter<double>("depth");
170  info.minE = cut.getParameter<double>("minEnergy");
171  info.maxE = cut.getParameter<double>("maxEnergy");
172  info.minTime = cut.getParameter<double>("minTime");
173  info.maxTime = cut.getParameter<double>("maxTime");
174  info.endcap = cut.getParameter<bool>("endcap");
175  cutInfo_.push_back(info);
176  }
177 
178  produces<reco::PFClusterCollection>();
179  produces<reco::PFClusterCollection>("OOT");
180 }
181 
184  iEvent.getByToken(clusters_, clusters);
185  auto out = std::make_unique<reco::PFClusterCollection>();
186  auto outOOT = std::make_unique<reco::PFClusterCollection>();
187 
188  for (const auto& cluster : *clusters) {
189  const double energy = cluster.energy();
190  const double time = cluster.time();
191  const double depth = cluster.depth();
192  const PFLayer::Layer layer = cluster.layer();
193  for (const auto& info : cutInfo_) {
194  if (energy < info.minE || energy > info.maxE)
195  continue;
196  if (depth < 0.9 * info.depth || depth > 1.1 * info.depth)
197  continue;
198  if ((info.endcap &&
200  (((!info.endcap) && (layer == PFLayer::ECAL_ENDCAP || layer == PFLayer::HCAL_ENDCAP))))
201  continue;
202 
203  if (time > info.minTime && time < info.maxTime)
204  out->push_back(cluster);
205  else
206  outOOT->push_back(cluster);
207 
208  break;
209  }
210  }
211 
212  iEvent.put(std::move(out));
213  iEvent.put(std::move(outOOT), "OOT");
214 }
215 
217 
218 // ------------ method called once each job just before starting event loop ------------
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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
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
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Layer
layer definition
Definition: PFLayer.h:29
std::vector< l1t::PFCluster > PFClusterCollection
Definition: PFCluster.h:87
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::PFClusterCollection > clusters_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45