]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/PlotTrackToPointResidualsSDD.C
Update master to aliroot
[u/mrichter/AliRoot.git] / ITS / macrosSDD / PlotTrackToPointResidualsSDD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TCanvas.h>
3 #include <TFile.h>
4 #include <TStyle.h>
5 #include <TProfile.h>
6 #include <TH1F.h>
7 #include <TH1D.h>
8 #include <TH2F.h>
9 #include <TPaveStats.h>
10 #include <TGraphErrors.h>
11 #include <TLatex.h>
12 #include <TString.h>
13 #endif
14
15
16 // Macro to plot the output of the ITS aligment task for SDD detector
17 // For the moment just plot quantities summed over modules
18 // Orgin: F. Prino
19
20
21 void PlotTrackToPointResidualsSDD(){
22   TFile *fil=new TFile("AnalysisResults.root");
23   TDirectoryFile* df=(TDirectoryFile*)fil->Get("ITSAlignQA");
24   TList* l=(TList*)df->Get("clistITSAlignQA");
25   TH2F* hSDDResidXvsX=0x0;
26   TH2F* hSDDResidX=0x0;
27   for(Int_t iMod=240; iMod<500; iMod++){
28     TString hname=Form("hSDDResidXvsX%d",iMod);
29     TH2F* h=(TH2F*)l->FindObject(hname);
30     if(h){
31       if(hSDDResidXvsX==0x0){
32         hSDDResidXvsX=(TH2F*)h->Clone("hSDDResidXvsXAll");
33       }else{
34         hSDDResidXvsX->Add(h);
35       }
36     }
37     TString hname2=Form("hSDDResidX%d",iMod);
38     TH2F* h2=(TH2F*)l->FindObject(hname2);
39     if(h2){
40       if(hSDDResidX==0x0){
41         hSDDResidX=(TH2F*)h2->Clone("hSDDResidXvsPtAll");
42       }else{
43         hSDDResidX->Add(h2);
44       }
45     }
46   }
47
48   gStyle->SetPalette(1);
49   gStyle->SetOptTitle(0);
50   TCanvas* c1=new TCanvas("c1","",700,1000);
51   c1->Divide(1,2);
52   c1->cd(1);
53   hSDDResidXvsX->Draw("colz");
54   hSDDResidXvsX->GetYaxis()->SetTitle("Track-to-point residual (cm)");
55   hSDDResidXvsX->GetXaxis()->SetTitle("X local (cm)");
56   TProfile* hpfx=hSDDResidXvsX->ProfileX();
57   hpfx->SetLineWidth(2);
58   hpfx->Draw("same");
59   c1->cd(2);
60   Int_t midbin=hSDDResidXvsX->GetXaxis()->FindBin(0.);
61   TH1D* hleft=hSDDResidXvsX->ProjectionY("hleft",1,midbin-2);
62   TH1D* hright=hSDDResidXvsX->ProjectionY("hright",midbin+2,hSDDResidXvsX->GetNbinsX());
63   hleft->Draw();
64   hleft->GetYaxis()->SetTitle("Track-to-point residual (cm)");
65   c1->Update();
66   TPaveStats *st1=(TPaveStats*)hleft->GetListOfFunctions()->FindObject("stats");
67   st1->SetY1NDC(0.71);
68   st1->SetY2NDC(0.9);
69   hright->SetLineColor(2);
70   hright->Draw("sames");
71   c1->Update();
72   TPaveStats *st2=(TPaveStats*)hright->GetListOfFunctions()->FindObject("stats");
73   st2->SetY1NDC(0.51);
74   st2->SetY2NDC(0.7);
75   st2->SetTextColor(2);
76   c1->Modified();
77   c1->Update();
78   TLatex *t2=new TLatex(0.14,0.7,"Left sides");
79   t2->SetTextSize(0.048);
80   t2->SetNDC();
81   TLatex *t3=new TLatex(0.14,0.5,"Right sides");
82   t3->SetNDC();
83   t3->SetTextSize(0.048);
84   t3->SetTextColor(2);
85   t2->Draw();
86   t3->Draw();
87
88   TCanvas* c2=new TCanvas("c2","",700,1000);
89   c2->Divide(1,2);
90   c2->cd(1);
91   hSDDResidX->GetXaxis()->SetRangeUser(0.,10.);
92   hSDDResidX->Draw("colz");
93   hSDDResidX->GetXaxis()->SetTitle("p_{t} (GeV/c)");
94   hSDDResidX->GetYaxis()->SetTitle("Track-to-point residual (cm)");
95   c2->cd(2);
96   TGraphErrors* grms=new TGraphErrors(0);
97   for(Int_t iBin=1; iBin<=hSDDResidX->GetNbinsX(); iBin++){
98     TH1D* htmp=hSDDResidX->ProjectionY("htmp",iBin,iBin);
99     if(htmp->GetEntries()>10.){
100       Int_t npt=grms->GetN();
101       grms->SetPoint(npt,hSDDResidX->GetXaxis()->GetBinCenter(iBin),htmp->GetRMS());
102       grms->SetPointError(npt,0.5*hSDDResidX->GetXaxis()->GetBinWidth(iBin),htmp->GetRMSError());
103     }
104   }
105   grms->SetMarkerStyle(20);
106   grms->GetXaxis()->SetLimits(0.,10.);
107   grms->Draw("AP");
108   grms->GetXaxis()->SetTitle("p_{t} (GeV/c)");
109   grms->GetYaxis()->SetTitle("RMS of Track-to-point residual (cm)");
110 }