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