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