]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDigitizer.cxx
Minor changes and cosmetics according to Anders.
[u/mrichter/AliRoot.git] / TPC / AliTPCDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 */
19
20
21
22
23 #include <stdlib.h>
24 #include <TTree.h> 
25 #include <TObjArray.h>
26 #include <TFile.h>
27 #include <TDirectory.h>
28 #include <iostream.h>
29
30 #include "AliTPCDigitizer.h"
31
32 #include "AliTPC.h"
33 #include "AliTPCParam.h"
34 #include "AliTPCParamSR.h" 
35 #include "AliRun.h"
36 #include "AliPDG.h"
37 #include "AliRunDigitizer.h"
38 #include "AliSimDigits.h"
39
40
41 ClassImp(AliTPCDigitizer)
42
43 //___________________________________________
44 AliTPCDigitizer::AliTPCDigitizer() :AliDigitizer()
45 {
46 // Default ctor - don't use it
47   fDebug =0;
48 }
49
50 //___________________________________________
51 AliTPCDigitizer::AliTPCDigitizer(AliRunDigitizer* manager) 
52     :AliDigitizer(manager)
53 {
54 // ctor which should be used
55   fDebug =0;
56   if (GetDebug()>2) 
57     cerr<<"AliTPCDigitizer::AliTPCDigitizer"
58         <<"(AliRunDigitizer* manager) was processed"<<endl;
59 }
60
61 //------------------------------------------------------------------------
62 AliTPCDigitizer::~AliTPCDigitizer()
63 {
64 // Destructor
65 }
66
67
68
69 //------------------------------------------------------------------------
70 Bool_t AliTPCDigitizer::Init()
71 {
72 // Initialization 
73     
74  return kTRUE;
75 }
76
77
78 //------------------------------------------------------------------------
79 void AliTPCDigitizer::Exec(Option_t* option)
80 {
81   ExecFast(option);  
82 }
83 //------------------------------------------------------------------------
84 void AliTPCDigitizer::ExecFast(Option_t* option)
85 {
86   
87   // merge input tree's with summable digits
88   //output stored in TreeTPCD
89   char s[100]; 
90   char ss[100];
91   TString optionString = option;
92   if (optionString.Data() == "deb") {
93     cout<<"AliTPCDigitizer::Exec: called with option deb "<<endl;
94     fDebug = 3;
95   }
96   //get detector and geometry
97   AliTPC *pTPC  = (AliTPC *) gAlice->GetModule("TPC");
98   AliTPCParam * param = pTPC->GetParam();
99   
100   sprintf(s,param->GetTitle());
101   sprintf(ss,"75x40_100x60");
102   if(strcmp(s,ss)==0){
103     printf("2 pad-length geom hits with 3 pad-lenght geom digits...\n");
104     delete param;
105     param=new AliTPCParamSR();
106   }
107   else{
108    sprintf(ss,"75x40_100x60_150x60");
109    if(strcmp(s,ss)!=0) {
110      printf("No TPC parameters found...\n");
111      exit(2); 
112    }
113   }
114   
115   pTPC->GenerNoise(500000); //create teble with noise
116   //
117   Int_t nInputs = fManager->GetNinputs();
118   Int_t * masks = new Int_t[nInputs];
119   for (Int_t i=0; i<nInputs;i++)
120     masks[i]= fManager->GetMask(i);
121   Short_t **pdig= new Short_t*[nInputs];   //pointers to the expanded digits array
122   Int_t **ptr=  new Int_t*[nInputs];       //pointers to teh expanded tracks array
123
124   //create digits array for given sectors
125   // make indexes
126   AliSimDigits ** digarr = new AliSimDigits*[nInputs]; 
127   for (Int_t i1=0;i1<nInputs; i1++){
128     digarr[i1]=0;
129     //    intree[i1]
130     TTree * treear =  fManager->GetInputTreeTPCS(i1);
131     if (!treear) {
132       cerr<<"AliTPCDigitizer: Input tree with SDigits not found in"
133           <<" input "<< i1<<endl;
134       return;
135     }
136     if (treear->GetIndex()==0) 
137       treear->BuildIndex("fSegmentID","fSegmentID");
138     if (!treear) {      
139       cerr<<" TPC -  not existing input = \n"<<i1<<" ";      
140     }
141     treear->GetBranch("Segment")->SetAddress(&digarr[i1]);
142   }
143   Stat_t nentries = fManager->GetInputTreeTPCS(0)->GetEntries();
144   
145
146   //create branch's in TPC treeD
147   AliSimDigits * digrow = new AliSimDigits;
148   TTree * tree  = fManager->GetTreeDTPC();
149   tree->Branch("Segment","AliSimDigits",&digrow);
150   //
151
152   param->SetZeroSup(2);
153
154   Int_t zerosup = param->GetZeroSup(); 
155   //
156   //Loop over segments of the TPC
157     
158   for (Int_t n=0; n<nentries; n++) {
159   //    for (Int_t n=0; n<300; n++) {
160     fManager->GetInputTreeTPCS(0)->GetEvent(n);      
161     digarr[0]->ExpandBuffer();
162     digarr[0]->ExpandTrackBuffer();
163            
164     for (Int_t i=1;i<nInputs; i++){ 
165       fManager->GetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());
166       if ((digarr[0]->GetID()-digarr[i]->GetID())>0) 
167         printf("problem - not corresponding segment in background event\n");
168       
169       digarr[i]->ExpandBuffer();
170       digarr[i]->ExpandTrackBuffer();
171     }   
172     Int_t sec, row;
173     if (!param->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
174       cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
175       continue;
176     }
177
178     digrow->SetID(digarr[0]->GetID());
179
180     Int_t nrows = digarr[0]->GetNRows();
181     Int_t ncols = digarr[0]->GetNCols();
182     digrow->Allocate(nrows,ncols);
183     digrow->AllocateTrack(3);
184
185     Float_t q=0;
186     Int_t label[1000]; //stack for 300 events 
187     Int_t labptr = 0;
188
189     Int_t nElems = nrows*ncols;     
190  
191     for (Int_t i=0;i<nInputs; i++){ 
192       pdig[i] = digarr[i]->GetDigits();
193       ptr[i]  = digarr[i]->GetTracks();
194     }
195     Short_t *pdig1= digrow->GetDigits();
196     Int_t   *ptr1= digrow->GetTracks() ;
197
198     
199
200     //    for (Int_t rows=0;rows<nrows; rows++){
201     //  for (Int_t col=0;col<ncols; col++){
202     for (Int_t elem=0;elem<nElems; elem++){    
203       //for (Int_t elem=nElems;elem<nElems; elem++){
204
205         q=0;
206         labptr=0;
207         // looop over digits 
208         for (Int_t i=0;i<nInputs; i++){ 
209           //          q  += digarr[i]->GetDigitFast(rows,col);
210           q  += *(pdig[i]);
211           
212           for (Int_t tr=0;tr<3;tr++) {
213             //             Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
214             Int_t lab = ptr[i][tr*nElems];
215             if ( (lab > 1) && *(pdig[i])>zerosup) {
216               label[labptr]=lab+masks[i];
217               labptr++;
218             }      
219           }
220           pdig[i]++;
221           ptr[i]++;
222           
223         }
224         q/=16.;  //conversion factor
225         //      Float_t noise  = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());  
226         Float_t noise  = pTPC->GetNoise();
227         q+=noise;
228         q=TMath::Nint(q);
229         if (q > zerosup){ 
230           
231           if(q > param->GetADCSat()) q = (Short_t)(param->GetADCSat());
232           //digrow->SetDigitFast((Short_t)q,rows,col);  
233           *pdig1 =Short_t(q);
234           for (Int_t tr=0;tr<3;tr++){
235             if (tr<labptr) 
236               // ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
237               ptr1[tr*nElems] = label[tr];
238             //else
239               //            ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);          
240             //  ptr1[tr*nElems] = 1;
241           }
242         }
243         pdig1++;
244         ptr1++;
245     }
246     
247     digrow->CompresBuffer(1,zerosup);
248     digrow->CompresTrackBuffer(1);
249     tree->Fill();
250     if (fDebug>0) cerr<<sec<<"\t"<<row<<"\n";  
251   } 
252   fManager->GetTreeDTPC()->Write(0,TObject::kOverwrite);
253   delete digrow;     
254   for (Int_t i1=0;i1<nInputs; i1++) delete digarr[i1];
255   delete []masks;
256   delete digarr;  
257 }
258
259
260
261 //------------------------------------------------------------------------
262 void AliTPCDigitizer::ExecSave(Option_t* option)
263 {
264   
265   // merge input tree's with summable digits
266   //output stored in TreeTPCD
267
268   TString optionString = option;
269   if (optionString.Data() == "deb") {
270     cout<<"AliTPCDigitizer::Exec: called with option deb "<<endl;
271     fDebug = 3;
272   }
273   //get detector and geometry 
274   printf("TPC merging -1  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
275   AliTPC *pTPC  = (AliTPC *) gAlice->GetModule("TPC");
276   AliTPCParam * param = pTPC->GetParam();
277   pTPC->GenerNoise(500000); //create teble with noise
278   printf("noise %f \n",  param->GetNoise()*param->GetNoiseNormFac());
279   //
280   Int_t nInputs = fManager->GetNinputs();
281   Int_t * masks = new Int_t[nInputs];
282   for (Int_t i=0; i<nInputs;i++)
283     masks[i]= fManager->GetMask(i);
284   //  Short_t **pdig= new Short_t*[nInputs];   //pointers to the expanded digits array
285   //Int_t **ptr=  new Int_t*[nInputs];       //pointers to teh expanded tracks array
286
287   //create digits array for given sectors
288   // make indexes
289    printf("TPC merging -2  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
290   AliSimDigits ** digarr = new AliSimDigits*[nInputs]; 
291   for (Int_t i1=0;i1<nInputs; i1++){
292     digarr[i1]=0;
293     //    intree[i1]
294     TTree * treear =  fManager->GetInputTreeTPCS(i1);
295     printf("TPC merging -2.7  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
296     TBranch * br = treear->GetBranch("fSegmentID");
297     if (br) br->GetFile()->cd();
298     printf("TPC merging -2.75  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
299     //treear->BuildIndex("fSegmentID","fSegmentID");
300     if (!treear) {      
301       cerr<<" TPC -  not existing input = \n"<<i1<<" ";      
302     } 
303     printf("TPC merging -2.8  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
304     treear->GetBranch("Segment")->SetAddress(&digarr[i1]);
305      printf("TPC merging -2.9  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
306   }
307   Stat_t nentries = fManager->GetInputTreeTPCS(0)->GetEntries();
308   
309
310   //create branch's in TPC treeD
311   AliSimDigits * digrow = new AliSimDigits;
312   TTree * tree  = fManager->GetTreeDTPC();
313   //if (tree->GetBranch("Segment") ) tree->GetBranch("Segment")->SetAddress(&digrow);
314   //else
315   tree->Branch("Segment","AliSimDigits",&digrow);
316   //
317  printf("TPC merging -3  -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
318   param->SetZeroSup(2);
319
320   Int_t zerosup = param->GetZeroSup();
321   //Loop over segments of the TPC
322     
323   for (Int_t n=0; n<nentries; n++) {
324   //    for (Int_t n=0; n<300; n++) {
325     fManager->GetInputTreeTPCS(0)->GetEvent(n);      
326     digarr[0]->ExpandBuffer();
327     digarr[0]->ExpandTrackBuffer();
328            
329     for (Int_t i=1;i<nInputs; i++){ 
330       fManager->GetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());      
331       //fManager->GetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),1);      
332       digarr[i]->ExpandBuffer();
333       digarr[i]->ExpandTrackBuffer();
334       if ((digarr[0]->GetID()-digarr[i]->GetID())>0) 
335         printf("problem\n");
336     
337     }   
338     
339     Int_t sec, row;
340     if (!param->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
341       cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
342       continue;
343     }
344
345     digrow->SetID(digarr[0]->GetID());
346
347     Int_t nrows = digarr[0]->GetNRows();
348     Int_t ncols = digarr[0]->GetNCols();
349     digrow->Allocate(nrows,ncols);
350     digrow->AllocateTrack(3);
351
352     Float_t q=0;
353     Int_t label[1000]; //stack for 300 events 
354     Int_t labptr = 0;
355
356     
357
358     for (Int_t rows=0;rows<nrows; rows++){
359       for (Int_t col=0;col<ncols; col++){
360     
361         q=0;
362         labptr=0;
363         // looop over digits 
364         for (Int_t i=0;i<nInputs; i++){ 
365           q  += digarr[i]->GetDigitFast(rows,col);
366           //q  += *(pdig[i]);
367           
368           for (Int_t tr=0;tr<3;tr++) {
369             Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
370             //Int_t lab = ptr[i][tr*nElems];
371             if ( (lab > 1) ) {
372               label[labptr]=lab+masks[i];
373               labptr++;
374             }      
375           }
376           // pdig[i]++;
377           //ptr[i]++;
378           
379         }
380         q/=16.;  //conversion factor
381         //      Float_t noise  = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());  
382         Float_t noise  = pTPC->GetNoise();
383         q+=noise;
384         q=TMath::Nint(q);
385         if (q > zerosup){ 
386           
387           if(q > param->GetADCSat()) q = (Short_t)(param->GetADCSat());
388           digrow->SetDigitFast((Short_t)q,rows,col);  
389           // *pdig1 =Short_t(q);
390           for (Int_t tr=0;tr<3;tr++){
391             if (tr<labptr) 
392               ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
393             //ptr1[tr*nElems] = label[tr];
394             //else
395               //            ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);          
396             //  ptr1[tr*nElems] = 1;
397           }
398         }
399         //pdig1++;
400         //ptr1++;
401     }
402     }
403     
404     digrow->CompresBuffer(1,zerosup);
405     digrow->CompresTrackBuffer(1);
406     tree->Fill();
407     if (fDebug>0) cerr<<sec<<"\t"<<row<<"\n";  
408   } 
409   printf("end TPC merging - end -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
410   fManager->GetTreeDTPC()->Write(0,TObject::kOverwrite);
411   delete digrow;     
412   for (Int_t i1=0;i1<nInputs; i1++) delete digarr[i1];
413   delete []masks;
414   delete digarr;  
415 }