test
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
15 #include "SiPixelClusterProducer.h"
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 
46  //---------------------------------------------------------------------------
48  //---------------------------------------------------------------------------
50  :
51  conf_(conf),
52  theSiPixelGainCalibration_(0),
53  clusterMode_("None"), // bogus
54  clusterizer_(0), // the default, in case we fail to make one
55  readyToCluster_(false), // since we obviously aren't
56  src_( conf.getParameter<edm::InputTag>( "src" ) ),
57  maxTotalClusters_( conf.getParameter<int32_t>( "maxNumberOfClusters" ) )
58  {
59  tPixelDigi = consumes<edm::DetSetVector<PixelDigi>>(src_);
60  //--- Declare to the EDM what kind of collections we will be making.
61  produces<SiPixelClusterCollectionNew>();
62 
63  std::string payloadType = conf.getParameter<std::string>( "payloadType" );
64 
65  if (strcmp(payloadType.c_str(), "HLT") == 0)
67  else if (strcmp(payloadType.c_str(), "Offline") == 0)
69  else if (strcmp(payloadType.c_str(), "Full") == 0)
71 
72  //--- Make the algorithm(s) according to what the user specified
73  //--- in the ParameterSet.
75 
76  }
77 
78  // Destructor
80  delete clusterizer_;
82  }
83 
84 
85  //---------------------------------------------------------------------------
87  //---------------------------------------------------------------------------
89  {
90 
91  //Setup gain calibration service
93 
94  // Step A.1: get input data
95  //edm::Handle<PixelDigiCollection> pixDigis;
97  e.getByToken(tPixelDigi, input);
98 
99  // Step A.2: get event setup
101  es.get<TrackerDigiGeometryRecord>().get( geom );
102 
103  // Step B: create the final output collection
104  std::auto_ptr<SiPixelClusterCollectionNew> output( new SiPixelClusterCollectionNew() );
105  //FIXME: put a reserve() here
106 
107  // Step C: Iterate over DetIds and invoke the pixel clusterizer algorithm
108  // on each DetUnit
109  run(*input, geom, *output );
110 
111  // Step D: write output to file
112  output->shrink_to_fit();
113  e.put( output );
114 
115  }
116 
117  //---------------------------------------------------------------------------
121  //---------------------------------------------------------------------------
123  clusterMode_ =
124  conf_.getUntrackedParameter<std::string>("ClusterMode","PixelThresholdClusterizer");
125 
126  if ( clusterMode_ == "PixelThresholdClusterizer" ) {
129  readyToCluster_ = true;
130  }
131  else {
132  edm::LogError("SiPixelClusterProducer") << "[SiPixelClusterProducer]:"
133  <<" choice " << clusterMode_ << " is invalid.\n"
134  << "Possible choices:\n"
135  << " PixelThresholdClusterizer";
136  readyToCluster_ = false;
137  }
138  }
139 
140  //---------------------------------------------------------------------------
142  //---------------------------------------------------------------------------
146  if ( ! readyToCluster_ ) {
147  edm::LogError("SiPixelClusterProducer")
148  <<" at least one clusterizer is not ready -- can't run!" ;
149  // TO DO: throw an exception here? The user may want to know...
150  return; // clusterizer is invalid, bail out
151  }
152 
153  int numberOfDetUnits = 0;
154  int numberOfClusters = 0;
155 
156  // Iterate on detector units
158  for( ; DSViter != input.end(); DSViter++) {
159  ++numberOfDetUnits;
160 
161  // LogDebug takes very long time, get rid off.
162  //LogDebug("SiStripClusterizer") << "[SiPixelClusterProducer::run] DetID" << DSViter->id;
163 
164  std::vector<short> badChannels;
165  DetId detIdObject(DSViter->detId());
166 
167  // Comment: At the moment the clusterizer depends on geometry
168  // to access information as the pixel topology (number of columns
169  // and rows in a detector module).
170  // In the future the geometry service will be replaced with
171  // a ES service.
172  const GeomDetUnit * geoUnit = geom->idToDetUnit( detIdObject );
173  const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
174  if (! pixDet) {
175  // Fatal error! TO DO: throw an exception!
176  assert(0);
177  }
178  {
179  // Produce clusters for this DetUnit and store them in
180  // a DetSet
181  edmNew::DetSetVector<SiPixelCluster>::FastFiller spc(output, DSViter->detId());
182  clusterizer_->clusterizeDetUnit(*DSViter, pixDet, badChannels, spc);
183  if ( spc.empty() ) {
184  spc.abort();
185  } else {
186  numberOfClusters += spc.size();
187  }
188  } // spc is not deleted and detsetvector updated
189  if ((maxTotalClusters_ >= 0) && (numberOfClusters > maxTotalClusters_)) {
190  edm::LogError("TooManyClusters") << "Limit on the number of clusters exceeded. An empty cluster collection will be produced instead.\n";
192  empty.swap(output);
193  break;
194  }
195  } // end of DetUnit loop
196 
197  //LogDebug ("SiPixelClusterProducer") << " Executing "
198  // << clusterMode_ << " resulted in " << numberOfClusters
199  // << " SiPixelClusters in " << numberOfDetUnits << " DetUnits.";
200  }
201 
202 
203 
204 
207 
209 
PixelClusterizerBase * clusterizer_
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual void setESObjects(const edm::EventSetup &es)=0
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
assert(m_qm.get())
void swap(DetSetVector &rh)
static std::string const input
Definition: EdmProvDump.cc:44
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
SiPixelClusterProducer(const edm::ParameterSet &conf)
Constructor: set the ParameterSet and defer all thinking to setupClusterizer().
int32_t maxTotalClusters_
Optional limit on the total number of clusters.
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:361
EDProducer to cluster PixelDigis into SiPixelClusters.
Definition: DetId.h:18
virtual void produce(edm::Event &e, const edm::EventSetup &c) override
The &quot;Event&quot; entrypoint: gets called by framework for every event.
const T & get() const
Definition: EventSetup.h:56
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > tPixelDigi
edmNew::DetSetVector< SiPixelCluster > SiPixelClusterCollectionNew
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:346
SiPixelGainCalibrationServiceBase * theSiPixelGainCalibration_
volatile std::atomic< bool > shutdown_flag false
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
void run(const edm::DetSetVector< PixelDigi > &input, edm::ESHandle< TrackerGeometry > &geom, edmNew::DetSetVector< SiPixelCluster > &output)
Iterate over DetUnits, and invoke the PixelClusterizer on each.
void setSiPixelGainCalibrationService(SiPixelGainCalibrationServiceBase *in)
A specific threshold-based pixel clustering algorithm.
virtual void clusterizeDetUnit(const edm::DetSet< PixelDigi > &input, const PixelGeomDetUnit *pixDet, const std::vector< short > &badChannels, edmNew::DetSetVector< SiPixelCluster >::FastFiller &output)=0