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