CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelClusterProducer.cc
Go to the documentation of this file.
1 
14 // Our own stuff
15 #include "SiPixelClusterProducer.h"
18 
19 // Geometry
21 
22 // Data Formats
26 
27 // Database payloads
31 
32 // Framework
36 
37 // STL
38 #include <vector>
39 #include <memory>
40 #include <string>
41 #include <iostream>
42 
43 // MessageLogger
45 
46 //---------------------------------------------------------------------------
48 //---------------------------------------------------------------------------
50  : tPutPixelClusters(produces<SiPixelClusterCollectionNew>()),
51  clusterMode_(conf.getParameter<std::string>("ClusterMode")),
52  maxTotalClusters_(conf.getParameter<int32_t>("maxNumberOfClusters")) {
53  if (clusterMode_ == "PixelThresholdReclusterizer")
54  tPixelClusters = consumes<SiPixelClusterCollectionNew>(conf.getParameter<edm::InputTag>("src"));
55  else
56  tPixelDigi = consumes<edm::DetSetVector<PixelDigi>>(conf.getParameter<edm::InputTag>("src"));
57 
58  trackerTopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
59  trackerGeomToken_ = esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>();
60 
61  const auto& payloadType = conf.getParameter<std::string>("payloadType");
62  if (payloadType == "HLT")
63  theSiPixelGainCalibration_ = std::make_unique<SiPixelGainCalibrationForHLTService>(conf, consumesCollector());
64  else if (payloadType == "Offline")
65  theSiPixelGainCalibration_ = std::make_unique<SiPixelGainCalibrationOfflineService>(conf, consumesCollector());
66  else if (payloadType == "Full")
67  theSiPixelGainCalibration_ = std::make_unique<SiPixelGainCalibrationService>(conf, consumesCollector());
68 
69  //--- Make the algorithm(s) according to what the user specified
70  //--- in the ParameterSet.
71  setupClusterizer(conf);
72 }
73 
74 // Destructor
76 
79 
80  desc.add<edm::InputTag>("src", edm::InputTag("siPixelDigis"));
81  desc.add<std::string>("ClusterMode", "PixelThresholdClusterizer");
82  desc.add<int>("maxNumberOfClusters", -1)->setComment("-1 means no limit");
83  desc.add<std::string>("payloadType", "Offline")
84  ->setComment("Options: HLT - column granularity, Offline - gain:col/ped:pix");
85 
87  SiPixelGainCalibrationServiceBase::fillPSetDescription(desc); // no-op, but in principle the structures are there...
88 
89  descriptions.add("SiPixelClusterizerDefault", desc);
90 }
91 
92 //---------------------------------------------------------------------------
94 //---------------------------------------------------------------------------
96  //Setup gain calibration service
97  theSiPixelGainCalibration_->setESObjects(es);
98 
99  // Step A.1: get input data
102  if (clusterMode_ == "PixelThresholdReclusterizer")
103  e.getByToken(tPixelClusters, inputClusters);
104  else
105  e.getByToken(tPixelDigi, inputDigi);
106 
107  // Step A.2: get event setup
109 
110  edm::ESHandle<TrackerTopology> trackerTopologyHandle = es.getHandle(trackerTopoToken_);
111  tTopo_ = trackerTopologyHandle.product();
112 
113  // Step B: create the final output collection
114  auto output = std::make_unique<SiPixelClusterCollectionNew>();
115  //FIXME: put a reserve() here
116 
117  // Step C: Iterate over DetIds and invoke the pixel clusterizer algorithm
118  // on each DetUnit
119  if (clusterMode_ == "PixelThresholdReclusterizer")
120  run(*inputClusters, geom, *output);
121  else
122  run(*inputDigi, geom, *output);
123 
124  // Step D: write output to file
125  output->shrink_to_fit();
126 
127  // set sequential identifier
128  for (auto& clusters : *output) {
129  uint16_t id = 0;
130  for (auto& cluster : clusters) {
131  cluster.setOriginalId(id++);
132  }
133  }
134  e.put(tPutPixelClusters, std::move(output));
135 }
136 
137 //---------------------------------------------------------------------------
141 //---------------------------------------------------------------------------
143  if (clusterMode_ == "PixelThresholdReclusterizer" || clusterMode_ == "PixelThresholdClusterizer") {
144  clusterizer_ = std::make_unique<PixelThresholdClusterizer>(conf);
145  clusterizer_->setSiPixelGainCalibrationService(theSiPixelGainCalibration_.get());
146  } else if (clusterMode_ == "PixelThresholdClusterizerForBricked") {
147  clusterizer_ = std::make_unique<PixelThresholdClusterizerForBricked>(conf);
148  clusterizer_->setSiPixelGainCalibrationService(theSiPixelGainCalibration_.get());
149  } else {
150  throw cms::Exception("Configuration") << "[SiPixelClusterProducer]:"
151  << " choice " << clusterMode_ << " is invalid.\n"
152  << "Possible choices:\n"
153  << " PixelThresholdClusterizer";
154  }
155 }
156 
157 //---------------------------------------------------------------------------
159 //---------------------------------------------------------------------------
160 template <typename T>
164  int numberOfDetUnits = 0;
165  int numberOfClusters = 0;
166 
167  // Iterate on detector units
168  for (auto const& dsv : input) {
169  ++numberOfDetUnits;
170 
171  // LogDebug takes very long time, get rid off.
172  //LogDebug("SiStripClusterizer") << "[SiPixelClusterProducer::run] DetID" << dsv.id;
173 
174  std::vector<short> badChannels;
175  DetId detIdObject(dsv.detId());
176 
177  // Comment: At the moment the clusterizer depends on geometry
178  // to access information as the pixel topology (number of columns
179  // and rows in a detector module).
180  // In the future the geometry service will be replaced with
181  // a ES service.
182  const GeomDetUnit* geoUnit = geom->idToDetUnit(detIdObject);
183  const PixelGeomDetUnit* pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
184  if (!pixDet) {
185  // Fatal error! TO DO: throw an exception!
186  assert(0);
187  }
188  {
189  // Produce clusters for this DetUnit and store them in
190  // a DetSet
191  edmNew::DetSetVector<SiPixelCluster>::FastFiller spc(output, dsv.detId());
192  clusterizer_->clusterizeDetUnit(dsv, pixDet, tTopo_, badChannels, spc);
193  if (spc.empty()) {
194  spc.abort();
195  } else {
196  numberOfClusters += spc.size();
197  }
198  } // spc is not deleted and detsetvector updated
199  if ((maxTotalClusters_ >= 0) && (numberOfClusters > maxTotalClusters_)) {
200  edm::LogError("TooManyClusters")
201  << "Limit on the number of clusters exceeded. An empty cluster collection will be produced instead.\n";
203  empty.swap(output);
204  break;
205  }
206  } // end of DetUnit loop
207 
208  //LogDebug ("SiPixelClusterProducer") << " Executing "
209  // << clusterMode_ << " resulted in " << numberOfClusters
210  // << " SiPixelClusters in " << numberOfDetUnits << " DetUnits.";
211 }
212 
215 
static void fillPSetDescription(edm::ParameterSetDescription &desc)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > trackerTopoToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void swap(DetSetVector &rh)
Log< level::Error, false > LogError
assert(be >=bs)
static std::string const input
Definition: EdmProvDump.cc:47
std::unique_ptr< PixelClusterizerBase > clusterizer_
def move
Definition: eostools.py:511
SiPixelClusterProducer(const edm::ParameterSet &conf)
Constructor: set the ParameterSet and defer all thinking to setupClusterizer().
std::unique_ptr< SiPixelGainCalibrationServiceBase > theSiPixelGainCalibration_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDPutTokenT< SiPixelClusterCollectionNew > tPutPixelClusters
const int32_t maxTotalClusters_
Optional limit on the total number of clusters.
EDProducer to cluster PixelDigis into SiPixelClusters.
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Definition: DetId.h:17
void produce(edm::Event &e, const edm::EventSetup &c) override
The &quot;Event&quot; entrypoint: gets called by framework for every event.
void setupClusterizer(const edm::ParameterSet &conf)
T const * product() const
Definition: ESHandle.h:86
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > tPixelDigi
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< SiPixelClusterCollectionNew > tPixelClusters
~SiPixelClusterProducer() override
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
void run(const T &input, const edm::ESHandle< TrackerGeometry > &geom, edmNew::DetSetVector< SiPixelCluster > &output)
Iterate over DetUnits, and invoke the PixelClusterizer on each.
long double T
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeomToken_
const TrackerTopology * tTopo_
const std::string clusterMode_