]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/MakeITSdEdxAnalysis.C
Fixed striong (Svein)
[u/mrichter/AliRoot.git] / ITS / MakeITSdEdxAnalysis.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TCanvas.h>
3 #include <TClassTable.h>
4 #include <TGraph.h>
5 #include <TTree.h>
6 #include <TH1.h>
7 #include <TFile.h>
8 #include <TH2.h>
9 #include <TInterpreter.h>
10 #include <TStyle.h>
11 #include "AliHeader.h"
12 #include "AliITSdEdxAnalyzer.h"
13 #include "AliRunLoader.h"
14 #include "AliESDEvent.h"
15 #endif
16
17 Bool_t MakeITSdEdxAnalysis(TString path=".",Int_t nmombins=16){
18
19   if (gClassTable->GetID("AliRun") < 0) {
20     gInterpreter->ExecuteMacro("loadlibs.C");
21   }
22  
23   TString galicefile = path+"/galice.root";
24   TString esdfile = path+"/AliESDs.root";
25
26   AliRunLoader* runLoader = AliRunLoader::Open(galicefile.Data());
27   if (!runLoader) {
28     printf("Error in getting run loader");
29     return kFALSE;
30   }
31   runLoader->LoadHeader();
32   runLoader->LoadKinematics();
33
34   TFile* esdFile = TFile::Open(esdfile.Data());
35   if (!esdFile || !esdFile->IsOpen()) {
36     printf("Error in opening ESD file");
37     return kFALSE;
38   }
39
40   AliESDEvent * esd = new AliESDEvent;
41   TTree* tree = (TTree*) esdFile->Get("esdTree");
42   if (!tree) {
43     printf("Error: no ESD tree found");
44     return kFALSE;
45   }
46   esd->ReadFromTree(tree);
47   AliITSdEdxAnalyzer* enan=new AliITSdEdxAnalyzer(nmombins,0.1,3.1);
48   enan->SetMIPdEdx(79.);
49
50   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
51     runLoader->GetEvent(iEvent);
52     tree->GetEvent(iEvent);
53     if (!esd) {
54       printf("Error: no ESD object found for event %d", iEvent);
55       return kFALSE;
56     }
57     AliStack* stack = runLoader->Stack();
58     enan->ReadEvent(esd,stack);
59   }
60   enan->WriteHistos();
61
62
63
64   TCanvas* c1=new TCanvas("c1","dEdx - pions - layer 4");
65   c1->Divide(4,4);
66   for(Int_t i=0; i<nmombins;i++){
67     TH1F* h=enan->GetSingleLayerdEdxHisto(4,211,i);
68     c1->cd(i+1);
69     h->Draw();
70   }
71
72
73   TCanvas* c1b=new TCanvas("c1b","dEdx - pions - Truncated mean");
74   c1b->Divide(4,4);
75   for(Int_t i=0; i<nmombins;i++){
76     TH1F* h=enan->GetTruncatedMeandEdxHisto(211,i);
77     c1b->cd(i+1);
78     h->Draw();
79   }
80
81   TCanvas* c2=new TCanvas("c2","Delta dEdx - pions - layer 4");
82   c2->Divide(4,4);
83   for(Int_t i=0; i<nmombins;i++){
84     TH1F* h=enan->GetSingleLayerDeltadEdxHisto(4,211,i);
85     c2->cd(i+1);
86     h->Draw();
87   }
88
89   TCanvas* c2b=new TCanvas("c2b","Delta dEdx - pions - Truncated Mean");
90   c2b->Divide(4,4);
91   for(Int_t i=0; i<nmombins;i++){
92     TH1F* h=enan->GetTruncatedMeanDeltadEdxHisto(211,i);
93     c2b->cd(i+1);
94     h->Draw();
95   }
96
97   gStyle->SetPalette(1);
98
99   TCanvas* c3=0;
100   c3=new TCanvas("c3","dEdx vs. p - pions - layer 4");
101   TH2F* h2=enan->GetSingleLayerdEdxVsPHisto(4,211);
102   h2->Draw("colz");
103   enan->SetUseBBFromAliExternalTrackParam();
104   TGraph* gpion1=enan->GetBetheBlochGraph(211);
105   gpion1->Draw("LSAME");
106   enan->SetUseBBFromAliITSpidESD();
107   TGraph* gpion2=enan->GetBetheBlochGraph(211);
108   gpion2->SetLineColor(2);
109   gpion2->Draw("LSAME");
110
111   TCanvas* c3b=0;
112   c3b=new TCanvas("c3b","dEdx vs. p - pions - Truncated Mean");
113   TH2F* h2b=enan->GetTruncatedMeandEdxVsPHisto(211);
114   h2b->Draw("colz");
115   gpion1->Draw("LSAME");
116   gpion2->Draw("LSAME");
117
118   return kTRUE;
119 }