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