]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMuonForwardTrackAnalysis.cxx
Updates to the cluster finder and the track finder
[u/mrichter/AliRoot.git] / MFT / AliMuonForwardTrackAnalysis.cxx
1 //====================================================================================================================================================
2 //
3 //      Class for the analysis of the ALICE muon forward tracks (MUON + MFT)
4 //
5 //      Contact author: antonio.uras@cern.ch
6 //
7 //====================================================================================================================================================
8
9 #include "TObject.h"
10 #include "TClonesArray.h"
11 #include "AliMuonForwardTrack.h"
12 #include "AliMuonForwardTrackPair.h"
13 #include "TMatrixD.h"
14 #include "TTree.h"
15 #include "TH1D.h"
16 #include "AliLog.h"
17 #include "TFile.h"
18 #include "TParticle.h"
19 #include "AliMUONTrackParam.h"
20 #include "AliMUONTrackExtrap.h"
21 #include "TGeoManager.h"
22 #include "TRandom.h"
23 #include "TLorentzVector.h"
24 #include "TDatabasePDG.h"
25 #include "TGraph.h"
26 #include "AliMuonForwardTrackAnalysis.h"
27
28 ClassImp(AliMuonForwardTrackAnalysis)
29
30 //====================================================================================================================================================
31
32 AliMuonForwardTrackAnalysis::AliMuonForwardTrackAnalysis():
33   TObject(),
34   fInputDir(0),
35   fOutputDir(0),
36   fInputTree(0),
37   fMuonForwardTracks(0),
38   fMuonForwardTrackPairs(0),
39   fMFTTrack(0),
40   fMFTTrackPair(0),
41   fMCRefTrack(0),
42   fEv(0),
43   fFirstEvent(-1),
44   fLastEvent(-1),
45   fNTracksOfEvent(0),
46   fNTracksAnalyzedOfEvent(0),
47   fNTracksAnalyzed(0),
48   fNPairsOfEvent(0),
49   fNPairsAnalyzedOfEvent(0),
50   fHistOffsetSingleMuonsX(0x0),
51   fHistOffsetSingleMuonsY(0x0),
52   fHistOffsetSingleMuons(0x0),
53   fHistWOffsetSingleMuons(0x0),
54   fHistErrorSingleMuonsX(0x0),
55   fHistErrorSingleMuonsY(0x0),
56   fHistOffsetSingleMuonsX_vsPtRapidity(0x0),
57   fHistOffsetSingleMuonsY_vsPtRapidity(0x0),
58   fHistSingleMuonsPtRapidity(0x0),
59   fHistSingleMuonsOffsetChi2(0x0),
60   fHistWOffsetMuonPairs(0x0),
61   fHistMassMuonPairs(0x0),
62   fHistMassMuonPairsWithoutMFT(0x0),
63   fHistMassMuonPairsMC(0x0),
64   fHistRapidityPtMuonPairsMC(0x0),
65   fGraphSingleMuonsOffsetChi2(0x0),
66   fNMassBins(1000),
67   fMassMin(0),
68   fMassMax(10),
69   fSingleMuonAnalysis(1),
70   fMuonPairAnalysis(1),
71   fMatchTrigger(0),
72   fOption(0),
73   fXVertResMC(50.e-4),
74   fYVertResMC(50.e-4),
75   fZVertResMC(50.e-4),
76   fMaxNWrongClustersMC(999),
77   fPtMinSingleMuons(0)
78 {
79
80   // default constructor
81
82   for (Int_t rapBin=0; rapBin<fNRapBinsOffsetSingleMuons; rapBin++) {
83     for (Int_t ptBin=0; ptBin<fNPtBinsOffsetSingleMuons; ptBin++) {
84       fHistOffsetSingleMuonsX_tmp[rapBin][ptBin] = 0x0;
85       fHistOffsetSingleMuonsY_tmp[rapBin][ptBin] = 0x0;
86     }
87   }
88
89 }
90
91 //====================================================================================================================================================
92
93 Bool_t AliMuonForwardTrackAnalysis::Init(Char_t *inputFileName) {
94
95   BookHistos();
96
97   TFile *inputFile = new TFile(Form("%s/%s",fInputDir.Data(),inputFileName));
98   if (!inputFile || !inputFile->IsOpen()) {
99     AliError(Form("Error opening file %s", inputFileName));
100     return kFALSE;
101   }
102   fInputTree = (TTree*) inputFile->Get("AliMuonForwardTracks");
103   if (!fInputTree) {
104     AliError("Error reading input tree");
105     return kFALSE;
106   }
107
108   if (fFirstEvent<0 || fLastEvent<0 || fFirstEvent>fLastEvent || fFirstEvent>=fInputTree->GetEntries()) {
109     fFirstEvent = 0;
110     fLastEvent  = fInputTree->GetEntries()-1;
111   }
112   else {
113     fLastEvent = TMath::Min(fLastEvent, Int_t(fInputTree->GetEntries()-1));
114   }
115
116   AliInfo(Form("Analysing events %d to %d", fFirstEvent, fLastEvent));
117
118   fMuonForwardTracks = new TClonesArray("AliMuonForwardTrack");
119   fInputTree->SetBranchAddress("tracks", &fMuonForwardTracks);  
120
121   TGeoManager::Import(Form("%s/geometry.root",fInputDir.Data()));
122
123   AliMUONTrackExtrap::SetField();
124
125   fMuonForwardTrackPairs = new TClonesArray("AliMuonForwardTrackPair");
126
127   return kTRUE;
128
129 }
130
131 //====================================================================================================================================================
132
133 Bool_t AliMuonForwardTrackAnalysis::LoadNextEvent() {
134
135   if (fEv>fLastEvent) return kFALSE;
136   if (fEv<fFirstEvent) { fEv++; return kTRUE; }
137   fMuonForwardTracks -> Clear();
138   fInputTree->GetEvent(fEv);
139   AliInfo(Form("**** analyzing event # %4d (%3d tracks) ****", fEv, fMuonForwardTracks->GetEntries()));
140
141   if (fSingleMuonAnalysis) {
142     fNTracksAnalyzedOfEvent = 0;
143     fNTracksOfEvent = fMuonForwardTracks->GetEntries();
144     while (AnalyzeSingleMuon()) continue;
145   }
146   
147   if (fMuonPairAnalysis) {
148     if (fMuonForwardTrackPairs) {
149       fMuonForwardTrackPairs->Clear();
150     }
151     BuildMuonPairs();
152     fNPairsAnalyzedOfEvent = 0;
153     fNPairsOfEvent = fMuonForwardTrackPairs->GetEntries();
154     while (AnalyzeMuonPair()) continue;
155   }
156
157   fEv++;
158   
159   return kTRUE;
160
161 }
162
163 //====================================================================================================================================================
164
165 Bool_t AliMuonForwardTrackAnalysis::AnalyzeSingleMuon() {
166
167   if (fNTracksAnalyzedOfEvent>=fNTracksOfEvent) return kFALSE;
168
169   fMFTTrack = (AliMuonForwardTrack*) fMuonForwardTracks->At(fNTracksAnalyzedOfEvent);
170   fNTracksAnalyzedOfEvent++;
171   if (fMatchTrigger && !fMFTTrack->GetMatchTrigger()) return kTRUE;
172   fMCRefTrack = fMFTTrack->GetMCTrackRef();
173
174   if (!fMCRefTrack) return kTRUE;
175   if (fMFTTrack->GetNWrongClustersMC()>fMaxNWrongClustersMC) return kTRUE;
176
177   Double_t xOrig=gRandom->Gaus(0., fXVertResMC);
178   Double_t yOrig=gRandom->Gaus(0., fYVertResMC);
179 //   Double_t xOrig = 0.;
180 //   Double_t yOrig = 0.;
181   Double_t zOrig=gRandom->Gaus(0., fZVertResMC);
182
183   AliMUONTrackParam *param = fMFTTrack->GetTrackParamAtMFTCluster(0);
184   AliMUONTrackExtrap::ExtrapToZCov(param, zOrig);
185
186   TLorentzVector pMu;
187   Double_t mMu = TDatabasePDG::Instance()->GetParticle("mu-")->Mass();
188   Double_t energy = TMath::Sqrt(param->P()*param->P() + mMu*mMu);
189   pMu.SetPxPyPzE(param->Px(), param->Py(), param->Pz(), energy);
190
191   if (fMFTTrack->Pt()<fPtMinSingleMuons) return kTRUE;
192
193   TMatrixD cov(5,5);
194   cov = param->GetCovariances();
195
196   fHistErrorSingleMuonsX -> Fill(1.e4*TMath::Sqrt(cov(0,0)));
197   fHistErrorSingleMuonsY -> Fill(1.e4*TMath::Sqrt(cov(2,2)));
198
199   Double_t dX = fMFTTrack->GetOffsetX(xOrig, zOrig);
200   Double_t dY = fMFTTrack->GetOffsetY(yOrig, zOrig);
201   
202   Double_t offset = fMFTTrack->GetOffset(xOrig, yOrig, zOrig);
203   Double_t weightedOffset = fMFTTrack->GetWeightedOffset(xOrig, yOrig, zOrig);
204
205   //  AliDebug(2, Form("pdg code = %d\n", fMCRefTrack->GetPdgCode()));
206
207   fHistOffsetSingleMuonsX -> Fill(1.e4*dX);
208   fHistOffsetSingleMuonsY -> Fill(1.e4*dY);
209   Int_t rapBin = fHistOffsetSingleMuonsX_vsPtRapidity->GetXaxis()->FindBin(pMu.Rapidity());
210   Int_t ptBin  = fHistOffsetSingleMuonsX_vsPtRapidity->GetYaxis()->FindBin(pMu.Pt());
211   //  AliDebug(1, Form("pt = %f (%d), rap = %f (%d)\n", pMu.Pt(), pMu.Rapidity(), ptBin, rapBin));
212   if (0<rapBin && rapBin<=fNRapBinsOffsetSingleMuons && 0<ptBin && ptBin<=fNPtBinsOffsetSingleMuons) {
213     fHistOffsetSingleMuonsX_tmp[rapBin-1][ptBin-1]->Fill(1.e4*dX);
214     fHistOffsetSingleMuonsY_tmp[rapBin-1][ptBin-1]->Fill(1.e4*dY);
215   }
216   fHistSingleMuonsPtRapidity  -> Fill(pMu.Rapidity(), pMu.Pt());
217   fHistOffsetSingleMuons      -> Fill(1.e4*offset);
218   fHistWOffsetSingleMuons     -> Fill(weightedOffset);
219   Double_t chi2OverNdf = fMFTTrack->GetGlobalChi2()/Double_t(fMFTTrack->GetNMFTClusters()+fMFTTrack->GetNMUONClusters());
220   fHistSingleMuonsOffsetChi2  -> Fill(1.e4*offset, chi2OverNdf);
221   fGraphSingleMuonsOffsetChi2 -> SetPoint(fGraphSingleMuonsOffsetChi2->GetN(),1.e4*offset, chi2OverNdf);
222
223   fNTracksAnalyzed++;
224
225   return kTRUE;
226
227 }
228
229 //====================================================================================================================================================
230
231 Bool_t AliMuonForwardTrackAnalysis::AnalyzeMuonPair() {
232
233   if (fNPairsAnalyzedOfEvent>=fNPairsOfEvent) return kFALSE;
234
235   fMFTTrackPair = (AliMuonForwardTrackPair*) fMuonForwardTrackPairs->At(fNPairsAnalyzedOfEvent);
236
237   if (fOption==kResonanceOnly && !fMFTTrackPair->IsResonance()) {
238     fNPairsAnalyzedOfEvent++;
239     return kTRUE;
240   }
241
242   Double_t xOrig=gRandom->Gaus(0., fXVertResMC);
243   Double_t yOrig=gRandom->Gaus(0., fYVertResMC);
244   Double_t zOrig=gRandom->Gaus(0., fZVertResMC);
245   AliDebug(1, Form("origin = (%f, %f, %f)", xOrig, yOrig, zOrig));
246
247   fHistMassMuonPairs           -> Fill(fMFTTrackPair->GetMass(zOrig));
248   fHistWOffsetMuonPairs        -> Fill(fMFTTrackPair->GetWeightedOffset(xOrig, yOrig, zOrig));
249   fHistMassMuonPairsWithoutMFT -> Fill(fMFTTrackPair->GetMassWithoutMFT(xOrig, yOrig, zOrig));
250   fHistMassMuonPairsMC         -> Fill(fMFTTrackPair->GetMassMC());
251   fHistRapidityPtMuonPairsMC   -> Fill(fMFTTrackPair->GetRapidityMC(), fMFTTrackPair->GetPtMC());
252
253   AliDebug(1, Form("mass = %f   MC = %f", fMFTTrackPair->GetMass(zOrig), fMFTTrackPair->GetMassMC()));
254
255   fNPairsAnalyzedOfEvent++;
256
257   return kTRUE;
258
259 }
260
261 //====================================================================================================================================================
262
263 void AliMuonForwardTrackAnalysis::BuildMuonPairs() {
264
265   for (Int_t iTrack=0; iTrack<fMuonForwardTracks->GetEntries(); iTrack++) {
266     for (Int_t jTrack=0; jTrack<iTrack; jTrack++) {
267     
268       AliMuonForwardTrack *track0 = (AliMuonForwardTrack*) fMuonForwardTracks->At(iTrack);
269       AliMuonForwardTrack *track1 = (AliMuonForwardTrack*) fMuonForwardTracks->At(jTrack);
270
271       if (fMatchTrigger) if (!track0->GetMatchTrigger() || !track1->GetMatchTrigger()) continue;
272       if (!track0->GetMCTrackRef() || !track1->GetMCTrackRef()) continue;
273       if (track0->GetNWrongClustersMC()>fMaxNWrongClustersMC || track1->GetNWrongClustersMC()>fMaxNWrongClustersMC) continue;
274       if (track0->Pt()<fPtMinSingleMuons || track1->Pt()<fPtMinSingleMuons) continue;
275
276       AliMuonForwardTrackPair *trackPair = new AliMuonForwardTrackPair(track0, track1);
277       if (fOption==kResonanceOnly && !trackPair->IsResonance()) {
278         delete trackPair;
279         continue;
280       }
281       new ((*fMuonForwardTrackPairs)[fMuonForwardTrackPairs->GetEntries()]) AliMuonForwardTrackPair(*trackPair);
282
283     }
284   }
285
286 }
287
288 //====================================================================================================================================================
289
290 void AliMuonForwardTrackAnalysis::Terminate(Char_t *outputFileName) {
291
292   for (Int_t rapBin=0; rapBin<fNRapBinsOffsetSingleMuons; rapBin++) {
293     for (Int_t ptBin=0; ptBin<fNPtBinsOffsetSingleMuons; ptBin++) {
294       Int_t binMin_x = fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->FindBin(-3*fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->GetRMS());
295       Int_t binMax_x = fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->FindBin(+3*fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->GetRMS());
296       for (Int_t bin=1; bin<=fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->GetNbinsX(); bin++) {
297         if (bin<binMin_x || bin>binMax_x) fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->SetBinContent(bin, 0);
298       }
299       Int_t binMin_y = fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->FindBin(-3*fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->GetRMS());
300       Int_t binMax_y = fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->FindBin(+3*fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->GetRMS());
301       for (Int_t bin=1; bin<=fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->GetNbinsX(); bin++) {
302         if (bin<binMin_y || bin>binMax_y) fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->SetBinContent(bin, 0);
303       }
304       fHistOffsetSingleMuonsX_vsPtRapidity->SetBinContent(rapBin+1, ptBin+1, fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->GetRMS());
305       fHistOffsetSingleMuonsX_vsPtRapidity->SetBinError(rapBin+1, ptBin+1, fHistOffsetSingleMuonsX_tmp[rapBin][ptBin]->GetRMSError());
306       fHistOffsetSingleMuonsY_vsPtRapidity->SetBinContent(rapBin+1, ptBin+1, fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->GetRMS());
307       fHistOffsetSingleMuonsY_vsPtRapidity->SetBinError(rapBin+1, ptBin+1, fHistOffsetSingleMuonsY_tmp[rapBin][ptBin]->GetRMSError());
308     }
309   }
310
311   TFile *fileOut = new TFile(Form("%s/%s",fOutputDir.Data(),outputFileName), "recreate");
312
313   printf("Writing output objects to file %s\n", fileOut->GetName());
314
315   fHistOffsetSingleMuonsX -> Write();
316   fHistOffsetSingleMuonsY -> Write();
317   fHistOffsetSingleMuons  -> Write();
318   fHistWOffsetSingleMuons -> Write();
319   fHistErrorSingleMuonsX  -> Write();
320   fHistErrorSingleMuonsY  -> Write();
321
322   fHistOffsetSingleMuonsX_vsPtRapidity -> Write();
323   fHistOffsetSingleMuonsY_vsPtRapidity -> Write();
324
325 //   for (Int_t rapBin=0; rapBin<fNRapBinsOffsetSingleMuons; rapBin++) {
326 //     for (Int_t ptBin=0; ptBin<fNPtBinsOffsetSingleMuons; ptBin++) {
327 //       fHistOffsetSingleMuonsX_tmp[rapBin][ptBin] -> Write();
328 //       fHistOffsetSingleMuonsY_tmp[rapBin][ptBin] -> Write();
329 //     }
330 //   }
331
332   fHistSingleMuonsPtRapidity -> Write();
333   fHistSingleMuonsOffsetChi2 -> Write();
334
335   fGraphSingleMuonsOffsetChi2 -> Write();
336
337   fHistWOffsetMuonPairs        -> Write();
338   fHistMassMuonPairs           -> Write();
339   fHistMassMuonPairsWithoutMFT -> Write();
340   fHistMassMuonPairsMC         -> Write();
341   fHistRapidityPtMuonPairsMC   -> Write();
342
343   fileOut -> Close();
344
345 }
346
347 //====================================================================================================================================================
348
349 void AliMuonForwardTrackAnalysis::BookHistos() {
350
351   fHistOffsetSingleMuonsX = new TH1D("fHistOffsetSingleMuonsX", "Offset for single muons along X",  200, -1000, 1000);
352   fHistOffsetSingleMuonsY = new TH1D("fHistOffsetSingleMuonsY", "Offset for single muons along Y",  200, -1000, 1000);
353   fHistErrorSingleMuonsX  = new TH1D("fHistErrorSingleMuonsX",  "Coordinate Error for single muons along X",  200, 0, 1000);
354   fHistErrorSingleMuonsY  = new TH1D("fHistErrorSingleMuonsY",  "Coordinate Error for single muons along Y",  200, 0, 1000);
355   fHistOffsetSingleMuons  = new TH1D("fHistOffsetSingleMuons",  "Offset for single muons",          200, 0, 2000);
356   fHistWOffsetSingleMuons = new TH1D("fHistWOffsetSingleMuons", "Weighted Offset for single muons", 300, 0, 15);  
357
358   fHistOffsetSingleMuonsX_vsPtRapidity = new TH2D("fHistOffsetSingleMuonsX_vsPtRapidity", "Offset for single muons along X", 
359                                                   10, -4, -2.5, 10, 0.5, 5.5);
360   fHistOffsetSingleMuonsY_vsPtRapidity = new TH2D("fHistOffsetSingleMuonsY_vsPtRapidity", "Offset for single muons along Y", 
361                                                   10, -4, -2.5, 10, 0.5, 5.5);
362
363   for (Int_t rapBin=0; rapBin<fNRapBinsOffsetSingleMuons; rapBin++) {
364     for (Int_t ptBin=0; ptBin<fNPtBinsOffsetSingleMuons; ptBin++) {
365       fHistOffsetSingleMuonsX_tmp[rapBin][ptBin] = new TH1D(Form("fHistOffsetSingleMuonsX_tmp_%02d_%02d",rapBin,ptBin), "",  1000, -1000, 1000);
366       fHistOffsetSingleMuonsY_tmp[rapBin][ptBin] = new TH1D(Form("fHistOffsetSingleMuonsY_tmp_%02d_%02d",rapBin,ptBin), "",  1000, -1000, 1000);
367     }
368   }
369
370   fHistSingleMuonsPtRapidity = new TH2D("fHistSingleMuonsPtRapidity", "Phase Space for single muons", 10, -4, -2.5, 10, 0.5, 5.5);
371   fHistSingleMuonsOffsetChi2 = new TH2D("fHistSingleMuonsOffsetChi2", "Offset vs #chi^{2}/ndf for single muons", 400, 0, 4000, 100, 0, 20);
372
373   fHistOffsetSingleMuonsX -> SetXTitle("Offset(X)  [#mum]");
374   fHistOffsetSingleMuonsY -> SetXTitle("Offset(Y)  [#mum]");
375   fHistErrorSingleMuonsX  -> SetXTitle("Err. on track position at z_{vtx} (X)  [#mum]");
376   fHistErrorSingleMuonsY  -> SetXTitle("Err. on track position at z_{vtx} (Y)  [#mum]");
377   fHistOffsetSingleMuons  -> SetXTitle("Offset  [#mum]");
378   fHistWOffsetSingleMuons -> SetXTitle("Weighted Offset");
379
380   fHistOffsetSingleMuonsX_vsPtRapidity -> SetXTitle("y^{#mu}");
381   fHistOffsetSingleMuonsY_vsPtRapidity -> SetXTitle("y^{#mu}");
382   fHistOffsetSingleMuonsX_vsPtRapidity -> SetYTitle("p_{T}^{#mu}  [GeV/c]");
383   fHistOffsetSingleMuonsY_vsPtRapidity -> SetYTitle("p_{T}^{#mu}  [GeV/c]");
384
385   fHistSingleMuonsPtRapidity -> SetXTitle("y^{#mu}");
386   fHistSingleMuonsPtRapidity -> SetYTitle("p_{T}^{#mu}  [GeV/c]");
387   fHistSingleMuonsOffsetChi2 -> SetXTitle("Offset  [#mum]");
388   fHistSingleMuonsOffsetChi2 -> SetYTitle("#chi^{2}/ndf");
389
390   fHistOffsetSingleMuonsX -> Sumw2();
391   fHistOffsetSingleMuonsY -> Sumw2();
392   fHistErrorSingleMuonsX  -> Sumw2();
393   fHistErrorSingleMuonsY  -> Sumw2();
394   fHistOffsetSingleMuons  -> Sumw2();
395   fHistWOffsetSingleMuons -> Sumw2();
396
397   fHistOffsetSingleMuonsX_vsPtRapidity -> Sumw2();
398   fHistOffsetSingleMuonsY_vsPtRapidity -> Sumw2();
399   fHistSingleMuonsPtRapidity           -> Sumw2();
400   fHistSingleMuonsOffsetChi2           -> Sumw2();
401     
402   //--------------------------------------------
403
404   fGraphSingleMuonsOffsetChi2 = new TGraph("fGraphSingleMuonsOffsetChi2");
405   fGraphSingleMuonsOffsetChi2 -> SetName("fGraphSingleMuonsOffsetChi2");
406
407   //--------------------------------------------
408
409   fHistWOffsetMuonPairs        = new TH1D("fHistWOffsetMuonPairs",        "Weighted Offset for Muon Pairs", 300, 0, 60);
410   fHistMassMuonPairs           = new TH1D("fHistMassMuonPairs",           "Dimuon Mass (MUON+MFT)",  fNMassBins, fMassMin, fMassMax);
411   fHistMassMuonPairsWithoutMFT = new TH1D("fHistMassMuonPairsWithoutMFT", "Dimuon Mass (MUON only)", fNMassBins, fMassMin, fMassMax);
412   fHistMassMuonPairsMC         = new TH1D("fHistMassMuonPairsMC",         "Dimuon Mass (MC)",        fNMassBins, fMassMin, fMassMax);
413   fHistRapidityPtMuonPairsMC   = new TH2D("fHistRapidityPtMuonPairsMC",   "Dimuon Phase Space (MC)", 100, -4.5, -2., 100, 0., 10.); 
414
415   fHistWOffsetMuonPairs        -> SetXTitle("Weighted Offset");
416   fHistMassMuonPairs           -> SetXTitle("Mass  [GeV/c^{2}]");
417   fHistMassMuonPairsWithoutMFT -> SetXTitle("Mass  [GeV/c^{2}]");
418   fHistMassMuonPairsMC         -> SetXTitle("Mass  [GeV/c^{2}]");
419   fHistRapidityPtMuonPairsMC   -> SetXTitle("Rapidity");
420   fHistRapidityPtMuonPairsMC   -> SetYTitle("p_{T}  [GeV/c]");
421
422   fHistWOffsetMuonPairs        -> Sumw2();
423   fHistMassMuonPairs           -> Sumw2();
424   fHistMassMuonPairsWithoutMFT -> Sumw2();
425   fHistMassMuonPairsMC         -> Sumw2();
426   fHistRapidityPtMuonPairsMC   -> Sumw2();
427
428 }
429
430 //====================================================================================================================================================
431