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