ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/L1RpcTriggerAnalysis/src/SynchroSelectorTrack.cc
Revision: 1.2
Committed: Thu May 17 20:39:45 2012 UTC (12 years, 11 months ago) by konec
Content type: text/plain
Branch: MAIN
CVS Tags: Artur_11_07_2013_B, Artur_11_07_2013_A, Artur_11_07_2013, Artur_28_06_2013, Mikolaj_cmssw533, Mikolaj_cmssw52x, HEAD
Changes since 1.1: +1 -1 lines
Error occurred while calculating annotation data.
Log Message:
*** empty log message ***

File Contents

# Content
1 #include "UserCode/L1RpcTriggerAnalysis/interface/SynchroSelectorTrack.h"
2
3
4 #include "FWCore/Framework/interface/Event.h"
5 #include "FWCore/Framework/interface/EventSetup.h"
6 #include "FWCore/Framework/interface/ESHandle.h"
7 #include "FWCore/Utilities/interface/InputTag.h"
8
9 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuRegionalCand.h"
10 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h"
11 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTExtendedCand.h"
12 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTReadoutCollection.h"
13
14 #include "DataFormats/TrackReco/interface/TrackFwd.h"
15 #include "DataFormats/TrackReco/interface/Track.h"
16
17 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
18 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
19
20 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
21 #include "MagneticField/Engine/interface/MagneticField.h"
22 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
23 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
24 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
25
26 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
27 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
28
29
30 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
31 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
32
33 #include "UserCode/L1RpcTriggerAnalysis/interface/TrackToL1ObjMatcher.h"
34 #include "TObjArray.h"
35 #include "TH1F.h"
36
37 #include "TH1.h"
38 #include <cmath>
39 using namespace edm;
40 using namespace std;
41
42 SynchroSelectorTrack::SynchroSelectorTrack(const edm::ParameterSet& cfg, TObjArray& histos)
43 : SynchroSelector(cfg)
44 {
45 hDxy = new TH1F("hDxy_SSTrack","hDxy_SSTrack",100.,0.,1.); histos.Add(hDxy);
46 SynchroSelector::hPullX = new TH1F("hPullX_SSTrack","hPullX_SSTrack",100,-10.,10.); histos.Add(hPullX);
47 SynchroSelector::hDistX = new TH1F("hDistX_SSTrack","hDistX_SSTrack",100,-100.,100.); histos.Add(hDistX);
48 hDeltaPhi = new TH1F("hDeltaPhi_SSTrack","hDeltaPhi_SSTrack",100, 0., 1.); histos.Add(hDeltaPhi);
49 hDeltaEta = new TH1F("hDeltaEta_SSTrack","hDeltaEta_SSTrack",100, 0., 1.); histos.Add(hDeltaEta);
50 hNum = new TH1F("hNumTracks","hNumTracks",50,0.,250.); histos.Add(hNum);
51 }
52
53
54
55 bool SynchroSelectorTrack::takeIt(const RPCDetId & det, const edm::Event&ev, const edm::EventSetup& es)
56 {
57 bool result = false;
58 float theMaxTIP = theConfig.getParameter<double>("maxTIP");
59 float theMinPt = theConfig.getParameter<double>("minPt");
60
61 edm::Handle<reco::TrackCollection> trackCollection;
62 ev.getByLabel(InputTag( theConfig.getParameter<std::string>("collection")),trackCollection);
63
64 edm::ESHandle<GlobalTrackingGeometry> globalGeometry;
65 es.get<GlobalTrackingGeometryRecord>().get(globalGeometry);
66
67 edm::ESHandle<MagneticField> magField;
68 es.get<IdealMagneticFieldRecord>().get(magField);
69
70 math::XYZPoint reference(0.,0.,0.);
71 if (theConfig.exists("beamSpot")) {
72 edm::InputTag beamSpot = theConfig.getParameter<edm::InputTag>("beamSpot");
73 edm::Handle<reco::BeamSpot> bsHandle;
74 ev.getByLabel( beamSpot, bsHandle);
75 if (bsHandle.isValid()) reference = math::XYZPoint(bsHandle->x0(), bsHandle->y0(), bsHandle->z0());
76 }
77
78 reco::TrackCollection tracks = *(trackCollection.product());
79 if (hNum) hNum->Fill(tracks.size());
80
81 mindetaForStudy = 100.;
82 mindphiForStudy = 100.;
83
84 typedef reco::TrackCollection::const_iterator IT;
85 for (IT it = tracks.begin(); it !=tracks.end(); ++it) {
86 const reco::Track & track = *it;
87 if (track.pt() < theMinPt ) continue;
88 double dxy = track.dxy(reference);
89 if (fabs(dxy) > theMaxTIP) continue;
90
91 TrajectoryStateOnSurface aTSOS = trajectoryStateTransform::outerStateOnSurface(track, *globalGeometry, magField.product());
92 if (!checkL1RpcMatching(aTSOS, ev,es) ) continue;
93 if (!checkRpcDetMatching(aTSOS,det,ev,es)) continue;
94 if (theConfig.getParameter<bool>("checkUniqueRecHitMatching") && !checkUniqueRecHitMatching(aTSOS,det,ev,es)) continue;
95
96 if (hDxy) hDxy->Fill(fabs(dxy));
97 result = true;
98 }
99
100 if (mindetaForStudy < 99.) hDeltaEta->Fill(mindetaForStudy);
101 if (mindphiForStudy < 99.) hDeltaPhi->Fill(mindphiForStudy);
102
103 return result;
104 }
105
106 bool SynchroSelectorTrack::checkL1RpcMatching( const TrajectoryStateOnSurface & tsos, const edm::Event&ev, const edm::EventSetup& es)
107 {
108 edm::Handle<L1MuGMTReadoutCollection> pCollection;
109 InputTag l1MuReadout = theConfig.getParameter<edm::InputTag>("l1MuReadout");
110 ev.getByLabel(l1MuReadout,pCollection);
111 L1MuGMTReadoutCollection const* gmtrc = pCollection.product();
112 if (!gmtrc) return false;
113
114 TrackToL1ObjMatcher matcher(theConfig.getParameter<edm::ParameterSet>("rpcMatcherPSet"));
115
116 vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
117 vector<L1MuGMTReadoutRecord>::const_iterator RRItr;
118 typedef vector<L1MuRegionalCand>::const_iterator ITC;
119 for( RRItr = gmt_records.begin() ; RRItr != gmt_records.end() ; RRItr++ ) {
120 for (int i=1; i<=2; ++i) {
121 vector<L1MuRegionalCand> Cands = (i==1) ? RRItr->getBrlRPCCands() : RRItr->getFwdRPCCands();
122 for(ITC it = Cands.begin() ; it != Cands.end() ; ++it ) {
123 if (it->empty()) continue;
124 //propagate and check matching for candidate
125 float rpcEta = it->etaValue();
126 float rpcPhi = it->phiValue();
127 bool matched = matcher(rpcEta, rpcPhi, tsos, ev, es);
128 TrackToL1ObjMatcher::LastResult matchingResult = matcher.lastResult();
129 if (matchingResult.isValid) {
130 if (matchingResult.deltaEta < mindetaForStudy) mindetaForStudy = matchingResult.deltaEta;
131 if (matchingResult.deltaPhi < mindphiForStudy) mindphiForStudy = matchingResult.deltaPhi;
132 }
133 if (matched) return true;
134 }
135 }
136 }
137 return false;
138 }
139
140
141
142