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