1 |
dgele |
1.1 |
#include "PhysicsTools/RecoAlgos/interface/TrackSelector.h"
|
2 |
|
|
|
3 |
|
|
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
|
4 |
|
|
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
|
5 |
|
|
#include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
|
6 |
|
|
#include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
|
7 |
|
|
|
8 |
|
|
#include <DataFormats/SiStripDetId/interface/SiStripDetId.h>
|
9 |
|
|
|
10 |
|
|
template<typename RecHitType, typename ClusterRefType>
|
11 |
|
|
void
|
12 |
|
|
helper::TrackCollectionStoreManager::ClusterHitRecord<RecHitType,ClusterRefType>::
|
13 |
|
|
rekey(TrackingRecHitCollection &hits, const ClusterRefType &newRef) const
|
14 |
|
|
{
|
15 |
|
|
//std::cout << "Rekeying hit with index " << index_ << ", detid = " << detid_ << std::endl;
|
16 |
|
|
TrackingRecHit & genericHit = hits[index_];
|
17 |
|
|
RecHitType * hit = 0;
|
18 |
|
|
if (genericHit.geographicalId().rawId() == detid_) { // a hit on this det, so it's simple
|
19 |
|
|
hit = static_cast<RecHitType *>(&genericHit);
|
20 |
|
|
} else { // projected or matched rechit
|
21 |
|
|
assert ( typeid(RecHitType) == typeid(SiStripRecHit2D) ); // must not happen for pixel
|
22 |
|
|
if (typeid(genericHit) == typeid(SiStripMatchedRecHit2D)) {
|
23 |
|
|
SiStripMatchedRecHit2D & mhit = static_cast<SiStripMatchedRecHit2D &>(genericHit);
|
24 |
|
|
// I need the reinterpret_cast because the compiler sees this code even if RecHitType = PixelRecHit
|
25 |
|
|
hit = reinterpret_cast<RecHitType *>(SiStripDetId(detid_).stereo() ? mhit.stereoHit() : mhit.monoHit());
|
26 |
|
|
} else {
|
27 |
|
|
assert (typeid(genericHit) == typeid(ProjectedSiStripRecHit2D)) ; // no option left, so this shoud be true
|
28 |
|
|
ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D &>(genericHit);
|
29 |
|
|
hit = reinterpret_cast<RecHitType *>(& phit.originalHit());
|
30 |
|
|
}
|
31 |
|
|
}
|
32 |
|
|
assert (hit != 0);
|
33 |
|
|
assert ( hit->cluster() == ref_ ); // otherwise something went wrong
|
34 |
|
|
hit->setClusterRef(newRef);
|
35 |
|
|
}
|
36 |
|
|
|
37 |
|
|
using namespace reco;
|
38 |
|
|
|
39 |
|
|
namespace helper
|
40 |
|
|
{
|
41 |
|
|
|
42 |
|
|
TrackCollectionStoreManager::
|
43 |
|
|
TrackCollectionStoreManager(const edm::Handle<reco::TrackCollection> & )
|
44 |
|
|
:
|
45 |
|
|
selTracks_( new reco::TrackCollection ),
|
46 |
|
|
selTrackExtras_( new reco::TrackExtraCollection ),
|
47 |
|
|
selHits_( new TrackingRecHitCollection ),
|
48 |
|
|
selStripClusters_( new edmNew::DetSetVector<SiStripCluster> ),
|
49 |
|
|
selPixelClusters_( new edmNew::DetSetVector<SiPixelCluster> ),
|
50 |
|
|
rTracks_(), rTrackExtras_(), rHits_(), rStripClusters_(), rPixelClusters_(),
|
51 |
|
|
idx_(0), hidx_(0),
|
52 |
|
|
cloneClusters_ (true)
|
53 |
|
|
{
|
54 |
|
|
}
|
55 |
|
|
|
56 |
|
|
|
57 |
|
|
|
58 |
|
|
//------------------------------------------------------------------
|
59 |
|
|
//! Process a single track.
|
60 |
|
|
//------------------------------------------------------------------
|
61 |
|
|
void
|
62 |
|
|
TrackCollectionStoreManager::
|
63 |
|
|
processTrack( const Track & trk )
|
64 |
|
|
{
|
65 |
|
|
//std::cout << "Process track " << idx_ << ", pt = " << trk.pt() << ", eta = " << trk.eta() << ", phi = " << trk.phi() << std::endl;
|
66 |
|
|
selTracks_->push_back( Track( trk ) );
|
67 |
|
|
selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, idx_ ++ ) );
|
68 |
|
|
selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
|
69 |
|
|
trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
|
70 |
|
|
trk.outerStateCovariance(), trk.outerDetId(),
|
71 |
|
|
trk.innerStateCovariance(), trk.innerDetId(),
|
72 |
|
|
trk.seedDirection() ) );
|
73 |
|
|
TrackExtra & tx = selTrackExtras_->back();
|
74 |
|
|
for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
|
75 |
|
|
//std::cout << "\\-- Process hit " << hidx_ << ", detid = " << (*hit)->geographicalId()() << std::endl;
|
76 |
|
|
selHits_->push_back( (*hit)->clone() );
|
77 |
|
|
TrackingRecHit * newHit = & (selHits_->back());
|
78 |
|
|
tx.add( TrackingRecHitRef( rHits_, hidx_ ++ ) );
|
79 |
|
|
|
80 |
|
|
//--- Skip the rest for this hit if we don't want to clone the cluster.
|
81 |
|
|
//--- The copy constructer in the rec hit will copy the link properly.
|
82 |
|
|
//
|
83 |
|
|
if (cloneClusters() == false)
|
84 |
|
|
continue; // go to the next hit on the track
|
85 |
|
|
|
86 |
|
|
//std::cout << "| I'm cloing clusters" << std::endl;
|
87 |
|
|
|
88 |
|
|
const DetId detId( (*hit)->geographicalId() );
|
89 |
|
|
if (newHit->isValid() && (detId.det() == DetId::Tracker)) {
|
90 |
|
|
//std::cout << "| It is a tracker hit" << std::endl;
|
91 |
|
|
|
92 |
|
|
const std::type_info & hit_type = typeid(*newHit);
|
93 |
|
|
if (hit_type == typeid(SiPixelRecHit)) {
|
94 |
|
|
//std::cout << "| It is a Pixel hit !!" << std::endl;
|
95 |
|
|
pixelClusterRecords_.push_back( PixelClusterHitRecord( static_cast<SiPixelRecHit &>(*newHit), hidx_ - 1) );
|
96 |
|
|
} else if (hit_type == typeid(SiStripRecHit2D)) {
|
97 |
|
|
//std::cout << "| It is a SiStripRecHit2D hit !!" << std::endl;
|
98 |
|
|
stripClusterRecords_.push_back( StripClusterHitRecord( static_cast<SiStripRecHit2D &>(*newHit), hidx_ - 1) );
|
99 |
|
|
} else if (hit_type == typeid(SiStripMatchedRecHit2D)) {
|
100 |
|
|
//std::cout << "| It is a SiStripMatchedRecHit2D hit !!" << std::endl;
|
101 |
|
|
SiStripMatchedRecHit2D & mhit = static_cast<SiStripMatchedRecHit2D &>(*newHit);
|
102 |
|
|
stripClusterRecords_.push_back( StripClusterHitRecord( *mhit.monoHit() , hidx_ - 1) );
|
103 |
|
|
stripClusterRecords_.push_back( StripClusterHitRecord( *mhit.stereoHit(), hidx_ - 1) );
|
104 |
|
|
} else if (hit_type == typeid(ProjectedSiStripRecHit2D)) {
|
105 |
|
|
//std::cout << "| It is a ProjectedSiStripRecHit2D hit !!" << std::endl;
|
106 |
|
|
ProjectedSiStripRecHit2D & phit = static_cast<ProjectedSiStripRecHit2D &>(*newHit);
|
107 |
|
|
stripClusterRecords_.push_back( StripClusterHitRecord( phit.originalHit(), hidx_ - 1) );
|
108 |
|
|
} else {
|
109 |
|
|
//std::cout << "| It is a " << hit_type.name() << " hit !?" << std::endl;
|
110 |
|
|
// do nothing. We might end up here for FastSim hits.
|
111 |
|
|
} // end 'switch' on hit type
|
112 |
|
|
} // end if it was a tracker hit
|
113 |
|
|
} // end of for loop over tracking rec hits on this track
|
114 |
|
|
} // end of track, and function
|
115 |
|
|
|
116 |
|
|
void
|
117 |
|
|
TrackCollectionStoreManager::
|
118 |
|
|
processAllClusters()
|
119 |
|
|
{
|
120 |
|
|
if (!pixelClusterRecords_.empty()) {
|
121 |
|
|
processClusters<SiPixelRecHit, SiPixelCluster>(pixelClusterRecords_, *selPixelClusters_, rPixelClusters_ );
|
122 |
|
|
}
|
123 |
|
|
if (!stripClusterRecords_.empty()) {
|
124 |
|
|
processClusters<SiStripRecHit2D,SiStripCluster>(stripClusterRecords_, *selStripClusters_, rStripClusters_ );
|
125 |
|
|
}
|
126 |
|
|
}
|
127 |
|
|
|
128 |
|
|
template<typename HitType, typename ClusterType>
|
129 |
|
|
void
|
130 |
|
|
TrackCollectionStoreManager::
|
131 |
|
|
processClusters( std::vector<ClusterHitRecord<HitType> > & clusterRecords,
|
132 |
|
|
edmNew::DetSetVector<ClusterType> & dsv,
|
133 |
|
|
edm::RefProd< edmNew::DetSetVector<ClusterType> > & refprod )
|
134 |
|
|
{
|
135 |
|
|
std::sort(clusterRecords.begin(), clusterRecords.end()); // this sorts them by detid
|
136 |
|
|
typedef typename std::vector<ClusterHitRecord<HitType> >::const_iterator RIT;
|
137 |
|
|
RIT it = clusterRecords.begin(), end = clusterRecords.end();
|
138 |
|
|
size_t clusters = 0;
|
139 |
|
|
while (it != end) {
|
140 |
|
|
RIT it2 = it;
|
141 |
|
|
uint32_t detid = it->detid();
|
142 |
|
|
|
143 |
|
|
// first isolate all clusters on the same detid
|
144 |
|
|
while ( (it2 != end) && (it2->detid() == detid)) { ++it2; }
|
145 |
|
|
// now [it, it2] bracket one detid
|
146 |
|
|
|
147 |
|
|
// then prepare to copy the clusters
|
148 |
|
|
typename edmNew::DetSetVector<ClusterType>::FastFiller filler(dsv, detid);
|
149 |
|
|
typename HitType::ClusterRef lastRef, newRef;
|
150 |
|
|
for ( ; it != it2; ++it) { // loop on the detid
|
151 |
|
|
// first check if we need to clone the hit
|
152 |
|
|
if (it->clusterRef() != lastRef) {
|
153 |
|
|
lastRef = it->clusterRef();
|
154 |
|
|
// clone cluster
|
155 |
|
|
filler.push_back( *lastRef );
|
156 |
|
|
// make new ref
|
157 |
|
|
newRef = typename HitType::ClusterRef( refprod, clusters++ );
|
158 |
|
|
}
|
159 |
|
|
// then fixup the reference
|
160 |
|
|
it->rekey( *selHits_, newRef );
|
161 |
|
|
|
162 |
|
|
} // end of the loop on a single detid
|
163 |
|
|
|
164 |
|
|
} // end of the loop on all clusters
|
165 |
|
|
|
166 |
|
|
clusterRecords.clear();
|
167 |
|
|
} // end of the function
|
168 |
|
|
|
169 |
|
|
//------------------------------------------------------------------
|
170 |
|
|
//! Put tracks, track extras and hits+clusters into the event.
|
171 |
|
|
//------------------------------------------------------------------
|
172 |
|
|
edm::OrphanHandle<reco::TrackCollection>
|
173 |
|
|
TrackCollectionStoreManager::
|
174 |
|
|
put( edm::Event & evt ) {
|
175 |
|
|
edm::OrphanHandle<reco::TrackCollection>
|
176 |
|
|
h = evt.put( selTracks_ );
|
177 |
|
|
evt.put( selTrackExtras_ );
|
178 |
|
|
evt.put( selHits_ );
|
179 |
|
|
evt.put( selStripClusters_ );
|
180 |
|
|
evt.put( selPixelClusters_ );
|
181 |
|
|
return h;
|
182 |
|
|
}
|
183 |
|
|
|
184 |
|
|
|
185 |
|
|
} // end of namespace helper
|
186 |
|
|
|