CMS 3D CMS Logo

HGCalClusterTestProducer.cc
Go to the documentation of this file.
1 #ifndef __RecoLocalCalo_HGCRecProducers_PFClusterProducer_H__
2 #define __RecoLocalCalo_HGCRecProducers_PFClusterProducer_H__
3 
4 // user include files
7 
10 
17 
18 #include <memory>
19 
22 
26 
28 
30  public:
33 
34  virtual void produce(edm::Event&, const edm::EventSetup&);
35 
36  private:
37 
41 
43 
44  std::unique_ptr<HGCalImagingAlgo> algo;
45  std::unique_ptr<HGCalDepthPreClusterer> multicluster_algo;
46  bool doSharing;
48 
50  // edm::EDGetTokenT<std::vector<reco::PFCluster> > hydraTokens[2];
51 };
52 
54 
56  algoId(reco::CaloCluster::undefined),
57  doSharing(ps.getParameter<bool>("doSharing")),
58  detector(ps.getParameter<std::string >("detector")), //one of EE, EF or "both"
59  verbosity((HGCalImagingAlgo::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity",3)){
60  double ecut = ps.getParameter<double>("ecut");
61  double delta_c = ps.getParameter<double>("deltac");
62  double kappa = ps.getParameter<double>("kappa");
63  double multicluster_radius = ps.getParameter<double>("multiclusterRadius");
64  double minClusters = ps.getParameter<unsigned>("minClusters");
65 
66 
67  if(detector=="all") {
68  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
69  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
70  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
72  }else if(detector=="EE") {
73  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
75  }else if(detector=="FH") {
76  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
78  } else {
79  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
81  }
82 
83 
84  if(doSharing){
85  double showerSigma = ps.getParameter<double>("showerSigma");
86  algo = std::make_unique<HGCalImagingAlgo>(delta_c, kappa, ecut, showerSigma, algoId, verbosity);
87  }else{
88  algo = std::make_unique<HGCalImagingAlgo>(delta_c, kappa, ecut, algoId, verbosity);
89  }
90 
91  auto sumes = consumesCollector();
92 
93  multicluster_algo = std::make_unique<HGCalDepthPreClusterer>(ps, sumes, multicluster_radius, minClusters);
94 
95  // hydraTokens[0] = consumes<std::vector<reco::PFCluster> >( edm::InputTag("FakeClusterGen") );
96  // hydraTokens[1] = consumes<std::vector<reco::PFCluster> >( edm::InputTag("FakeClusterCaloFace") );
97 
98  //std::cout << "Constructing HGCalClusterTestProducer" << std::endl;
99 
100  produces<std::vector<reco::BasicCluster> >();
101  produces<std::vector<reco::BasicCluster> >("sharing");
102 
103  produces<std::vector<reco::HGCalMultiCluster> >();
104  produces<std::vector<reco::HGCalMultiCluster> >("sharing");
105 }
106 
108  const edm::EventSetup& es) {
109 
113 
114 
115  std::unique_ptr<std::vector<reco::BasicCluster> > clusters( new std::vector<reco::BasicCluster> ),
116  clusters_sharing( new std::vector<reco::BasicCluster> );
117 
118  algo->reset();
119 
120  algo->getEventSetup(es);
121 
122  multicluster_algo->getEvent(evt);
123  multicluster_algo->getEventSetup(es);
124 
125  switch(algoId){
127  evt.getByToken(hits_ee_token,ee_hits);
128  algo->makeClusters(*ee_hits);
129  break;
131  evt.getByToken(hits_fh_token,fh_hits);
132  evt.getByToken(hits_bh_token,bh_hits);
133  if( fh_hits.isValid() ) {
134  algo->makeClusters(*fh_hits);
135  } else if ( bh_hits.isValid() ) {
136  algo->makeClusters(*bh_hits);
137  }
138  break;
140  evt.getByToken(hits_ee_token,ee_hits);
141  algo->makeClusters(*ee_hits);
142  evt.getByToken(hits_fh_token,fh_hits);
143  algo->makeClusters(*fh_hits);
144  evt.getByToken(hits_bh_token,bh_hits);
145  algo->makeClusters(*bh_hits);
146  break;
147  default:
148  break;
149  }
150  *clusters = algo->getClusters(false);
151  if(doSharing)
152  *clusters_sharing = algo->getClusters(true);
153 
154  //std::cout << "Density based cluster size: " << clusters->size() << std::endl;
155  //if(doSharing)
156  //std::cout << "Sharing clusters size : " << clusters_sharing->size() << std::endl;
157 
158  // edm::Handle<std::vector<reco::PFCluster> > hydra[2];
159  std::vector<std::string> names;
160  names.push_back(std::string("gen"));
161  names.push_back(std::string("calo_face"));
162  // for( unsigned i = 0 ; i < 2; ++i ) {
163  // evt.getByToken(hydraTokens[i],hydra[i]);
164  // std::cout << "hydra " << names[i] << " size : " << hydra[i]->size() << std::endl;
165  // }
166 
167 
168 
169  auto clusterHandle = evt.put(std::move(clusters));
170  auto clusterHandleSharing = evt.put(std::move(clusters_sharing),"sharing");
171 
172  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
173  for( unsigned i = 0; i < clusterHandle->size(); ++i ) {
174  edm::Ptr<reco::BasicCluster> ptr(clusterHandle,i);
175  clusterPtrs.push_back(ptr);
176  }
177 
178  for( unsigned i = 0; i < clusterHandleSharing->size(); ++i ) {
179  edm::Ptr<reco::BasicCluster> ptr(clusterHandleSharing,i);
180  clusterPtrsSharing.push_back(ptr);
181  }
182 
183  std::unique_ptr<std::vector<reco::HGCalMultiCluster> >
184  multiclusters( new std::vector<reco::HGCalMultiCluster> ),
185  multiclusters_sharing( new std::vector<reco::HGCalMultiCluster> );
186 
187  *multiclusters = std::move(multicluster_algo->makePreClusters(clusterPtrs));
188  *multiclusters_sharing = std::move(multicluster_algo->makePreClusters(clusterPtrsSharing));
189 
190  evt.put(std::move(multiclusters));
191  evt.put(std::move(multiclusters_sharing),"sharing");
192 
193 }
194 
195 #endif
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
HGCalClusterTestProducer(const edm::ParameterSet &)
static const HistoName names[]
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
reco::CaloCluster::AlgoId algoId
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:141
edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token
std::unique_ptr< HGCalImagingAlgo > algo
edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
HGCalImagingAlgo::VerbosityLevel verbosity
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token
fixed size matrix
virtual void produce(edm::Event &, const edm::EventSetup &)
static const G4double kappa
std::unique_ptr< HGCalDepthPreClusterer > multicluster_algo
def move(src, dest)
Definition: eostools.py:510