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