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