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