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