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