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