CMS 3D CMS Logo

HGCalClusterProducer.cc
Go to the documentation of this file.
1 #ifndef __RecoLocalCalo_HGCRecProducers_HGCalClusterProducer_H__
2 #define __RecoLocalCalo_HGCRecProducers_HGCalClusterProducer_H__
3 
4 // user include files
8 
11 
18 
19 #include <memory>
20 #include <chrono>
21 
25 
29 
31 
33  public:
36 
37  virtual void produce(edm::Event&, const edm::EventSetup&);
38 
39  private:
40 
44 
46 
47  std::unique_ptr<HGCalImagingAlgo> algo;
48  std::unique_ptr<HGCal3DClustering> multicluster_algo;
49  bool doSharing;
51 
53 };
54 
56 
58  algoId(reco::CaloCluster::undefined),
59  doSharing(ps.getParameter<bool>("doSharing")),
60  detector(ps.getParameter<std::string >("detector")), // one of EE, FH, BH or "all"
61  verbosity((HGCalImagingAlgo::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity",3)){
62  double ecut = ps.getParameter<double>("ecut");
63  std::vector<double> vecDeltas = ps.getParameter<std::vector<double> >("deltac");
64  double kappa = ps.getParameter<double>("kappa");
65  std::vector<double> multicluster_radii = ps.getParameter<std::vector<double> >("multiclusterRadii");
66  double minClusters = ps.getParameter<unsigned>("minClusters");
67  std::vector<double> dEdXweights = ps.getParameter<std::vector<double> >("dEdXweights");
68  std::vector<double> thicknessCorrection = ps.getParameter<std::vector<double> >("thicknessCorrection");
69  std::vector<double> fcPerMip = ps.getParameter<std::vector<double> >("fcPerMip");
70  double fcPerEle = ps.getParameter<double>("fcPerEle");
71  std::vector<double> nonAgedNoises = ps.getParameter<std::vector<double> >("nonAgedNoises");
72  double noiseMip = ps.getParameter<double>("noiseMip");
73  bool dependSensor = ps.getParameter<bool>("dependSensor");
74 
75 
76  if(detector=="all") {
77  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
78  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
79  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
81  }else if(detector=="EE") {
82  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
84  }else if(detector=="FH") {
85  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
87  } else {
88  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
90  }
91 
92 
93  if(doSharing){
94  double showerSigma = ps.getParameter<double>("showerSigma");
95  algo = std::make_unique<HGCalImagingAlgo>(vecDeltas, kappa, ecut, showerSigma, algoId, dependSensor, dEdXweights, thicknessCorrection, fcPerMip, fcPerEle, nonAgedNoises, noiseMip, verbosity);
96  }else{
97  algo = std::make_unique<HGCalImagingAlgo>(vecDeltas, kappa, ecut, algoId, dependSensor, dEdXweights, thicknessCorrection, fcPerMip, fcPerEle, nonAgedNoises, noiseMip, verbosity);
98  }
99 
100  auto sumes = consumesCollector();
101 
102  multicluster_algo = std::make_unique<HGCal3DClustering>(ps, sumes, multicluster_radii, minClusters);
103 
104  produces<std::vector<reco::BasicCluster> >();
105  produces<std::vector<reco::BasicCluster> >("sharing");
106 
107  produces<std::vector<reco::HGCalMultiCluster> >();
108  produces<std::vector<reco::HGCalMultiCluster> >("sharing");
109 }
110 
112  const edm::EventSetup& es) {
113 
117 
118 
119  std::unique_ptr<std::vector<reco::BasicCluster> > clusters( new std::vector<reco::BasicCluster> ),
120  clusters_sharing( new std::vector<reco::BasicCluster> );
121 
122  algo->reset();
123 
124  algo->getEventSetup(es);
125 
126  multicluster_algo->getEvent(evt);
127  multicluster_algo->getEventSetup(es);
128 
129  switch(algoId){
131  evt.getByToken(hits_ee_token,ee_hits);
132  algo->populate(*ee_hits);
133  break;
135  evt.getByToken(hits_fh_token,fh_hits);
136  evt.getByToken(hits_bh_token,bh_hits);
137  if( fh_hits.isValid() ) {
138  algo->populate(*fh_hits);
139  } else if ( bh_hits.isValid() ) {
140  algo->populate(*bh_hits);
141  }
142  break;
144  evt.getByToken(hits_ee_token,ee_hits);
145  algo->populate(*ee_hits);
146  evt.getByToken(hits_fh_token,fh_hits);
147  algo->populate(*fh_hits);
148  evt.getByToken(hits_bh_token,bh_hits);
149  algo->populate(*bh_hits);
150  break;
151  default:
152  break;
153  }
154  algo->makeClusters();
155  *clusters = algo->getClusters(false);
156  if(doSharing)
157  *clusters_sharing = algo->getClusters(true);
158 
159  std::vector<std::string> names;
160  names.push_back(std::string("gen"));
161  names.push_back(std::string("calo_face"));
162 
163  auto clusterHandle = evt.put(std::move(clusters));
164  auto clusterHandleSharing = evt.put(std::move(clusters_sharing),"sharing");
165 
166  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
167  for( unsigned i = 0; i < clusterHandle->size(); ++i ) {
168  edm::Ptr<reco::BasicCluster> ptr(clusterHandle,i);
169  clusterPtrs.push_back(ptr);
170  }
171 
172  if(doSharing){
173  for( unsigned i = 0; i < clusterHandleSharing->size(); ++i ) {
174  edm::Ptr<reco::BasicCluster> ptr(clusterHandleSharing,i);
175  clusterPtrsSharing.push_back(ptr);
176  }
177  }
178 
179  std::unique_ptr<std::vector<reco::HGCalMultiCluster> >
180  multiclusters( new std::vector<reco::HGCalMultiCluster> ),
181  multiclusters_sharing( new std::vector<reco::HGCalMultiCluster> );
182 
183  std::chrono::high_resolution_clock::time_point then = std::chrono::high_resolution_clock::now();
184  *multiclusters = multicluster_algo->makeClusters(clusterPtrs);
185  if(doSharing)
186  *multiclusters_sharing = multicluster_algo->makeClusters(clusterPtrsSharing);
187  evt.put(std::move(multiclusters));
188  if(doSharing)
189  evt.put(std::move(multiclusters_sharing),"sharing");
190  std::chrono::high_resolution_clock::time_point now = std::chrono::high_resolution_clock::now();
191  std::chrono::duration<float> time_span = std::chrono::duration_cast<std::chrono::milliseconds>(now - then);
192  // delta += float (now.tv_usec - then.tv_usec)/1000.;
193  edm::LogInfo ("HGCalClusterProducer") << "Time taken by multiclustering " << time_span.count() << " ms";
194 }
195 
196 #endif
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
static const HistoName names[]
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
std::unique_ptr< HGCalImagingAlgo > algo
HGCalImagingAlgo::VerbosityLevel verbosity
bool isValid() const
Definition: HandleBase.h:74
std::unique_ptr< HGCal3DClustering > multicluster_algo
edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token
edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token
HGCalClusterProducer(const edm::ParameterSet &)
fixed size matrix
virtual void produce(edm::Event &, const edm::EventSetup &)
reco::CaloCluster::AlgoId algoId
static const G4double kappa
def move(src, dest)
Definition: eostools.py:510