CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
cms::SiPixelRecHitConverter Class Reference

#include <SiPixelRecHitConverter.h>

Inheritance diagram for cms::SiPixelRecHitConverter:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

virtual void beginJob ()
 
virtual void produce (edm::Event &e, const edm::EventSetup &c)
 The "Event" entrypoint: gets called by framework for every event. More...
 
void run (const edmNew::DetSetVector< SiPixelCluster > &input, SiPixelRecHitCollectionNew &output, edm::ESHandle< TrackerGeometry > &geom)
 
void run (edm::Handle< edmNew::DetSetVector< SiPixelCluster > > inputhandle, SiPixelRecHitCollectionNew &output, edm::ESHandle< TrackerGeometry > &geom)
 
 SiPixelRecHitConverter (const edm::ParameterSet &conf)
 Constructor: set the ParameterSet and defer all thinking to setupCPE(). More...
 
virtual ~SiPixelRecHitConverter ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Attributes

edm::ParameterSet conf_
 
const
PixelClusterParameterEstimator
cpe_
 
std::string cpeName_
 
bool m_newCont
 
bool ready_
 
edm::InputTag src_
 
int theVerboseLevel
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 62 of file SiPixelRecHitConverter.h.

Constructor & Destructor Documentation

SiPixelRecHitConverter::SiPixelRecHitConverter ( const edm::ParameterSet conf)
explicit

Constructor: set the ParameterSet and defer all thinking to setupCPE().

Definition at line 48 of file SiPixelRecHitConverter.cc.

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  }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const PixelClusterParameterEstimator * cpe_
SiPixelRecHitConverter::~SiPixelRecHitConverter ( )
virtual

Definition at line 63 of file SiPixelRecHitConverter.cc.

64  {
65  }

Member Function Documentation

void SiPixelRecHitConverter::beginJob ( void  )
virtual

Reimplemented from edm::EDProducer.

Definition at line 71 of file SiPixelRecHitConverter.cc.

72  {
73  }
void SiPixelRecHitConverter::produce ( edm::Event e,
const edm::EventSetup c 
)
virtual

The "Event" entrypoint: gets called by framework for every event.

Implements edm::EDProducer.

Definition at line 78 of file SiPixelRecHitConverter.cc.

References conf_, cpe_, cpeName_, relativeConstraints::geom, edm::EventSetup::get(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), LaserDQM_cfg::input, convertSQLitetoXML_cfg::output, edm::Event::put(), ready_, run(), and src_.

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  }
T getParameter(std::string const &) const
void run(const edmNew::DetSetVector< SiPixelCluster > &input, SiPixelRecHitCollectionNew &output, edm::ESHandle< TrackerGeometry > &geom)
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
const PixelClusterParameterEstimator * cpe_
void cms::SiPixelRecHitConverter::run ( const edmNew::DetSetVector< SiPixelCluster > &  input,
SiPixelRecHitCollectionNew output,
edm::ESHandle< TrackerGeometry > &  geom 
)

Referenced by produce().

void SiPixelRecHitConverter::run ( edm::Handle< edmNew::DetSetVector< SiPixelCluster > >  inputhandle,
SiPixelRecHitCollectionNew output,
edm::ESHandle< TrackerGeometry > &  geom 
)

Iterate over DetUnits, then over Clusters and invoke the CPE on each, and make a RecHit to store the result. New interface reading DetSetVector by V.Chiochia (May 30th, 2006)

Definition at line 117 of file SiPixelRecHitConverter.cc.

References edmNew::DetSetVector< T >::begin(), cpe_, cpeName_, cond::rpcobgas::detid, edmNew::DetSetVector< T >::end(), LaserDQM_cfg::input, asciidump::le, ClusterParameterEstimator< T >::localParameters(), LogDebug, edmNew::makeRefTo(), edmNew::DetSetVector< T >::FastFiller::push_back(), PixelCPEBase::rawQualityWord(), ready_, edmNew::DetSetVector< T >::FastFiller::size(), and theVerboseLevel.

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  }
#define LogDebug(id)
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)
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator begin() const
SiPixelRecHitQuality::QualWordType rawQualityWord() const
data_type const * const_iterator
Definition: DetSetNew.h:25
virtual LocalValues localParameters(const T &, const GeomDetUnit &) const =0
const_iterator end() const
Definition: DetId.h:20
Pixel cluster – collection of neighboring pixels above threshold.
const PixelClusterParameterEstimator * cpe_
Pixel Reconstructed Hit.

Member Data Documentation

edm::ParameterSet cms::SiPixelRecHitConverter::conf_
private

Definition at line 95 of file SiPixelRecHitConverter.h.

Referenced by produce().

const PixelClusterParameterEstimator* cms::SiPixelRecHitConverter::cpe_
private

Definition at line 98 of file SiPixelRecHitConverter.h.

Referenced by produce(), and run().

std::string cms::SiPixelRecHitConverter::cpeName_
private

Definition at line 97 of file SiPixelRecHitConverter.h.

Referenced by produce(), and run().

bool cms::SiPixelRecHitConverter::m_newCont
private

Definition at line 103 of file SiPixelRecHitConverter.h.

bool cms::SiPixelRecHitConverter::ready_
private

Definition at line 100 of file SiPixelRecHitConverter.h.

Referenced by produce(), and run().

edm::InputTag cms::SiPixelRecHitConverter::src_
private

Definition at line 101 of file SiPixelRecHitConverter.h.

Referenced by produce().

int cms::SiPixelRecHitConverter::theVerboseLevel
private

Definition at line 102 of file SiPixelRecHitConverter.h.

Referenced by run().