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