]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/Cal/AliTRDplotNoiseBaseline.C
end-of-line normalization
[u/mrichter/AliRoot.git] / TRD / Cal / AliTRDplotNoiseBaseline.C
CommitLineData
a65a7e70 1#if !defined( __CINT__) || defined(__MAKECINT__)
2
3
4#include <Riostream.h>
5#include <TSystem.h>
6#include <TProfile2D.h>
7#include <TCanvas.h>
8#include <TH1F.h>
9#include <TH2I.h>
10#include <TStyle.h>
11#include <TFile.h>
12#include <TGrid.h>
13
14
15#include "AliCDBManager.h"
16#include "AliCDBStorage.h"
17#include "AliCDBEntry.h"
18
19
20#include "../TRD/AliTRDarrayF.h"
21#include "../TRD/AliTRDCalibPadStatus.h"
22#include "../TRD/Cal/AliTRDCalPadStatus.h"
23#include "../TRD/Cal/AliTRDCalDet.h"
24#include "../TRD/Cal/AliTRDCalPad.h"
25#include "../TRD/Cal/AliTRDCalROC.h"
26#include "../TRD/AliTRDcalibDB.h"
27
28
29#endif
30
31
32//void PlotNoiseBaseline(Int_t run, Int_t sm, Int_t det, const char * pathdatabase="local:///d/alice12/bailhache/TestShuttle/database/", const char * pathreferencefile="local:///d/alice12/bailhache/TestShuttle/reference")
33//void PlotNoiseBaseline(Int_t run=34529, Int_t sm=0, Int_t det=0, const char * pathdatabase="alien://Folder=/alice/data/2008/LHC08b/OCDB/", const char * pathreferencedatabase="alien://Folder=/alice/data/2008/LHC08b/Reference/")
34//void PlotNoiseBaseline(Int_t run=1, Int_t sm=0, Int_t det=0, const char * pathdatabase="local:///d/alice12/bailhache/AliAnalysisTask/v4-13-Head/SHUTTLE/TestShuttle/TestCDB/", const char * pathreferencedatabase="local:///d/alice12/bailhache/AliAnalysisTask/v4-13-Head/SHUTTLE/TestShuttle/TestReference/")
35void AliTRDplotNoiseBaseline(Int_t run=34529, Int_t sm=0, Int_t det=0, const char * pathdatabase="alien://Folder=/alice/data/2008/LHC08b/OCDB/", const char * pathreferencedatabase="alien://Folder=/alice/data/2008/LHC08b/Reference/")
36{
37
38 //TGrid::Connect("alien://",0,0,"t");
39
40 AliCDBManager *CDB = AliCDBManager::Instance();
41 CDB->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
42 CDB->SetSpecificStorage("TRD/Calib/PadNoise",pathdatabase);
43 CDB->SetSpecificStorage("TRD/Calib/DetNoise",pathdatabase);
44 CDB->SetSpecificStorage("TRD/Calib/PadStatus",pathdatabase);
45 CDB->SetRun(run);
46
47 AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
48
49 //const AliTRDCalDet *u = cal->GetNoiseDet();
50
51 AliTRDCalDet *u = new AliTRDCalDet("u","u");
52 for(Int_t k = 0; k < 540; k++){
53 u->SetValue(k,10.0);
54 }
55
56 //Style
57 //************************
58 gStyle->SetPalette(1);
59 gStyle->SetOptStat(1111);
60 gStyle->SetPadBorderMode(0);
61 gStyle->SetCanvasColor(10);
62 gStyle->SetPadLeftMargin(0.13);
63 gStyle->SetPadRightMargin(0.13);
64
65 //Build the Cal Pad
66 //********************************
67 Int_t smi = sm*30;
68 AliTRDCalPad *ki = new AliTRDCalPad("testnoise","testnoise");
69 for(Int_t k = 0; k < 540; k++){
70 ki->SetCalROC(k,(AliTRDCalROC *) cal->GetNoiseROC(k));
71 }
72
73 // padstatus 2D
74 Int_t smn = (Int_t) det/30;
75 if((smn==0) || (smn==1) || (smn==2) || (smn==9) || (smn==10) || (smn==11)) smn = 1;
76 if((smn==3) || (smn==4) || (smn==5) || (smn==12) || (smn==13) || (smn==14)) smn = 2;
77 if((smn==6) || (smn==7) || (smn==8) || (smn==15) || (smn==16) || (smn==17)) smn = 3;
78 TString name("TRD/DAQData/PadStatus");
79 name += smn;
80 //name += 3;
81 AliCDBEntry *entrypadstatus = AliCDBManager::Instance()->Get("TRD/Calib/PadStatus",run);
82 if(!entrypadstatus) return;
83 AliTRDCalPadStatus *lo = (AliTRDCalPadStatus *)entrypadstatus->GetObject();
84 AliCDBEntry *entryo = AliCDBManager::Instance()->GetStorage(pathreferencedatabase)->Get(name, run);
85 if(!entryo) return;
86 AliTRDCalibPadStatus *calpad = (AliTRDCalibPadStatus *) entryo->GetObject();
87 if(!calpad) return;
88
89
90 // Plot
91 //***********
92
93
94 // noise 2D
95 TCanvas *cnoise = new TCanvas((const char*)"noise1",(const char*)"noise1",50,50,600,800);
96 cnoise->Divide(3,2);
97 cnoise->cd(1);
98 ((TH2F *)ki->MakeHisto2DSmPl(sm,0,u,0,0.0,3.5,-1))->Draw("colz");
99 cnoise->cd(2);
100 ((TH2F *)ki->MakeHisto2DSmPl(sm,1,u,0,0.0,3.5,-1))->Draw("colz");
101 cnoise->cd(3);
102 ((TH2F *)ki->MakeHisto2DSmPl(sm,2,u,0,0.0,3.5,-1))->Draw("colz");
103 cnoise->cd(4);
104 ((TH2F *)ki->MakeHisto2DSmPl(sm,3,u,0,0.0,3.5,-1))->Draw("colz");
105 cnoise->cd(5);
106 ((TH2F *)ki->MakeHisto2DSmPl(sm,4,u,0,0.0,3.5,-1))->Draw("colz");
107 cnoise->cd(6);
108 ((TH2F *)ki->MakeHisto2DSmPl(sm,5,u,0,0.0,3.5,-1))->Draw("colz");
109
110
111 // Pad Status
112 TCanvas *cpadstatus = new TCanvas((const char*)"padstatus",(const char*)"padstatus",50,50,600,800);
113 cpadstatus->Divide(3,2);
114 cpadstatus->cd(1);
115 ((TH2F *)lo->MakeHisto2DSmPl(sm,0))->Draw("colz");
116 cpadstatus->cd(2);
117 ((TH2F *)lo->MakeHisto2DSmPl(sm,1))->Draw("colz");
118 cpadstatus->cd(3);
119 ((TH2F *)lo->MakeHisto2DSmPl(sm,2))->Draw("colz");
120 cpadstatus->cd(4);
121 ((TH2F *)lo->MakeHisto2DSmPl(sm,3))->Draw("colz");
122 cpadstatus->cd(5);
123 ((TH2F *)lo->MakeHisto2DSmPl(sm,4))->Draw("colz");
124 cpadstatus->cd(6);
125 ((TH2F *)lo->MakeHisto2DSmPl(sm,5))->Draw("colz");
126
127
128
129 // reference data
130
131 TCanvas *cpoui = new TCanvas((const char*)"cpoui",(const char*)"cpoui",50,50,600,800);
132 cpoui->cd();
133 ((TH2F *)calpad->GetHisto(det))->Draw("lego");
134
135
136 AliTRDCalROC *ouip = calpad->GetCalRocMean(det);
137 TCanvas *cpouilo = new TCanvas((const char*)"cpouilo",(const char*)"cpouilo",50,50,600,800);
138 cpouilo->Divide(2,1);
139 cpouilo->cd(1);
140 ((TH1F *)ouip->MakeHisto1D(8.5,10.5,-1,10.0))->Draw();
141 //((TH1F *)ouip->MakeHisto1D(0.85,1.05,-1))->Draw();
142 cpouilo->cd(2);
143 ((TH2F *)ouip->MakeHisto2D(8.5,10.5,-1,10.0))->Draw("colz");
144 //((TH2F *)ouip->MakeHisto2D(0.85,1.05,-1))->Draw("colz");
145
146 AliTRDCalROC *ouiphy = calpad->GetCalRocRMS(det);
147 TCanvas *cpouiloh = new TCanvas((const char*)"cpouiloh",(const char*)"cpouiloh",50,50,600,800);
148 cpouiloh->Divide(2,1);
149 cpouiloh->cd(1);
150 ((TH1F *)ouiphy->MakeHisto1D(0.1,4.5,-1,10.0))->Draw();
151 //((TH1F *)ouiphy->MakeHisto1D(0.01,0.45,-1))->Draw();
152 cpouiloh->cd(2);
153 ((TH2F *)ouiphy->MakeHisto2D(0.1,4.5,-1,10.0))->Draw("colz");
154 //((TH2F *)ouiphy->MakeHisto2D(0.01,0.45,-1))->Draw("colz");
155
156
157
158
159}