#include <L1TkClusterBuilder.h>
//////////////////////////////////////// Stacked Tracker Simulations /// / Nicola Pozzobon, UNIPD /// / 2011, June /// ////////////////////////////////////// ************************ DECLARATION OF CLASS ************************
Definition at line 34 of file L1TkClusterBuilder.h.
L1TkClusterBuilder< T >::L1TkClusterBuilder | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Constructors.
Close class.
NOTE since pattern hit correlation must be performed within a stacked module, one must store Clusters in a proper way, providing easy access to them in a detector/member-wise way
***************************** IMPLEMENTATION OF METHODS ***************************** Constructors Default is for PixelDigis
Definition at line 74 of file L1TkClusterBuilder.h.
References edm::ParameterSet::getParameter().
{ rawHitInputTags = iConfig.getParameter< std::vector< edm::InputTag > >("rawHits"); simTrackInputTag = iConfig.getParameter< edm::InputTag >("simTrackHits"); ADCThreshold = iConfig.getParameter< unsigned int >("ADCThreshold"); produces< std::vector< L1TkCluster< T > > >(); }
L1TkClusterBuilder< T >::~L1TkClusterBuilder | ( | ) |
L1TkClusterBuilder< Ref_PSimHit_ >::L1TkClusterBuilder | ( | const edm::ParameterSet & | iConfig | ) |
Constructors Specialize for PSimHits
Definition at line 85 of file L1TkClusterBuilder.h.
References edm::ParameterSet::getParameter().
{ rawHitInputTags = iConfig.getParameter< std::vector< edm::InputTag > >("rawHits"); simTrackInputTag = iConfig.getParameter< edm::InputTag >("simTrackHits"); produces< std::vector< L1TkCluster< Ref_PSimHit_ > > >(); }
void L1TkClusterBuilder< T >::beginRun | ( | edm::Run & | run, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Mandatory methods.
Begin run.
Get the geometry
Get the clustering algorithm
Print some information when loaded
Reimplemented from edm::EDProducer.
Definition at line 98 of file L1TkClusterBuilder.h.
References gather_cfg::cout, edm::EventSetup::get(), and edm::ESHandle< T >::product().
{ edm::ESHandle< StackedTrackerGeometry > StackedTrackerGeomHandle; iSetup.get< StackedTrackerGeometryRecord >().get( StackedTrackerGeomHandle ); theStackedTrackers = StackedTrackerGeomHandle.product(); iSetup.get< ClusteringAlgorithmRecord >().get( ClusteringAlgoHandle ); std::cout << std::endl; std::cout << "L1TkClusterBuilder<" << templateNameFinder<T>() << "> loaded modules:" << "\n\tClusteringAlgorithm:\t" << ClusteringAlgoHandle->AlgorithmName() << std::endl; std::cout << std::endl; }
void L1TkClusterBuilder< T >::endRun | ( | edm::Run & | run, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
End run.
Reimplemented from edm::EDProducer.
Definition at line 118 of file L1TkClusterBuilder.h.
{}
void L1TkClusterBuilder< T >::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implement the producer.
Prepare output
Get the PixelDigiSimLink
Get the SimTracks
This is a map containing hits: a vector of type T is mapped wrt the DetId
Loop over the detector elements
Temp vectors containing the vectors of the hits used to build each cluster
Find the hits in each stack member
If there are hits, cluster them It is the ClusteringAlgorithm::Cluster method which calls the constructor to the Cluster class!
Create L1TkCluster objects and store them
If MC, check also fake or non/fake
End of loop over detector elements
Put output in the event
Implements edm::EDProducer.
Definition at line 122 of file L1TkClusterBuilder.h.
References funct::false, edm::Event::getByLabel(), i, StackedTrackerDetUnit::Id(), edm::EventBase::isRealData(), edm::Event::put(), StackedTrackerDetUnit::stackMember(), and groupFilesInBlocks::temp.
{ std::auto_ptr< std::vector< L1TkCluster< T > > > ClustersForOutput( new std::vector< L1TkCluster< T > > ); iEvent.getByLabel( "simSiPixelDigis", PixelDigiSimLinkHandle ); iEvent.getByLabel( simTrackInputTag, SimTrackHandle ); std::map< DetId, std::vector< T > > rawHits; this->RetrieveRawHits( rawHits, iEvent ); StackedTrackerGeometry::StackContainerIterator StackedTrackerIterator; for ( StackedTrackerIterator = theStackedTrackers->stacks().begin(); StackedTrackerIterator != theStackedTrackers->stacks().end(); ++StackedTrackerIterator ) { StackedTrackerDetUnit* Unit = *StackedTrackerIterator; StackedTrackerDetId Id = Unit->Id(); assert(Unit == theStackedTrackers->idToStack(Id)); std::vector< std::vector< T > > innerHits, outerHits; typename std::map< DetId, std::vector< T > >::const_iterator innerHitFind = rawHits.find( Unit->stackMember(0) ); typename std::map< DetId, std::vector< T > >::const_iterator outerHitFind = rawHits.find( Unit->stackMember(1) ); if ( innerHitFind != rawHits.end() ) ClusteringAlgoHandle->Cluster( innerHits, innerHitFind->second ); if ( outerHitFind != rawHits.end() ) ClusteringAlgoHandle->Cluster( outerHits, outerHitFind->second ); for ( unsigned int i = 0; i < innerHits.size(); i++ ) { L1TkCluster< T > temp( innerHits.at(i), Id, 0 ); if ( iEvent.isRealData() == false ) theStackedTrackers->checkSimTrack(&temp, PixelDigiSimLinkHandle, SimTrackHandle ); ClustersForOutput->push_back( temp ); } for ( unsigned int i = 0; i < outerHits.size(); i++ ) { L1TkCluster< T > temp( outerHits.at(i), Id, 1 ); if ( iEvent.isRealData() == false ) theStackedTrackers->checkSimTrack( &temp, PixelDigiSimLinkHandle, SimTrackHandle ); ClustersForOutput->push_back( temp ); } } iEvent.put( ClustersForOutput ); }
void L1TkClusterBuilder< T >::RetrieveRawHits | ( | std::map< DetId, std::vector< T > > & | mRawHits, |
const edm::Event & | iEvent | ||
) | [private] |
Get hits.
void L1TkClusterBuilder< Ref_PixelDigi_ >::RetrieveRawHits | ( | std::map< DetId, std::vector< Ref_PixelDigi_ > > & | mRawHits, |
const edm::Event & | iEvent | ||
) | [private] |
Retrieve hits from the event Specialize template for PixelDigis
Loop over the tags used to identify hits in the cfg file
For each tag, get the corresponding handle
Loop over detector elements identifying PixelDigis
Is it Pixel?
Loop over Digis in this specific detector element
If the Digi is over threshold, accept it as a raw hit and put into map
End of threshold selection
End of loop over digis
End of "is Pixel"
End of loop over detector elements
End of loop over tags
Definition at line 190 of file L1TkClusterBuilder.h.
References edm::DetSetVector< T >::begin(), edm::DetSetVector< T >::end(), edm::Event::getByLabel(), and edm::makeRefTo().
{ mRawHits.clear(); std::vector< edm::InputTag >::iterator it; for ( it = rawHitInputTags.begin(); it != rawHitInputTags.end(); ++it ) { edm::Handle< edm::DetSetVector< PixelDigi > > HitHandle; iEvent.getByLabel( *it, HitHandle ); edm::DetSetVector<PixelDigi>::const_iterator detsIter; edm::DetSet<PixelDigi>::const_iterator hitsIter; for ( detsIter = HitHandle->begin(); detsIter != HitHandle->end(); detsIter++ ) { DetId id = detsIter->id; if ( id.subdetId()==1 || id.subdetId()==2 ) { for ( hitsIter = detsIter->data.begin(); hitsIter != detsIter->data.end(); hitsIter++ ) { if ( hitsIter->adc() >= ADCThreshold ) { mRawHits[id].push_back( makeRefTo( HitHandle, id , hitsIter ) ); } } } } } }
void L1TkClusterBuilder< Ref_PSimHit_ >::RetrieveRawHits | ( | std::map< DetId, std::vector< Ref_PSimHit_ > > & | mRawHits, |
const edm::Event & | iEvent | ||
) | [private] |
Retrieve hits from the event Specialize template for SimHits
Loop over the tags used to identify hits in the cfg file
For each tag, get the corresponding handle
Loop over the hits
If Pixel, directly map it
End of loop over hits
End of loop over tags
Definition at line 237 of file L1TkClusterBuilder.h.
References edm::Event::getByLabel(), and i.
{ mRawHits.clear(); for ( std::vector<edm::InputTag>::iterator it = rawHitInputTags.begin(); it != rawHitInputTags.end(); ++it ) { edm::Handle< edm::PSimHitContainer > HitHandle; iEvent.getByLabel( *it, HitHandle ); for ( unsigned int i=0 ; i != HitHandle->size() ; ++i ) { DetId id = DetId(HitHandle->at(i).detUnitId()); if ( id.subdetId()==1 || id.subdetId()==2 ) mRawHits[id].push_back( edm::Ref<edm::PSimHitContainer>( HitHandle, i ) ); } } }
unsigned int L1TkClusterBuilder< T >::ADCThreshold [private] |
Definition at line 53 of file L1TkClusterBuilder.h.
edm::ESHandle< ClusteringAlgorithm< T > > L1TkClusterBuilder< T >::ClusteringAlgoHandle [private] |
Definition at line 48 of file L1TkClusterBuilder.h.
edm::Handle< edm::DetSetVector< PixelDigiSimLink > > L1TkClusterBuilder< T >::PixelDigiSimLinkHandle [private] |
Definition at line 49 of file L1TkClusterBuilder.h.
std::vector< edm::InputTag > L1TkClusterBuilder< T >::rawHitInputTags [private] |
Definition at line 51 of file L1TkClusterBuilder.h.
edm::Handle< edm::SimTrackContainer > L1TkClusterBuilder< T >::SimTrackHandle [private] |
Definition at line 50 of file L1TkClusterBuilder.h.
edm::InputTag L1TkClusterBuilder< T >::simTrackInputTag [private] |
Definition at line 52 of file L1TkClusterBuilder.h.
const StackedTrackerGeometry* L1TkClusterBuilder< T >::theStackedTrackers [private] |
Data members.
Definition at line 47 of file L1TkClusterBuilder.h.