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