]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHMerger.cxx
removed a cout
[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.2  2001/03/14 18:16:08  jbarbosa
19 Corrected bug (more to correct).
20 File "points.dat" is no longer created.
21
22 Revision 1.1  2001/02/27 22:13:34  jbarbosa
23 Implementing merger class.
24
25 */
26
27 #include <TTree.h> 
28 #include <TVector.h>
29 #include <TObjArray.h>
30 #include <TFile.h>
31 #include <TDirectory.h>
32 #include <TParticle.h>
33
34
35 // #include "AliMerger.h"
36 // #include "AliMergable.h"
37 #include "AliRICHMerger.h"
38 #include "AliRICHChamber.h"
39 #include "AliHitMap.h"
40 #include "AliRICHHitMapA1.h"
41 #include "AliRICH.h"
42 #include "AliRICHHit.h"
43 #include "AliRICHSDigit.h"
44 #include "AliRICHDigit.h"
45 #include "AliRICHTransientDigit.h"
46 #include "AliRun.h"
47 #include "AliPDG.h"
48
49 ClassImp(AliRICHMerger)
50
51 //___________________________________________
52 AliRICHMerger::AliRICHMerger()
53 {
54 // Default constructor    
55     fEvNrSig = 0;
56     fEvNrBgr = 0;
57     fMerge   = kDigitize;
58     fFnBgr   = 0;
59     fHitMap  = 0;
60     fList    = 0;
61     fTrH1    = 0;
62     fHitsBgr = 0;
63     fSDigitsBgr = 0;
64     fHitMap     = 0;
65     fList       = 0;
66     fTrList     = 0;
67     fAddress    = 0; 
68 }
69
70 //------------------------------------------------------------------------
71 AliRICHMerger::~AliRICHMerger()
72 {
73 // Destructor
74     if (fTrH1)       delete fTrH1;
75     if (fHitsBgr)    delete fHitsBgr;
76     if (fSDigitsBgr) delete fSDigitsBgr;
77     if (fHitMap)     delete fHitMap;
78     if (fList)       delete fList;
79     if (fTrList)     delete fTrList;
80     if (fAddress)    delete fAddress; 
81 }
82
83 //------------------------------------------------------------------------
84 Bool_t AliRICHMerger::Exists(const AliRICHSDigit *mergable)
85 {
86     AliRICHSDigit *padhit = (AliRICHSDigit*) mergable;
87     return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
88 }
89
90 //------------------------------------------------------------------------
91 void AliRICHMerger::Update(AliRICHSDigit *mergable)
92 {
93     AliRICHSDigit *padhit = (AliRICHSDigit*) mergable;    
94     AliRICHTransientDigit* pdigit;
95     Int_t ipx      = padhit->PadX();        // pad number on X
96     Int_t ipy      = padhit->PadY();        // pad number on Y
97     Int_t iqpad    = Int_t(padhit->QPad()); // charge per pad
98
99     pdigit = (AliRICHTransientDigit*) fHitMap[fNch]->GetHit(ipx, ipy);
100     // update charge
101     //
102     (*pdigit).AddSignal(iqpad);
103     (*pdigit).AddPhysicsSignal(iqpad);          
104     // update list of tracks
105     //
106     TObjArray* fTrList = (TObjArray*)pdigit->TrackList();
107     Int_t lastEntry = fTrList->GetLast();
108     TVector *pTrack = (TVector*)fTrList->At(lastEntry);
109     TVector &ptrk   = *pTrack;
110     TVector &trinfo = *((TVector*) (*fAddress)[fCountadr-1]);
111     Int_t lastTrack = Int_t(ptrk(0));
112
113     if (trinfo(0) != kBgTag) {
114         if (lastTrack == fTrack) {
115             Int_t lastCharge = Int_t(ptrk(1));
116             lastCharge += iqpad;
117             fTrList->RemoveAt(lastEntry);
118             trinfo(1) = lastCharge;
119             fTrList->AddAt(&trinfo,lastEntry);
120         } else {
121             fTrList->Add(&trinfo);
122         }
123     } else {
124         if (lastTrack != -1) fTrList->Add(&trinfo);
125     }
126 }
127
128 //------------------------------------------------------------------------
129 void AliRICHMerger::CreateNew(AliRICHSDigit *mergable)
130 {
131     AliRICHSDigit *padhit = (AliRICHSDigit*) mergable;    
132     AliRICHTransientDigit* pdigit;
133
134     Int_t ipx      = padhit->PadX();       // pad number on X
135     Int_t ipy      = padhit->PadY();       // pad number on Y
136     fList->AddAtAndExpand(
137         new AliRICHTransientDigit(fNch,fDigits),fCounter);
138     fHitMap[fNch]->SetHit(ipx, ipy, fCounter);
139     fCounter++;
140     pdigit = (AliRICHTransientDigit*)fList->At(fList->GetLast());
141     // list of tracks
142     TObjArray *fTrList = (TObjArray*)pdigit->TrackList();
143     TVector &trinfo    = *((TVector*) (*fAddress)[fCountadr-1]);
144     fTrList->Add(&trinfo);
145 }
146
147
148 //------------------------------------------------------------------------
149 void AliRICHMerger::Init()
150 {
151 // Initialisation
152     
153     if (fMerge) 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   //FILE* points; //these will be the digits...
180   
181   //points=fopen("points.dat","w");
182   
183   AliRICHChamber*       iChamber;
184   AliSegmentation*  segmentation;
185   
186     Int_t digitise=0;
187     Int_t trk[50];
188     Int_t chtrk[50];  
189     TObjArray *list=new TObjArray;
190     static TClonesArray *pAddress=0;
191     if(!pAddress) pAddress=new TClonesArray("TVector",1000);
192     Int_t digits[5]; 
193     
194     AliRICH *pRICH = (AliRICH *) gAlice->GetDetector("RICH");
195     AliHitMap* pHitMap[10];
196     Int_t i;
197     for (i=0; i<10; i++) {pHitMap[i]=0;}
198     
199     if (fMerge ) {
200       fBgrFile->cd();
201       //fBgrFile->ls();
202       //
203       // Get Hits Tree header from file
204       if(fHitsBgr) fHitsBgr->Clear();
205       if(fSDigitsBgr) fSDigitsBgr->Clear();
206       if(fTrH1) delete fTrH1;
207       fTrH1 = 0;
208       
209       char treeName[20];
210       sprintf(treeName,"TreeH%d",fEvNrBgr);
211       fTrH1 = (TTree*)gDirectory->Get(treeName);
212       if (!fTrH1) {
213         printf("ERROR: cannot find Hits Tree for event:%d\n",fEvNrBgr);
214       }
215       //
216       // Set branch addresses
217       TBranch *branch;
218       char branchname[20];
219       sprintf(branchname,"%s",pRICH->GetName());
220       if (fTrH1 && fHitsBgr) {
221         branch = fTrH1->GetBranch(branchname);
222         if (branch) branch->SetAddress(&fHitsBgr);
223       }
224       if (fTrH1 && fSDigitsBgr) {
225         branch = fTrH1->GetBranch("RICHSDigits");
226         if (branch) branch->SetAddress(&fSDigitsBgr);
227       }
228     }
229     
230     AliHitMap* hm;
231     Int_t countadr=0;
232     Int_t counter=0;
233     for (i =0; i<kNCH; i++) {
234       iChamber= &(pRICH->Chamber(i));
235       segmentation=iChamber->GetSegmentationModel(1);
236       pHitMap[i] = new AliRICHHitMapA1(segmentation, list);
237     }
238     //
239     //   Loop over tracks
240     //
241     
242     TTree *treeH = gAlice->TreeH();
243     Int_t ntracks =(Int_t) treeH->GetEntries();
244     for (Int_t track=0; track<ntracks; track++) {
245       gAlice->ResetHits();
246       treeH->GetEvent(track);
247       //
248       //   Loop over hits
249       for(AliRICHHit* mHit=(AliRICHHit*)pRICH->FirstHit(-1); 
250           mHit;
251           mHit=(AliRICHHit*)pRICH->NextHit()) 
252         {
253           
254           Int_t   nch   = mHit->fChamber-1;  // chamber number
255           Int_t   index = mHit->Track();
256           if (nch >kNCH) continue;
257           iChamber = &(pRICH->Chamber(nch));
258           
259           TParticle *current = (TParticle*)gAlice->Particle(index);
260           
261           if (current->GetPdgCode() >= 50000050)
262             {
263               TParticle *motherofcurrent = (TParticle*)gAlice->Particle(current->GetFirstMother());
264               particle = motherofcurrent->GetPdgCode();
265             }
266           else
267             {
268               particle = current->GetPdgCode();
269             }
270
271           
272           //printf("Flag:%d\n",flag);
273           //printf("Track:%d\n",track);
274           //printf("Particle:%d\n",particle);
275           
276           digitise=1;
277           
278           if (flag == 1) 
279             if(TMath::Abs(particle)==211 || TMath::Abs(particle)==111)
280               digitise=0;
281           
282           if (flag == 2)
283             if(TMath::Abs(particle)==321 || TMath::Abs(particle)==130 || TMath::Abs(particle)==310 
284                || TMath::Abs(particle)==311)
285               digitise=0;
286           
287           if (flag == 3 && TMath::Abs(particle)==2212)
288             digitise=0;
289           
290           if (flag == 4 && TMath::Abs(particle)==13)
291             digitise=0;
292           
293           if (flag == 5 && TMath::Abs(particle)==11)
294             digitise=0;
295           
296           if (flag == 6 && TMath::Abs(particle)==2112)
297             digitise=0;
298           
299           
300           //printf ("Particle: %d, Flag: %d, Digitise: %d\n",particle,flag,digitise); 
301           
302           
303           if (digitise)
304             {
305               
306               //
307               // Loop over pad hits
308               for (AliRICHSDigit* mPad=
309                      (AliRICHSDigit*)pRICH->FirstPad(mHit,pRICH->SDigits());
310                    mPad;
311                    mPad=(AliRICHSDigit*)pRICH->NextPad(pRICH->SDigits()))
312                 {
313                   Int_t ipx      = mPad->fPadX;       // pad number on X
314                   Int_t ipy      = mPad->fPadY;       // pad number on Y
315                   Int_t iqpad    = mPad->fQpad;       // charge per pad
316                   //
317                   //
318                   //printf("X:%d, Y:%d, Q:%d\n",ipx,ipy,iqpad);
319                   
320                   Float_t thex, they, thez;
321                   segmentation=iChamber->GetSegmentationModel(0);
322                   segmentation->GetPadC(ipx,ipy,thex,they,thez);
323                   new((*pAddress)[countadr++]) TVector(2);
324                   TVector &trinfo=*((TVector*) (*pAddress)[countadr-1]);
325                   trinfo(0)=(Float_t)track;
326                   trinfo(1)=(Float_t)iqpad;
327                   
328                   digits[0]=ipx;
329                   digits[1]=ipy;
330                   digits[2]=iqpad;
331                   
332                   AliRICHTransientDigit* pdigit;
333                   // build the list of fired pads and update the info
334                   if (!pHitMap[nch]->TestHit(ipx, ipy)) {
335                     list->AddAtAndExpand(new AliRICHTransientDigit(nch,digits),counter);
336                     pHitMap[nch]->SetHit(ipx, ipy, counter);
337                     counter++;
338                     pdigit=(AliRICHTransientDigit*)list->At(list->GetLast());
339                     // list of tracks
340                     TObjArray *trlist=(TObjArray*)pdigit->TrackList();
341                     trlist->Add(&trinfo);
342                   } else {
343                     pdigit=(AliRICHTransientDigit*) pHitMap[nch]->GetHit(ipx, ipy);
344                     // update charge
345                     (*pdigit).fSignal+=iqpad;
346                     // update list of tracks
347                     TObjArray* trlist=(TObjArray*)pdigit->TrackList();
348                     Int_t lastEntry=trlist->GetLast();
349                     TVector *ptrkP=(TVector*)trlist->At(lastEntry);
350                     TVector &ptrk=*ptrkP;
351                     Int_t lastTrack=Int_t(ptrk(0));
352                     Int_t lastCharge=Int_t(ptrk(1));
353                     if (lastTrack==track) {
354                       lastCharge+=iqpad;
355                       trlist->RemoveAt(lastEntry);
356                       trinfo(0)=lastTrack;
357                       trinfo(1)=lastCharge;
358                       trlist->AddAt(&trinfo,lastEntry);
359                     } else {
360                       trlist->Add(&trinfo);
361                     }
362                     // check the track list
363                     Int_t nptracks=trlist->GetEntriesFast();
364                     if (nptracks > 2) {
365                       printf("Attention - tracks:  %d (>2)\n",nptracks);
366                       //printf("cat,nch,ix,iy %d %d %d %d  \n",icat+1,nch,ipx,ipy);
367                       for (Int_t tr=0;tr<nptracks;tr++) {
368                         TVector *pptrkP=(TVector*)trlist->At(tr);
369                         TVector &pptrk=*pptrkP;
370                         trk[tr]=Int_t(pptrk(0));
371                         chtrk[tr]=Int_t(pptrk(1));
372                       }
373                     } // end if nptracks
374                   } //  end if pdigit
375                 } //end loop over clusters
376             }// track type condition
377         } // hit loop
378     } // track loop
379     
380     // open the file with background
381     
382     if (fMerge) {
383       ntracks = (Int_t)fTrH1->GetEntries();
384       //
385       //   Loop over tracks
386       //
387       for (fTrack = 0; fTrack < ntracks; fTrack++) {
388         
389         if (fHitsBgr)       fHitsBgr->Clear();
390         if (fSDigitsBgr)    fSDigitsBgr->Clear();
391         
392         fTrH1->GetEvent(fTrack);
393         //
394         //   Loop over hits
395         AliRICHHit* mHit;
396         for(int i = 0; i < fHitsBgr->GetEntriesFast(); ++i) 
397           {     
398             mHit   = (AliRICHHit*) (*fHitsBgr)[i];
399             fNch   = mHit->Chamber()-1;  // chamber number
400             iChamber = &(pRICH->Chamber(fNch));
401             //Float_t xbgr = mHit->X();
402             //Float_t ybgr = mHit->Y();
403             Bool_t cond  = kFALSE;
404             cond  = kTRUE;
405             //
406             // Loop over pad hits
407             for (AliRICHSDigit* mPad =
408                    (AliRICHSDigit*)pRICH->FirstPad(mHit,fSDigitsBgr);
409                  mPad;
410                  mPad = (AliRICHSDigit*)pRICH->NextPad(fSDigitsBgr))
411               {
412                 Int_t ipx      = mPad->PadX();        // pad number on X
413                 Int_t ipy      = mPad->PadY();        // pad number on Y
414                 Int_t iqpad    = Int_t(mPad->QPad()); // charge per pad
415                 
416                 new((*fAddress)[fCountadr++]) TVector(2);
417                 TVector &trinfo = *((TVector*) (*fAddress)[fCountadr-1]);
418                 trinfo(0) = kBgTag;  // tag background
419                 trinfo(1) = kBgTag;
420                 
421                 fDigits[0] = ipx;
422                 fDigits[1] = ipy;
423                 fDigits[3] = iqpad;
424                 
425                 // build the list of fired pads and update the info
426                 if (!Exists(mPad)) {
427                   CreateNew(mPad);
428                 } else {
429                   Update(mPad);
430                 } //  end if !Exists
431               } //end loop over clusters
432           } // hit loop
433       } // track loop
434       
435       TTree *fAli = gAlice->TreeK();
436       TFile *file = NULL;
437       
438       if (fAli) file = fAli->GetCurrentFile();
439       file->cd();
440     } // if fMerge
441     
442     Int_t tracks[10];
443     Int_t charges[10];
444     //cout<<"Start filling digits \n "<<endl;
445     Int_t nentries=list->GetEntriesFast();
446     //printf(" \n \n nentries %d \n",nentries);
447     
448     // start filling the digits
449     
450     for (Int_t nent=0;nent<nentries;nent++) {
451       AliRICHTransientDigit *address=(AliRICHTransientDigit*)list->At(nent);
452       if (address==0) continue; 
453       
454       Int_t ich=address->fChamber;
455       Int_t q=address->fSignal; 
456       iChamber=&(pRICH->Chamber(ich));
457       AliRICHResponse * response=iChamber->GetResponseModel();
458       Int_t adcmax= (Int_t) response->MaxAdc();
459       
460       
461       // add white noise and do zero-suppression and signal truncation (new electronics,old electronics gaus 1.2,0.2)
462       //printf("Treshold: %d\n",iChamber->fTresh->GetHitIndex(address->fPadX,address->fPadY));
463       Int_t pedestal = iChamber->fTresh->GetHitIndex(address->fPadX,address->fPadY);
464
465       //printf("Pedestal:%d\n",pedestal);
466       //Int_t pedestal=0;
467       Float_t treshold = (pedestal + 4*2.2);
468       
469       Float_t meanNoise = gRandom->Gaus(2.2, 0.3);
470       Float_t noise     = gRandom->Gaus(0, meanNoise);
471       q+=(Int_t)(noise + pedestal);
472       //q+=(Int_t)(noise);
473       //          magic number to be parametrised !!! 
474       if ( q <= treshold) 
475         {
476           q = q - pedestal;
477           continue;
478         }
479       q = q - pedestal;
480       if ( q >= adcmax) q=adcmax;
481       digits[0]=address->fPadX;
482       digits[1]=address->fPadY;
483       digits[2]=q;
484       
485       TObjArray* trlist=(TObjArray*)address->TrackList();
486       Int_t nptracks=trlist->GetEntriesFast();
487       
488       // this was changed to accomodate the real number of tracks
489       if (nptracks > 10) {
490         printf("Attention - tracks > 10 %d\n",nptracks);
491         nptracks=10;
492       }
493       if (nptracks > 2) {
494         printf("Attention - tracks > 2  %d \n",nptracks);
495         //printf("cat,ich,ix,iy,q %d %d %d %d %d \n",
496         //icat,ich,digits[0],digits[1],q);
497       }
498       for (Int_t tr=0;tr<nptracks;tr++) {
499         TVector *ppP=(TVector*)trlist->At(tr);
500         TVector &pp  =*ppP;
501         tracks[tr]=Int_t(pp(0));
502         charges[tr]=Int_t(pp(1));
503       }      //end loop over list of tracks for one pad
504       if (nptracks < 10 ) {
505         for (Int_t t=nptracks; t<10; t++) {
506           tracks[t]=0;
507           charges[t]=0;
508         }
509       }
510       //write file
511       //if (ich==2)
512         //fprintf(points,"%4d,      %4d,      %4d\n",digits[0],digits[1],digits[2]);
513       
514       // fill digits
515       pRICH->AddDigits(ich,tracks,charges,digits);
516     }   
517     gAlice->TreeD()->Fill();
518     
519     list->Delete();
520     for(Int_t ii=0;ii<kNCH;++ii) {
521       if (pHitMap[ii]) {
522         hm=pHitMap[ii];
523         delete hm;
524         pHitMap[ii]=0;
525       }
526     }
527     
528     //TTree *TD=gAlice->TreeD();
529     //Stat_t ndig=TD->GetEntries();
530     //cout<<"number of digits  "<<ndig<<endl;
531     TClonesArray *fDch;
532     for (int k=0;k<kNCH;k++) {
533       fDch= pRICH->DigitsAddress(k);
534       int ndigit=fDch->GetEntriesFast();
535       printf ("Chamber %d digits %d \n",k,ndigit);
536     }
537     pRICH->ResetDigits();
538     // char hname[30];
539     // sprintf(hname,"TreeD%d",nev);
540     // gAlice->TreeD()->Write(hname);
541     gAlice->TreeD()->Write(0,TObject::kOverwrite);
542     // reset tree
543     //    gAlice->TreeD()->Reset();
544     delete list;
545     pAddress->Clear();
546     // gObjectTable->Print();
547
548 }
549
550
551
552
553