1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
|
#include "CommonTools/RecoAlgos/interface/TrackSelector.h"
using namespace reco;
namespace helper {
TrackCollectionStoreManager::TrackCollectionStoreManager(const edm::Handle<reco::TrackCollection>&)
: selTracks_(new reco::TrackCollection),
selTrackExtras_(new reco::TrackExtraCollection),
selHits_(new TrackingRecHitCollection),
selStripClusters_(new edmNew::DetSetVector<SiStripCluster>),
selPixelClusters_(new edmNew::DetSetVector<SiPixelCluster>),
selPhase2OTClusters_(new edmNew::DetSetVector<Phase2TrackerCluster1D>),
rTracks_(),
rTrackExtras_(),
rHits_(),
clusterStorer_(),
idx_(0),
hidx_(0),
cloneClusters_(true) {}
//------------------------------------------------------------------
//! Process a single track.
//------------------------------------------------------------------
void TrackCollectionStoreManager::processTrack(const Track& trk) {
selTracks_->push_back(Track(trk));
selTracks_->back().setExtra(TrackExtraRef(rTrackExtras_, idx_++));
selTrackExtras_->push_back(TrackExtra(trk.outerPosition(),
trk.outerMomentum(),
trk.outerOk(),
trk.innerPosition(),
trk.innerMomentum(),
trk.innerOk(),
trk.outerStateCovariance(),
trk.outerDetId(),
trk.innerStateCovariance(),
trk.innerDetId(),
trk.seedDirection()));
TrackExtra& tx = selTrackExtras_->back();
auto const firstHitIndex = hidx_;
unsigned int nHitsAdded = 0;
for (trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++hit, ++hidx_) {
selHits_->push_back((*hit)->clone());
TrackingRecHit* newHit = &(selHits_->back());
++nHitsAdded;
//--- Skip the rest for this hit if we don't want to clone the cluster.
//--- The copy constructer in the rec hit will copy the link properly.
if (cloneClusters() && newHit->isValid() && ((*hit)->geographicalId().det() == DetId::Tracker)) {
clusterStorer_.addCluster(*selHits_, hidx_);
}
} // end of for loop over tracking rec hits on this track
tx.setHits(rHits_, firstHitIndex, nHitsAdded);
} // end of track, and function
//------------------------------------------------------------------
//! Put tracks, track extras and hits+clusters into the event.
//------------------------------------------------------------------
edm::OrphanHandle<reco::TrackCollection> TrackCollectionStoreManager::put(edm::Event& evt) {
edm::OrphanHandle<reco::TrackCollection> h = evt.put(std::move(selTracks_));
evt.put(std::move(selTrackExtras_));
evt.put(std::move(selHits_));
evt.put(std::move(selStripClusters_));
evt.put(std::move(selPixelClusters_));
evt.put(std::move(selPhase2OTClusters_));
return h;
}
} // end of namespace helper
|