CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelClusterProducer.cc
Go to the documentation of this file.
1 
14 // Our own stuff
17 
18 // Geometry
21 
22 // Data Formats
26 
27 // Database payloads
31 
32 // Framework
35 
36 // STL
37 #include <vector>
38 #include <memory>
39 #include <string>
40 #include <iostream>
41 
42 // MessageLogger
44 
45 namespace cms
46 {
47 
48  //---------------------------------------------------------------------------
50  //---------------------------------------------------------------------------
52  :
53  conf_(conf),
54  theSiPixelGainCalibration_(0),
55  clusterMode_("None"), // bogus
56  clusterizer_(0), // the default, in case we fail to make one
57  readyToCluster_(false), // since we obviously aren't
58  src_( conf.getParameter<edm::InputTag>( "src" ) ),
59  maxTotalClusters_( conf.getParameter<int32_t>( "maxNumberOfClusters" ) )
60  {
61  //--- Declare to the EDM what kind of collections we will be making.
62  produces<SiPixelClusterCollectionNew>();
63 
64  std::string payloadType = conf.getParameter<std::string>( "payloadType" );
65 
66  if (strcmp(payloadType.c_str(), "HLT") == 0)
68  else if (strcmp(payloadType.c_str(), "Offline") == 0)
70  else if (strcmp(payloadType.c_str(), "Full") == 0)
72 
73  //--- Make the algorithm(s) according to what the user specified
74  //--- in the ParameterSet.
76 
77  }
78 
79  // Destructor
81  delete clusterizer_;
83  }
84 
85  //void SiPixelClusterProducer::beginJob( const edm::EventSetup& es )
87  {
88  edm::LogInfo("SiPixelClusterizer") << "[SiPixelClusterizer::beginJob]";
90  }
91 
92  //---------------------------------------------------------------------------
94  //---------------------------------------------------------------------------
96  {
97 
98  //Setup gain calibration service
100 
101  // Step A.1: get input data
102  //edm::Handle<PixelDigiCollection> pixDigis;
104  e.getByLabel( src_, input);
105 
106  // Step A.2: get event setup
108  es.get<TrackerDigiGeometryRecord>().get( geom );
109 
110  // Step B: create the final output collection
111  std::auto_ptr<SiPixelClusterCollectionNew> output( new SiPixelClusterCollectionNew() );
112  //FIXME: put a reserve() here
113 
114  // Step C: Iterate over DetIds and invoke the pixel clusterizer algorithm
115  // on each DetUnit
116  run(*input, geom, *output );
117 
118  // Step D: write output to file
119  e.put( output );
120 
121  }
122 
123  //---------------------------------------------------------------------------
127  //---------------------------------------------------------------------------
129  clusterMode_ =
130  conf_.getUntrackedParameter<std::string>("ClusterMode","PixelThresholdClusterizer");
131 
132  if ( clusterMode_ == "PixelThresholdClusterizer" ) {
134  readyToCluster_ = true;
135  }
136  else {
137  edm::LogError("SiPixelClusterProducer") << "[SiPixelClusterProducer]:"
138  <<" choice " << clusterMode_ << " is invalid.\n"
139  << "Possible choices:\n"
140  << " PixelThresholdClusterizer";
141  readyToCluster_ = false;
142  }
143  }
144 
145  //---------------------------------------------------------------------------
147  //---------------------------------------------------------------------------
151  if ( ! readyToCluster_ ) {
152  edm::LogError("SiPixelClusterProducer")
153  <<" at least one clusterizer is not ready -- can't run!" ;
154  // TO DO: throw an exception here? The user may want to know...
155  return; // clusterizer is invalid, bail out
156  }
157 
158  int numberOfDetUnits = 0;
159  int numberOfClusters = 0;
160 
161  // Iterate on detector units
163  for( ; DSViter != input.end(); DSViter++) {
164  ++numberOfDetUnits;
165 
166  // LogDebug takes very long time, get rid off.
167  //LogDebug("SiStripClusterizer") << "[SiPixelClusterProducer::run] DetID" << DSViter->id;
168 
169  std::vector<short> badChannels;
170  DetId detIdObject(DSViter->detId());
171 
172  // Comment: At the moment the clusterizer depends on geometry
173  // to access information as the pixel topology (number of columns
174  // and rows in a detector module).
175  // In the future the geometry service will be replaced with
176  // a ES service.
177  const GeomDetUnit * geoUnit = geom->idToDetUnit( detIdObject );
178  const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
179  if (! pixDet) {
180  // Fatal error! TO DO: throw an exception!
181  assert(0);
182  }
183  // Produce clusters for this DetUnit and store them in
184  // a DetSet
185  edmNew::DetSetVector<SiPixelCluster>::FastFiller spc(output, DSViter->detId());
186  clusterizer_->clusterizeDetUnit(*DSViter, pixDet, badChannels, spc);
187  if ( spc.empty() ) {
188  spc.abort();
189  } else {
190  numberOfClusters += spc.size();
191  }
192 
193  if ((maxTotalClusters_ >= 0) && (numberOfClusters > maxTotalClusters_)) {
194  edm::LogError("TooManyClusters") << "Limit on the number of clusters exceeded. An empty cluster collection will be produced instead.\n";
196  empty.swap(output);
197  break;
198  }
199  } // end of DetUnit loop
200 
201  //LogDebug ("SiPixelClusterProducer") << " Executing "
202  // << clusterMode_ << " resulted in " << numberOfClusters
203  // << " SiPixelClusters in " << numberOfDetUnits << " DetUnits.";
204  }
205 
206 } // end of namespace cms
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
PixelClusterizerBase * clusterizer_
virtual void setESObjects(const edm::EventSetup &es)=0
void swap(DetSetVector &rh)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:356
tuple conf
Definition: dbtoconf.py:185
Definition: DetId.h:20
const T & get() const
Definition: EventSetup.h:55
SiPixelClusterProducer(const edm::ParameterSet &conf)
Constructor: set the ParameterSet and defer all thinking to setupClusterizer().
edmNew::DetSetVector< SiPixelCluster > SiPixelClusterCollectionNew
void run(const edm::DetSetVector< PixelDigi > &input, edm::ESHandle< TrackerGeometry > &geom, edmNew::DetSetVector< SiPixelCluster > &output)
Iterate over DetUnits, and invoke the PixelClusterizer on each.
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:341
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
virtual void produce(edm::Event &e, const edm::EventSetup &c)
The &quot;Event&quot; entrypoint: gets called by framework for every event.
int32_t maxTotalClusters_
Optional limit on the total number of clusters.
SiPixelGainCalibrationServiceBase * theSiPixelGainCalibration_
void setSiPixelGainCalibrationService(SiPixelGainCalibrationServiceBase *in)
An explicit threshold-based clustering algorithm.
virtual void clusterizeDetUnit(const edm::DetSet< PixelDigi > &input, const PixelGeomDetUnit *pixDet, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)=0