CMS 3D CMS Logo

MTDClusterProducer.cc
Go to the documentation of this file.
1 //---------------------------------------------------------------------------
6 //---------------------------------------------------------------------------
7 // Our own stuff
10 
11 // Data Formats
13 
14 // STL
15 #include <vector>
16 #include <memory>
17 #include <string>
18 #include <iostream>
19 
20 // MessageLogger
29 
38 
40 public:
41  //--- Constructor, virtual destructor (just in case)
42  explicit MTDClusterProducer(const edm::ParameterSet& conf);
43  ~MTDClusterProducer() override = default;
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
46  //--- The top-level event method.
47  void produce(edm::Event& e, const edm::EventSetup& c) override;
48 
49  //--- Execute the algorithm(s).
50  template <typename T>
51  void run(const T& input, FTLClusterCollection& output);
52 
53 private:
56 
57  const std::string ftlbInstance_; // instance name of barrel clusters
58  const std::string ftleInstance_; // instance name of endcap clusters
59 
60  const std::string clusterMode_; // user's choice of the clusterizer
61  std::unique_ptr<MTDClusterizerBase> clusterizer_; // what we got (for now, one ptr to base class)
62 
67 };
68 
69 //---------------------------------------------------------------------------
71 //---------------------------------------------------------------------------
73  : btlHits_(consumes<FTLRecHitCollection>(conf.getParameter<edm::InputTag>("srcBarrel"))),
74  etlHits_(consumes<FTLRecHitCollection>(conf.getParameter<edm::InputTag>("srcEndcap"))),
75  ftlbInstance_(conf.getParameter<std::string>("BarrelClusterName")),
76  ftleInstance_(conf.getParameter<std::string>("EndcapClusterName")),
77  clusterMode_(conf.getParameter<std::string>("ClusterMode")) {
78  //--- Declare to the EDM what kind of collections we will be making.
79  produces<FTLClusterCollection>(ftlbInstance_);
80  produces<FTLClusterCollection>(ftleInstance_);
81 
82  mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
83  mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
84 
85  //--- Make the algorithm(s) according to what the user specified
86  //--- in the ParameterSet.
87  if (clusterMode_ == "MTDThresholdClusterizer") {
88  clusterizer_ = std::make_unique<MTDThresholdClusterizer>(conf);
89  } else {
90  throw cms::Exception("MTDClusterProducer") << "[MTDClusterProducer]:"
91  << " choice " << clusterMode_ << " is invalid.\n"
92  << "Possible choices:\n"
93  << " MTDThresholdClusterizer";
94  }
95 }
96 
97 // Configuration descriptions
100  desc.add<edm::InputTag>("srcBarrel", edm::InputTag("mtdRecHits:FTLBarrel"));
101  desc.add<edm::InputTag>("srcEndcap", edm::InputTag("mtdRecHits:FTLEndcap"));
102  desc.add<std::string>("BarrelClusterName", "FTLBarrel");
103  desc.add<std::string>("EndcapClusterName", "FTLEndcap");
104  desc.add<std::string>("ClusterMode", "MTDThresholdClusterizer");
106  descriptions.add("mtdClusterProducer", desc);
107 }
108 
109 //---------------------------------------------------------------------------
111 //---------------------------------------------------------------------------
113  // Step A.1: get input data
116  e.getByToken(btlHits_, inputBarrel);
117  e.getByToken(etlHits_, inputEndcap);
118 
119  // Step A.2: get event setup
121  geom_ = geom.product();
122 
123  auto mtdTopo = es.getTransientHandle(mtdtopoToken_);
124  topo_ = mtdTopo.product();
125 
126  // Step B: create the final output collection
127  auto outputBarrel = std::make_unique<FTLClusterCollection>();
128  auto outputEndcap = std::make_unique<FTLClusterCollection>();
129 
130  if (inputBarrel.isValid())
131  run(*inputBarrel, *outputBarrel);
132  else
133  edm::LogWarning("MTDReco") << "MTDClusterProducer:: Missing Barrel Digis";
134  if (inputEndcap.isValid())
135  run(*inputEndcap, *outputEndcap);
136  else
137  edm::LogWarning("MTDReco") << "MTDClusterProducer:: Missing Endcap Digis";
138 
139  e.put(std::move(outputBarrel), ftlbInstance_);
140  e.put(std::move(outputEndcap), ftleInstance_);
141 }
142 
143 //---------------------------------------------------------------------------
145 //---------------------------------------------------------------------------
146 template <typename T>
148  if (!clusterizer_) {
149  throw cms::Exception("MTDClusterProducer") << " at least one clusterizer is not ready -- can't run!";
150  }
151 
152  clusterizer_->clusterize(input, geom_, topo_, output);
153 
154  LogDebug("MTDClusterProducer") << " Executing " << clusterMode_ << " resulted in " << output.size()
155  << " MTDClusters for " << input.size() << " Hits.";
156 }
157 
160 
static void fillPSetDescription(edm::ParameterSetDescription &desc)
edm::ESGetToken< MTDTopology, MTDTopologyRcd > mtdtopoToken_
std::unique_ptr< MTDClusterizerBase > clusterizer_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &e, const edm::EventSetup &c) override
The "Event" entrypoint: gets called by framework for every event.
const edm::EDGetTokenT< FTLRecHitCollection > btlHits_
const edm::EDGetTokenT< FTLRecHitCollection > etlHits_
static std::string const input
Definition: EdmProvDump.cc:50
MTDClusterProducer(const edm::ParameterSet &conf)
Constructor: set the ParameterSet and defer all thinking to setupClusterizer().
const MTDTopology * topo_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
EDProducer to cluster FTLRecHits into FTLClusters.
void run(const T &input, FTLClusterCollection &output)
Iterate over DetUnits, and invoke the PixelClusterizer on each.
const std::string clusterMode_
const std::string ftlbInstance_
~MTDClusterProducer() override=default
edm::ESGetToken< MTDGeometry, MTDDigiGeometryRecord > mtdgeoToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
const std::string ftleInstance_
HLT enums.
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:141
Definition: output.py:1
Log< level::Warning, false > LogWarning
long double T
def move(src, dest)
Definition: eostools.py:511
const MTDGeometry * geom_
#define LogDebug(id)