ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/Utils/interface/MuonTools.hh
Revision: 1.1
Committed: Mon Nov 3 13:19:17 2008 UTC (16 years, 6 months ago) by pharris
Content type: text/plain
Branch: MAIN
Log Message:
 Added the Muon Tools

File Contents

# User Rev Content
1 pharris 1.1 #ifndef MITWLNU_MUONTOOLS_HH
2     #define MITWLNU_MUONTOOLS_HH
3     #include <vector>
4     #include <iostream>
5    
6     #include "TH2D.h"
7     #include "TFile.h"
8     #include "MitAna/DataTree/interface/Collections.h"
9     #include "MitCommon/MathTools/interface/MathUtil.h"
10    
11     namespace mithep {
12     class MuonTools {
13     public:
14     MuonTools();
15     ~MuonTools();
16     double getCaloCompatability(mithep::Muon* iMuon,bool iEMSpecial, bool iCorrectedHCAL);
17     enum selection {
18     AllArbitrated,
19     PromptTight,
20     TMLastStationLoose,
21     TMLastStationTight,
22     TMOneStationLoose,
23     TMOneStationTight,
24     TM2DCompatibilityLoose,
25     TM2DCompatibilityTight
26     };
27     static inline double sigWeight(double iVal0,double iVal1) {
28     if(iVal1 < 1.) return 1;
29     if(iVal0 < 3. && iVal1 > 3.) {
30     double lVal = TMath::Max(iVal0,1.);
31     return 1/TMath::Power(lVal,0.25);
32     }
33     double lVal = TMath::Max(iVal1,1.);
34     return 1/TMath::Power(lVal,0.25);
35     }
36     static inline int lastHit(mithep::Muon *iMuon) {
37     int lId = -1;
38     for(int i0 = 0; i0 < 8; i0++) {
39     if(iMuon->GetDX(i0) < 99999 || iMuon->GetDY(i0) < 99999) lId = i0;
40     }
41     return lId;
42     }
43     static inline int maxChamberId(mithep::Muon* iMuon,double iMaxD,double iMaxP) {
44     int lId = -1;
45     for(int i0 = 0; i0 < 8; i0++) {
46     if(iMuon->GetTrackDist(i0) < iMaxD &&
47     iMuon->GetTrackDist(i0)/iMuon->GetTrackDistErr(i0) < iMaxP) {
48     if(iMuon->GetDX(i0) >= 999999) {lId = i0;} else {lId = -1;}
49     }
50     }
51     return lId;
52     }
53     static inline int lastStation(mithep::Muon* iMuon,double iMaxD,double iMaxP) {
54     int lId = -1;
55     for(int i0 = 0; i0 < 8; i0++) {
56     if((lId % 4) > (i0 % 4)) continue;
57     if(iMuon->GetTrackDist(i0) < iMaxD &&
58     iMuon->GetTrackDist(i0)/iMuon->GetTrackDistErr(i0) < iMaxP) lId = i0;
59     }
60     return lId;
61     }
62     static inline int lastStation(mithep::Muon* iMuon,int iMax=8) {
63     int lId = -1; if(iMax > 8) iMax = 8;
64     for(int i0 = 0; i0 < iMax; i0++) {if(iMuon->StationBit(i0) && ((lId % 4) < (i0 % 4)))lId = i0;}
65     return lId;
66     }
67     static inline bool overflow(TH2D *iHist,double lVal0,double lVal1) {
68     if(iHist == 0) return true;
69     if(iHist ->GetXaxis()->FindBin(lVal0) == 0 ||
70     iHist ->GetXaxis()->FindBin(lVal0) > iHist->GetNbinsX() ||
71     iHist ->GetYaxis()->FindBin(lVal0) == 0 ||
72     iHist ->GetYaxis()->FindBin(lVal0) > iHist->GetNbinsY()) return true;
73     return false;
74     }
75     static inline MCParticle* etaPhiMatch(MCParticleCol* iMCParts,Muon *iMuon) {
76     mithep::MCParticle *lMC = 0; double lDR = 1.; double lPt = -1.;
77     for(unsigned int i0 = 0; i0 < iMCParts->GetEntries(); i0++) {
78     mithep::MCParticle* pMC = iMCParts->At(i0);
79     if(pMC->Status() != 1) continue;
80     double pDR = mithep::MathUtils::DeltaR(pMC->Mom(),iMuon->Mom());
81     if(pDR > lDR && pDR > 0.01) continue;
82     if(fabs(pMC->Pt()) < lPt) continue;
83     lDR = pDR;
84     lMC = pMC;
85     lPt = pMC->Pt();
86     }
87     return lMC;
88     }
89     static inline Muon* match(MuonCol* iMuons,Track *iTrack) {
90     mithep::Muon *lMuon = 0; double lDR = 1.; double lPt = -1.;
91     for(unsigned int i0 = 0; i0 < iMuons->GetEntries(); i0++) {
92     mithep::Muon* pMuon = iMuons->At(i0);
93     double pDR = mithep::MathUtils::DeltaR(pMuon->Mom(),iTrack->Mom4(0.109));
94     if(pDR > lDR && pDR > 0.01) continue;
95     if(fabs(pMuon->Pt()) < lPt) continue;
96     lDR = pDR;
97     lMuon = pMuon;
98     lPt = pMuon->Pt();
99     }
100     return lMuon;
101     }
102     static inline int NSegments(mithep::Muon *iMuon) {
103     int lNSegs = 0;
104     for(int i0 = 0; i0 < 8; i0++) {
105     if(iMuon->StationBit(i0)) lNSegs++;
106     }
107     return lNSegs;
108     }
109     static inline bool promptTight(int iId,Muon *iMuon) {
110     const mithep::Track *lTrack = 0;
111     if(iId == 0 ) lTrack = iMuon->GlobalTrk();
112     if(iId == 1 || (iId == -1 && lTrack == 0)) lTrack = iMuon->TrackerTrk();
113     if(iId == 2 || (iId == -1 && lTrack == 0)) lTrack = iMuon->StandaloneTrk();
114     if(lTrack == 0) return 0;
115     if(lTrack->NHits() < 11) return false;
116     if(lTrack->Chi2()/lTrack->Ndof() > 10.) return false;
117     if(lTrack->D0() > 0.2) return false;
118     if((lastHit(iMuon)% 4) == 0) return false;
119     return true;
120     }
121     static inline bool TMLastStation(Muon *iMuon,double iDYMin = 3., double iPYMin = 3.,double iDXMin = 3., double iPXMin = 3.,int iN = 2) {
122     if(iMuon->NSegments() < iN) return false; //2Last 1One
123     int lLast = lastStation(iMuon,-3.,-3.); //Last Required Station
124     if(lLast < 0) return false;
125     if(iMuon->GetDX(lLast) > 9999.) return false;
126     lLast = lastStation(iMuon); //No Requirements Imply StationMask (Segment and Track Arbitrated)
127     if(lLast < 0) return false;
128     if(!(fabs(iMuon->GetDX(lLast)) < iDXMin ||
129     fabs(iMuon->GetPullX(lLast)) < iPXMin)) return false;
130     if(lLast == 3) lLast = lastStation(iMuon,3);
131     if(lLast < 0) return false;
132     if(!(fabs(iMuon->GetDY(lLast)) < iDYMin || //Old Code it was 9999
133     fabs(iMuon->GetPullY(lLast)) < iPYMin)) return false;
134     return true;
135     }
136     static inline bool TMOneStation(Muon *iMuon,double iDYMin = 3., double iPYMin = 3.,double iDXMin = 3., double iPXMin = 3.,int iN = 1) {
137     if(iMuon->NSegments() < iN) return false; //2Last 1One
138     bool pGoodX = false; bool pBadY = false;
139     for(int i0 = 0; i0 < 8; i0++) {
140     if((fabs(iMuon->GetDX(i0)) < iDXMin ||
141     fabs(iMuon->GetPullX(i0)) < iPXMin)) pGoodX = true;
142     if(pGoodX &&
143     (fabs(iMuon->GetDY(i0)) < iDYMin ||
144     fabs(iMuon->GetPullY(i0)) < iPYMin)) return true;
145     if(fabs(iMuon->GetDY(i0)) < 999999) pBadY = true;
146     if(i0 == 3 && pGoodX && !pBadY) return true;
147     }
148     return false;
149     }
150     static inline double getCaloCompatabilitySlowly(mithep::Muon* iMuon,bool iEMSpecial, bool iCorrectedHCAL) {
151     //if(fpion_em_etaEmi == 0) {
152     std::cout << "Warning Loading Compatability Slowly" << std::endl;
153     TFile* lPion_templates = new TFile("MitCommon/MuonTools/data/PionCaloTemplate.root","READ");
154     TFile* lMuon_templates = new TFile("MitCommon/MuonTools/data/MuonCaloTemplate.root","READ");
155     TH2D* lpion_em_etaEmi = (TH2D*) lPion_templates->Get("em_etaEmi");
156     TH2D* lpion_had_etaEmi = (TH2D*) lPion_templates->Get("had_etaEmi");
157     TH2D* lpion_had_etaTmi = (TH2D*) lPion_templates->Get("had_etaTmi");
158     TH2D* lpion_em_etaB = (TH2D*) lPion_templates->Get("em_etaB") ;
159     TH2D* lpion_had_etaB = (TH2D*) lPion_templates->Get("had_etaB") ;
160     TH2D* lpion_ho_etaB = (TH2D*) lPion_templates->Get("ho_etaB") ;
161     TH2D* lpion_had_etaTpl = (TH2D*) lPion_templates->Get("had_etaTpl");
162     TH2D* lpion_em_etaEpl = (TH2D*) lPion_templates->Get("em_etaEpl") ;
163     TH2D* lpion_had_etaEpl = (TH2D*) lPion_templates->Get("had_etaEpl");
164     TH2D* lmuon_em_etaEmi = (TH2D*) lMuon_templates->Get("em_etaEmi") ;
165     TH2D* lmuon_had_etaEmi = (TH2D*) lMuon_templates->Get("had_etaEmi");
166     TH2D* lmuon_had_etaTmi = (TH2D*) lMuon_templates->Get("had_etaTmi");
167     TH2D* lmuon_em_etaB = (TH2D*) lMuon_templates->Get("em_etaB") ;
168     TH2D* lmuon_had_etaB = (TH2D*) lMuon_templates->Get("had_etaB") ;
169     TH2D* lmuon_ho_etaB = (TH2D*) lMuon_templates->Get("ho_etaB") ;
170     TH2D* lmuon_had_etaTpl = (TH2D*) lMuon_templates->Get("had_etaTpl");
171     TH2D* lmuon_em_etaEpl = (TH2D*) lMuon_templates->Get("em_etaEpl");
172     TH2D* lmuon_had_etaEpl = (TH2D*) lMuon_templates->Get("had_etaEpl");
173     //}
174     double lEta = -1.; double lP = -1;
175     double lEM = -5.; double lHad = 0; double lHO = 0;
176     lEta = iMuon->Eta();
177     lP = iMuon->P();
178     if(lP >= 2000.) lP = 1999.9;
179     if(!iEMSpecial || iMuon->EmEnergy() != 0.) lEM = iMuon->EmEnergy();
180     lHad = iMuon->HadEnergy();
181     lHO = iMuon->HoEnergy();
182     if(lP < 0. ) return 0.5;
183     if(fabs(lEta) > 2.5 ) return 0.5;
184     TH2D* lTMuonHad = NULL;
185     TH2D* lTPionHad = NULL;
186     TH2D* lTMuonHo = NULL;
187     TH2D* lTPionHo = NULL;
188     TH2D* lTMuonEm = NULL;
189     TH2D* lTPionEm = NULL;
190    
191     if(fabs(lEta) >= 1.27) {
192     if(iCorrectedHCAL) lHad *= 1.8/2.2;
193     if(lEta > 0) {
194     lTPionHad = lpion_had_etaEpl;
195     lTMuonHad = lmuon_had_etaEpl;
196     } else {
197     lTPionHad = lpion_had_etaEmi;
198     lTMuonHad = lmuon_had_etaEmi;
199     }
200     }
201     if(fabs(lEta) < 1.27 && fabs(lEta) >= 1.1 ) {
202     if(iCorrectedHCAL) lHad *= (1.8/(-2.2*fabs(lEta)+5.5));
203     if(lEta > 0) {
204     lTPionHad = lpion_had_etaTpl;
205     lTMuonHad = lmuon_had_etaTpl;
206     } else {
207     lTPionHad = lpion_had_etaTmi;
208     lTMuonHad = lmuon_had_etaTmi;
209     }
210     }
211     if(fabs(lEta) < 1.1) {
212     if(iCorrectedHCAL) lHad *= sin(2*atan(exp(iMuon->Eta())));
213     lTPionHad = lpion_had_etaB;
214     lTMuonHad = lmuon_had_etaB;
215     }
216     if(lEta > 1.479 ) {
217     lTPionEm = lpion_em_etaEpl;
218     lTMuonEm = lmuon_em_etaEpl;
219     }
220     if(fabs(lEta) <= 1.479) {
221     lTPionEm = lpion_em_etaB;
222     lTMuonEm = lmuon_em_etaB;
223     }
224     if(lEta < -1.479 ) {
225     lTPionEm = lpion_em_etaEmi;
226     lTMuonEm = lmuon_em_etaEmi;
227     }
228     if(fabs(lEta) < 1.28) {
229     lTPionHo = lpion_ho_etaB;
230     lTMuonHo = lmuon_ho_etaB;
231     }
232    
233     double lPBX = 1.; double lPSX = 1.;
234     double lPBY = 1.; double lPSY = 1.;
235     double lPBZ = 1.; double lPSZ = 1.;
236     if(!overflow(lTPionEm, lP,lEM)) lPBX = lTPionEm ->GetBinContent(lTPionEm ->GetXaxis()->FindBin(lP),lTPionEm ->GetYaxis()->FindBin(lEM) );
237     if(!overflow(lTPionHad,lP,lHad)) lPBY = lTPionHad->GetBinContent(lTPionHad->GetXaxis()->FindBin(lP),lTPionHad->GetYaxis()->FindBin(lHad));
238     if(!overflow(lTPionHo, lP,lHO)) lPBZ = lTPionHo ->GetBinContent(lTPionHo ->GetXaxis()->FindBin(lP),lTPionHo ->GetYaxis()->FindBin(lHO) );
239     if(!overflow(lTMuonEm, lP,lEM )) lPSX = lTMuonEm ->GetBinContent(lTMuonEm ->GetXaxis()->FindBin(lP),lTMuonEm ->GetYaxis()->FindBin(lEM) );
240     if(!overflow(lTMuonHad,lP,lHad)) lPSY = lTMuonHad->GetBinContent(lTMuonHad->GetXaxis()->FindBin(lP),lTMuonHad->GetYaxis()->FindBin(lHad));
241     if(!overflow(lTMuonHo ,lP,lHO)) lPSZ = lTMuonHo ->GetBinContent(lTMuonHo ->GetXaxis()->FindBin(lP),lTMuonHo ->GetYaxis()->FindBin(lHO) );
242    
243     if(lPSX == 0. || lPBX == 0. || (lEM <= 0. && !iEMSpecial)) {lPSX = 1.; lPBX = 1.;}
244     if(lPSY == 0. || lPBY == 0. || lHad == 0.) {lPSY = 1.; lPBY = 1.;}
245     if(lPSZ == 0. || lPBZ == 0. || lHO == 0.) {lPSZ = 1.; lPBZ = 1.;}
246     lPion_templates->Close();
247     lMuon_templates->Close();
248     if((lPSX*lPSY*lPSZ+lPBX*lPBY*lPBZ) > 0.) return lPSX*lPSY*lPSZ/(lPSX*lPSY*lPSZ+lPBX*lPBY*lPBZ);
249     return 0.5;
250     }
251     static inline float getSegmentCompatability(const mithep::Muon* iMuon) {
252     int lNStationsCrossed = 0;
253     int lNStationsSegment = 0;
254    
255     std::vector<int> lStTrack(8);
256     std::vector<int> lStSegmentmatch(8);
257     std::vector<int> lStCrossed(8);
258     std::vector<float> lStBoundary(8);
259     std::vector<float> lStWeight(8);
260     float lWeight = 0.;
261     bool lAdjust = true;
262     for(int i0 = 0; i0 < 8; i0++) {
263     if(iMuon->GetTrackDist(i0) < 999999. ) {
264     lNStationsCrossed++;
265     lStCrossed[i0] = 1;
266     lStBoundary[i0] = 0.;
267     if(iMuon->GetTrackDist(i0) > -10. ) lStBoundary[i0] = iMuon->GetTrackDist(i0);
268     }
269     if(iMuon->GetDX(i0) < 999999.) { //Use iMuon->GetSegmentX--> CHECK
270     lNStationsSegment++;
271     lStSegmentmatch[i0] = 1;
272     }
273     if(iMuon->GetDY(i0) < 999999.) lAdjust = false;
274     }
275     int lPCross = -1;
276     const float lAtWeight = 0.5;
277     for(int i0 = 0; i0< 8; i0++) {
278     lStWeight[i0] = 0;
279     if(lStCrossed[i0] > 0 ) { lPCross++;
280     switch ( lNStationsCrossed ) {
281     case 1 :
282     lStWeight[i0] = 1.;
283     break;
284     case 2 :
285     if ( lPCross == 0 ) lStWeight[i0] = 0.33;
286     else lStWeight[i0] = 0.67;
287     break;
288     case 3 :
289     if ( lPCross == 0 ) lStWeight[i0] = 0.23;
290     else if( lPCross == 1 ) lStWeight[i0] = 0.33;
291     else lStWeight[i0] = 0.44;
292     break;
293     case 4 :
294     if ( lPCross == 0 ) lStWeight[i0] = 0.10;
295     else if( lPCross == 1 ) lStWeight[i0] = 0.20;
296     else if( lPCross == 2 ) lStWeight[i0] = 0.30;
297     else lStWeight[i0] = 0.40;
298     break;
299    
300     default :
301     lStWeight[i0] = 1./lNStationsCrossed;
302     }
303    
304     if(lStSegmentmatch[i0] <= 0 && lStBoundary[i0] != 0. ) {
305     lStWeight[i0] *= lAtWeight*0.5*(TMath::Erf(lStBoundary[i0]/6.)+1.);
306     } else if(lStSegmentmatch[i0] <= 0 && lStBoundary[i0] == 0) {lStWeight[i0] = 0.;}
307    
308     if( lStSegmentmatch[i0] > 0) {
309     double lP2X = TMath::Power(iMuon->GetPullX(i0),2.);
310     double lP2Y = TMath::Power(iMuon->GetPullY(i0),2.);
311     double lD2X = TMath::Power(iMuon->GetDX(i0),2.);
312     double lD2Y = TMath::Power(iMuon->GetDY(i0),2.);
313     if( iMuon->GetDY(i0) < 999999 && iMuon->GetDX(i0) < 999999) {
314     lStWeight[i0] *= sigWeight(sqrt(lD2X+lD2Y),sqrt(lP2X+lP2Y));
315     } else if (iMuon->GetDY(i0) >= 999999 && i0 < 4) {
316     lStWeight[i0] *= sigWeight(iMuon->GetDX(i0),iMuon->GetPullX(i0));
317     } else if(i0 < 4) {
318     lStWeight[i0] *= sigWeight(iMuon->GetDY(i0),iMuon->GetPullY(i0));
319     }
320     }
321     }
322     lWeight += lStWeight[i0];
323     }
324     if(lNStationsCrossed == 0) lWeight = 0.5;
325     return lWeight;
326     }
327     static inline bool isGoodMuon(mithep::Muon *iMuon,selection iSelection) {
328     double lVal = 0;
329     switch(iSelection) {
330     case AllArbitrated:
331     if(iMuon->StandaloneTrk() != 0 || iMuon->GlobalTrk()!= 0) return true;
332     if(iMuon->NSegments() > 0) return true;
333     return false;
334     break;
335     case PromptTight:
336     return promptTight(-1,iMuon);
337     break;
338     case TMOneStationLoose:
339     return TMOneStation(iMuon,99999,999999);
340     case TMOneStationTight:
341     return TMOneStation(iMuon);
342     case TMLastStationLoose:
343     return TMLastStation(iMuon,999999,999999);
344     break;
345     case TMLastStationTight:
346     return TMLastStation(iMuon);
347     break;
348     case TM2DCompatibilityLoose:
349     lVal = 1.2*getSegmentCompatability(iMuon);
350     if(lVal/1.2 == 0.5) return false;
351     //if(fClass ) lVal += 0.8*getCaloCompatability(iMuon,true,true);
352     lVal += 0.8*getCaloCompatabilitySlowly(iMuon,true,true);
353     if(lVal > 0.7) return true;
354     return false;
355     break;
356     case TM2DCompatibilityTight:
357     lVal = 1.2*getSegmentCompatability(iMuon);
358     if(lVal/1.2 == 0.5) return false;
359     //if(fClass ) lVal += 0.8*getCaloCompatability(iMuon,true,true);
360     lVal += 0.8*getCaloCompatabilitySlowly(iMuon,true,true);
361     if(lVal > 1.0) return true;
362     return false;
363     break;
364     default:
365     return false;
366     }
367     return false;
368     }
369     bool isGood(mithep::Muon *iMuon,selection iSelection);
370     private:
371     TFile *fMuon_templates;
372     TFile *fPion_templates;
373     TH2D* fpion_em_etaEmi;
374     TH2D* fpion_had_etaEmi;
375     TH2D* fpion_had_etaTmi;
376     TH2D* fpion_em_etaB;
377     TH2D* fpion_had_etaB;
378     TH2D* fpion_ho_etaB;
379     TH2D* fpion_had_etaTpl;
380     TH2D* fpion_em_etaEpl;
381     TH2D* fpion_had_etaEpl;
382     TH2D* fmuon_em_etaEmi;
383     TH2D* fmuon_had_etaEmi;
384     TH2D* fmuon_had_etaTmi;
385     TH2D* fmuon_em_etaB;
386     TH2D* fmuon_had_etaB;
387     TH2D* fmuon_ho_etaB;
388     TH2D* fmuon_had_etaTpl;
389     TH2D* fmuon_em_etaEpl;
390     TH2D* fmuon_had_etaEpl;
391     };
392     }
393    
394     #endif