]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONMerger.cxx
8c4dc4f6c8097e7fdb94397bf37b87bd4d3e8df9
[u/mrichter/AliRoot.git] / MUON / AliMUONMerger.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.4  2001/03/20 13:36:11  egangler
19 TFile memory leak and "too many files open" problem solved (?):
20 the backgroud file is open only once forever, and not once for each event.
21
22 Revision 1.3  2001/03/05 23:57:44  morsch
23 Writing of digit tree moved to macro.
24
25 Revision 1.2  2001/03/05 08:40:25  morsch
26 Method SortTracks(..) imported from AliMUON.
27
28 Revision 1.1  2001/02/02 14:11:53  morsch
29 AliMUONMerger prototype to be called by the merge manager.
30
31 */
32
33 #include <TTree.h> 
34 #include <TVector.h>
35 #include <TObjArray.h>
36 #include <TFile.h>
37 #include <TDirectory.h>
38
39
40 // #include "AliMerger.h"
41 // #include "AliMergable.h"
42 #include "AliMUONMerger.h"
43 #include "AliMUONConstants.h"
44 #include "AliMUONChamber.h"
45 #include "AliHitMap.h"
46 #include "AliMUONHitMapA1.h"
47 #include "AliMUON.h"
48 #include "AliMUONHit.h"
49 #include "AliMUONPadHit.h"
50 #include "AliMUONDigit.h"
51 #include "AliMUONTransientDigit.h"
52 #include "AliRun.h"
53 #include "AliPDG.h"
54
55 ClassImp(AliMUONMerger)
56
57 //___________________________________________
58 AliMUONMerger::AliMUONMerger()
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     fPadHitsBgr = 0;
70     fHitMap     = 0;
71     fList       = 0;
72     fBgrFile    = 0;
73 }
74
75 //------------------------------------------------------------------------
76 AliMUONMerger::~AliMUONMerger()
77 {
78 // Destructor
79     if (fTrH1)       delete fTrH1;
80     if (fHitsBgr)    delete fHitsBgr;
81     if (fPadHitsBgr) delete fPadHitsBgr;
82     if (fHitMap)     delete fHitMap;
83     if (fList)       delete fList;
84     if (fBgrFile)    delete fBgrFile;
85 }
86
87 //------------------------------------------------------------------------
88 Bool_t AliMUONMerger::Exists(const AliMUONPadHit *mergable)
89 {
90     AliMUONPadHit *padhit = (AliMUONPadHit*) mergable;
91     return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
92 }
93
94 //------------------------------------------------------------------------
95 void AliMUONMerger::Update(AliMUONPadHit *padhit)
96 {
97     AliMUONTransientDigit *pdigit = 
98       static_cast<AliMUONTransientDigit*>(
99       fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
100
101     // update charge
102     //
103     Int_t iqpad    = padhit->QPad();        // charge per pad
104     pdigit->AddSignal(iqpad);
105     pdigit->AddPhysicsSignal(iqpad);            
106
107     // update list of tracks
108     //
109     Int_t track, charge;    
110     if (fSignal) {
111       track = fTrack;
112       charge = iqpad;
113     } else {
114       track = kBgTag;
115       charge = kBgTag;
116     }
117     pdigit->UpdateTrackList(track,charge);
118 }
119
120 //------------------------------------------------------------------------
121 void AliMUONMerger::CreateNew(AliMUONPadHit *padhit)
122 {
123     fList->AddAtAndExpand(
124         new AliMUONTransientDigit(fNch,fDigits),fCounter);
125     fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
126     AliMUONTransientDigit* pdigit = 
127       (AliMUONTransientDigit*)fList->At(fList->GetLast());
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 AliMUONMerger::Init()
144 {
145 // Initialisation
146     // open only once the background file !!
147     if (fMerge && !fBgrFile) fBgrFile = InitBgr();
148 }
149
150
151
152 //------------------------------------------------------------------------
153 TFile* AliMUONMerger::InitBgr()
154 {
155 // Initialise background event
156     TFile *file = new TFile(fFnBgr);
157 // add error checking later
158     printf("\n AliMUONMerger has opened %s file with background event \n", fFnBgr);
159     fHitsBgr     = new TClonesArray("AliMUONHit",1000);
160     fPadHitsBgr  = new TClonesArray("AliMUONPadHit",1000);
161     return file;
162 }
163
164 //------------------------------------------------------------------------
165 void AliMUONMerger::Digitise()
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     AliMUONChamber*   iChamber;
172     AliSegmentation*  segmentation;
173
174     fList = new TObjArray;
175
176     AliMUON *pMUON  = (AliMUON *) gAlice->GetModule("MUON");
177     fHitMap= new AliHitMap* [AliMUONConstants::NCh()];
178     if (fMerge ) {
179         fBgrFile->cd();
180         //
181         // Get Hits Tree header from file
182         //if(fHitsBgr) fHitsBgr->Clear();     // Useless because line 327
183         //if(fPadHitsBgr) fPadHitsBgr->Clear(); // Useless because line 328
184         if(fTrH1) delete fTrH1;
185         fTrH1 = 0;
186         
187         char treeName[20];
188         sprintf(treeName,"TreeH%d",fEvNrBgr);
189         fTrH1 = (TTree*)gDirectory->Get(treeName);
190         if (!fTrH1) {
191             printf("ERROR: cannot find Hits Tree for event:%d\n",fEvNrBgr);
192         }
193         //
194         // Set branch addresses
195         TBranch *branch;
196         char branchname[20];
197         sprintf(branchname,"%s",pMUON->GetName());
198         if (fTrH1 && fHitsBgr) {
199             branch = fTrH1->GetBranch(branchname);
200             if (branch) branch->SetAddress(&fHitsBgr);
201         }
202         if (fTrH1 && fPadHitsBgr) {
203             branch = fTrH1->GetBranch("MUONCluster");
204             if (branch) branch->SetAddress(&fPadHitsBgr);
205         }
206     }
207     //
208     // loop over cathodes
209     //
210     AliHitMap* hm;
211     fSignal = kTRUE;
212     for (int icat = 0; icat < 2; icat++) { 
213         fCounter = 0;
214         Int_t * nmuon = new Int_t [AliMUONConstants::NCh()];
215         for (Int_t i = 0; i < AliMUONConstants::NCh(); i++) {
216             iChamber = &(pMUON->Chamber(i));
217             if (iChamber->Nsec() == 1 && icat == 1) {
218                 continue;
219             } else {
220                 segmentation = iChamber->SegmentationModel(icat+1);
221             }
222             fHitMap[i] = new AliMUONHitMapA1(segmentation, fList);
223             nmuon[i] = 0;
224         }
225
226 //
227 //   Loop over tracks
228 //
229
230         TTree *treeH  = gAlice->TreeH();
231         Int_t ntracks = (Int_t) treeH->GetEntries();
232
233         for (fTrack = 0; fTrack < ntracks; fTrack++) {
234             gAlice->ResetHits();
235             treeH->GetEvent(fTrack);
236 //
237 //   Loop over hits
238             for(AliMUONHit* mHit = (AliMUONHit*)pMUON->FirstHit(-1); 
239                 mHit;
240                 mHit = (AliMUONHit*)pMUON->NextHit()) 
241             {
242                 fNch = mHit->Chamber()-1;  // chamber number
243                 if (fNch > AliMUONConstants::NCh()-1) continue;
244                 iChamber = &(pMUON->Chamber(fNch));
245                 
246 //
247 // Loop over pad hits
248                 for (AliMUONPadHit* mPad =
249                          (AliMUONPadHit*)pMUON->FirstPad(mHit,pMUON->PadHits());
250                      mPad;
251                      mPad = (AliMUONPadHit*)pMUON->NextPad(pMUON->PadHits()))
252                 {
253                     Int_t cathode  = mPad->Cathode();      // cathode number
254                     Int_t ipx      = mPad->PadX();         // pad number on X
255                     Int_t ipy      = mPad->PadY();         // pad number on Y
256                     Int_t iqpad    = Int_t(mPad->QPad());  // charge per pad
257                     if (cathode != (icat+1)) continue;
258
259                     segmentation = iChamber->SegmentationModel(cathode);
260
261                     fDigits[0] = ipx;
262                     fDigits[1] = ipy;
263                     fDigits[2] = icat;
264                     fDigits[3] = iqpad;
265                     fDigits[4] = iqpad;
266                     if (mHit->Particle() == kMuonPlus ||
267                         mHit->Particle() == kMuonMinus) {
268                         fDigits[5] = mPad->HitNumber();
269                     } else fDigits[5] = -1;
270
271                     // build the list of fired pads and update the info
272
273                     if (!Exists(mPad)) {
274                         CreateNew(mPad);
275                     } else {
276                         Update(mPad);
277                     } //  end if pdigit
278                 } //end loop over clusters
279             } // hit loop
280         } // track loop
281
282         // open the file with background
283        
284         if (fMerge) {
285             fSignal = kFALSE;
286             ntracks = (Int_t)fTrH1->GetEntries();
287 //
288 //   Loop over tracks
289 //
290             for (fTrack = 0; fTrack < ntracks; fTrack++) {
291
292                 if (fHitsBgr)       fHitsBgr->Clear();
293                 if (fPadHitsBgr)    fPadHitsBgr->Clear();
294
295                 fTrH1->GetEvent(fTrack);
296 //
297 //   Loop over hits
298                 AliMUONHit* mHit;
299                 for(int i = 0; i < fHitsBgr->GetEntriesFast(); ++i) 
300                 {       
301                     mHit   = (AliMUONHit*) (*fHitsBgr)[i];
302                     fNch   = mHit->Chamber()-1;  // chamber number
303                     iChamber = &(pMUON->Chamber(fNch));
304 //
305 // Loop over pad hits
306                     for (AliMUONPadHit* mPad =
307                              (AliMUONPadHit*)pMUON->FirstPad(mHit,fPadHitsBgr);
308                          mPad;
309                          mPad = (AliMUONPadHit*)pMUON->NextPad(fPadHitsBgr))
310                     {
311                         Int_t cathode  = mPad->Cathode();     // cathode number
312                         Int_t ipx      = mPad->PadX();        // pad number on X
313                         Int_t ipy      = mPad->PadY();        // pad number on Y
314                         Int_t iqpad    = Int_t(mPad->QPad()); // charge per pad
315
316                         if (cathode != (icat+1)) continue;
317                         fDigits[0] = ipx;
318                         fDigits[1] = ipy;
319                         fDigits[2] = icat;
320                         fDigits[3] = iqpad;
321                         fDigits[4] = 0;
322                         fDigits[5] = -1;
323                         
324                         // build the list of fired pads and update the info
325                         if (!Exists(mPad)) {
326                             CreateNew(mPad);
327                         } else {
328                             Update(mPad);
329                         } //  end if !Exists
330                     } //end loop over clusters
331                 } // hit loop
332             } // track loop
333
334             TTree *fAli = gAlice->TreeK();
335             TFile *file = NULL;
336             
337             if (fAli) file = fAli->GetCurrentFile();
338             file->cd();
339         } // if fMerge
340
341         Int_t tracks[kMAXTRACKS];
342         Int_t charges[kMAXTRACKS];
343         Int_t nentries = fList->GetEntriesFast();
344         
345         for (Int_t nent = 0; nent < nentries; nent++) {
346             AliMUONTransientDigit *address = (AliMUONTransientDigit*)fList->At(nent);
347             if (address == 0) continue; 
348             Int_t ich = address->Chamber();
349             Int_t   q = address->Signal(); 
350             iChamber = &(pMUON->Chamber(ich));
351 //
352 //  Digit Response (noise, threshold, saturation, ...)
353             AliMUONResponse * response = iChamber->ResponseModel();
354             q = response->DigitResponse(q);
355             
356             if (!q) continue;
357             
358             fDigits[0] = address->PadX();
359             fDigits[1] = address->PadY();
360             fDigits[2] = address->Cathode();
361             fDigits[3] = q;
362             fDigits[4] = address->Physics();
363             fDigits[5] = address->Hit();
364             
365             Int_t nptracks = address->GetNTracks();
366
367             if (nptracks > kMAXTRACKS) {
368                 printf("\n Attention - nptracks > kMAXTRACKS %d \n", nptracks);
369                 nptracks = kMAXTRACKS;
370             }
371             if (nptracks > 2) {
372                 printf("Attention - nptracks > 2  %d \n",nptracks);
373                 printf("cat,ich,ix,iy,q %d %d %d %d %d \n",icat,ich,fDigits[0],fDigits[1],q);
374             }
375             for (Int_t tr = 0; tr < nptracks; tr++) {
376                 tracks[tr]   = address->GetTrack(tr);
377                 charges[tr]  = address->GetCharge(tr);
378             }      //end loop over list of tracks for one pad
379             // Sort list of tracks according to charge
380             if (nptracks > 1) {
381                 SortTracks(tracks,charges,nptracks);
382             }
383             if (nptracks < kMAXTRACKS ) {
384                 for (Int_t i = nptracks; i < kMAXTRACKS; i++) {
385                     tracks[i]  = 0;
386                     charges[i] = 0;
387                 }
388             }
389             
390             // fill digits
391             pMUON->AddDigits(ich,tracks,charges,fDigits);
392         }
393         gAlice->TreeD()->Fill();
394         pMUON->ResetDigits();
395         fList->Delete();
396
397         
398         for(Int_t ii = 0; ii < AliMUONConstants::NCh(); ++ii) {
399             if (fHitMap[ii]) {
400                 hm=fHitMap[ii];
401                 delete hm;
402                 fHitMap[ii] = 0;
403             }
404         }
405         delete [] nmuon;    
406     } //end loop over cathodes
407     delete [] fHitMap;
408     delete fList;
409     
410 //  no need to delete ... and it makes a crash also
411 //    if (fHitsBgr)    fHitsBgr->Delete();
412 //    if (fPadHitsBgr) fPadHitsBgr->Delete();
413     // gObjectTable->Print();
414 }
415
416
417
418 void AliMUONMerger::SortTracks(Int_t *tracks,Int_t *charges,Int_t ntr)
419 {
420   //
421   // Sort the list of tracks contributing to a given digit
422   // Only the 3 most significant tracks are acctually sorted
423   //
424   
425   //
426   //  Loop over signals, only 3 times
427   //
428   
429   Int_t qmax;
430   Int_t jmax;
431   Int_t idx[3] = {-2,-2,-2};
432   Int_t jch[3] = {-2,-2,-2};
433   Int_t jtr[3] = {-2,-2,-2};
434   Int_t i,j,imax;
435   
436   if (ntr<3) imax=ntr;
437   else imax=3;
438   for(i=0;i<imax;i++){
439     qmax=0;
440     jmax=0;
441     
442     for(j=0;j<ntr;j++){
443       
444       if((i == 1 && j == idx[i-1]) 
445          ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
446       
447       if(charges[j] > qmax) {
448         qmax = charges[j];
449         jmax=j;
450       }       
451     } 
452     
453     if(qmax > 0) {
454       idx[i]=jmax;
455       jch[i]=charges[jmax]; 
456       jtr[i]=tracks[jmax]; 
457     }
458     
459   } 
460   
461   for(i=0;i<3;i++){
462     if (jtr[i] == -2) {
463          charges[i]=0;
464          tracks[i]=0;
465     } else {
466          charges[i]=jch[i];
467          tracks[i]=jtr[i];
468     }
469   }
470 }
471
472
473