]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDigitizer.cxx
Use debug methods
[u/mrichter/AliRoot.git] / RICH / AliRICHDigitizer.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.1  2001/11/02 15:37:26  hristov
19    Digitizer class created. Code cleaning and bug fixes (J.Chudoba)
20 */
21 #include <iostream> 
22
23 #include <TTree.h> 
24 #include <TObjArray.h>
25 #include <TFile.h>
26 #include <TDirectory.h>
27 #include <TParticle.h>
28
29 #include "AliRICHDigitizer.h"
30 #include "AliRICHChamber.h"
31 #include "AliHitMap.h"
32 #include "AliRICHHitMapA1.h"
33 #include "AliRICH.h"
34 #include "AliRICHHit.h"
35 #include "AliRICHSDigit.h"
36 #include "AliRICHDigit.h"
37 #include "AliRICHTransientDigit.h"
38 #include "AliRun.h"
39 #include "AliPDG.h"
40 #include "AliRunDigitizer.h"
41
42 ClassImp(AliRICHDigitizer)
43
44 //___________________________________________
45   AliRICHDigitizer::AliRICHDigitizer()
46 {
47 // Default constructor - don't use it   
48   fHits = 0;
49   fSDigits = 0;
50   fHitMap = 0;
51   fTDList = 0;
52 }
53
54 ////////////////////////////////////////////////////////////////////////
55 AliRICHDigitizer::AliRICHDigitizer(AliRunDigitizer* manager) 
56     :AliDigitizer(manager)
57 {
58 // ctor which should be used
59   fHits = 0;
60   fSDigits = 0;
61   fHitMap = 0;
62   fTDList = 0;
63   fDebug   = 0; 
64   if (GetDebug()>2) 
65     cerr<<"AliRICHDigitizer::AliRICHDigitizer"
66         <<"(AliRunDigitizer* manager) was processed"<<endl;
67 }
68
69 ////////////////////////////////////////////////////////////////////////
70 AliRICHDigitizer::~AliRICHDigitizer()
71 {
72 // Destructor
73   if (fHits)       delete fHits;
74   if (fSDigits)    delete fSDigits;
75   if (fHitMap)     delete fHitMap;
76   if (fTDList)     delete fTDList;
77 }
78
79 //------------------------------------------------------------------------
80 Bool_t AliRICHDigitizer::Exists(const AliRICHSDigit *padhit)
81 {
82   return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
83 }
84
85 //------------------------------------------------------------------------
86 void AliRICHDigitizer::Update(AliRICHSDigit *padhit)
87 {
88   AliRICHTransientDigit *pdigit = 
89     static_cast<AliRICHTransientDigit*>(
90       fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
91
92   // update charge
93   //
94   Int_t iqpad    = Int_t(padhit->QPad()); // charge per pad
95   pdigit->AddSignal(iqpad);
96   pdigit->AddPhysicsSignal(iqpad);              
97
98   // update list of tracks
99   //
100   Int_t track, charge;
101   track = fTrack+fMask;
102   if (fSignal) {
103     charge = iqpad;
104   } else {
105     charge = kBgTag;
106   }
107   pdigit->UpdateTrackList(track,charge);
108 }
109
110 //------------------------------------------------------------------------
111 void AliRICHDigitizer::CreateNew(AliRICHSDigit *padhit)
112 {
113   fTDList->AddAtAndExpand(
114     new AliRICHTransientDigit(fNch,fDigits),fCounter);
115   fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
116
117   AliRICHTransientDigit* pdigit = 
118     static_cast<AliRICHTransientDigit*>(fTDList->Last());
119   // list of tracks
120   Int_t track, charge;
121   if (fSignal) {
122     track = fTrack;
123     charge = padhit->QPad();
124   } else {
125     track = kBgTag;
126     charge = kBgTag;
127   }
128   pdigit->AddToTrackList(track,charge);
129   fCounter++;
130 }
131
132
133 ////////////////////////////////////////////////////////////////////////
134 Bool_t AliRICHDigitizer::Init()
135 {
136 // Initialisation
137   fHits     = new TClonesArray("AliRICHHit",1000);
138   fSDigits  = new TClonesArray("AliRICHSDigit",1000);
139   return kTRUE;
140 }
141
142 ////////////////////////////////////////////////////////////////////////
143 //void AliRICHDigitizer::Digitise(Int_t nev, Int_t flag)
144 void AliRICHDigitizer::Exec(Option_t* option)
145 {
146   TString optionString = option;
147   if (optionString.Data() == "deb") {
148     cout<<"AliMUONDigitizer::Exec: called with option deb "<<endl;
149     fDebug = 3;
150   }
151
152   AliRICHChamber*       iChamber;
153   AliSegmentation*  segmentation;
154   
155   fTDList = new TObjArray;
156     
157   AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
158   pRICH->MakeBranchInTreeD(fManager->GetTreeD());
159   fHitMap= new AliHitMap* [kNCH];
160         
161   for (Int_t i =0; i<kNCH; i++) {
162     iChamber= &(pRICH->Chamber(i));
163     segmentation=iChamber->GetSegmentationModel(1);
164     fHitMap[i] = new AliRICHHitMapA1(segmentation, fTDList);
165   }
166
167 // Loop over files to digitize
168   fSignal = kTRUE;
169   fCounter = 0;
170   for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
171        inputFile++) {
172
173 // Connect RICH branches
174
175     if (inputFile > 0 ) fSignal = kFALSE;
176     TBranch *branchHits = 0;
177     TBranch *branchSDigits = 0;
178     TTree *treeH = fManager->GetInputTreeH(inputFile);
179     if (GetDebug()>2) {
180       cerr<<"   inputFile  "<<inputFile<<endl;
181       cerr<<"   treeH, fHits "<<treeH<<" "<<fHits<<endl;
182     }
183     if (treeH && fHits) {
184       branchHits = treeH->GetBranch("RICH");
185       if (branchHits) {
186         fHits->Clear();
187         branchHits->SetAddress(&fHits);
188       }
189       else
190         Error("Exec","branch RICH was not found");
191     }
192     if (GetDebug()>2) cerr<<"   branchHits = "<<branchHits<<endl;
193
194     if (treeH && fSDigits) {
195       branchSDigits = treeH->GetBranch("RICHSDigits");
196       if (branchSDigits) 
197         branchSDigits->SetAddress(&fSDigits);
198       else
199         Error("exec","branch RICHSDigits was not found");
200     }
201     if (GetDebug()>2) cerr<<"   branchSDigits = "<<branchSDigits<<endl;
202
203
204     //
205     //   Loop over tracks
206     //
207     
208     Int_t ntracks =(Int_t) treeH->GetEntries();
209     for (fTrack=0; fTrack<ntracks; fTrack++) {
210       fHits->Clear();
211       fSDigits->Clear();  
212       branchHits->GetEntry(fTrack);
213       branchSDigits->GetEntry(fTrack);
214
215
216       //
217       //   Loop over hits
218       for(Int_t i = 0; i < fHits->GetEntriesFast(); ++i) {
219         AliRICHHit* mHit = static_cast<AliRICHHit*>(fHits->At(i));
220         fNch = mHit->Chamber()-1;  // chamber number
221         if (fNch >= kNCH) {
222           cerr<<"AliRICHDigitizer: chamber nr. fNch out of range: "<<fNch<<endl;
223           cerr<<"               track: "<<fTrack<<endl;
224           continue;
225         }
226         iChamber = &(pRICH->Chamber(fNch));
227           
228         //
229         // Loop over pad hits
230         for (AliRICHSDigit* mPad=
231                (AliRICHSDigit*)pRICH->FirstPad(mHit,fSDigits);
232              mPad;
233              mPad=(AliRICHSDigit*)pRICH->NextPad(fSDigits))
234         {
235           Int_t iqpad    = mPad->QPad();       // charge per pad
236           fDigits[0]=mPad->PadX();
237           fDigits[1]=mPad->PadY();
238           fDigits[2]=iqpad;
239           fDigits[3]=iqpad;
240           fDigits[4]=mPad->HitNumber();
241                   
242           // build the list of fired pads and update the info
243           if (Exists(mPad)) {
244             Update(mPad);
245           } else {
246             CreateNew(mPad);
247           }
248         } //end loop over clusters
249       } // hit loop
250     } // track loop
251   } // end file loop
252   if (GetDebug()>2) cerr<<"END OF FILE LOOP"<<endl;
253
254   Int_t tracks[kMAXTRACKSPERRICHDIGIT];
255   Int_t charges[kMAXTRACKSPERRICHDIGIT];
256   Int_t nentries=fTDList->GetEntriesFast();
257     
258   for (Int_t nent=0;nent<nentries;nent++) {
259     AliRICHTransientDigit *transDigit=(AliRICHTransientDigit*)fTDList->At(nent);
260     if (transDigit==0) continue; 
261     Int_t ich=transDigit->GetChamber();
262     Int_t q=transDigit->Signal(); 
263     iChamber=&(pRICH->Chamber(ich));
264     AliRICHResponse * response=iChamber->GetResponseModel();
265     Int_t adcmax= (Int_t) response->MaxAdc();
266       
267       
268     // add white noise and do zero-suppression and signal truncation (new electronics,old electronics gaus 1.2,0.2)
269     //printf("Treshold: %d\n",iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY()));
270
271 // tmp change
272 //    Int_t pedestal = iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY());
273     Int_t pedestal = 0;
274
275
276     //printf("Pedestal:%d\n",pedestal);
277     //Int_t pedestal=0;
278     Float_t treshold = (pedestal + 4*2.2);
279       
280     Float_t meanNoise = gRandom->Gaus(2.2, 0.3);
281     Float_t noise     = gRandom->Gaus(0, meanNoise);
282     q+=(Int_t)(noise + pedestal);
283     //q+=(Int_t)(noise);
284     //          magic number to be parametrised !!! 
285     if ( q <= treshold) 
286     {
287       q = q - pedestal;
288       continue;
289     }
290     q = q - pedestal;
291     if ( q >= adcmax) q=adcmax;
292     fDigits[0]=transDigit->PadX();
293     fDigits[1]=transDigit->PadY();
294     fDigits[2]=q;
295     fDigits[3]=transDigit->Physics();
296     fDigits[4]=transDigit->Hit();
297
298     Int_t nptracks = transDigit->GetNTracks();  
299
300     // this was changed to accomodate the real number of tracks
301     if (nptracks > kMAXTRACKSPERRICHDIGIT) {
302       printf("Attention - tracks > 10 %d\n",nptracks);
303       nptracks=kMAXTRACKSPERRICHDIGIT;
304     }
305     if (nptracks > 2) {
306       printf("Attention - tracks > 2  %d \n",nptracks);
307     }
308     for (Int_t tr=0;tr<nptracks;tr++) {
309       tracks[tr]=transDigit->GetTrack(tr);
310       charges[tr]=transDigit->GetCharge(tr);
311     }      //end loop over list of tracks for one pad
312     if (nptracks < kMAXTRACKSPERRICHDIGIT ) {
313       for (Int_t t=nptracks; t<kMAXTRACKSPERRICHDIGIT; t++) {
314         tracks[t]=0;
315         charges[t]=0;
316       }
317     }
318     //write file
319     //if (ich==2)
320     //fprintf(points,"%4d,      %4d,      %4d\n",digits[0],digits[1],digits[2]);
321       
322     // fill digits
323     pRICH->AddDigits(ich,tracks,charges,fDigits);
324   }     
325   fManager->GetTreeD()->Fill();
326
327   pRICH->ResetDigits();
328   fTDList->Delete();    // or fTDList->Clear(); ???
329   for(Int_t ii=0;ii<kNCH;++ii) {
330     if (fHitMap[ii]) {
331       delete fHitMap[ii];
332       fHitMap[ii]=0;
333     }
334   }
335     
336   TClonesArray *richDigits;
337   for (Int_t k=0;k<kNCH;k++) {
338     richDigits = pRICH->DigitsAddress(k);
339     Int_t ndigit=richDigits->GetEntriesFast();
340     printf ("Chamber %d digits %d \n",k,ndigit);
341   }
342   pRICH->ResetDigits(); /// ??? should it be here???
343   fManager->GetTreeD()->Write(0,TObject::kOverwrite);
344   delete [] fHitMap;
345   delete fTDList;
346
347   if (fHits)    fHits->Delete();
348   if (fSDigits) fSDigits->Delete();
349
350 }