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