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