]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHMerger.cxx
Template for single ring run configuration. (J. Barbosa)
[u/mrichter/AliRoot.git] / RICH / AliRICHMerger.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.9  2002/10/14 14:57:32  hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
20
21 Revision 1.7.6.1  2002/07/24 10:07:52  alibrary
22 Updating VirtualMC
23
24 Revision 1.7  2001/11/02 15:37:26  hristov
25 Digitizer class created. Code cleaning and bug fixes (J.Chudoba)
26
27 Revision 1.5  2001/10/23 13:03:35  hristov
28 The access to several data members was changed from public to protected. The digitisation was adapted to the multi-event case (J.Chudoba)
29
30 Revision 1.4  2001/10/21 18:31:24  hristov
31 Several pointers were set to zero in the default constructors to avoid memory management problems
32
33 Revision 1.3  2001/05/16 14:57:20  alibrary
34 New files for folders and Stack
35
36 Revision 1.2  2001/03/14 18:16:08  jbarbosa
37 Corrected bug (more to correct).
38 File "points.dat" is no longer created.
39
40 Revision 1.1  2001/02/27 22:13:34  jbarbosa
41 Implementing merger class.
42
43 */
44 #include <Riostream.h> 
45
46 #include <TTree.h> 
47 #include <TObjArray.h>
48 #include <TFile.h>
49 #include <TDirectory.h>
50 #include <TParticle.h>
51
52 #include "AliRICHMerger.h"
53 #include "AliRICHChamber.h"
54 #include "AliHitMap.h"
55 #include "AliRICHHitMapA1.h"
56 #include "AliRICH.h"
57 #include "AliRICHHit.h"
58 #include "AliRICHSDigit.h"
59 #include "AliRICHDigit.h"
60 #include "AliRICHTransientDigit.h"
61 #include "AliRun.h"
62 #include "AliPDG.h"
63
64 ClassImp(AliRICHMerger)
65
66 //___________________________________________
67   AliRICHMerger::AliRICHMerger()
68 {
69 // Default constructor    
70   fEvNrSig = 0;
71   fEvNrBgr = 0;
72   fMerge   = kDigitize;
73   fFnBgr   = 0;
74   fHitMap  = 0;
75   fList    = 0;
76   fTrH1    = 0;
77   fHitsBgr = 0;
78   fSDigitsBgr = 0;
79   fHitMap     = 0;
80   fList       = 0;
81   fBgrFile    = 0;
82 }
83
84 //------------------------------------------------------------------------
85 AliRICHMerger::~AliRICHMerger()
86 {
87 // Destructor
88   if (fTrH1)       delete fTrH1;
89   if (fHitsBgr)    delete fHitsBgr;
90   if (fSDigitsBgr) delete fSDigitsBgr;
91   if (fHitMap)     delete fHitMap;
92   if (fList)       delete fList;
93   if (fBgrFile)    delete fBgrFile;
94 }
95
96 //------------------------------------------------------------------------
97 Bool_t AliRICHMerger::Exists(const AliRICHSDigit *padhit)
98 {
99   return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
100 }
101
102 //------------------------------------------------------------------------
103 void AliRICHMerger::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   if (fSignal) {
119     track = fTrack;
120     charge = iqpad;
121   } else {
122     track = kBgTag;
123     charge = kBgTag;
124   }
125   pdigit->UpdateTrackList(track,charge);
126 }
127
128 //------------------------------------------------------------------------
129 void AliRICHMerger::CreateNew(AliRICHSDigit *padhit)
130 {
131   fList->AddAtAndExpand(
132     new AliRICHTransientDigit(fNch,fDigits),fCounter);
133   fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
134
135   AliRICHTransientDigit* pdigit = 
136     static_cast<AliRICHTransientDigit*>(fList->Last());
137   // list of tracks
138   Int_t track, charge;
139   if (fSignal) {
140     track = fTrack;
141     charge = padhit->QPad();
142   } else {
143     track = kBgTag;
144     charge = kBgTag;
145   }
146   pdigit->AddToTrackList(track,charge);
147   fCounter++;
148 }
149
150
151 //------------------------------------------------------------------------
152 void AliRICHMerger::Init()
153 {
154 // Initialisation
155     
156   if (fMerge && !fBgrFile) fBgrFile = InitBgr();
157 }
158
159
160
161 //------------------------------------------------------------------------
162 TFile* AliRICHMerger::InitBgr()
163 {
164 // Initialise background event
165   TFile *file = new TFile(fFnBgr);
166 // add error checking later
167   printf("\n AliRICHMerger has opened %s file with background event \n", fFnBgr);
168   fHitsBgr     = new TClonesArray("AliRICHHit",1000);
169   fSDigitsBgr  = new TClonesArray("AliRICHSDigit",1000);
170   return file;
171 }
172
173 //------------------------------------------------------------------------
174 void AliRICHMerger::Digitise(Int_t nev, Int_t flag)
175 {
176
177   // keep galice.root for signal and name differently the file for 
178   // background when add! otherwise the track info for signal will be lost !
179
180   Int_t particle;
181
182   AliRICHChamber*       iChamber;
183   AliSegmentation*  segmentation;
184   
185   fList = new TObjArray;
186     
187   AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
188   fHitMap= new AliHitMap* [kNCH];
189     
190   if (fMerge ) {
191     fBgrFile->cd();
192     //
193     // Get Hits Tree header from file
194     if(fHitsBgr) fHitsBgr->Clear();
195     if(fSDigitsBgr) fSDigitsBgr->Clear();
196     if(fTrH1) delete fTrH1;
197     fTrH1 = 0;
198       
199     char treeName[20];
200     sprintf(treeName,"TreeH%d",fEvNrBgr);
201     fTrH1 = (TTree*)gDirectory->Get(treeName);
202     if (!fTrH1) {
203       printf("ERROR: cannot find Hits Tree for event:%d\n",fEvNrBgr);
204     }
205     //
206     // Set branch addresses
207     TBranch *branch;
208     char branchname[20];
209     sprintf(branchname,"%s",pRICH->GetName());
210     if (fTrH1 && fHitsBgr) {
211       branch = fTrH1->GetBranch(branchname);
212       if (branch) branch->SetAddress(&fHitsBgr);
213     }
214     if (fTrH1 && fSDigitsBgr) {
215       branch = fTrH1->GetBranch("RICHSDigits");
216       if (branch) branch->SetAddress(&fSDigitsBgr);
217     }
218   }
219     
220   for (Int_t i =0; i<kNCH; i++) {
221     iChamber= &(pRICH->Chamber(i));
222     segmentation=iChamber->GetSegmentationModel(1);
223     fHitMap[i] = new AliRICHHitMapA1(segmentation, fList);
224   }
225   //
226   //   Loop over tracks
227   //
228     
229   fSignal = kTRUE;
230   fCounter = 0;
231   TTree *treeH = gAlice->TreeH();
232   Int_t ntracks =(Int_t) treeH->GetEntries();
233   for (fTrack=0; fTrack<ntracks; fTrack++) {
234     gAlice->ResetHits();
235     treeH->GetEvent(fTrack);
236     //
237     //   Loop over hits
238     for(AliRICHHit* mHit=(AliRICHHit*)pRICH->FirstHit(-1); 
239         mHit;
240         mHit=(AliRICHHit*)pRICH->NextHit()) 
241     {
242           
243       fNch = mHit->Chamber()-1;  // chamber number
244       Int_t trackIndex = mHit->Track();
245       if (fNch >= kNCH) {
246         cerr<<"AliRICHMerger: chamber nr. fNch out of range: "<<fNch<<endl;
247         cerr<<"               track: "<<fTrack<<endl;
248         continue;
249       }
250       iChamber = &(pRICH->Chamber(fNch));
251           
252
253 // 
254 // If flag is set, create digits only for some particles
255 //
256       TParticle *current = (TParticle*)gAlice->Particle(trackIndex);
257       if (current->GetPdgCode() >= 50000050)
258       {
259         TParticle *motherofcurrent = (TParticle*)gAlice->Particle(current->GetFirstMother());
260         particle = motherofcurrent->GetPdgCode();
261       }
262       else
263       {
264         particle = current->GetPdgCode();
265       }
266
267           
268       //printf("Flag:%d\n",flag);
269       //printf("Track:%d\n",track);
270       //printf("Particle:%d\n",particle);
271           
272       Int_t digitise=1;
273           
274       if (flag == 1) 
275         if(TMath::Abs(particle)==211 || TMath::Abs(particle)==111)
276           digitise=0;
277           
278       if (flag == 2)
279         if(TMath::Abs(particle)==321 || TMath::Abs(particle)==130 || TMath::Abs(particle)==310 
280            || TMath::Abs(particle)==311)
281           digitise=0;
282           
283       if (flag == 3 && TMath::Abs(particle)==2212)
284         digitise=0;
285           
286       if (flag == 4 && TMath::Abs(particle)==13)
287         digitise=0;
288           
289       if (flag == 5 && TMath::Abs(particle)==11)
290         digitise=0;
291           
292       if (flag == 6 && TMath::Abs(particle)==2112)
293         digitise=0;
294           
295           
296       //printf ("Particle: %d, Flag: %d, Digitise: %d\n",particle,flag,digitise); 
297           
298           
299       if (digitise)
300       {
301               
302         //
303         // Loop over pad hits
304         for (AliRICHSDigit* mPad=
305                (AliRICHSDigit*)pRICH->FirstPad(mHit,pRICH->SDigits());
306              mPad;
307              mPad=(AliRICHSDigit*)pRICH->NextPad(pRICH->SDigits()))
308         {
309           Int_t iqpad    = mPad->QPad();       // charge per pad
310           fDigits[0]=mPad->PadX();
311           fDigits[1]=mPad->PadY();
312           fDigits[2]=iqpad;
313           fDigits[3]=iqpad;
314           fDigits[4]=mPad->HitNumber();
315                   
316           // build the list of fired pads and update the info
317           if (Exists(mPad)) {
318             Update(mPad);
319           } else {
320             CreateNew(mPad);
321           }
322         } //end loop over clusters
323       }// track type condition
324     } // hit loop
325   } // track loop
326     
327   // open the file with background
328     
329   if (fMerge) {
330     fSignal = kFALSE;
331     ntracks = (Int_t)fTrH1->GetEntries();
332     //
333     //   Loop over tracks
334     //
335     for (fTrack = 0; fTrack < ntracks; fTrack++) {
336         
337       if (fHitsBgr)       fHitsBgr->Clear();
338       if (fSDigitsBgr)    fSDigitsBgr->Clear();
339         
340       fTrH1->GetEvent(fTrack);
341       //
342       //   Loop over hits
343       AliRICHHit* mHit;
344       for(Int_t i = 0; i < fHitsBgr->GetEntriesFast(); ++i) 
345       { 
346         mHit   = (AliRICHHit*) (*fHitsBgr)[i];
347         fNch   = mHit->Chamber()-1;  // chamber number
348         iChamber = &(pRICH->Chamber(fNch));
349
350         //
351         // Loop over pad hits
352         for (AliRICHSDigit* mPad =
353                (AliRICHSDigit*)pRICH->FirstPad(mHit,fSDigitsBgr);
354              mPad;
355              mPad = (AliRICHSDigit*)pRICH->NextPad(fSDigitsBgr))
356         {
357           fDigits[0] = mPad->PadX(); 
358           fDigits[1] = mPad->PadY();  
359           fDigits[2] = mPad->QPad();
360           fDigits[3] = 0;
361           fDigits[4] = -1; // set hit number to -1 for bgr
362                 
363           // build the list of fired pads and update the info
364           if (!Exists(mPad)) {
365             CreateNew(mPad);
366           } else {
367             Update(mPad);
368           } //  end if !Exists
369         } //end loop over clusters
370       } // hit loop
371     } // track loop
372       
373     TTree *treeK = gAlice->TreeK();
374     TFile *file = NULL;
375       
376     if (treeK) file = treeK->GetCurrentFile();
377     file->cd();
378   } // if fMerge
379     
380   Int_t tracks[kMAXTRACKSPERRICHDIGIT];
381   Int_t charges[kMAXTRACKSPERRICHDIGIT];
382   Int_t nentries=fList->GetEntriesFast();
383     
384   for (Int_t nent=0;nent<nentries;nent++) {
385     AliRICHTransientDigit *transDigit=(AliRICHTransientDigit*)fList->At(nent);
386     if (transDigit==0) continue; 
387     Int_t ich=transDigit->GetChamber();
388     Int_t q=transDigit->Signal(); 
389     iChamber=&(pRICH->Chamber(ich));
390     AliRICHResponse * response=iChamber->GetResponseModel();
391     Int_t adcmax= (Int_t) response->MaxAdc();
392       
393       
394     // add white noise and do zero-suppression and signal truncation (new electronics,old electronics gaus 1.2,0.2)
395     //printf("Treshold: %d\n",iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY()));
396     Int_t pedestal = iChamber->fTresh->GetHitIndex(transDigit->PadX(),transDigit->PadY());
397
398     //printf("Pedestal:%d\n",pedestal);
399     //Int_t pedestal=0;
400     Float_t treshold = (pedestal + 4*2.2);
401       
402     Float_t meanNoise = gRandom->Gaus(2.2, 0.3);
403     Float_t noise     = gRandom->Gaus(0, meanNoise);
404     q+=(Int_t)(noise + pedestal);
405     //q+=(Int_t)(noise);
406     //          magic number to be parametrised !!! 
407     if ( q <= treshold) 
408     {
409       q = q - pedestal;
410       continue;
411     }
412     q = q - pedestal;
413     if ( q >= adcmax) q=adcmax;
414     fDigits[0]=transDigit->PadX();
415     fDigits[1]=transDigit->PadY();
416     fDigits[2]=q;
417     fDigits[3]=transDigit->Physics();
418     fDigits[4]=transDigit->Hit();
419
420     Int_t nptracks = transDigit->GetNTracks();  
421
422     // this was changed to accomodate the real number of tracks
423     if (nptracks > kMAXTRACKSPERRICHDIGIT) {
424       printf("Attention - tracks > 10 %d\n",nptracks);
425       nptracks=kMAXTRACKSPERRICHDIGIT;
426     }
427     if (nptracks > 2) {
428       printf("Attention - tracks > 2  %d \n",nptracks);
429       //printf("cat,ich,ix,iy,q %d %d %d %d %d \n",
430       //icat,ich,digits[0],digits[1],q);
431     }
432     for (Int_t tr=0;tr<nptracks;tr++) {
433       tracks[tr]=transDigit->GetTrack(tr);
434       charges[tr]=transDigit->GetCharge(tr);
435     }      //end loop over list of tracks for one pad
436     if (nptracks < kMAXTRACKSPERRICHDIGIT ) {
437       for (Int_t t=nptracks; t<kMAXTRACKSPERRICHDIGIT; t++) {
438         tracks[t]=0;
439         charges[t]=0;
440       }
441     }
442     //write file
443     //if (ich==2)
444     //fprintf(points,"%4d,      %4d,      %4d\n",digits[0],digits[1],digits[2]);
445       
446     // fill digits
447     pRICH->AddDigits(ich,tracks,charges,fDigits);
448   }     
449   gAlice->TreeD()->Fill();
450
451   fList->Delete();
452   for(Int_t ii=0;ii<kNCH;++ii) {
453     if (fHitMap[ii]) {
454       delete fHitMap[ii];
455       fHitMap[ii]=0;
456     }
457   }
458     
459   TClonesArray *richDigits;
460   for (Int_t k=0;k<kNCH;k++) {
461     richDigits = pRICH->DigitsAddress(k);
462     Int_t ndigit=richDigits->GetEntriesFast();
463     printf ("Chamber %d digits %d \n",k,ndigit);
464   }
465   pRICH->ResetDigits(); /// ??? should it be here???
466   gAlice->TreeD()->Write(0,TObject::kOverwrite);
467   // reset tree
468   //    gAlice->TreeD()->Reset();
469   //    delete fList; // deleted in dtor
470   // gObjectTable->Print();
471
472 }