]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCHits2Digits.C
Coding convention rules obeyed
[u/mrichter/AliRoot.git] / TPC / AliTPCHits2Digits.C
1 Int_t AliTPCHits2Digits()
2 {
3
4   // new version by J.Belikov
5
6   // Connect the Root Galice file containing Geometry, Kine and Hits
7
8   const char * inFile = "galice.root";  
9   TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(inFile);
10   if (file) {file->Close(); delete file;}
11   file = new TFile(inFile,"UPDATE");
12   if (!file->IsOpen()) {
13     cerr<<"Can't open "<<inFile<<" !\n";
14     return 1;
15   }
16
17   // Get AliRun object from file or create it if not on file
18   if (gAlice) delete gAlice;
19   gAlice = (AliRun*)file->Get("gAlice");
20   if (!gAlice) {
21     cerr<<"AliTPCHits2Digits.C : AliRun object not found on file\n";
22     return 2;
23   }
24
25   gAlice->GetEvent(0);
26   AliTPC *TPC = (AliTPC*)gAlice->GetDetector("TPC");      
27
28 //Set response functions
29
30   AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
31   AliTPCPRF2D    * prfinner   = new AliTPCPRF2D;
32   AliTPCPRF2D    * prfouter   = new AliTPCPRF2D;
33   AliTPCRF1D     * rf    = new AliTPCRF1D(kTRUE);
34   rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
35   rf->SetOffset(3*param->GetZSigma());
36   rf->Update();
37
38   TDirectory *savedir=gDirectory;
39   TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
40   if (!f->IsOpen()) { 
41      cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n"
42      return 3;
43   }
44   prfinner->Read("prf_07504_Gati_056068_d02");
45   prfouter->Read("prf_10006_Gati_047051_d03");
46   f->Close();
47   savedir->cd();
48
49   param->SetInnerPRF(prfinner);
50   param->SetOuterPRF(prfouter); 
51   param->SetTimeRF(rf);
52   TPC->SetParam(param);
53    
54   cerr<<"Digitizing TPC...\n";
55
56   //setup TPCDigitsArray 
57   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
58   arr->SetClass("AliSimDigits");
59   arr->Setup(param);
60   TPC->SetParam(param);
61   arr->MakeTree();
62
63   TPC->SetDigitsArray(arr);
64   TPC->Hits2DigitsSector(1);             
65   TPC->Hits2DigitsSector(2);             
66   TPC->Hits2DigitsSector(3);             
67   TPC->Hits2DigitsSector(1+18);             
68   TPC->Hits2DigitsSector(2+18);             
69   TPC->Hits2DigitsSector(3+18);             
70
71   TPC->Hits2DigitsSector(36+1);             
72   TPC->Hits2DigitsSector(36+2);             
73   TPC->Hits2DigitsSector(36+3);             
74   TPC->Hits2DigitsSector(36+1+18);             
75   TPC->Hits2DigitsSector(36+2+18);             
76   TPC->Hits2DigitsSector(36+3+18); 
77             
78   //write results
79
80   char treeName[100];
81   sprintf(treeName,"TreeD_%s",param->GetTitle());
82   TPC->GetDigitsArray()->GetTree()->Write(treeName);
83
84   delete gAlice; gAlice=0;
85   file->Close(); delete file;
86   return 0;
87 };
88