]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCDDL.C
preliminary update of the rec-hlt-tpc.C macro: adding the CA tracker and merger
[u/mrichter/AliRoot.git] / TPC / AliTPCDDL.C
CommitLineData
b62e2a95 1#if !defined(__CINT__)
2#include <TFile.h>
3#include <TTree.h>
4#include "AliTPCParamSR.h"
5#include "AliTPCDigitsArray.h"
6#include "AliSimDigits.h"
2e9f335b 7#include "AliTPCBuffer.h"
8#endif
9
10
04fa961a 11void AliTPCDDL(Int_t eventNumber=0, Int_t eth=0){
2e9f335b 12 //eth is a threshold.
a79660fb 13 //Digits stored into a file have an amplitude value greater than "eth"
04fa961a 14
15 const char * inFile_new = "galice.root";
16 AliRunLoader *rl = AliRunLoader::Open(inFile_new,"Event","read");
17
18 Int_t nevents=rl->GetNumberOfEvents();
19 cout<<"Number of Events:"<<nevents<<endl;
20 while (eventNumber<=0 || eventNumber>nevents){
21 cout<<"Insert the event number:";
22 cin>>eventNumber;
23 cout<<endl;
24 }
25 rl->GetEvent(eventNumber-1);
26 AliLoader *tpcloader=rl->GetLoader("TPCLoader");
27 tpcloader->LoadDigits();
28 TTree *digitsTree=tpcloader->TreeD();
29
30 AliSimDigits digrows, *dummy=&digrows;
31 digitsTree->GetBranch("Segment")->SetAddress(&dummy);
32 Stat_t nrows = digitsTree->GetEntries();
33 cout<<"Number of entries (rows):"<<nrows<<endl;
34 // get the TPC parameters
35 rl->CdGAFile();
36 AliTPCParamSR* param = AliTPC::LoadTPCParam(gFile);
37 if (!param)
38 cout<<"No TPC parameter"<<endl;
2e9f335b 39 AliTPCDigitsArray *digarr=new AliTPCDigitsArray;
40 digarr->Setup(param);
04fa961a 41 digarr->ConnectTree(digitsTree);
42
2e9f335b 43 AliTPCBuffer *b=new AliTPCBuffer("AliTPCDDL.dat");
9f992f70 44
45 //Verbose level
46 // 0: Silent
47 // 1: cout messages
48 // 2: txt files with digits
49 //BE CAREFUL, verbose level 2 MUST be used only for debugging and
50 //it is highly suggested to use this mode only for debugging digits files
51 //reasonably small, because otherwise the size of the txt files can reach
52 //quickly several MB wasting time and disk space.
53 b->SetVerbose(0);
54
55
04fa961a 56 nrows=Int_t(digarr->GetTree()->GetEntries());
2e9f335b 57 cout<<"Number of entries "<<nrows<<endl;
58 Int_t PSector=-1;
59 Int_t SubSec=0;
60 for (Int_t n=0; n<nrows; n++) {
61 AliSimDigits *digrow=(AliSimDigits*)digarr->LoadEntry(n);
62
63 Int_t sec,row; // sector and row number (in the TPC)
64 param->AdjustSectorRow(digrow->GetID(),sec,row);
65 // cout<<sec<<" row "<<row<<endl;
66 if(PSector!=sec){
67 SubSec=0;
68 PSector=sec;
69 }//end if
70
71 if(sec<36){
72 //inner sector [0;35]
73 if(row!=30)
74 //the whole row is written into the output file
75 b->WriteRowBinary(eth,digrow,0,0,0,sec,SubSec,row);
76 else{
77 //only the pads in the range [37;48] are written into the output file
78 b->WriteRowBinary(eth,digrow,37,48,1,sec,SubSec,row);
79 SubSec=1;
80 //only the pads outside the range [37;48] are written into the output file
81 b->WriteRowBinary(eth,digrow,37,48,2,sec,SubSec,row);
82 }//end else
83 }//end if
84 else{
85 //outer sector [36;71]
86 if(row==54)SubSec=2;
87 if((row!=27)&&(row!=76))
88 b->WriteRowBinary(eth,digrow,0,0,0,sec,SubSec,row);
89 else{
90 if(row==27){
91 //only the pads outside the range [43;46] are written into the output file
92 b->WriteRowBinary(eth,digrow,43,46,2,sec,SubSec,row);
93 SubSec=1;
94 //only the pads in the range [43;46] are written into the output file
95 b->WriteRowBinary(eth,digrow,43,46,1,sec,SubSec,row);
96 }
97 if(row==76){
98 //only the pads outside the range [33;88] are written into the output file
99 b->WriteRowBinary(eth,digrow,33,88,2,sec,SubSec,row);
100 SubSec=3;
101 //only the pads in the range [33;88] are written into the output file
102 b->WriteRowBinary(eth,digrow,33,88,1,sec,SubSec,row);
103 }
104 }
105 }//end else
106 }//end for
2e9f335b 107 cout<<"File created !"<<endl;
108 cout<<"Total number of digits: "<<b->GetDigNumber()<<endl;
109 delete b;
110 return;
111}//end AliTPCDataChallenge