#include <SiPixelRecHitConverter.h>
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. | |
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(). | |
virtual | ~SiPixelRecHitConverter () |
Private Attributes | |
edm::ParameterSet | conf_ |
const PixelClusterParameterEstimator * | cpe_ |
std::string | cpeName_ |
bool | m_newCont |
bool | ready_ |
edm::InputTag | src_ |
int | theVerboseLevel |
Definition at line 62 of file SiPixelRecHitConverter.h.
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.
: conf_(conf), cpeName_("None"), // bogus cpe_(0), // the default, in case we fail to make one ready_(false), // since we obviously aren't src_( conf.getParameter<edm::InputTag>( "src" ) ), theVerboseLevel(conf.getUntrackedParameter<int>("VerboseLevel",0)) { //--- Declare to the EDM what kind of collections we will be making. produces<SiPixelRecHitCollection>(); }
SiPixelRecHitConverter::~SiPixelRecHitConverter | ( | ) | [virtual] |
Definition at line 63 of file SiPixelRecHitConverter.cc.
{ }
void SiPixelRecHitConverter::beginJob | ( | void | ) | [virtual] |
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(), collect_tpl::input, convertSQLitetoXML_cfg::output, edm::Event::put(), ready_, run(), and src_.
{ // Step A.1: get input data edm::Handle< edmNew::DetSetVector<SiPixelCluster> > input; e.getByLabel( src_, input); // Step A.2: get event setup edm::ESHandle<TrackerGeometry> geom; es.get<TrackerDigiGeometryRecord>().get( geom ); // Step B: create empty output collection std::auto_ptr<SiPixelRecHitCollectionNew> output(new SiPixelRecHitCollectionNew); // Step B*: create CPE edm::ESHandle<PixelClusterParameterEstimator> hCPE; std::string cpeName_ = conf_.getParameter<std::string>("CPE"); es.get<TkPixelCPERecord>().get(cpeName_,hCPE); const PixelClusterParameterEstimator &cpe(*hCPE); cpe_ = &cpe; if(cpe_) ready_ = true; // Step C: Iterate over DetIds and invoke the strip CPE algorithm // on each DetUnit run( input, *output, geom ); e.put(output); }
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(), collect_tpl::input, asciidump::le, ClusterParameterEstimator< T >::localParameters(), LogDebug, edm::makeRefTo(), edmNew::DetSetVector< T >::FastFiller::push_back(), PixelCPEBase::rawQualityWord(), ready_, edmNew::DetSetVector< T >::FastFiller::size(), and theVerboseLevel.
{ if ( ! ready_ ) { edm::LogError("SiPixelRecHitConverter") << " at least one CPE is not ready -- can't run!"; // TO DO: throw an exception here? The user may want to know... assert(0); return; // clusterizer is invalid, bail out } int numberOfDetUnits = 0; int numberOfClusters = 0; const edmNew::DetSetVector<SiPixelCluster>& input = *inputhandle; edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=input.begin(); for ( ; DSViter != input.end() ; DSViter++) { numberOfDetUnits++; unsigned int detid = DSViter->detId(); DetId detIdObject( detid ); const GeomDetUnit * genericDet = geom->idToDetUnit( detIdObject ); const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet); assert(pixDet); SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(output,detid); edmNew::DetSet<SiPixelCluster>::const_iterator clustIt = DSViter->begin(), clustEnd = DSViter->end(); for ( ; clustIt != clustEnd; clustIt++) { numberOfClusters++; std::pair<LocalPoint, LocalError> lv = cpe_->localParameters( *clustIt, *genericDet ); LocalPoint lp( lv.first ); LocalError le( lv.second ); // Create a persistent edm::Ref to the cluster edm::Ref< edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster > cluster = edmNew::makeRefTo( inputhandle, clustIt); // Make a RecHit and add it to the DetSet // old : recHitsOnDetUnit.push_back( new SiPixelRecHit( lp, le, detIdObject, &*clustIt) ); SiPixelRecHit hit( lp, le, detIdObject, cluster); // Copy the extra stuff; unfortunately, only the derivatives of PixelCPEBase // are capable of doing that. So until we get rid of CPEFromDetPosition // we'll have to dynamic_cast :( // &&& This cast can be moved to the setupCPE, so that it's done once per job. const PixelCPEBase * cpeBase = dynamic_cast< const PixelCPEBase* >( cpe_ ); if (cpeBase) { hit.setRawQualityWord( cpeBase->rawQualityWord() ); // hit.setProbabilityX( cpeBase->probabilityX() ); // hit.setProbabilityY( cpeBase->probabilityY() ); // hit.setQBin( cpeBase->qBin() ); // hit.setCotAlphaFromCluster( cpeBase->cotAlphaFromCluster() ); // hit.setCotBetaFromCluster ( cpeBase->cotBetaFromCluster() ); } // // Now save it ================= recHitsOnDetUnit.push_back(hit); // ============================= } // <-- End loop on Clusters if ( recHitsOnDetUnit.size()>0 ) { if (theVerboseLevel > 2) LogDebug("SiPixelRecHitConverter") << " Found " << recHitsOnDetUnit.size() << " RecHits on " << detid; } } // <-- End loop on DetUnits if ( theVerboseLevel > 2 ) LogDebug ("SiPixelRecHitConverter") << cpeName_ << " converted " << numberOfClusters << " SiPixelClusters into SiPixelRecHits, in " << numberOfDetUnits << " DetUnits."; }
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.
std::string cms::SiPixelRecHitConverter::cpeName_ [private] |
Definition at line 97 of file SiPixelRecHitConverter.h.
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.
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().