ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerTracks.cc
Revision: 1.16
Committed: Fri Aug 29 02:50:02 2008 UTC (16 years, 8 months ago) by loizides
Content type: text/plain
Branch: MAIN
Changes since 1.15: +8 -19 lines
Log Message:
Try to reuse as much of track filler as possible.

File Contents

# User Rev Content
1 loizides 1.16 // $Id: FillerTracks.cc,v 1.15 2008/08/28 22:21:01 loizides Exp $
2 loizides 1.1
3 bendavid 1.11 #include "MitProd/TreeFiller/interface/FillerTracks.h"
4 loizides 1.1 #include "FWCore/MessageLogger/interface/MessageLogger.h"
5     #include "DataFormats/Common/interface/Handle.h"
6 bendavid 1.14 #include "DataFormats/TrackReco/interface/HitPattern.h"
7 loizides 1.1 #include "DataFormats/TrackReco/interface/Track.h"
8     #include "DataFormats/TrackReco/interface/TrackFwd.h"
9     #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
10     #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
11     #include "DataFormats/RecoCandidate/interface/TrackAssociation.h"
12    
13     using namespace std;
14     using namespace edm;
15     using namespace mithep;
16    
17 loizides 1.4 //--------------------------------------------------------------------------------------------------
18 loizides 1.16 FillerTracks::FillerTracks(const ParameterSet &cfg, const char *name, bool active,
19     const char *edmName, const char *mitName) :
20 loizides 1.7 BaseFiller(cfg,name,active),
21 loizides 1.16 edmName_(Conf().getUntrackedParameter<string>("edmName",edmName)),
22     mitName_(Conf().getUntrackedParameter<string>("mitName",mitName)),
23 loizides 1.1 edmSimAssociationName_(Conf().getUntrackedParameter<string>("edmSimAssociationName","")),
24 loizides 1.13 simMapName_(Conf().getUntrackedParameter<string>("simMapName","SimMap")),
25     trackMapName_(Conf().getUntrackedParameter<string>("trackMapName",
26     Form("%sMapName",mitName_.c_str()))),
27     simMap_(0),
28 loizides 1.7 tracks_(new mithep::TrackArr(250)),
29 loizides 1.1 trackMap_(new mithep::TrackMap)
30     {
31 loizides 1.13 // Constructor.
32 loizides 1.16
33     InitLayerMap();
34 loizides 1.1 }
35    
36 loizides 1.4 //--------------------------------------------------------------------------------------------------
37 loizides 1.1 FillerTracks::~FillerTracks()
38     {
39 loizides 1.13 // Destructor.
40    
41 loizides 1.5 delete tracks_;
42 loizides 1.1 delete trackMap_;
43     }
44    
45 loizides 1.4 //--------------------------------------------------------------------------------------------------
46 loizides 1.1 void FillerTracks::BookDataBlock(TreeWriter &tws)
47     {
48 loizides 1.13 // Add tracks branch to tree, publish and get our objects.
49    
50 loizides 1.1 tws.AddBranch(mitName_.c_str(),&tracks_);
51 loizides 1.12
52 loizides 1.13 simMap_ = OS()->get<SimParticleMap>(simMapName_.c_str());
53     OS()->add<TrackMap>(trackMap_,trackMapName_.c_str());
54     OS()->add<TrackArr>(tracks_,mitName_.c_str());
55 loizides 1.1 }
56    
57 loizides 1.4 //--------------------------------------------------------------------------------------------------
58 loizides 1.1 void FillerTracks::FillDataBlock(const edm::Event &event,
59     const edm::EventSetup &setup)
60     {
61 loizides 1.13 // Fill tracks from edm collection into our collection.
62    
63 paus 1.10 tracks_ ->Reset();
64 loizides 1.1 trackMap_->Reset();
65 loizides 1.6
66     Handle<reco::TrackCollection> hTrackProduct;
67     GetProduct(edmName_, hTrackProduct, event);
68 loizides 1.1
69 loizides 1.6 trackMap_->SetEdmProductId(hTrackProduct.id().id());
70     const reco::TrackCollection inTracks = *(hTrackProduct.product());
71 loizides 1.1
72 paus 1.10 // for MC SimParticle association (reco->sim mappings)
73 loizides 1.1 reco::RecoToSimCollection simAssociation;
74     if (simMap_ && !edmSimAssociationName_.empty()) {
75     Handle<reco::RecoToSimCollection> simAssociationProduct;
76 loizides 1.6 GetProduct(edmSimAssociationName_, simAssociationProduct, event);
77 loizides 1.1 simAssociation = *(simAssociationProduct.product());
78     }
79    
80 paus 1.10 // loop through all tracks and fill the information
81     for (reco::TrackCollection::const_iterator it = inTracks.begin();
82     it != inTracks.end(); ++it) {
83 loizides 1.5 mithep::Track *outTrack = tracks_->Allocate();
84 paus 1.10 // create track and set the core parameters
85 bendavid 1.14 new (outTrack) mithep::Track(it->qoverp(),it->lambda(),it->phi(),it->dxy(),it->dsz());
86     outTrack->SetErrors(it->qoverpError(),it->lambdaError(),it->phiError(),it->dxyError(),it->dszError());
87 loizides 1.1
88 bendavid 1.14 //Fill track quality information
89     outTrack->SetChi2(it->chi2());
90 loizides 1.16 outTrack->SetNdof(static_cast<Int_t>(it->ndof()));
91 bendavid 1.14
92     //Fill hit layer map
93     const reco::HitPattern &hits = it->hitPattern();
94     for (Int_t i=0; i<hits.numberOfHits(); i++) {
95     uint32_t hit = hits.getHitPattern(i);
96     if ( hits.validHitFilter(hit) )
97     if ( hits.trackerHitFilter(hit) )
98     outTrack->SetHit(layerMap_[hit]);
99    
100     // if ( hits.muonDTHitFilter(hit) )
101     // printf("Muon DT Layer = %i\n", hits.getLayer(hit));
102     // if ( hits.muonCSCHitFilter(hit) )
103     // printf("Muon CSC Layer = %i\n", hits.getLayer(hit));
104     // if ( hits.muonRPCHitFilter(hit) )
105     // printf("Muon RPC Layer = %i\n", hits.getLayer(hit));
106     }
107    
108    
109 loizides 1.1 // add reference between mithep and edm object
110 paus 1.10 reco::TrackRef theRef(hTrackProduct, it - inTracks.begin());
111 loizides 1.1 trackMap_->Add(theRef, outTrack);
112    
113     if (simMap_ && !edmSimAssociationName_.empty()) {
114     reco::TrackBaseRef theBaseRef(theRef);
115     vector<pair<TrackingParticleRef, double> > simRefs;
116 paus 1.10 Bool_t noSimParticle = 0;
117 loizides 1.1 try {
118     simRefs = simAssociation[theBaseRef]; //try to get the sim references if existing
119     }
120 loizides 1.5 catch (edm::Exception &ex) {
121 paus 1.10 noSimParticle = 1;
122 loizides 1.1 }
123 paus 1.10
124 loizides 1.1 if (!noSimParticle) { //loop through sim match candidates
125     for (vector<pair<TrackingParticleRef, double> >::const_iterator simRefPair=simRefs.begin();
126     simRefPair != simRefs.end(); ++simRefPair)
127    
128 paus 1.10 if (simRefPair->second > 0.5) // require more than 50% shared hits between reco and sim
129 bendavid 1.9 outTrack->SetMCPart(simMap_->GetMit(simRefPair->first)); //add reco->sim reference
130 loizides 1.1 }
131     }
132     }
133     tracks_->Trim();
134     }
135 bendavid 1.14
136     //--------------------------------------------------------------------------------------------------
137     void FillerTracks::InitLayerMap()
138     {
139     // Initialize mapping between hit layer format in reco::HitPattern and the one used in
140 loizides 1.15 // mithep::Track. Note in 21x stereo layers are treated separatelely.
141 bendavid 1.14
142 loizides 1.15 layerMap_[1160] = mithep::Track::PXB1;
143     layerMap_[1168] = mithep::Track::PXB2;
144     layerMap_[1176] = mithep::Track::PXB3;
145     layerMap_[1288] = mithep::Track::PXF1;
146     layerMap_[1296] = mithep::Track::PXF2;
147     layerMap_[1416] = mithep::Track::TIB1;
148     layerMap_[1420] = mithep::Track::TIB1S;
149     layerMap_[1424] = mithep::Track::TIB2;
150     layerMap_[1428] = mithep::Track::TIB2S;
151     layerMap_[1432] = mithep::Track::TIB3;
152     layerMap_[1440] = mithep::Track::TIB4;
153     layerMap_[1544] = mithep::Track::TID1;
154     layerMap_[1548] = mithep::Track::TID1S;
155     layerMap_[1552] = mithep::Track::TID2;
156     layerMap_[1556] = mithep::Track::TID2S;
157     layerMap_[1560] = mithep::Track::TID3;
158     layerMap_[1672] = mithep::Track::TOB1;
159     layerMap_[1676] = mithep::Track::TOB1S;
160     layerMap_[1680] = mithep::Track::TOB2;
161     layerMap_[1684] = mithep::Track::TOB2S;
162     layerMap_[1688] = mithep::Track::TOB3;
163     layerMap_[1696] = mithep::Track::TOB4;
164     layerMap_[1704] = mithep::Track::TOB5;
165     layerMap_[1712] = mithep::Track::TOB6;
166     layerMap_[1800] = mithep::Track::TEC1;
167     layerMap_[1804] = mithep::Track::TEC1S;
168     layerMap_[1808] = mithep::Track::TEC2;
169     layerMap_[1812] = mithep::Track::TEC2S;
170     layerMap_[1816] = mithep::Track::TEC3;
171     layerMap_[1824] = mithep::Track::TEC4;
172     layerMap_[1832] = mithep::Track::TEC5;
173     layerMap_[1836] = mithep::Track::TEC5S;
174     layerMap_[1840] = mithep::Track::TEC6;
175     layerMap_[1848] = mithep::Track::TEC7;
176     layerMap_[1856] = mithep::Track::TEC8;
177     layerMap_[1864] = mithep::Track::TEC9;
178 bendavid 1.14 }