]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/scripts/CompareMatDep.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / FMD / scripts / CompareMatDep.C
1 Int_t 
2 DetIndex(Int_t d, Char_t r)
3 {
4   return (d == 1 ? 0 : (d-1) * 2 - 1) + (r == 'I' ? 0 : 1);
5 }
6 Int_t
7 RingColor(Int_t d, Char_t r)
8 {
9   return (d == 1 ? kRed : d == 2 ? kGreen : kBlue) + (r == 'I' ? 1 : 3);
10 }
11
12 void CompareMatDep() {
13   gStyle->SetOptTitle(kFALSE);
14
15   Int_t spread = 10;
16   TMultiGraph* mg = new TMultiGraph;
17
18   for(Int_t d = 1; d<=3;d++) {
19     Int_t nr = (d == 1 ? 1 : 2);
20     for(Int_t q = 0; q < nr; q++) {
21       Char_t        r = (q == 0 ? 'I' : 'O');
22       TGraphErrors* g = new TGraphErrors();
23       g->SetName(Form("fmd%d%c", d, r));
24       g->SetTitle(Form("FMD%d%c (+%d)", d, r, spread*(DetIndex(d,r)+1)));
25       g->SetLineColor(RingColor(d, r));
26       g->SetMarkerColor(RingColor(d, r));
27       g->SetFillColor(RingColor(d, r));
28       g->SetFillStyle(3002);
29       g->SetMarkerStyle(20);
30       g->SetMarkerSize(1.3);
31       mg->Add(g);
32     }
33   }
34   
35   TGraphErrors* g = new TGraphErrors();
36   g->SetName("all");
37   g->SetTitle("All");
38   g->SetLineColor(kBlack);
39   g->SetMarkerColor(kBlack);
40   g->SetFillColor(kGray);
41   g->SetFillStyle(3002);
42   g->SetMarkerStyle(20);
43   g->SetMarkerSize(1.3);
44   mg->Add(g);
45
46   const char* base = "/mnt/samsung/oldhome/ALICE/FMDanalysis/BackgroundCorrection/BackgroundStudy/TestOfMaterialBudget/";
47
48   Int_t  samples[] = { -20, -15, -10, -5, 0, 5, 10 , 15, 20, 999 };
49   Int_t* sample    = samples;
50   Int_t  i         = 0;
51   while (*sample < 999) { 
52     TFile* fin = 0;
53     if(*sample < 0) 
54       fin = TFile::Open(Form("%s/minus%dpercent/fmd_dNdeta.root",
55                              base,-*sample));
56     else if(*sample > 0) 
57       fin = TFile::Open(Form("%s/plus%dpercent/fmd_dNdeta.root",base,*sample));
58     if(*sample == 0) 
59       fin = TFile::Open(Form("%s/refdist/fmd_dNdeta.root", base));
60     if (!fin) { 
61       Error("CompareMatDep", "Failed to open input file for %d", *sample);
62       sample++;
63       i++;
64       continue;
65     }
66     Float_t twa = 0;
67     Float_t tsw = 0;
68
69     for(Int_t d = 1; d<=3;d++) {
70       Int_t nr = (d == 1 ? 1 : 2);
71       for(Int_t q = 0; q < nr; q++) {
72         Char_t        r = (q == 0 ? 'I' : 'O');
73         TH1F* hRingMult = (TH1F*)fin->Get(Form("hRingMult_FMD%d%c",
74                                               d,r));
75         Float_t wav    = 0;
76         Float_t sumofw = 0;
77
78         // Weighted average 
79         for(Int_t n = 1;n<=hRingMult->GetNbinsX();n++) {
80           if(hRingMult->GetBinContent(n) > -1 && 
81              hRingMult->GetBinError(n) > 0 ) {
82             Float_t weight = 1/TMath::Power(hRingMult->GetBinError(n),2);
83             wav            += weight*hRingMult->GetBinContent(n);
84             sumofw         += weight;
85             twa            += weight*hRingMult->GetBinContent(n);
86             tsw             += weight;
87           }
88         }
89         
90         TGraphErrors* g = 
91           static_cast<TGraphErrors*>(mg->GetListOfGraphs()->At(DetIndex(d, r)));
92         g->SetPoint(i, *sample, 100*wav/sumofw+spread*(DetIndex(d,r)+1));
93         g->SetPointError(i, 0, 100*1/TMath::Sqrt(sumofw)); 
94         Info("CompareMatDep", "Setting point %d @ %d of %s to %f +/- %f", 
95              i, *sample, g->GetName(), 100*wav/sumofw, 
96              100*1/TMath::Sqrt(sumofw));
97       } // for q 
98     } // for d
99
100     TGraphErrors* g = 
101       static_cast<TGraphErrors*>(mg->GetListOfGraphs()->At(5));
102     g->SetPoint(i, *sample, 100*twa/tsw);
103     g->SetPointError(i, 0, 100*1/TMath::Sqrt(tsw)); 
104     Info("CompareMatDep", "Setting point %d @ %d of %s to %f +/- %f", 
105          i, *sample, g->GetName(), 100*twa/tsw, 
106          100*1/TMath::Sqrt(tsw));
107
108     sample++;
109     i++;
110   }
111   
112   TCanvas* c  = new TCanvas("c", "C", 800, 800);
113   c->SetTopMargin(0.05);
114   c->SetRightMargin(0.05);
115   c->SetLeftMargin(0.13);
116   c->SetFillColor(0);
117   c->SetBorderSize(0);
118   c->SetBorderMode(0);
119   c->SetGridx();
120   c->SetGridy();
121
122   mg->Draw("a pl");
123   mg->GetHistogram()->SetXTitle("#Delta#rho [%]");
124   mg->GetHistogram()->SetYTitle("Average relative difference [%]");
125   mg->GetHistogram()->GetXaxis()->SetTitleFont(132);
126   mg->GetHistogram()->GetXaxis()->SetLabelFont(132);
127   mg->GetHistogram()->GetXaxis()->SetNdivisions(10);
128   mg->GetHistogram()->GetYaxis()->SetTitleFont(132);
129   mg->GetHistogram()->GetYaxis()->SetLabelFont(132);
130   mg->GetHistogram()->GetYaxis()->SetNdivisions(1010);
131   mg->GetHistogram()->GetYaxis()->SetTitleOffset(1.2);
132   // c->Clear();
133   // mg->Draw("a e3");
134
135   TLegend* leg = gPad->BuildLegend(0.15,0.7,0.4,0.945);
136   leg->SetFillColor(0);
137   leg->SetBorderSize(0);
138   leg->SetTextFont(132);
139   leg->SetFillStyle(0);
140
141   TIter next(mg->GetListOfGraphs());
142   TObject* o = 0;
143
144   c->Print("materialDependence.png");
145 }