CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFClusterSelector.cc
Go to the documentation of this file.
2 #include <memory>
3 
7 
8 
9 using namespace std;
10 using namespace edm;
11 
12 
14  clusters_(consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("src"))),
15  energyRanges_(iConfig.getParameter<std::vector<double> >("energyRanges")),
16  timingCutsLowBarrel_(iConfig.getParameter<std::vector<double> >("timingCutsLowBarrel")),
17  timingCutsHighBarrel_(iConfig.getParameter<std::vector<double> >("timingCutsHighBarrel")),
18  timingCutsLowEndcap_(iConfig.getParameter<std::vector<double> >("timingCutsLowEndcap")),
19  timingCutsHighEndcap_(iConfig.getParameter<std::vector<double> >("timingCutsHighEndcap"))
20 {
21  produces<reco::PFClusterCollection>();
22 }
23 
24 
26  const edm::EventSetup& iSetup) {
27 
29  iEvent.getByToken(clusters_,clusters);
30  std::auto_ptr<reco::PFClusterCollection> out(new reco::PFClusterCollection);
31 
32  const std::vector<double>* timingCutsLow = NULL;
33  const std::vector<double>* timingCutsHigh = NULL;
34  for(const auto& cluster : *clusters ) {
35  const double energy = cluster.energy();
36  const double time = cluster.time();
37 
38  switch(cluster.layer()) {
42  timingCutsLow = &timingCutsLowBarrel_;
43  timingCutsHigh = &timingCutsHighBarrel_;
44  break;
47  timingCutsLow = &timingCutsLowEndcap_;
48  timingCutsHigh = &timingCutsHighEndcap_;
49  break;
50  default:
51  continue;
52  }
53  auto e_bin = std::lower_bound(energyRanges_.begin(),energyRanges_.end(),energy);
54  // Returns iterator to the first value that is *greater* than the value
55  // we are comparing to.
56  // The timing cuts are indexed by the bin high edge so we just need the
57  // distance between this bin timing cut vectors are padded to avoid
58  // overflows.
59  const unsigned idx = std::distance(energyRanges_.begin(),e_bin);
60  if( time > (*timingCutsLow)[idx] && time < (*timingCutsHigh)[idx] ) {
61  out->push_back(cluster);
62  }
63  }
64  iEvent.put( out);
65 
66 }
67 
69 
70 // ------------ method called once each job just before starting event loop ------------
71 void
73  const EventSetup& es) {
74 
75 
76 }
77 
78 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::vector< double > timingCutsHighBarrel_
#define NULL
Definition: scimark2.h:8
std::vector< double > timingCutsLowBarrel_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
edm::EDGetTokenT< reco::PFClusterCollection > clusters_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
tuple out
Definition: dbtoconf.py:99
std::vector< double > timingCutsHighEndcap_
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
std::vector< double > energyRanges_
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
std::vector< double > timingCutsLowEndcap_
virtual void beginRun(const edm::Run &run, const edm::EventSetup &es)
PFClusterSelector(const edm::ParameterSet &)
Definition: Run.h:41