]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFSDigits2Digits.C
8de37e776daaa6439179719c1859e8eb8280a4c9
[u/mrichter/AliRoot.git] / TOF / AliTOFSDigits2Digits.C
1 Int_t AliTOFSDigits2Digits(Int_t nev=5) {
2
3   // Adapted to the NewIO by I.Belikov (Jouri.Belikov@cern.ch)
4
5   // number 3 is a legacy from AliDigit object
6     
7    if (gAlice) {
8       delete gAlice->GetRunLoader();
9       delete gAlice;//if everything was OK here it is already NULL
10       gAlice = 0x0;
11    }
12
13    AliRunLoader* rl = AliRunLoader::Open("galice.root");
14    if (rl == 0x0) {
15       cerr<<"Can not open session"<<endl;
16       return 1;
17    }
18    
19    if (rl->LoadgAlice()) {
20       cerr<<"Error occured while loading gAlice"<<endl;
21       return 1;
22    }
23
24    AliLoader *tofl = rl->GetLoader("TOFLoader");
25    if (tofl == 0x0) {
26       cerr<<"Can not get the TOF Loader"<<endl;
27       return 1;
28    }
29
30    gAlice=rl->GetAliRun();
31    if (!gAlice) {
32       cerr<<"Can't get gAlice !\n";
33       return 1;
34    }
35
36    tofl->LoadSDigits("read");
37    tofl->LoadDigits("recreate");
38
39    Int_t totndig=0;   // total number of digits
40    Int_t tottracks=0; // total number of tracks contributing to totndig
41
42    if (nev>rl->GetNumberOfEvents()) nev=rl->GetNumberOfEvents();
43    
44    TClonesArray *fSDigits=new TClonesArray("AliTOFSDigit",  1000); 
45    TClonesArray *fDigits =new TClonesArray("AliTOFdigit",  1000);
46    TClonesArray &da=*fDigits; 
47
48    for (Int_t ievent = 0; ievent < nev; ievent++) {
49      rl->GetEvent(ievent);
50
51      TTree *sTree=tofl->TreeS();
52      if (sTree == 0) {
53        cerr<<"Can't get the sdigit tree !\n";
54        return 1;
55      }
56      TBranch *branch=sTree->GetBranch("TOF");
57      if (!branch) {
58        cerr<<"Cant' get the branch !\n";
59        return 1;
60      }
61      branch->SetAddress(&fSDigits);           
62     
63     TTree *dTree=tofl->TreeD();
64     if (dTree == 0) {
65        tofl->MakeTree("D");
66        dTree=tofl->TreeD();
67     }
68     branch=dTree->GetBranch("TOF");
69     if (!branch) dTree->Branch("TOF",&fDigits);    
70     else branch->SetAddress(&fDigits);
71
72     Int_t nEntries = sTree->GetEntries();                                
73     for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
74       sTree->GetEvent(iEntry);
75       
76       Int_t ndig = fSDigits->GetEntriesFast();
77       cout << "----------------<AliTOFSDigits2Digits>----------------" << endl;
78       cout << "Found " << ndig << " TOF SDigits for event " << ievent << endl;
79       cout << "------------------------------------------------------" << endl;
80
81       for (Int_t k = 0; k < ndig; k++) {
82         Int_t    vol[5];       // location for a digit
83         
84         // Get the information for this digit
85         AliTOFSDigit *tofsdigit = (AliTOFSDigit *)fSDigits->UncheckedAt(k);
86
87         Int_t nslot=tofsdigit->GetNDigits(); // get the number of slots
88                                              // for current sdigit
89
90         // TOF sdigit volumes (always the same for all slots)
91         Int_t sector    = tofsdigit->GetSector(); // range [1-18]
92         Int_t plate     = tofsdigit->GetPlate();  // range [1- 5]
93         Int_t strip     = tofsdigit->GetStrip();  // range [1-20]
94         Int_t padz      = tofsdigit->GetPadz();   // range [1- 2]
95         Int_t padx      = tofsdigit->GetPadx();   // range [1-48]
96
97         vol[0] = sector;
98         vol[1] = plate;
99         vol[2] = strip;
100         vol[3] = padx;
101         vol[4] = padz;
102
103         //--------------------- QA section ----------------------
104         // in the while, I perform QA
105         Bool_t isSDigitBad = (sector<1 || sector>18 || 
106                                plate<1 || plate >5  || 
107                                 padz<1 || padz>2    || 
108                                 padx<1 || padx>48);
109         
110         if (isSDigitBad) {
111           cout << "<AliTOFSDigits2Digits>  strange sdigit found" << endl;
112           return 3;
113         }
114         //-------------------------------------------------------
115
116         // start loop on number of slots for current sdigit
117         for (Int_t islot = 0; islot < nslot; islot++) {
118           Float_t  digit[2];     // TOF digit variables
119           const Int_t kMAXDIGITS = 3;
120           Int_t tracknum[kMAXDIGITS];//contributing tracks for the current slot
121         
122           Float_t tdc=tofsdigit->GetTdc(islot); digit[0]=tdc;
123           Float_t adc=tofsdigit->GetAdc(islot); digit[1]=adc;
124
125           tracknum[0]=tofsdigit->GetTrack(islot,0);
126           tracknum[1]=tofsdigit->GetTrack(islot,1);
127           tracknum[2]=tofsdigit->GetTrack(islot,2);
128
129           for (Int_t i = 0; i < kMAXDIGITS; i++) {
130             tottracks++;
131             // search for the first empty location
132             if(tracknum[i]==-1){
133               tottracks--;
134               break;
135             }
136           }
137
138           // adding a TOF digit for each slot
139           {
140           Int_t ndigits=da.GetEntriesFast();
141           //cerr<<ndigits<<" ndig \n";
142           new (da[ndigits]) AliTOFdigit(tracknum, vol, digit);
143           }          
144           totndig++;
145         }
146
147
148       } // end loop on sdigits
149       fSDigits->Clear();
150       
151     } // end loop on entries
152
153     dTree->Fill();
154     tofl->WriteDigits("OVERWRITE");
155
156     // free used memory
157     fDigits->Clear();
158
159   } // end loop on events 
160
161   delete fSDigits;
162   delete fDigits;
163
164   cout << "----------------------------------------------------------" << endl;
165   cout << "<AliTOFSDigits2Digits> Summary" << endl;
166   cout << "contributing tracks to " << totndig << " digits: " << tottracks << endl; 
167   cout << "----------------------------------------------------------" << endl;
168
169   return 0;
170   
171
172 }