CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/SLHCUpgradeSimulations/L1TrackTrigger/interface/L1TkClusterBuilder.h

Go to the documentation of this file.
00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 #ifndef CLUSTER_BUILDER_H
00010 #define CLUSTER_BUILDER_H
00011 
00012 #include "FWCore/Framework/interface/Frameworkfwd.h"
00013 #include "FWCore/Framework/interface/EDProducer.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 
00020 #include "SLHCUpgradeSimulations/L1TrackTrigger/interface/ClusteringAlgorithm.h"
00021 #include "SLHCUpgradeSimulations/L1TrackTrigger/interface/ClusteringAlgorithmRecord.h"
00022 #include "Geometry/TrackerGeometryBuilder/interface/StackedTrackerDetUnit.h"
00023 #include <memory>
00024 #include <map>
00025 #include <vector>
00026 
00033 template< typename T >
00034 class L1TkClusterBuilder : public edm::EDProducer
00035 {
00038   public:
00040     explicit L1TkClusterBuilder( const edm::ParameterSet& iConfig );
00041 
00043     ~L1TkClusterBuilder();
00044 
00045   private:
00047     const StackedTrackerGeometry              *theStackedTrackers;
00048     edm::ESHandle< ClusteringAlgorithm< T > > ClusteringAlgoHandle; // Handles are needed in ::produce()
00049     edm::Handle< edm::DetSetVector< PixelDigiSimLink > >   PixelDigiSimLinkHandle;
00050     edm::Handle< edm::SimTrackContainer >                  SimTrackHandle;
00051     std::vector< edm::InputTag >                           rawHitInputTags;
00052     edm::InputTag                                          simTrackInputTag;
00053     unsigned int                                           ADCThreshold;  
00054 
00056     virtual void beginRun( edm::Run& run, const edm::EventSetup& iSetup );
00057     virtual void endRun( edm::Run& run, const edm::EventSetup& iSetup );
00058     virtual void produce( edm::Event& iEvent, const edm::EventSetup& iSetup );
00059 
00061     void RetrieveRawHits( std::map< DetId, std::vector< T > > &mRawHits, const edm::Event& iEvent );
00062 
00063 }; 
00064 
00071 
00072 
00073 template< typename T >
00074 L1TkClusterBuilder< T >::L1TkClusterBuilder( const edm::ParameterSet& iConfig )
00075 {
00076   rawHitInputTags  = iConfig.getParameter< std::vector< edm::InputTag > >("rawHits");
00077   simTrackInputTag = iConfig.getParameter< edm::InputTag >("simTrackHits");
00078   ADCThreshold     = iConfig.getParameter< unsigned int >("ADCThreshold");
00079   produces< std::vector< L1TkCluster< T > > >();
00080 }
00081 
00084 template<>
00085 L1TkClusterBuilder< Ref_PSimHit_ >::L1TkClusterBuilder( const edm::ParameterSet& iConfig )
00086 {
00087   rawHitInputTags  = iConfig.getParameter< std::vector< edm::InputTag > >("rawHits");
00088   simTrackInputTag = iConfig.getParameter< edm::InputTag >("simTrackHits");
00089   produces< std::vector< L1TkCluster< Ref_PSimHit_ > > >();
00090 }
00091 
00093 template< typename T >
00094 L1TkClusterBuilder< T >::~L1TkClusterBuilder(){}
00095 
00097 template< typename T >
00098 void L1TkClusterBuilder< T >::beginRun( edm::Run& run, const edm::EventSetup& iSetup )
00099 {
00101   edm::ESHandle< StackedTrackerGeometry > StackedTrackerGeomHandle;
00102   iSetup.get< StackedTrackerGeometryRecord >().get( StackedTrackerGeomHandle );
00103   theStackedTrackers = StackedTrackerGeomHandle.product();
00104 
00106   iSetup.get< ClusteringAlgorithmRecord >().get( ClusteringAlgoHandle );
00107 
00109   std::cout << std::endl;
00110   std::cout << "L1TkClusterBuilder<" << templateNameFinder<T>() << "> loaded modules:"
00111             << "\n\tClusteringAlgorithm:\t" << ClusteringAlgoHandle->AlgorithmName()
00112             << std::endl;
00113   std::cout << std::endl;
00114 }
00115 
00117 template< typename T >
00118 void L1TkClusterBuilder< T >::endRun( edm::Run& run, const edm::EventSetup& iSetup ){}
00119 
00121 template< typename T >
00122 void L1TkClusterBuilder< T >::produce( edm::Event& iEvent, const edm::EventSetup& iSetup )
00123 {
00125   std::auto_ptr< std::vector< L1TkCluster< T > > > ClustersForOutput( new std::vector< L1TkCluster< T > > );
00126 
00128   iEvent.getByLabel( "simSiPixelDigis", PixelDigiSimLinkHandle );
00129 
00131   iEvent.getByLabel( simTrackInputTag, SimTrackHandle );
00132 
00133   std::map< DetId, std::vector< T > > rawHits; 
00134 
00135 
00136   this->RetrieveRawHits( rawHits, iEvent );
00137 
00139   StackedTrackerGeometry::StackContainerIterator StackedTrackerIterator;
00140   for ( StackedTrackerIterator = theStackedTrackers->stacks().begin();
00141         StackedTrackerIterator != theStackedTrackers->stacks().end();
00142         ++StackedTrackerIterator )
00143   {
00144     StackedTrackerDetUnit* Unit = *StackedTrackerIterator;
00145     StackedTrackerDetId Id = Unit->Id();
00146     assert(Unit == theStackedTrackers->idToStack(Id));
00147 
00150     std::vector< std::vector< T > > innerHits, outerHits;
00151 
00153     typename std::map< DetId, std::vector< T > >::const_iterator innerHitFind = rawHits.find( Unit->stackMember(0) );
00154     typename std::map< DetId, std::vector< T > >::const_iterator outerHitFind = rawHits.find( Unit->stackMember(1) );
00155 
00159     if ( innerHitFind != rawHits.end() ) ClusteringAlgoHandle->Cluster( innerHits, innerHitFind->second );
00160     if ( outerHitFind != rawHits.end() ) ClusteringAlgoHandle->Cluster( outerHits, outerHitFind->second );
00161 
00163     for ( unsigned int i = 0; i < innerHits.size(); i++ )
00164     {
00165       L1TkCluster< T > temp( innerHits.at(i), Id, 0 );
00167       if ( iEvent.isRealData() == false )
00168         theStackedTrackers->checkSimTrack(&temp, PixelDigiSimLinkHandle, SimTrackHandle ); 
00169 
00170       ClustersForOutput->push_back( temp );
00171     }
00172     for ( unsigned int i = 0; i < outerHits.size(); i++ )
00173     {
00174       L1TkCluster< T > temp( outerHits.at(i), Id, 1 );
00175       if ( iEvent.isRealData() == false )
00176         theStackedTrackers->checkSimTrack( &temp, PixelDigiSimLinkHandle, SimTrackHandle ); 
00177 
00178       ClustersForOutput->push_back( temp );
00179     }
00180 
00181   } 
00182 
00184   iEvent.put( ClustersForOutput );
00185 }
00186 
00189 template<>
00190 void L1TkClusterBuilder< Ref_PixelDigi_ >::RetrieveRawHits( std::map< DetId, std::vector< Ref_PixelDigi_ > > &mRawHits,
00191                                                                          const edm::Event& iEvent )
00192 {
00193   mRawHits.clear();
00195   std::vector< edm::InputTag >::iterator it;
00196   for ( it = rawHitInputTags.begin();
00197         it != rawHitInputTags.end();
00198         ++it )
00199   {
00201     edm::Handle< edm::DetSetVector< PixelDigi > > HitHandle;
00202     iEvent.getByLabel( *it, HitHandle );
00203 
00204     edm::DetSetVector<PixelDigi>::const_iterator detsIter;
00205     edm::DetSet<PixelDigi>::const_iterator       hitsIter;
00206 
00208     for ( detsIter = HitHandle->begin();
00209           detsIter != HitHandle->end();
00210           detsIter++ )
00211     {
00212       DetId id = detsIter->id;
00213 
00215       if ( id.subdetId()==1 || id.subdetId()==2 )
00216       {
00218         for ( hitsIter = detsIter->data.begin();
00219               hitsIter != detsIter->data.end();
00220               hitsIter++ )
00221         {
00222           if ( hitsIter->adc() >= ADCThreshold )
00223           {
00226             mRawHits[id].push_back( makeRefTo( HitHandle, id , hitsIter ) );
00227           } 
00228         } 
00229       } 
00230     } 
00231   } 
00232 }
00233 
00236 template<>
00237 void L1TkClusterBuilder< Ref_PSimHit_ >::RetrieveRawHits( std::map< DetId, std::vector< Ref_PSimHit_ > > &mRawHits,
00238                                                                        const edm::Event& iEvent )
00239 {
00240   mRawHits.clear();
00242   for ( std::vector<edm::InputTag>::iterator it = rawHitInputTags.begin();
00243         it != rawHitInputTags.end();
00244         ++it )
00245   {
00247     edm::Handle< edm::PSimHitContainer > HitHandle;
00248     iEvent.getByLabel( *it, HitHandle );
00249 
00251     for ( unsigned int i=0 ; i != HitHandle->size() ; ++i )
00252     {
00253       DetId id = DetId(HitHandle->at(i).detUnitId());
00255       if ( id.subdetId()==1 || id.subdetId()==2 )
00256         mRawHits[id].push_back( edm::Ref<edm::PSimHitContainer>( HitHandle, i ) );
00257     } 
00258   } 
00259 }
00260 
00261 #endif
00262