ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerMuons.cc
Revision: 1.44
Committed: Sat May 5 16:49:59 2012 UTC (13 years ago) by paus
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_032, Mit_031, Mit_030, Mit_029c, Mit_029b, Mit_030_pre1, Mit_029a, Mit_029, Mit_029_pre1, Mit_028a, Mit_028, Mit_027a, Mit_027, HEAD
Changes since 1.43: +62 -56 lines
Log Message:
Version 027 - complete version for ICHEP 2012.

File Contents

# User Rev Content
1 paus 1.44 // $Id: FillerMuons.cc,v 1.43 2012/03/29 23:41:59 paus Exp $
2 loizides 1.1
3     #include "MitProd/TreeFiller/interface/FillerMuons.h"
4     #include "DataFormats/MuonReco/interface/Muon.h"
5 sixie 1.26 #include "DataFormats/MuonReco/interface/MuonQuality.h"
6 loizides 1.1 #include "DataFormats/MuonReco/interface/MuonFwd.h"
7 sixie 1.26 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
8 bendavid 1.15 #include "DataFormats/Common/interface/RefToPtr.h"
9 bendavid 1.27 #include "DataFormats/VertexReco/interface/VertexFwd.h"
10     #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
11     #include "TrackingTools/TransientTrack/plugins/TransientTrackBuilderESProducer.h"
12     #include "TrackingTools/IPTools/interface/IPTools.h"
13 paus 1.44 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexTrackCompatibilityEstimator.h"
14 loizides 1.22 #include "MitAna/DataTree/interface/Names.h"
15     #include "MitAna/DataTree/interface/MuonCol.h"
16 loizides 1.1 #include "MitAna/DataTree/interface/Track.h"
17 loizides 1.22 #include "MitProd/ObjectService/interface/ObjectService.h"
18 bendavid 1.37 #include "RecoVertex/VertexTools/interface/LinearizedTrackStateFactory.h"
19     #include "RecoVertex/VertexTools/interface/VertexTrackFactory.h"
20     #include "RecoVertex/VertexPrimitives/interface/VertexTrack.h"
21     #include "RecoVertex/VertexPrimitives/interface/CachingVertex.h"
22     #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexUpdator.h"
23 ksung 1.40 #include "MitEdm/Tools/interface/VertexReProducer.h"
24 loizides 1.1
25     using namespace std;
26     using namespace edm;
27     using namespace mithep;
28    
29 loizides 1.4 //--------------------------------------------------------------------------------------------------
30 loizides 1.18 FillerMuons::FillerMuons(const edm::ParameterSet &cfg, const char *name, bool active) :
31 paus 1.43 BaseFiller (cfg,name,active),
32     edmName_ (Conf().getUntrackedParameter<string>("edmName","muons")),
33     expectedHitsName_ (Conf().getUntrackedParameter<string>("expectedHitsName","")),
34     mitName_ (Conf().getUntrackedParameter<string>("mitName",Names::gkMuonBrn)),
35     globalTrackMapName_ (Conf().getUntrackedParameter<string>("globalTrackMapName","")),
36     staTrackMapName_ (Conf().getUntrackedParameter<string>("staTrackMapName","")),
37     staVtxTrackMapName_ (Conf().getUntrackedParameter<string>("staVtxTrackMapName","")),
38     trackerTrackMapName_ (Conf().getUntrackedParameter<string>("trackerTrackMapName","")),
39     muonMapName_ (Conf().getUntrackedParameter<string>("muonMapName","")),
40     pvEdmName_ (Conf().getUntrackedParameter<string>("pvEdmName","offlinePrimaryVertices")),
41     pvBSEdmName_ (Conf().getUntrackedParameter<string>("pvBSEdmName","offlinePrimaryVerticesWithBS")),
42     fitUnbiasedVertex_ (Conf().getUntrackedParameter<bool>("fitUnbiasedVertex",true)),
43     globalTrackMap_ (0),
44     standaloneTrackMap_ (0),
45 bendavid 1.19 standaloneVtxTrackMap_(0),
46 paus 1.43 trackerTrackMap_ (0),
47     muonMap_ (new mithep::MuonMap),
48     muons_ (new mithep::MuonArr(16))
49 loizides 1.1 {
50     // Constructor.
51     }
52    
53 loizides 1.4 //--------------------------------------------------------------------------------------------------
54 paus 1.44 FillerMuons::~FillerMuons()
55 loizides 1.20 {
56     // Destructor.
57    
58 loizides 1.5 delete muons_;
59 bendavid 1.19 delete muonMap_;
60 loizides 1.1 }
61    
62 loizides 1.4 //--------------------------------------------------------------------------------------------------
63 paus 1.44 void FillerMuons::BookDataBlock(TreeWriter &tws)
64     {
65 loizides 1.9 // Add muons branch to tree and get pointers to maps.
66 loizides 1.1
67 loizides 1.20 tws.AddBranch(mitName_,&muons_);
68     OS()->add<mithep::MuonArr>(muons_,mitName_);
69 bendavid 1.19
70 loizides 1.20 if (!globalTrackMapName_.empty()) {
71     globalTrackMap_ = OS()->get<TrackMap>(globalTrackMapName_);
72     if (globalTrackMap_)
73     AddBranchDep(mitName_,globalTrackMap_->GetBrName());
74     }
75     if (!staTrackMapName_.empty()) {
76     standaloneTrackMap_ = OS()->get<TrackMap>(staTrackMapName_);
77     if (standaloneTrackMap_)
78     AddBranchDep(mitName_,standaloneTrackMap_->GetBrName());
79     }
80     if (!staVtxTrackMapName_.empty()) {
81     standaloneVtxTrackMap_ = OS()->get<TrackMap>(staVtxTrackMapName_);
82     if (standaloneVtxTrackMap_)
83     AddBranchDep(mitName_,standaloneVtxTrackMap_->GetBrName());
84     }
85     if (!trackerTrackMapName_.empty()) {
86     trackerTrackMap_ = OS()->get<TrackMap>(trackerTrackMapName_);
87     if (trackerTrackMap_)
88     AddBranchDep(mitName_,trackerTrackMap_->GetBrName());
89     }
90     if (!muonMapName_.empty()) {
91     muonMap_->SetBrName(mitName_);
92     OS()->add<MuonMap>(muonMap_,muonMapName_);
93     }
94 loizides 1.1 }
95    
96 loizides 1.4 //--------------------------------------------------------------------------------------------------
97 paus 1.44 void FillerMuons::FillDataBlock(const edm::Event &event,
98 loizides 1.4 const edm::EventSetup &setup)
99 loizides 1.1 {
100 paus 1.44 // Fill muon information.
101 loizides 1.1
102 paus 1.43 muons_ ->Delete();
103 bendavid 1.19 muonMap_->Reset();
104 paus 1.44
105 loizides 1.6 Handle<reco::MuonCollection> hMuonProduct;
106 paus 1.44 GetProduct(edmName_, hMuonProduct, event);
107     const reco::MuonCollection inMuons = *(hMuonProduct.product());
108 loizides 1.1
109 bendavid 1.27 edm::Handle<reco::VertexCollection> hVertex;
110     event.getByLabel(pvEdmName_, hVertex);
111     const reco::VertexCollection *pvCol = hVertex.product();
112    
113     edm::Handle<reco::VertexCollection> hVertexBS;
114     event.getByLabel(pvBSEdmName_, hVertexBS);
115     const reco::VertexCollection *pvBSCol = hVertexBS.product();
116 paus 1.44
117 bendavid 1.27 edm::ESHandle<TransientTrackBuilder> hTransientTrackBuilder;
118     setup.get<TransientTrackRecord>().get("TransientTrackBuilder",hTransientTrackBuilder);
119     const TransientTrackBuilder *transientTrackBuilder = hTransientTrackBuilder.product();
120    
121 bendavid 1.28 KalmanVertexTrackCompatibilityEstimator<5> kalmanEstimator;
122    
123 bendavid 1.37 LinearizedTrackStateFactory lTrackFactory;
124     VertexTrackFactory<5> vTrackFactory;
125     KalmanVertexUpdator<5> updator;
126    
127 paus 1.44 for (reco::MuonCollection::const_iterator iM = inMuons.begin(); iM != inMuons.end(); ++iM) {
128 bendavid 1.25 mithep::Muon* outMuon = muons_->AddNew();
129    
130 paus 1.43 outMuon->SetPtEtaPhi (iM->pt(),iM->eta(),iM->phi());
131     outMuon->SetCharge (iM->charge());
132     outMuon->SetIsoR03SumPt (iM->isolationR03().sumPt);
133     outMuon->SetIsoR03EmEt (iM->isolationR03().emEt);
134     outMuon->SetIsoR03HadEt (iM->isolationR03().hadEt);
135     outMuon->SetIsoR03HoEt (iM->isolationR03().hoEt);
136     outMuon->SetIsoR03NTracks (iM->isolationR03().nTracks);
137     outMuon->SetIsoR03NJets (iM->isolationR03().nJets);
138     outMuon->SetIsoR05SumPt (iM->isolationR05().sumPt);
139     outMuon->SetIsoR05EmEt (iM->isolationR05().emEt);
140     outMuon->SetIsoR05HadEt (iM->isolationR05().hadEt);
141     outMuon->SetIsoR05HoEt (iM->isolationR05().hoEt);
142     outMuon->SetIsoR05NTracks (iM->isolationR05().nTracks);
143     outMuon->SetIsoR05NJets (iM->isolationR05().nJets);
144     outMuon->SetEmEnergy (iM->calEnergy().em);
145     outMuon->SetHadEnergy (iM->calEnergy().had);
146     outMuon->SetHoEnergy (iM->calEnergy().ho);
147     outMuon->SetEmS9Energy (iM->calEnergy().emS9);
148     outMuon->SetHadS9Energy (iM->calEnergy().hadS9);
149     outMuon->SetHoS9Energy (iM->calEnergy().hoS9);
150     outMuon->SetIsGlobalMuon (iM->isGlobalMuon());
151     outMuon->SetIsTrackerMuon (iM->isTrackerMuon());
152 bendavid 1.17 outMuon->SetIsStandaloneMuon(iM->isStandAloneMuon());
153 paus 1.44 outMuon->SetIsPFMuon (iM->isPFMuon());
154 paus 1.43 outMuon->SetIsCaloMuon (iM->isCaloMuon());
155 bendavid 1.41 //kink algorithm
156 paus 1.43 outMuon->SetTrkKink (iM->combinedQuality().trkKink);
157     outMuon->SetGlbKink (iM->combinedQuality().glbKink);
158 sixie 1.26 //save results from the muon selector in the MuonQuality bitmask
159     outMuon->Quality().SetQuality(MuonQuality::All);
160     if (muon::isGoodMuon(*iM,muon::AllGlobalMuons))
161     outMuon->Quality().SetQuality(MuonQuality::AllGlobalMuons);
162     if (muon::isGoodMuon(*iM,muon::AllStandAloneMuons))
163     outMuon->Quality().SetQuality(MuonQuality::AllStandAloneMuons);
164     if (muon::isGoodMuon(*iM,muon::AllTrackerMuons))
165     outMuon->Quality().SetQuality(MuonQuality::AllTrackerMuons);
166     if (muon::isGoodMuon(*iM,muon::TrackerMuonArbitrated))
167     outMuon->Quality().SetQuality(MuonQuality::TrackerMuonArbitrated);
168     if (muon::isGoodMuon(*iM,muon::AllArbitrated))
169     outMuon->Quality().SetQuality(MuonQuality::AllArbitrated);
170     if (muon::isGoodMuon(*iM,muon::GlobalMuonPromptTight))
171     outMuon->Quality().SetQuality(MuonQuality::GlobalMuonPromptTight);
172     if (muon::isGoodMuon(*iM,muon::TMLastStationLoose))
173     outMuon->Quality().SetQuality(MuonQuality::TMLastStationLoose);
174     if (muon::isGoodMuon(*iM,muon::TMLastStationTight))
175     outMuon->Quality().SetQuality(MuonQuality::TMLastStationTight);
176     if (muon::isGoodMuon(*iM,muon::TM2DCompatibilityLoose))
177     outMuon->Quality().SetQuality(MuonQuality::TM2DCompatibilityLoose);
178     if (muon::isGoodMuon(*iM,muon::TM2DCompatibilityTight))
179     outMuon->Quality().SetQuality(MuonQuality::TM2DCompatibilityTight);
180     if (muon::isGoodMuon(*iM,muon::TMOneStationLoose))
181     outMuon->Quality().SetQuality(MuonQuality::TMOneStationLoose);
182     if (muon::isGoodMuon(*iM,muon::TMOneStationTight))
183     outMuon->Quality().SetQuality(MuonQuality::TMOneStationTight);
184     if (muon::isGoodMuon(*iM,muon::TMLastStationOptimizedLowPtLoose))
185     outMuon->Quality().SetQuality(MuonQuality::TMLastStationOptimizedLowPtLoose);
186     if (muon::isGoodMuon(*iM,muon::TMLastStationOptimizedLowPtTight))
187     outMuon->Quality().SetQuality(MuonQuality::TMLastStationOptimizedLowPtTight);
188     if (muon::isGoodMuon(*iM,muon::GMTkChiCompatibility))
189     outMuon->Quality().SetQuality(MuonQuality::GMTkChiCompatibility);
190     if (muon::isGoodMuon(*iM,muon::GMStaChiCompatibility))
191     outMuon->Quality().SetQuality(MuonQuality::GMStaChiCompatibility);
192     if (muon::isGoodMuon(*iM,muon::GMTkKinkTight))
193     outMuon->Quality().SetQuality(MuonQuality::GMTkKinkTight);
194     if (muon::isGoodMuon(*iM,muon::TMLastStationAngLoose))
195     outMuon->Quality().SetQuality(MuonQuality::TMLastStationAngLoose);
196     if (muon::isGoodMuon(*iM,muon::TMLastStationAngTight))
197     outMuon->Quality().SetQuality(MuonQuality::TMLastStationAngTight);
198     if (muon::isGoodMuon(*iM,muon::TMOneStationAngLoose))
199     outMuon->Quality().SetQuality(MuonQuality::TMOneStationAngLoose);
200     if (muon::isGoodMuon(*iM,muon::TMOneStationAngTight))
201     outMuon->Quality().SetQuality(MuonQuality::TMOneStationAngTight);
202     if (muon::isGoodMuon(*iM,muon::TMLastStationOptimizedBarrelLowPtLoose))
203     outMuon->Quality().SetQuality(MuonQuality::TMLastStationOptimizedBarrelLowPtLoose);
204     if (muon::isGoodMuon(*iM,muon::TMLastStationOptimizedBarrelLowPtTight))
205     outMuon->Quality().SetQuality(MuonQuality::TMLastStationOptimizedBarrelLowPtTight);
206    
207 bendavid 1.30 if (globalTrackMap_ && iM->combinedMuon().isNonnull()) {
208 paus 1.44 outMuon->SetGlobalTrk (globalTrackMap_->GetMit(refToPtr(iM->combinedMuon())));
209 paus 1.43 outMuon->SetNValidHits(iM->globalTrack()->hitPattern().numberOfValidMuonHits());
210 bendavid 1.30 }
211 paus 1.44 if (standaloneTrackMap_ && standaloneVtxTrackMap_ && iM->standAloneMuon().isNonnull()) {
212 loizides 1.1 Int_t refProductId = iM->standAloneMuon().id().id();
213     if ( refProductId == standaloneVtxTrackMap_->GetEdmProductId())
214 paus 1.44 outMuon->SetStandaloneTrk(standaloneVtxTrackMap_->GetMit(refToPtr(iM->standAloneMuon())));
215 loizides 1.1 else if ( refProductId == standaloneTrackMap_->GetEdmProductId())
216 paus 1.44 outMuon->SetStandaloneTrk(standaloneTrackMap_->GetMit(refToPtr(iM->standAloneMuon())));
217 loizides 1.1 else throw edm::Exception(edm::errors::Configuration, "FillerMuons:FillDataBlock()\n")
218 loizides 1.14 << "Error! Track reference in unmapped collection " << edmName_ << endl;
219 loizides 1.1 }
220 bendavid 1.30 if (trackerTrackMap_ && iM->track().isNonnull()) {
221 bendavid 1.15 outMuon->SetTrackerTrk(trackerTrackMap_->GetMit(refToPtr(iM->track())));
222 bendavid 1.30 }
223 loizides 1.14
224 bendavid 1.27 //compute impact parameter with respect to PV
225     if (iM->track().isNonnull()) {
226 paus 1.44 const reco::TransientTrack &tt = transientTrackBuilder->build(iM->track());
227 bendavid 1.27
228 paus 1.43 reco::Vertex thevtx = pvCol ->at(0);
229     reco::Vertex thevtxbs = pvBSCol->at(0);
230 bendavid 1.37
231 paus 1.43 reco::Vertex thevtxub = pvCol ->at(0);
232 bendavid 1.41 reco::Vertex thevtxubbs = pvBSCol->at(0);
233 bendavid 1.37
234 ksung 1.40 reco::TrackCollection newTkCollection;
235     bool foundMatch = false;
236     for(reco::Vertex::trackRef_iterator itk = thevtx.tracks_begin(); itk!=thevtx.tracks_end(); itk++) {
237 paus 1.44 if(itk->get() == &*(iM->innerTrack())) {
238 ksung 1.40 foundMatch = true;
239     continue;
240     }
241     newTkCollection.push_back(*itk->get());
242     }
243    
244 paus 1.44 if (foundMatch && fitUnbiasedVertex_) {
245 ksung 1.40 edm::Handle<reco::BeamSpot> bs;
246     event.getByLabel(edm::InputTag("offlineBeamSpot"),bs);
247 paus 1.44
248 ksung 1.40 VertexReProducer revertex(hVertex,event);
249     edm::Handle<reco::BeamSpot> pvbeamspot;
250     event.getByLabel(revertex.inputBeamSpot(),pvbeamspot);
251     vector<TransientVertex> pvs = revertex.makeVertices(newTkCollection,*pvbeamspot,setup);
252 bendavid 1.41 if(pvs.size()>0) {
253     thevtxub = pvs.front(); // take the first in the list
254     }
255 paus 1.44
256 ksung 1.40 VertexReProducer revertexbs(hVertexBS,event);
257     edm::Handle<reco::BeamSpot> pvbsbeamspot;
258     event.getByLabel(revertexbs.inputBeamSpot(),pvbsbeamspot);
259     vector<TransientVertex> pvbss = revertexbs.makeVertices(newTkCollection,*pvbsbeamspot,setup);
260 bendavid 1.41 if(pvbss.size()>0) {
261     thevtxubbs = pvbss.front(); // take the first in the list
262     }
263 ksung 1.40 }
264    
265    
266 bendavid 1.37 //preserve sign of transverse impact parameter (cross-product definition from track, not lifetime-signing)
267     const double thesign = ( (-iM->track()->dxy(thevtx.position())) >=0 ) ? 1. : -1.;
268     const double thesignbs = ( (-iM->track()->dxy(thevtxbs.position())) >=0 ) ? 1. : -1.;
269    
270     const std::pair<bool,Measurement1D> &d0pv = IPTools::absoluteTransverseImpactParameter(tt,thevtx);
271 bendavid 1.27 if (d0pv.first) {
272 bendavid 1.37 outMuon->SetD0PV(thesign*d0pv.second.value());
273 bendavid 1.27 outMuon->SetD0PVErr(d0pv.second.error());
274     }
275 bendavid 1.28 else {
276     outMuon->SetD0PV(-99.0);
277     }
278 bendavid 1.27
279 bendavid 1.37 const std::pair<bool,Measurement1D> &ip3dpv = IPTools::absoluteImpactParameter3D(tt,thevtx);
280 bendavid 1.27 if (ip3dpv.first) {
281 bendavid 1.37 outMuon->SetIp3dPV(thesign*ip3dpv.second.value());
282 bendavid 1.27 outMuon->SetIp3dPVErr(ip3dpv.second.error());
283     }
284 bendavid 1.28 else {
285     outMuon->SetIp3dPV(-99.0);
286     }
287 bendavid 1.27
288 bendavid 1.37 const std::pair<bool,Measurement1D> &d0pvbs = IPTools::absoluteTransverseImpactParameter(tt,thevtxbs);
289 bendavid 1.27 if (d0pvbs.first) {
290 bendavid 1.37 outMuon->SetD0PVBS(thesignbs*d0pvbs.second.value());
291 bendavid 1.27 outMuon->SetD0PVBSErr(d0pvbs.second.error());
292     }
293 bendavid 1.28 else {
294     outMuon->SetD0PVBS(-99.0);
295     }
296 bendavid 1.27
297 bendavid 1.37 const std::pair<bool,Measurement1D> &ip3dpvbs = IPTools::absoluteImpactParameter3D(tt,thevtxbs);
298 bendavid 1.27 if (ip3dpvbs.first) {
299 bendavid 1.37 outMuon->SetIp3dPVBS(thesignbs*ip3dpvbs.second.value());
300 bendavid 1.27 outMuon->SetIp3dPVBSErr(ip3dpvbs.second.error());
301     }
302 bendavid 1.28 else {
303     outMuon->SetIp3dPVBS(-99.0);
304     }
305    
306 bendavid 1.41 const std::pair<bool,Measurement1D> &d0pvub = IPTools::absoluteTransverseImpactParameter(tt,thevtxub);
307     if (d0pvub.first) {
308     outMuon->SetD0PVUB(thesign*d0pvub.second.value());
309     outMuon->SetD0PVUBErr(d0pvub.second.error());
310     }
311     else {
312     outMuon->SetD0PVUB(-99.0);
313     }
314    
315     const std::pair<bool,Measurement1D> &ip3dpvub = IPTools::absoluteImpactParameter3D(tt,thevtxub);
316     if (ip3dpvub.first) {
317     outMuon->SetIp3dPVUB(thesign*ip3dpvub.second.value());
318     outMuon->SetIp3dPVUBErr(ip3dpvub.second.error());
319     }
320     else {
321     outMuon->SetIp3dPVUB(-99.0);
322     }
323    
324     const std::pair<bool,Measurement1D> &d0pvubbs = IPTools::absoluteTransverseImpactParameter(tt,thevtxubbs);
325     if (d0pvubbs.first) {
326     outMuon->SetD0PVUBBS(thesignbs*d0pvubbs.second.value());
327     outMuon->SetD0PVUBBSErr(d0pvubbs.second.error());
328     }
329     else {
330     outMuon->SetD0PVUBBS(-99.0);
331     }
332    
333     const std::pair<bool,Measurement1D> &ip3dpvubbs = IPTools::absoluteImpactParameter3D(tt,thevtxubbs);
334     if (ip3dpvubbs.first) {
335     outMuon->SetIp3dPVUBBS(thesignbs*ip3dpvubbs.second.value());
336     outMuon->SetIp3dPVUBBSErr(ip3dpvubbs.second.error());
337     }
338     else {
339     outMuon->SetIp3dPVUBBS(-99.0);
340     }
341    
342 bendavid 1.28 //compute compatibility with PV using taking into account also the case where muon track
343     //was included in the vertex fit
344     if (iM->track()->extra().isAvailable()) {
345 paus 1.44
346 bendavid 1.28 const std::pair<bool,double> &pvCompat = kalmanEstimator.estimate(pvCol->at(0),tt);
347     if (pvCompat.first) {
348     outMuon->SetPVCompatibility(pvCompat.second);
349     }
350     else {
351     outMuon->SetPVCompatibility(-99.0);
352     }
353 paus 1.44
354 bendavid 1.28 const std::pair<bool,double> &pvbsCompat = kalmanEstimator.estimate(pvBSCol->at(0),tt);
355     if (pvbsCompat.first) {
356     outMuon->SetPVBSCompatibility(pvbsCompat.second);
357     }
358     else {
359     outMuon->SetPVBSCompatibility(-99.0);
360     }
361     }
362 bendavid 1.27
363     }
364 bendavid 1.39
365 loizides 1.14 outMuon->SetNChambers (iM->numberOfChambers());
366 pharris 1.11 outMuon->SetStationMask(iM->stationMask(reco::Muon::SegmentAndTrackArbitration));
367 ceballos 1.33 outMuon->SetNMatches (iM->numberOfMatches());
368 loizides 1.14 for(int i0 = 0; i0 < 4; i0++) {
369 loizides 1.20 // DTs
370 pharris 1.11 outMuon->SetDX(i0, iM->dX(i0+1,1));
371     outMuon->SetDY(i0, iM->dY(i0+1,1));
372     outMuon->SetPullX(i0, iM->pullX(i0+1,1));
373     outMuon->SetPullY(i0, iM->pullY(i0+1,1));
374     outMuon->SetTrackDist(i0, iM->trackDist(i0+1,1));
375     outMuon->SetTrackDistErr(i0, iM->trackDistErr(i0+1,1));
376 pharris 1.34 outMuon->SetNSegments(i0, NumberOfSegments(&(*iM),i0+1,1));
377 loizides 1.20 // CSCs
378 pharris 1.11 outMuon->SetDX(4+i0, iM->dX (i0+1,2));
379     outMuon->SetDY(4+i0, iM->dY (i0+1,2));
380     outMuon->SetPullX(4+i0, iM->pullX (i0+1,2));
381     outMuon->SetPullY(4+i0, iM->pullY (i0+1,2));
382     outMuon->SetTrackDist(4+i0, iM->trackDist(i0+1,2));
383     outMuon->SetTrackDistErr(4+i0,iM->trackDistErr(i0+1,2));
384 pharris 1.34 outMuon->SetNSegments(4+i0, NumberOfSegments(&(*iM),i0+1,2));
385 pharris 1.11 }
386 bendavid 1.19
387 bendavid 1.32 reco::MuonRef theRef(hMuonProduct, iM - inMuons.begin());
388     // fill corrected expected inner hits
389 bendavid 1.36 if (iM->innerTrack().isNonnull()) {
390 paus 1.44 outMuon->SetCorrectedNExpectedHitsInner (iM->innerTrack()->trackerExpectedHitsInner().numberOfHits());
391     outMuon->SetValidFraction (iM->innerTrack()->validFraction());
392     outMuon->SetNTrkLayersHit (iM->innerTrack()->hitPattern().trackerLayersWithMeasurement());
393     outMuon->SetNTrkLayersNoHit (iM->innerTrack()->hitPattern().trackerLayersWithoutMeasurement());
394     outMuon->SetNPxlLayersHit (iM->innerTrack()->hitPattern().pixelLayersWithMeasurement());
395     outMuon->SetNTrkLostHitsIn (iM->innerTrack()->trackerExpectedHitsInner().numberOfLostTrackerHits());
396     outMuon->SetNTrkLostHitsOut (iM->innerTrack()->trackerExpectedHitsOuter().numberOfLostTrackerHits());
397 bendavid 1.32 }
398 bendavid 1.36
399 loizides 1.20 // add muon to map
400 bendavid 1.19 edm::Ptr<reco::Muon> thePtr(hMuonProduct, iM - inMuons.begin());
401     muonMap_->Add(thePtr, outMuon);
402 paus 1.44
403 bendavid 1.21 if (verbose_>1) {
404     if (!outMuon->HasGlobalTrk() && !outMuon->HasStandaloneTrk()) {
405     printf("mithep::Muon, pt=%5f, eta=%5f, phi=%5f, mass=%5f\n",outMuon->Pt(),outMuon->Eta(),outMuon->Phi(), outMuon->Mass());
406 paus 1.44 printf(" reco::Muon, pt=%5f, eta=%5f, phi=%5f, mass=%5f\n",iM->pt(),iM->eta(),iM->phi(),iM->mass());
407 bendavid 1.21 }
408     }
409 bendavid 1.19
410 pharris 1.11 }
411 loizides 1.1 muons_->Trim();
412     }
413 pharris 1.34
414     int FillerMuons::NumberOfSegments(const reco::Muon *iM, int station, int muonSubdetId, reco::Muon::ArbitrationType type ) {
415     if(!iM->isMatchesValid()) return 0;
416     int segments(0);
417 paus 1.44
418 pharris 1.34 for( std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = iM->matches().begin();
419     chamberMatch != iM->matches().end(); chamberMatch++ )
420     {
421     if(chamberMatch->segmentMatches.empty()) continue;
422     if(!(chamberMatch->station()==station && chamberMatch->detector()==muonSubdetId)) continue;
423     if(type == reco::Muon::NoArbitration) {
424 paus 1.44 segments += chamberMatch->segmentMatches.size();
425     continue;
426 pharris 1.34 }
427     for( std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin();
428 paus 1.44 segmentMatch != chamberMatch->segmentMatches.end(); segmentMatch++ )
429     {
430     if(type == reco::Muon::SegmentArbitration)
431     if(segmentMatch->isMask(reco::MuonSegmentMatch::BestInChamberByDR)) {
432     segments++;
433     break;
434     }
435     if(type == reco::Muon::SegmentAndTrackArbitration)
436     if(segmentMatch->isMask(reco::MuonSegmentMatch::BestInChamberByDR) &&
437     segmentMatch->isMask(reco::MuonSegmentMatch::BelongsToTrackByDR)) {
438     segments++;
439     break;
440     }
441     if(type == reco::Muon::SegmentAndTrackArbitrationCleaned)
442     if(segmentMatch->isMask(reco::MuonSegmentMatch::BestInChamberByDR) &&
443     segmentMatch->isMask(reco::MuonSegmentMatch::BelongsToTrackByDR) &&
444     segmentMatch->isMask(reco::MuonSegmentMatch::BelongsToTrackByCleaning)) {
445     segments++;
446     break;
447     }
448     if(type > 1<<7)
449     if(segmentMatch->isMask(type)) {
450     segments++;
451     break;
452     }
453     }
454 pharris 1.34 }
455     return segments;
456     }