ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitEdm/Producers/src/TrackToTrackAssociator.cc
Revision: 1.3
Committed: Wed Jul 15 20:38:24 2009 UTC (15 years, 9 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_025c_branch2, Mit_025c_branch1, Mit_030, Mit_029c, Mit_029b, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_025c_branch0, Mit_028, Mit_027a, Mit_027, Mit_026, Mit_025e, Mit_025d, Mit_025c, Mit_025b, Mit_025a, Mit_025, Mit_025pre2, Mit_024b, Mit_025pre1, Mit_024a, Mit_024, Mit_023, Mit_022a, Mit_022, Mit_020d, TMit_020d, Mit_020c, Mit_021, Mit_021pre2, Mit_021pre1, Mit_020b, Mit_020a, Mit_020, Mit_020pre1, Mit_018, Mit_017, Mit_017pre3, Mit_017pre2, Mit_017pre1, V07-05-00, Mit_016, Mit_015b, Mit_015a, Mit_015, Mit_014e, Mit_014d, Mit_014c, Mit_014b, ConvRejection-10-06-09, Mit_014a, Mit_014, Mit_014pre3, Mit_014pre2, Mit_014pre1, Mit_013d, Mit_013c, Mit_013b, Mit_013a, Mit_013, Mit_013pre1, Mit_012i, Mit_012h, Mit_012g, Mit_012f, Mit_012e, Mit_012d, Mit_012c, Mit_012b, Mit_012a, Mit_012, Mit_011a, Mit_011, Mit_010a, Mit_010, HEAD
Branch point for: Mit_025c_branch
Changes since 1.2: +1 -2 lines
Log Message:
Cleanup

File Contents

# User Rev Content
1 loizides 1.3 // $Id: TrackToTrackAssociator.cc,v 1.2 2009/03/20 18:01:48 loizides Exp $
2 bendavid 1.1
3     #include "MitEdm/Producers/interface/TrackToTrackAssociator.h"
4     #include <TSystem.h>
5     #include <TError.h>
6     #include "DataFormats/TrackReco/interface/TrackFwd.h"
7     #include "MitEdm/DataFormats/interface/Types.h"
8    
9     using namespace std;
10     using namespace edm;
11     using namespace reco;
12     using namespace mitedm;
13    
14     //--------------------------------------------------------------------------------------------------
15     TrackToTrackAssociator::TrackToTrackAssociator(const ParameterSet& cfg) :
16     fromTracksName_(cfg.getUntrackedParameter<string>("fromTracks")),
17     toTracksName_(cfg.getUntrackedParameter<string>("toTracks"))
18     {
19     // Constructor: Register your base product.
20 loizides 1.2
21 bendavid 1.1 produces<mitedm::TrackAssociation>();
22     }
23    
24     //--------------------------------------------------------------------------------------------------
25     void TrackToTrackAssociator::PrintErrorAndExit(const char *msg) const
26     {
27     // Print error message and then exit.
28    
29     Error("PrintErrorAndExit", msg);
30     gSystem->Exit(1);
31     }
32    
33     //--------------------------------------------------------------------------------------------------
34     void TrackToTrackAssociator::produce(Event &evt, const EventSetup &setup)
35     {
36 loizides 1.2 // Produce the track association.
37 bendavid 1.1
38     Handle<View<reco::Track> > hFromTrackProduct;
39     GetProduct(fromTracksName_, hFromTrackProduct, evt);
40     const View<reco::Track> fromTracks = *(hFromTrackProduct.product());
41    
42     Handle<View<reco::Track> > hToTrackProduct;
43     GetProduct(toTracksName_, hToTrackProduct, evt);
44     const View<reco::Track> toTracks = *(hToTrackProduct.product());
45    
46     auto_ptr<mitedm::TrackAssociation> association(new mitedm::TrackAssociation());
47    
48 loizides 1.2 //fill all for each fromTrack, fill an association for all toTracks which share some hits,
49     //using as the association quality the ratio number of matched hits/number of valid hits
50     //on toTrack
51 bendavid 1.1 for (View<reco::Track>::const_iterator tFrom = fromTracks.begin();
52     tFrom != fromTracks.end(); ++tFrom) {
53    
54     reco::TrackBaseRef refFrom = fromTracks.refAt(tFrom-fromTracks.begin());
55    
56     for (View<reco::Track>::const_iterator tTo = toTracks.begin();
57     tTo != toTracks.end(); ++tTo) {
58    
59 loizides 1.2 uint nShared = 0;
60     for (uint i=0; i<tTo->recHitsSize(); ++i) {
61     if (tTo->recHit(i)->isValid()) {
62     bool matchedHit = false;
63     for (uint j=0; j<tFrom->recHitsSize() && !matchedHit; ++j) {
64     if ( tTo->recHit(i)->sharesInput(tFrom->recHit(j).get(), TrackingRecHit::some) ) {
65     nShared++;
66     matchedHit=true;
67 bendavid 1.1 }
68 loizides 1.2 }
69     }
70     }
71    
72     if (nShared>0) {
73     double rShared = (double)nShared/(double)(tTo->numberOfValidHits());
74     reco::TrackBaseRef refTo = toTracks.refAt(tTo-toTracks.begin());
75     std::pair<TrackBaseRef, double> assocWQuality(refTo,rShared);
76     association->insert(refFrom, assocWQuality);
77     }
78 bendavid 1.1 }
79     }
80     evt.put(association);
81     }
82    
83 loizides 1.2 // define this as a plug-in
84 bendavid 1.1 DEFINE_FWK_MODULE(TrackToTrackAssociator);