]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDigitizer.cxx
Avoid warning.
[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 /* $Id$ */
17
18 //Piotr.Skowronski@cern.ch :
19 //Corrections applied in order to compile (only) with new I/O and folder structure
20 //To be implemented correctly by responsible
21
22 #include <Riostream.h> 
23
24 #include <TTree.h> 
25 #include <TObjArray.h>
26 #include <TFile.h>
27 #include <TDirectory.h>
28 #include <TParticle.h>
29
30 #include "AliRunLoader.h"
31 #include "AliLoader.h"
32
33 #include "AliRICHDigitizer.h"
34 #include "AliRICHChamber.h"
35 #include "AliHitMap.h"
36 #include "AliRICHHitMapA1.h"
37 #include "AliRICH.h"
38 #include "AliRICHHit.h"
39 #include "AliRICHSDigit.h"
40 #include "AliRICHDigit.h"
41 #include "AliRICHTransientDigit.h"
42 #include "AliRun.h"
43 #include "AliPDG.h"
44 #include "AliRunDigitizer.h"
45
46 ClassImp(AliRICHDigitizer)
47
48 //___________________________________________
49   AliRICHDigitizer::AliRICHDigitizer()
50 {
51 // Default constructor - don't use it   
52   fHits = 0;
53   fSDigits = 0;
54   fHitMap = 0;
55   fTDList = 0;
56 }
57
58 ////////////////////////////////////////////////////////////////////////
59 AliRICHDigitizer::AliRICHDigitizer(AliRunDigitizer* manager) 
60     :AliDigitizer(manager)
61 {
62 // ctor which should be used
63   fHits = 0;
64   fSDigits = 0;
65   fHitMap = 0;
66   fTDList = 0;
67   fDebug   = 0; 
68   if (GetDebug()>2) 
69     cerr<<"AliRICHDigitizer::AliRICHDigitizer"
70       <<"(AliRunDigitizer* manager) was processed"<<endl;
71 }
72
73 ////////////////////////////////////////////////////////////////////////
74 AliRICHDigitizer::~AliRICHDigitizer()
75 {
76 // Destructor
77     if (fHits) {
78         fHits->Delete();
79         delete fHits;
80     }
81     if (fSDigits) {
82         fSDigits->Delete();
83         delete fSDigits;
84     }
85     for (Int_t i=0; i<kNCH; i++ )
86         delete fHitMap[i];
87     delete [] fHitMap;
88     if (fTDList) {
89         fTDList->Delete();
90         delete fTDList;
91   }
92 }
93
94 //------------------------------------------------------------------------
95 Bool_t AliRICHDigitizer::Exists(const AliRICHSDigit *padhit)
96 {
97   return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
98 }
99
100 //------------------------------------------------------------------------
101 void AliRICHDigitizer::Update(AliRICHSDigit *padhit)
102 {
103   AliRICHTransientDigit *pdigit = 
104     static_cast<AliRICHTransientDigit*>(
105       fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
106
107   // update charge
108   //
109   Int_t iqpad    = Int_t(padhit->QPad()); // charge per pad
110   pdigit->AddSignal(iqpad);
111   pdigit->AddPhysicsSignal(iqpad);            
112
113   // update list of tracks
114   //
115   Int_t track, charge;
116   track = fTrack+fMask;
117   if (fSignal) {
118     charge = iqpad;
119   } else {
120     charge = kBgTag;
121   }
122   pdigit->UpdateTrackList(track,charge);
123 }
124
125 //------------------------------------------------------------------------
126 void AliRICHDigitizer::CreateNew(AliRICHSDigit *padhit)
127 {
128   fTDList->AddAtAndExpand(
129     new AliRICHTransientDigit(fNch,fDigits),fCounter);
130   fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
131
132   AliRICHTransientDigit* pdigit = 
133     static_cast<AliRICHTransientDigit*>(fTDList->Last());
134   // list of tracks
135   Int_t track, charge;
136   if (fSignal) {
137     track = fTrack;
138     charge = padhit->QPad();
139   } else {
140     track = kBgTag;
141     charge = kBgTag;
142   }
143   pdigit->AddToTrackList(track,charge);
144   fCounter++;
145 }
146
147
148 ////////////////////////////////////////////////////////////////////////
149 Bool_t AliRICHDigitizer::Init()
150 {
151 // Initialisation
152   fHits     = new TClonesArray("AliRICHHit",1000);
153   fSDigits  = new TClonesArray("AliRICHSDigit",1000);
154   return kTRUE;
155 }
156
157 ////////////////////////////////////////////////////////////////////////
158 //void AliRICHDigitizer::Digitise(Int_t nev, Int_t flag)
159 void AliRICHDigitizer::Exec(Option_t* option)
160 {
161   TString optionString = option;
162   if (optionString.Data() == "deb") {
163     cout<<"AliMUONDigitizer::Exec: called with option deb "<<endl;
164     fDebug = 3;
165   }
166
167   AliRICHChamber*       iChamber;
168   AliSegmentation*  segmentation;
169  
170   AliRunLoader *inRL, *outRL;//in and out Run Loaders
171   AliLoader *ingime, *outgime;// in and out ITSLoaders
172  
173   outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
174   outgime = outRL->GetLoader("RICHLoader");
175
176   fTDList = new TObjArray;
177   
178      
179   AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
180   
181   if (outgime->TreeD() == 0x0) outgime->MakeTree("D");
182   
183   pRICH->MakeBranchInTreeD(outgime->TreeD());
184   fHitMap= new AliHitMap* [kNCH];
185         
186   for (Int_t i =0; i<kNCH; i++) {
187     iChamber= &(pRICH->Chamber(i));
188     segmentation=iChamber->GetSegmentationModel(1);
189     fHitMap[i] = new AliRICHHitMapA1(segmentation, fTDList);
190   }
191
192 // Loop over files to digitize
193   fSignal = kTRUE;
194   fCounter = 0;
195   for (Int_t inputFile=0; inputFile<fManager->GetNinputs(); inputFile++) 
196    {
197     inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
198     ingime = inRL->GetLoader("RICHLoader");
199     
200 // Connect RICH branches
201
202     if (inputFile > 0 ) fSignal = kFALSE;
203     TBranch *branchHits = 0;
204     TBranch *branchSDigits = 0;
205
206     ingime->LoadHits("READ");
207     TTree *treeH = ingime->TreeH();
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        {
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   outgime->TreeD()->Fill();
359
360   //pRICH->ResetDigits();
361   fTDList->Delete();    // or fTDList->Clear(); ???
362   for(Int_t ii=0;ii<kNCH;++ii) {
363     if (fHitMap[ii]) {
364       delete fHitMap[ii];
365       fHitMap[ii]=0;
366     }
367   }
368     
369   TClonesArray *richDigits;
370   for (Int_t k=0;k<kNCH;k++) {
371     richDigits = pRICH->DigitsAddress(k);
372     Int_t ndigit=richDigits->GetEntriesFast();
373     printf ("Chamber %d digits %d \n",k,ndigit);
374   }
375   pRICH->ResetDigits(); /// ??? should it be here???
376   
377   outgime->WriteDigits("OVERWRITE");
378
379   delete [] fHitMap;
380   delete fTDList;
381
382   if (fHits)    fHits->Delete();
383   if (fSDigits) fSDigits->Delete();
384
385 }