CMS 3D CMS Logo

SiPixelRecHitConverter.cc

Go to the documentation of this file.
00001 
00012 // Our own stuff
00013 #include "RecoLocalTracker/SiPixelRecHits/interface/SiPixelRecHitConverter.h"
00014 #include "RecoLocalTracker/SiPixelRecHits/interface/CPEFromDetPosition.h"
00015 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEInitial.h"
00016 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEParmError.h"
00017 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPETemplateReco.h"
00018 #include "RecoLocalTracker/SiPixelRecHits/interface/PixelCPEGeneric.h"
00019 
00020 // Geometry
00021 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00022 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00023 
00024 // Data Formats
00025 #include "DataFormats/DetId/interface/DetId.h"
00026 #include "DataFormats/Common/interface/Ref.h"
00027 #include "DataFormats/Common/interface/DetSet2RangeMap.h"
00028 
00029 // Framework Already defined in the *.h file
00030 //#include "DataFormats/Common/interface/Handle.h"
00031 //#include "FWCore/Framework/interface/ESHandle.h"
00032 
00033 // STL
00034 #include <vector>
00035 #include <memory>
00036 #include <string>
00037 #include <iostream>
00038 
00039 // MessageLogger
00040 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00041 
00042 #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h"
00043 
00044 using namespace std;
00045 
00046 namespace cms
00047 {
00048   //---------------------------------------------------------------------------
00050   //---------------------------------------------------------------------------
00051   SiPixelRecHitConverter::SiPixelRecHitConverter(edm::ParameterSet const& conf) 
00052     : 
00053     conf_(conf),
00054     cpeName_("None"),     // bogus
00055     cpe_(0),              // the default, in case we fail to make one
00056     ready_(false),        // since we obviously aren't
00057     src_( conf.getParameter<edm::InputTag>( "src" ) ),
00058     theVerboseLevel(conf.getUntrackedParameter<int>("VerboseLevel",0)),
00059     m_newCont(conf.getUntrackedParameter<bool>("newContainer",false))
00060   {
00061     //--- Declare to the EDM what kind of collections we will be making.
00062     produces<SiPixelRecHitCollection>();
00063     if (m_newCont)
00064       produces<SiPixelRecHitCollectionNew>();
00065    
00066   }
00067   
00068   // Destructor
00069   SiPixelRecHitConverter::~SiPixelRecHitConverter() 
00070   { 
00071   }  
00072   
00073   //---------------------------------------------------------------------------
00074   // Begin job: get magnetic field
00075   //---------------------------------------------------------------------------
00076   void SiPixelRecHitConverter::beginJob(const edm::EventSetup& c) 
00077   {
00078   }
00079 
00080   //---------------------------------------------------------------------------
00082   //---------------------------------------------------------------------------
00083   void SiPixelRecHitConverter::produce(edm::Event& e, const edm::EventSetup& es)
00084   {
00085 
00086     // Step A.1: get input data
00087     edm::Handle< edmNew::DetSetVector<SiPixelCluster> > input;
00088     e.getByLabel( src_, input);
00089     
00090     // Step A.2: get event setup
00091     edm::ESHandle<TrackerGeometry> geom;
00092     es.get<TrackerDigiGeometryRecord>().get( geom );
00093 
00094                 // Step B: create empty output collection
00095     std::auto_ptr<SiPixelRecHitCollectionNew> output(new SiPixelRecHitCollectionNew);
00096 
00097                 // Step B*: create CPE
00098                 edm::ESHandle<PixelClusterParameterEstimator> hCPE;
00099                 std::string cpeName_ = conf_.getParameter<std::string>("CPE");
00100                 es.get<TkPixelCPERecord>().get(cpeName_,hCPE);
00101                 const PixelClusterParameterEstimator &cpe(*hCPE);
00102                 cpe_ = &cpe;
00103                 
00104                 if(cpe_) ready_ = true;
00105 
00106 
00107     // Step C: Iterate over DetIds and invoke the strip CPE algorithm
00108     // on each DetUnit
00109 
00110     run( input, *output, geom );
00111 
00112 
00113     // Step Z: temporary write also the old collection
00114     std::auto_ptr<SiPixelRecHitCollection> old(new SiPixelRecHitCollection);
00115     edmNew::copy(*output,*old);
00116     e.put(old);
00117     
00118     
00119     // Step D: write output to file
00120     if (m_newCont)
00121       e.put(output);
00122 
00123   }
00124 
00125   //---------------------------------------------------------------------------
00129   //---------------------------------------------------------------------------
00130         void SiPixelRecHitConverter::run(edm::Handle<edmNew::DetSetVector<SiPixelCluster> >  inputhandle,
00131                                         SiPixelRecHitCollectionNew &output,
00132                                         edm::ESHandle<TrackerGeometry> & geom) 
00133         {
00134                 if ( ! ready_ ) 
00135                 {
00136                         edm::LogError("SiPixelRecHitConverter") << " at least one CPE is not ready -- can't run!";
00137                         // TO DO: throw an exception here?  The user may want to know...
00138                         assert(0);
00139                         return;   // clusterizer is invalid, bail out
00140                 }
00141         
00142                 int numberOfDetUnits = 0;
00143                 int numberOfClusters = 0;
00144                 
00145                 const edmNew::DetSetVector<SiPixelCluster>& input = *inputhandle;
00146                 
00147                 edmNew::DetSetVector<SiPixelCluster>::const_iterator DSViter=input.begin();
00148 
00149                 for ( ; DSViter != input.end() ; DSViter++) 
00150                 {
00151                         numberOfDetUnits++;
00152                         unsigned int detid = DSViter->detId();
00153                         DetId detIdObject( detid );  
00154                         const GeomDetUnit * genericDet = geom->idToDetUnit( detIdObject );
00155                         const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(genericDet);
00156                         assert(pixDet); 
00157                         SiPixelRecHitCollectionNew::FastFiller recHitsOnDetUnit(output,detid);
00158                         
00159                         edmNew::DetSet<SiPixelCluster>::const_iterator clustIt = DSViter->begin(), clustEnd = DSViter->end();
00160 
00161                         for ( ; clustIt != clustEnd; clustIt++) 
00162                         {
00163                                 numberOfClusters++;
00164                                 std::pair<LocalPoint, LocalError> lv = cpe_->localParameters( *clustIt, *genericDet );
00165                                 LocalPoint lp( lv.first );
00166                                 LocalError le( lv.second );
00167                                 // Create a persistent edm::Ref to the cluster
00168                                 edm::Ref< edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster > cluster = edmNew::makeRefTo( inputhandle, clustIt);
00169                                 // Make a RecHit and add it to the DetSet
00170                                 // old : recHitsOnDetUnit.push_back( new SiPixelRecHit( lp, le, detIdObject, &*clustIt) );
00171                                 SiPixelRecHit hit( lp, le, detIdObject, cluster);
00172                                 #ifdef SIPIXELRECHIT_HAS_EXTRA_INFO
00173                                 // Copy the extra stuff; unfortunately, only the derivatives of PixelCPEBase
00174                                 // are capable of doing that.  So until we get rid of CPEFromDetPosition
00175                                 // we'll have to dynamic_cast :(
00176                                 // &&& This cast can be moved to the setupCPE, so that it's done once per job.
00177                                 PixelCPEBase * cpeBase = dynamic_cast< PixelCPEBase* >( cpe_ );
00178                                 if (cpeBase) {
00179                                         hit.setProbabilityX( cpeBase->probabilityX() );
00180                                         hit.setProbabilityY( cpeBase->probabilityY() );
00181                                         hit.setQBin( cpeBase->qBin() );
00182                                         hit.setCotAlphaFromCluster( cpeBase->cotAlphaFromCluster() );
00183                                         hit.setCotBetaFromCluster ( cpeBase->cotBetaFromCluster()  );
00184                                 }
00185                                 #endif
00186                                 // 
00187                                 // Now save it =================
00188                                 recHitsOnDetUnit.push_back(hit);
00189                                 // =============================
00190                         } //  <-- End loop on Clusters
00191                         
00192                         if ( recHitsOnDetUnit.size()>0 ) 
00193                         {
00194                                 if (theVerboseLevel > 2) 
00195                                 LogDebug("SiPixelRecHitConverter") << " Found " 
00196                                                                 << recHitsOnDetUnit.size() << " RecHits on " << detid;  
00197                         }
00198                         
00199                 } //    <-- End loop on DetUnits
00200                 
00201                 if ( theVerboseLevel > 2 ) LogDebug ("SiPixelRecHitConverter") 
00202                 << cpeName_ << " converted " << numberOfClusters 
00203                 << " SiPixelClusters into SiPixelRecHits, in " 
00204                 << numberOfDetUnits << " DetUnits.";    
00205         }
00206 }  // end of namespace cms

Generated on Tue Jun 9 17:43:59 2009 for CMSSW by  doxygen 1.5.4