CMS 3D CMS Logo

MultiClustersFromTrackstersProducer.cc
Go to the documentation of this file.
1 // user include files
4 
9 
12 
13 #include <algorithm>
14 #include <array>
15 #include <string>
16 #include <vector>
17 
19 public:
22  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
23 
24  void produce(edm::Event&, const edm::EventSetup&) override;
25 
26 private:
29 };
30 
32 
34  : layer_clusters_token_(consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("LayerClusters"))),
35  tracksters_token_(consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("Tracksters"))) {
36  produces<std::vector<reco::HGCalMultiCluster>>();
37 }
38 
40  // hgcalMultiClusters
42  desc.add<edm::InputTag>("Tracksters", edm::InputTag("Tracksters", "TrackstersByCA"));
43  desc.add<edm::InputTag>("LayerClusters", edm::InputTag("hgcalMergeLayerClusters"));
44  desc.addUntracked<unsigned int>("verbosity", 3);
45  descriptions.add("multiClustersFromTrackstersProducer", desc);
46 }
47 
49  auto multiclusters = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
51  evt.getByToken(tracksters_token_, tracksterHandle);
52 
53  edm::Handle<std::vector<reco::CaloCluster>> layer_clustersHandle;
54  evt.getByToken(layer_clusters_token_, layer_clustersHandle);
55 
56  auto const& tracksters = *tracksterHandle;
57  auto const& layerClusters = *layer_clustersHandle;
58 
60  for (unsigned i = 0; i < layerClusters.size(); ++i) {
61  edm::Ptr<reco::BasicCluster> ptr(layer_clustersHandle, i);
62  clusterPtrs.push_back(ptr);
63  }
64 
65  std::for_each(std::begin(tracksters), std::end(tracksters), [&](auto const& trackster) {
66  // Create an empty multicluster if the trackster has no layer clusters.
67  // This could happen when a seed leads to no trackster and a dummy one is produced.
68 
69  std::array<double, 3> baricenter{{0., 0., 0.}};
70  double total_weight = 0.;
72  int counter = 0;
73  if (!trackster.vertices().empty()) {
74  std::for_each(std::begin(trackster.vertices()), std::end(trackster.vertices()), [&](unsigned int idx) {
75  temp.push_back(clusterPtrs[idx]);
76  auto fraction = 1.f / trackster.vertex_multiplicity(counter++);
77  for (auto const& cell : clusterPtrs[idx]->hitsAndFractions()) {
78  temp.addHitAndFraction(cell.first, cell.second * fraction);
79  }
80  auto weight = clusterPtrs[idx]->energy() * fraction;
81  total_weight += weight;
82  baricenter[0] += clusterPtrs[idx]->x() * weight;
83  baricenter[1] += clusterPtrs[idx]->y() * weight;
84  baricenter[2] += clusterPtrs[idx]->z() * weight;
85  });
87  std::begin(baricenter), std::end(baricenter), std::begin(baricenter), [&total_weight](double val) -> double {
88  return val / total_weight;
89  });
90  }
91  temp.setEnergy(total_weight);
92  temp.setCorrectedEnergy(total_weight);
93  temp.setPosition(math::XYZPoint(baricenter[0], baricenter[1], baricenter[2]));
95  temp.setTime(trackster.time(), trackster.timeError());
96  multiclusters->push_back(temp);
97  });
98 
99  evt.put(std::move(multiclusters));
100 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
Definition: weight.py:1
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:536
edm::EDGetTokenT< std::vector< ticl::Trackster > > tracksters_token_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
Definition: Common.h:8
edm::EDGetTokenT< std::vector< reco::CaloCluster > > layer_clusters_token_
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned transform(const HcalDetId &id, unsigned transformCode)