CMS 3D CMS Logo

TICLLayerTileProducer.cc
Go to the documentation of this file.
1 // Author: Marco Rovere, marco.rovere@cern.ch
2 // Date: 05/2019
3 //
4 #include <memory> // unique_ptr
5 
11 
14 
16 
18 public:
19  explicit TICLLayerTileProducer(const edm::ParameterSet &ps);
20  ~TICLLayerTileProducer() override{};
21  void beginRun(edm::Run const &, edm::EventSetup const &) override;
22  void produce(edm::Event &, const edm::EventSetup &) override;
23  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
24 
25 private:
31  bool doNose_;
32 };
33 
35  : detector_(ps.getParameter<std::string>("detector")) {
36  geometry_token_ = esConsumes<CaloGeometry, CaloGeometryRecord, edm::Transition::BeginRun>();
37 
38  doNose_ = (detector_ == "HFNose");
39 
40  if (doNose_) {
42  consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layer_HFNose_clusters"));
43  produces<TICLLayerTilesHFNose>();
44  } else {
45  clusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layer_clusters"));
46  produces<TICLLayerTiles>();
47  }
48 }
49 
53 }
54 
56  std::unique_ptr<TICLLayerTilesHFNose> resultHFNose;
57  std::unique_ptr<TICLLayerTiles> result;
58  if (doNose_) {
59  resultHFNose = std::make_unique<TICLLayerTilesHFNose>();
60  } else {
61  result = std::make_unique<TICLLayerTiles>();
62  }
63 
65  if (doNose_)
66  evt.getByToken(clusters_HFNose_token_, cluster_h);
67  else
68  evt.getByToken(clusters_token_, cluster_h);
69 
70  const auto &layerClusters = *cluster_h;
71  int lcId = 0;
72  for (auto const &lc : layerClusters) {
73  const auto firstHitDetId = lc.hitsAndFractions()[0].first;
74  int layer = rhtools_.getLayerWithOffset(firstHitDetId) +
75  rhtools_.lastLayer(doNose_) * ((rhtools_.zside(firstHitDetId) + 1) >> 1) - 1;
76 
77  assert(layer >= 0);
78 
79  if (doNose_) {
80  resultHFNose->fill(layer, lc.eta(), lc.phi(), lcId);
81  } else {
82  result->fill(layer, lc.eta(), lc.phi(), lcId);
83  LogDebug("TICLLayerTileProducer") << "Adding layerClusterId: " << lcId << " into bin [eta,phi]: [ "
84  << (*result)[layer].etaBin(lc.eta()) << ", "
85  << (*result)[layer].phiBin(lc.phi()) << "] for layer: " << layer << std::endl;
86  }
87  lcId++;
88  }
89  if (doNose_)
90  evt.put(std::move(resultHFNose));
91  else
92  evt.put(std::move(result));
93 }
94 
97  desc.add<std::string>("detector", "HGCAL");
98  desc.add<edm::InputTag>("layer_clusters", edm::InputTag("hgcalMergeLayerClusters"));
99  desc.add<edm::InputTag>("layer_HFNose_clusters", edm::InputTag("hgcalLayerClustersHFNose"));
100  descriptions.add("ticlLayerTileProducer", desc);
101 }
102 
edm::EDGetTokenT< std::vector< reco::CaloCluster > > clusters_HFNose_token_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
hgcal::RecHitTools rhtools_
void produce(edm::Event &, const edm::EventSetup &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void beginRun(edm::Run const &, edm::EventSetup const &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:526
assert(be >=bs)
int zside(const DetId &id) const
Definition: RecHitTools.cc:174
TICLLayerTileProducer(const edm::ParameterSet &ps)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< std::vector< reco::CaloCluster > > clusters_token_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometry_token_
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
#define LogDebug(id)
unsigned int lastLayer(bool nose=false) const
Definition: RecHitTools.h:80
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381