CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelRecHitConverter.cc
Go to the documentation of this file.
1 
12 // Our own stuff
16 
17 // Geometry
20 
21 // Data Formats
25 
26 // Framework Already defined in the *.h file
27 //#include "DataFormats/Common/interface/Handle.h"
28 //#include "FWCore/Framework/interface/ESHandle.h"
29 
30 // STL
31 #include <vector>
32 #include <memory>
33 #include <string>
34 #include <iostream>
35 
36 // MessageLogger
38 
40 
41 using namespace std;
42 
43 namespace cms
44 {
45  //---------------------------------------------------------------------------
47  //---------------------------------------------------------------------------
48  SiPixelRecHitConverter::SiPixelRecHitConverter(edm::ParameterSet const& conf)
49  :
50  conf_(conf),
51  cpeName_("None"), // bogus
52  cpe_(0), // the default, in case we fail to make one
53  ready_(false), // since we obviously aren't
54  src_( conf.getParameter<edm::InputTag>( "src" ) ),
55  theVerboseLevel(conf.getUntrackedParameter<int>("VerboseLevel",0))
56  {
57  //--- Declare to the EDM what kind of collections we will be making.
58  produces<SiPixelRecHitCollection>();
59 
60  }
61 
62  // Destructor
64  {
65  }
66 
67  //---------------------------------------------------------------------------
68  // Begin job: get magnetic field
69  //---------------------------------------------------------------------------
70  //void SiPixelRecHitConverter::beginJob()
72  {
73  }
74 
75  //---------------------------------------------------------------------------
77  //---------------------------------------------------------------------------
79  {
80 
81  // Step A.1: get input data
83  e.getByLabel( src_, input);
84 
85  // Step A.2: get event setup
87  es.get<TrackerDigiGeometryRecord>().get( geom );
88 
89  // Step B: create empty output collection
90  std::auto_ptr<SiPixelRecHitCollectionNew> output(new SiPixelRecHitCollectionNew);
91 
92  // Step B*: create CPE
94  std::string cpeName_ = conf_.getParameter<std::string>("CPE");
95  es.get<TkPixelCPERecord>().get(cpeName_,hCPE);
96  const PixelClusterParameterEstimator &cpe(*hCPE);
97  cpe_ = &cpe;
98 
99  if(cpe_) ready_ = true;
100 
101 
102  // Step C: Iterate over DetIds and invoke the strip CPE algorithm
103  // on each DetUnit
104 
105  run( input, *output, geom );
106 
107 
108  e.put(output);
109 
110  }
111 
112  //---------------------------------------------------------------------------
116  //---------------------------------------------------------------------------
120  {
121  if ( ! ready_ )
122  {
123  edm::LogError("SiPixelRecHitConverter") << " at least one CPE is not ready -- can't run!";
124  // TO DO: throw an exception here? The user may want to know...
125  assert(0);
126  return; // clusterizer is invalid, bail out
127  }
128 
129  int numberOfDetUnits = 0;
130  int numberOfClusters = 0;
131 
132  const edmNew::DetSetVector<SiPixelCluster>& input = *inputhandle;
133 
135 
136  for ( ; DSViter != input.end() ; DSViter++)
137  {
138  numberOfDetUnits++;
139  unsigned int detid = DSViter->detId();
140  DetId detIdObject( detid );
141  const GeomDetUnit * genericDet = geom->idToDetUnit( detIdObject );
142  const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
143  assert(pixDet);
144  SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(output,detid);
145 
146  edmNew::DetSet<SiPixelCluster>::const_iterator clustIt = DSViter->begin(), clustEnd = DSViter->end();
147 
148  for ( ; clustIt != clustEnd; clustIt++)
149  {
150  numberOfClusters++;
151  std::pair<LocalPoint, LocalError> lv = cpe_->localParameters( *clustIt, *genericDet );
152  LocalPoint lp( lv.first );
153  LocalError le( lv.second );
154  // Create a persistent edm::Ref to the cluster
156  // Make a RecHit and add it to the DetSet
157  // old : recHitsOnDetUnit.push_back( new SiPixelRecHit( lp, le, detIdObject, &*clustIt) );
158  SiPixelRecHit hit( lp, le, detIdObject, cluster);
159  // Copy the extra stuff; unfortunately, only the derivatives of PixelCPEBase
160  // are capable of doing that. So until we get rid of CPEFromDetPosition
161  // we'll have to dynamic_cast :(
162  // &&& This cast can be moved to the setupCPE, so that it's done once per job.
163  const PixelCPEBase * cpeBase
164  = dynamic_cast< const PixelCPEBase* >( cpe_ );
165  if (cpeBase)
166  {
167  hit.setRawQualityWord( cpeBase->rawQualityWord() );
168  // hit.setProbabilityX( cpeBase->probabilityX() );
169  // hit.setProbabilityY( cpeBase->probabilityY() );
170  // hit.setQBin( cpeBase->qBin() );
171  // hit.setCotAlphaFromCluster( cpeBase->cotAlphaFromCluster() );
172  // hit.setCotBetaFromCluster ( cpeBase->cotBetaFromCluster() );
173  }
174  //
175  // Now save it =================
176  recHitsOnDetUnit.push_back(hit);
177  // =============================
178  } // <-- End loop on Clusters
179 
180  if ( recHitsOnDetUnit.size()>0 )
181  {
182  if (theVerboseLevel > 2)
183  LogDebug("SiPixelRecHitConverter") << " Found "
184  << recHitsOnDetUnit.size() << " RecHits on " << detid;
185  }
186 
187  } // <-- End loop on DetUnits
188 
189  if ( theVerboseLevel > 2 ) LogDebug ("SiPixelRecHitConverter")
190  << cpeName_ << " converted " << numberOfClusters
191  << " SiPixelClusters into SiPixelRecHits, in "
192  << numberOfDetUnits << " DetUnits.";
193  }
194 } // end of namespace cms
#define LogDebug(id)
T getParameter(std::string const &) const
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
void push_back(data_type const &d)
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator begin() const
void run(const edmNew::DetSetVector< SiPixelCluster > &input, SiPixelRecHitCollectionNew &output, edm::ESHandle< TrackerGeometry > &geom)
SiPixelRecHitQuality::QualWordType rawQualityWord() const
data_type const * const_iterator
Definition: DetSetNew.h:25
virtual LocalValues localParameters(const T &, const GeomDetUnit &) const =0
static int theVerboseLevel
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
const_iterator end() const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple conf
Definition: dbtoconf.py:185
Definition: DetId.h:20
virtual void produce(edm::Event &e, const edm::EventSetup &c)
The &quot;Event&quot; entrypoint: gets called by framework for every event.
const T & get() const
Definition: EventSetup.h:55
Pixel cluster – collection of neighboring pixels above threshold.
const PixelClusterParameterEstimator * cpe_
Pixel Reconstructed Hit.