]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/FindKrClustersRaw.C
Using custom streamer to guarantie back compatibility (Peter Christiansen)
[u/mrichter/AliRoot.git] / TPC / FindKrClustersRaw.C
1 //\r
2 \r
3 Int_t FindKrClusterCheck(const char *fileName="data.root");\r
4 \r
5 \r
6 Int_t FindKrClustersRaw(const char *fileName="data.root"){\r
7 \r
8   \r
9 \r
10   //\r
11   // remove Altro warnings\r
12   //\r
13   AliLog::SetClassDebugLevel("AliTPCRawStream",-5);\r
14   AliLog::SetClassDebugLevel("AliRawReaderDate",-5);\r
15   AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);\r
16   AliLog::SetModuleDebugLevel("RAW",-5);\r
17   //\r
18   // Get calibration\r
19   //\r
20   char *ocdbpath = gSystem->Getenv("OCDB_PATH");\r
21   //char *ocdbpath ="local:///afs/cern.ch/alice/tpctest/OCDB";\r
22   if (ocdbpath==0){\r
23     ocdbpath="alien://folder=/alice/data/2007/LHC07w/OCDB/";\r
24   }\r
25   printf("OCDB PATH = %s\n",ocdbpath); \r
26   AliCDBManager * man = AliCDBManager::Instance();\r
27   man->SetDefaultStorage(ocdbpath);\r
28   man->SetRun(100000000);\r
29 \r
30   AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();\r
31   AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();\r
32   //\r
33   //\r
34 \r
35 \r
36   //define tree\r
37   TFile *hfile=new TFile("adc.root","RECREATE","ADC file");\r
38   // Create a ROOT Tree\r
39   TTree *mytree = new TTree("Kr","Krypton cluster tree");\r
40 \r
41 \r
42   Int_t debugLevel=1;\r
43   if(debugLevel>0){\r
44     TH1F *histoRow   =new TH1F("histoRow","rows",100,0.,100.);\r
45     TH1F *histoPad   =new TH1F("histoPad","pads",150,0.,150.);\r
46     TH1F *histoTime  =new TH1F("histoTime","timebins",100,0.,1000.);\r
47     TH2F *histoRowPad=new TH2F("histoRowPad","pads-vs-rows" ,150,0.,150.,100,0.,100.);\r
48   }\r
49 \r
50 \r
51   //one general output\r
52   AliTPCclustererKr *clusters = new AliTPCclustererKr();\r
53   clusters->SetOutput(mytree);\r
54   clusters->SetRecoParam(0);\r
55 \r
56   if(debugLevel>0){\r
57     clusters->SetDebugLevel(debugLevel);\r
58     clusters->SetHistoRow(histoRow   );\r
59     clusters->SetHistoPad(histoPad   );\r
60     clusters->SetHistoTime(histoTime  );\r
61     clusters->SetHistoRowPad(histoRowPad);\r
62   }\r
63 \r
64 \r
65   AliTPCParamSR *param=new AliTPCParamSR();\r
66   //only for geometry parameters loading - temporarly\r
67 //  AliRunLoader* rl = AliRunLoader::Open("galice.root");\r
68 //  AliTPCParam *param=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");\r
69   //if (!param) {cerr<<"TPC parameters have not been found !\n"; return 4;}\r
70 \r
71   clusters->SetParam(param);\r
72 \r
73   //set cluster finder parameters (from data)\r
74   clusters->SetZeroSup(param->GetZeroSup());//zero suppression parameter\r
75   clusters->SetFirstBin(60);//first bin\r
76   clusters->SetLastBin(950);//last bin\r
77   clusters->SetMaxNoiseAbs(2);//maximal noise\r
78   clusters->SetMaxNoiseSigma(3);//maximal amount of sigma of noise\r
79 \r
80   //set cluster finder parameters (from MC)\r
81   clusters->SetMinAdc(3);//signal threshold (everything below is treated as 0)\r
82   clusters->SetMinTimeBins(2);//number of neighbouring timebins\r
83   clusters->SetMaxPadRangeCm(5.);//distance of the cluster center to the center of a pad (in cm)\r
84   clusters->SetMaxRowRangeCm(5.);//distance of the cluster center to the center of a padrow (in cm)\r
85   clusters->SetMaxTimeRange(5.);//distance of the cluster center to the max time bin on a pad (in tackts)\r
86   //ie. fabs(centerT - time)<7\r
87   clusters->SetValueToSize(7.);//cut reduce peak at 0\r
88   clusters->SetIsolCut(3);//set isolation cut threshold\r
89 \r
90   AliRawReader *reader = new AliRawReaderRoot(fileName);\r
91   reader->Reset();\r
92 \r
93   TStopwatch timer;\r
94   timer.Start();\r
95 \r
96   AliAltroRawStreamFast* stream = new AliAltroRawStreamFast(reader);\r
97   stream->SelectRawData("TPC");\r
98 \r
99   Int_t evtnr=0;\r
100   while (reader->NextEvent()) {\r
101     //output for each event\r
102 \r
103     //if(evtnr>4) break;\r
104     cout<<"Evt = "<<evtnr<<endl;\r
105     clusters->FinderIO(reader);\r
106     evtnr++;\r
107     AliSysInfo::AddStamp(Form("Event%d",evtnr),evtnr);\r
108   }\r
109 \r
110 \r
111   mytree->Print();//print rootuple summary \r
112   // Save all objects in this file\r
113   hfile->Write();\r
114   // Close the file\r
115   hfile->Close();\r
116 \r
117   timer.Stop();\r
118   timer.Print();\r
119   printf("Deleting clusterer\n");\r
120   delete clusters;\r
121   printf("Deleting stream\n");\r
122   delete stream;\r
123   printf("Deleting raw reader\n");\r
124   delete reader;\r
125 \r
126 //   TCanvas *c2=new TCanvas("c2","title",800,800);\r
127 //   c2->SetHighLightColor(2);\r
128 //   c2->Range(-458.9552,-2948.238,3296.642,26856.6);\r
129 //   c2->SetBorderSize(2);\r
130 //   c2->SetLeftMargin(0.15);\r
131 //   c2->SetRightMargin(0.06);\r
132 //   c2->SetFrameFillColor(0);\r
133 \r
134 //   gStyle->SetOptStat(111111);\r
135 //   histoRow->Draw();\r
136 //   c2->Print("rows.ps");\r
137 //   histoPad->Draw();\r
138 //   c2->Print("pads.ps");\r
139 //   histoTime->Draw();\r
140 //   c2->Print("timebins.ps");\r
141 //   histoRowPad->Draw();\r
142 //   c2->Print("row-pad.ps");\r
143 \r
144   return 0;\r
145 }\r
146 \r
147 \r
148 Int_t FindKrClusterCheck(const char *fileName){\r
149   //\r
150   //\r
151   gSystem->Load("$ROOTSYS/lib/libGui.so");\r
152   gSystem->Load("$ROOTSYS/lib/libTree.so");\r
153   gSystem->Load("$MEMSTAT/libMemStat.so");\r
154   {\r
155     TMemStat memstat(1000000,100000,kTRUE);\r
156     AliSysInfo::AddCallBack(TMemStatManager::GetInstance()->fStampCallBack);\r
157     FindKrClustersRaw(fileName); \r
158   }\r
159   // the output memstat.root file\r
160   TMemStat draw("memstat.root");\r
161   // Print some information\r
162   // code info\r
163   draw.MakeReport(0,0);\r
164 \r
165 }\r