]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/CalibMacros/CalibCosmic.C
forgot file in previous commit
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibCosmic.C
CommitLineData
91fd44c9 1/*
35a255d8 2 .x ~/NimStyle.C
91fd44c9 3 .x ~/UliStyle.C
4 gSystem->Load("libANALYSIS");
5 gSystem->Load("libTPCcalib");
6 .L $ALICE_ROOT/TPC/CalibMacros/CalibCosmic.C
7 // init
8 Init();
9 SetDefaultCut();
10
11*/
35a255d8 12
91fd44c9 13AliTPCcalibCosmic * cosmic =0;
14TObjArray fitArr;
15
16void Init(){
17 //
18 //
19 TFile fcalib("CalibObjects.root");
20 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
21 cosmic = ( AliTPCcalibCosmic *)array->FindObject("cosmicTPC");
22 TString axisName[9];
23 axisName[0] ="#Delta"; axisName[1] ="N_{cl}";
35a255d8 24 axisName[2] ="DCA_{r}(cm)";
91fd44c9 25 axisName[3] ="z (cm)"; axisName[4] ="sin(#phi)";
35a255d8 26 axisName[5] ="tan(#theta)"; axisName[6] ="1/p_{t} (1/GeV)";
27 axisName[7] ="p_{t} (GeV)"; axisName[8] ="alpha";
91fd44c9 28
29 {
30 for (Int_t ivar=0;ivar<6;ivar++){
31 for (Int_t ivar2=0;ivar2<9;ivar2++){
32 cosmic->fHistoDelta[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
33 cosmic->fHistoDelta[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
34 cosmic->fHistoPull[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
35 cosmic->fHistoPull[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
36 }
37 }
38 }
39}
40
35a255d8 41void SetRangeAll(Int_t axis, Float_t xmin, Float_t xmax){
42
91fd44c9 43 for (Int_t i=0;i<6;i++){
44 //
35a255d8 45 cosmic->fHistoDelta[i]->GetAxis(axis)->SetRangeUser(xmin,xmax);
46 cosmic->fHistoPull[i]->GetAxis(axis)->SetRangeUser(xmin,xmax);
47 }
48}
49
91fd44c9 50
35a255d8 51void SetDefaultCut(){
52 for (Int_t i=0;i<6;i++){
53 //
54 cosmic->fHistoDelta[i]->GetAxis(1)->SetRangeUser(120,200);
55 cosmic->fHistoPull[i]->GetAxis(1)->SetRangeUser(120,200);
56 cosmic->fHistoDelta[i]->GetAxis(3)->SetRangeUser(-0,250);
57 cosmic->fHistoPull[i]->GetAxis(3)->SetRangeUser(-0,250);
58 cosmic->fHistoDelta[i]->GetAxis(2)->SetRangeUser(0,60);
59 cosmic->fHistoPull[i]->GetAxis(2)->SetRangeUser(0,60);
91fd44c9 60 }
61}
62
63TH2 * GetDelta2D(Int_t type, Int_t var){
64 TH2 * his = cosmic->fHistoDelta[type]->Projection(0,var);
65 his->SetXTitle(cosmic->fHistoDelta[type]->GetAxis(var)->GetName());
66 his->SetYTitle(cosmic->fHistoDelta[type]->GetAxis(0)->GetName());
67 return his;
68}
69
70
71TH1* GetResol2DSigma(Int_t type, Int_t var){
72
73 TH2 * his = cosmic->fHistoDelta[type]->Projection(0,var);
74 his->SetXTitle(cosmic->fHistoDelta[type]->GetAxis(var)->GetName());
75 his->SetYTitle(cosmic->fHistoDelta[type]->GetAxis(0)->GetName());
76 his->FitSlicesY(0,0,-1,0,"QNR",&fitArr);
77 TH1 * hres = (TH1*)(fitArr.At(2)->Clone());
78 return hres;
79}
80
81
82TH1 * GetDelta(Int_t type){
83 TH1 * his = cosmic->fHistoDelta[type]->Projection(0);
84 his->SetXTitle(cosmic->fHistoDelta[type]->GetAxis(0)->GetName());
85 return his;
86}
87
88TH2 * GetPull2D(Int_t type, Int_t var){
89 TH2 * his = cosmic->fHistoPull[type]->Projection(0,var);
90 his->SetXTitle(cosmic->fHistoPull[type]->GetAxis(var)->GetName());
91 his->SetYTitle(cosmic->fHistoPull[type]->GetAxis(0)->GetName());
92 return his;
93}
94
95TH1* GetPull2DSigma(Int_t type, Int_t var){
96
97 TH2 * his = cosmic->fHistoPull[type]->Projection(0,var);
98 his->SetXTitle(cosmic->fHistoPull[type]->GetAxis(var)->GetName());
99 his->SetYTitle(cosmic->fHistoPull[type]->GetAxis(0)->GetName());
100 his->FitSlicesY(0,0,-1,0,"QNR",&fitArr);
101 TH1 * hres = (TH1*)(fitArr.At(2)->Clone());
102 return hres;
103}
104
105
106
107TH1 * GetPull(Int_t type){
108 TH1 * his = cosmic->fHistoPull[type]->Projection(0);
109 his->SetXTitle(cosmic->fHistoPull[type]->GetAxis(0)->GetName());
110 return his;
111}
35a255d8 112
113
114void DrawResoldEdx(AliTPCcalibCosmic *cosmic){
115 //
116 //
117 //
118 Int_t kmicolors[10]={1,2,3,6,7,8,9,10,11,12};
119 Int_t kmimarkers[10]={21,22,23,24,25,26,27,28,29,30};
120 TH2 *htemp;
121 TObjArray arr;
122 TH1 * hResolMax[4];
123 TH1 * hResolTot[4];
124 //
125 for (Int_t ipad=0;ipad<4;ipad++){
126 cosmic->fHistodEdxTot[ipad]->GetAxis(4)->SetRangeUser(-0.6,0.6);
127 cosmic->fHistodEdxMax[ipad]->GetAxis(4)->SetRangeUser(-0.6,0.6);
128 }
129 cosmic->fHistodEdxTot[0]->GetAxis(1)->SetRangeUser(30,62);
130 cosmic->fHistodEdxTot[1]->GetAxis(1)->SetRangeUser(30,62);
131 cosmic->fHistodEdxTot[2]->GetAxis(1)->SetRangeUser(10,35);
132 cosmic->fHistodEdxTot[3]->GetAxis(1)->SetRangeUser(10,150);
133 cosmic->fHistodEdxMax[0]->GetAxis(1)->SetRangeUser(30,62);
134 cosmic->fHistodEdxMax[1]->GetAxis(1)->SetRangeUser(30,62);
135 cosmic->fHistodEdxMax[2]->GetAxis(1)->SetRangeUser(10,35);
136 cosmic->fHistodEdxMax[3]->GetAxis(1)->SetRangeUser(10,150);
137 //
138
139 for (Int_t ipad=0;ipad<4;ipad++){
140 htemp = cosmic->fHistodEdxTot[ipad]->Projection(0,1);
141 htemp->FitSlicesY(0,0,-1,0,"QNR",&arr);
142 hResolTot[ipad] = (TH1*)(arr.At(2)->Clone());
143 delete htemp;
144 arr.SetOwner(kTRUE);
145 arr.Delete();
146 hResolTot[ipad]->Scale(1./TMath::Sqrt(2.));
147 //
148 htemp = cosmic->fHistodEdxMax[ipad]->Projection(0,1);
149 htemp->FitSlicesY(0,0,-1,0,"QNR",&arr);
150 hResolMax[ipad] = (TH1*)(arr.At(2)->Clone());
151 delete htemp;
152 arr.SetOwner(kTRUE);
153 arr.Delete();
154 hResolMax[ipad]->Scale(1./TMath::Sqrt(2.));
155 }
156 hResolTot[3]->GetXaxis()->SetRangeUser(0,160);
157 hResolTot[3]->SetXTitle("N_{cl}");
158 hResolTot[3]->SetYTitle("#sigma(dEdx/dEdx_{d})/#sqrt{2.}");
159 hResolTot[3]->SetTitle("Relative dEdx resolution");
160 for (Int_t ipad=3;ipad>=0;ipad--){
161 hResolTot[ipad]->SetMaximum(0.1);
162 hResolTot[ipad]->SetMinimum(0.);
163 hResolTot[ipad]->SetMarkerColor(kmicolors[ipad]+0);
164 hResolTot[ipad]->SetMarkerStyle(kmimarkers[ipad]+1);
165 if (ipad==3) hResolTot[ipad]->Draw();
166 hResolTot[ipad]->Draw("same");
167 //
168 hResolMax[ipad]->SetMaximum(0.1);
169 hResolMax[ipad]->SetMinimum(0.);
170 hResolMax[ipad]->SetMarkerColor(kmicolors[ipad]+0);
171 hResolMax[ipad]->SetMarkerStyle(kmimarkers[ipad]+4);
172 hResolMax[ipad]->Draw("same");
173 }
174
175}