]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONDigitizer.cxx
bug correction
[u/mrichter/AliRoot.git] / MUON / AliMUONDigitizer.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.6  2001/11/02 12:55:45  jchudoba
19 cleanup of the code, add const to Get methods
20
21 Revision 1.4  2001/10/18 14:44:09  jchudoba
22 Define constant MAXTRACKS for maximum number of tracks associated with 1 digit
23
24 Revision 1.3  2001/10/04 20:01:54  jchudoba
25 changes for TTask implementation, some other small editing
26
27 Revision 1.2  2001/07/28 10:46:04  hristov
28 AliRunDigitizer.h included; typos corrected
29
30 Revision 1.1  2001/07/27 15:41:01  jchudoba
31 merging/digitization classes
32
33 */
34 #include <TTree.h> 
35 #include <TObjArray.h>
36 #include <TFile.h>
37 #include <TDirectory.h>
38 #include <iostream.h>
39
40 #include "AliMUONDigitizer.h"
41 #include "AliMUONConstants.h"
42 #include "AliMUONChamber.h"
43 #include "AliMUONHitMapA1.h"
44 #include "AliMUON.h"
45 #include "AliMUONHit.h"
46 #include "AliMUONPadHit.h"
47 #include "AliMUONDigit.h"
48 #include "AliMUONTransientDigit.h"
49 #include "AliRun.h"
50 #include "AliPDG.h"
51 #include "AliRunDigitizer.h"
52
53 ClassImp(AliMUONDigitizer)
54
55 //___________________________________________
56 AliMUONDigitizer::AliMUONDigitizer() :AliDigitizer()
57 {
58 // Default ctor - don't use it
59   fHits = 0;
60   fPadHits = 0;
61   fHitMap = 0;
62   fTDList = 0;
63 }
64
65 //___________________________________________
66 AliMUONDigitizer::AliMUONDigitizer(AliRunDigitizer* manager) 
67     :AliDigitizer(manager)
68 {
69 // ctor which should be used
70   fHitMap  = 0;
71   fTDList  = 0;
72   fHits    = 0;
73   fPadHits = 0;
74   fDebug   = 0; 
75   if (GetDebug()>2) 
76     cerr<<"AliMUONDigitizer::AliMUONDigitizer"
77         <<"(AliRunDigitizer* manager) was processed"<<endl;
78 }
79
80 //------------------------------------------------------------------------
81 AliMUONDigitizer::~AliMUONDigitizer()
82 {
83 // Destructor
84   if (fHits)       delete fHits;
85   if (fPadHits)    delete fPadHits;
86   if (fHitMap)     delete [] fHitMap;
87   if (fTDList)     delete fTDList;
88 }
89
90 //------------------------------------------------------------------------
91 Bool_t AliMUONDigitizer::Exists(const AliMUONPadHit *padhit) const
92 {
93   return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
94 }
95
96 //------------------------------------------------------------------------
97 void AliMUONDigitizer::Update(AliMUONPadHit *padhit)
98 {
99     AliMUONTransientDigit *pdigit = 
100       static_cast<AliMUONTransientDigit*>(
101       fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
102
103     // update charge
104     //
105     Int_t iqpad    = padhit->QPad();        // charge per pad
106     pdigit->AddSignal(iqpad);
107     pdigit->AddPhysicsSignal(iqpad);            
108
109     // update list of tracks
110     //
111     Int_t track, charge;    
112     track = fTrack+fMask;
113     if (fSignal) {
114       charge = iqpad;
115     } else {
116       charge = kBgTag;
117     }
118     pdigit->UpdateTrackList(track,charge);
119 }
120
121 //------------------------------------------------------------------------
122 void AliMUONDigitizer::CreateNew(AliMUONPadHit *padhit)
123 {
124 // Create new AliMUONTransientDigit and add it to the fTDList
125
126   fTDList->AddAtAndExpand(
127     new AliMUONTransientDigit(fNch,fDigits),fCounter);
128   fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
129   AliMUONTransientDigit* pdigit = 
130     (AliMUONTransientDigit*)fTDList->Last();
131   // list of tracks
132   Int_t track, charge;    
133   track = fTrack+fMask;
134   if (fSignal) {
135     charge = padhit->QPad();
136   } else {
137     charge = kBgTag;
138   }
139   pdigit->AddToTrackList(track,charge);
140   fCounter++;
141 }
142
143
144 //------------------------------------------------------------------------
145 Bool_t AliMUONDigitizer::Init()
146 {
147 // Initialization 
148     
149   fHits     = new TClonesArray("AliMUONHit",1000);
150   fPadHits  = new TClonesArray("AliMUONPadHit",1000);
151   return kTRUE;
152 }
153
154
155 //------------------------------------------------------------------------
156 //void AliMUONDigitizer::Digitize()
157 void AliMUONDigitizer::Exec(Option_t* option)
158 {
159
160   TString optionString = option;
161   if (optionString.Data() == "deb") {
162     cout<<"AliMUONDigitizer::Exec: called with option deb "<<endl;
163     fDebug = 3;
164   }
165   AliMUONChamber*   iChamber;
166   AliSegmentation*  segmentation;
167   
168   if (GetDebug()>2) cerr<<"   AliMUONDigitizer::Digitize() starts"<<endl;
169   fTDList = new TObjArray;
170
171   AliMUON *pMUON  = (AliMUON *) gAlice->GetModule("MUON");
172   if (!pMUON) {
173     cerr<<"AliMUONDigitizer::Digitize Error:"
174         <<" module MUON not found in the input file"<<endl;
175     return;
176   }
177   pMUON->MakeBranchInTreeD(fManager->GetTreeD());
178   fHitMap= new AliMUONHitMapA1* [AliMUONConstants::NCh()];
179
180   //
181   // loop over cathodes
182   //
183
184   for (int icat = 0; icat < 2; icat++) { 
185     fCounter = 0;
186     for (Int_t i = 0; i < AliMUONConstants::NCh(); i++) {
187       iChamber = &(pMUON->Chamber(i));
188 //      if (!(iChamber->Nsec() == 1 && icat == 1)) {
189         segmentation = iChamber->SegmentationModel(icat+1);
190         fHitMap[i] = new AliMUONHitMapA1(segmentation, fTDList);
191 //      }
192     }
193
194
195 // Loop over files to digitize
196     fSignal = kTRUE;
197     for (Int_t inputFile=0; inputFile<fManager->GetNinputs();
198          inputFile++) {
199 // Connect MUON branches
200
201       if (inputFile > 0 ) fSignal = kFALSE;
202       TBranch *branchHits = 0;
203       TBranch *branchPadHits = 0;
204       TTree *treeH = fManager->GetInputTreeH(inputFile);
205       if (GetDebug()>2) {
206         cerr<<"   inputFile , cathode = "<<inputFile<<" "
207             <<icat<<endl;
208         cerr<<"   treeH, fHits "<<treeH<<" "<<fHits<<endl;
209       }
210       if (treeH && fHits) {
211         branchHits = treeH->GetBranch("MUON");
212         if (branchHits) {
213           fHits->Clear();
214           branchHits->SetAddress(&fHits);
215         }
216         else
217           Error("Exec","branch MUON was not found");
218       }
219       if (GetDebug()>2) cerr<<"   branchHits = "<<branchHits<<endl;
220
221       if (treeH && fPadHits) {
222         branchPadHits = treeH->GetBranch("MUONCluster");
223         if (branchPadHits) 
224           branchPadHits->SetAddress(&fPadHits);
225         else
226           Error("Exec","branch MUONCluster was not found");
227       }
228       if (GetDebug()>2) cerr<<"   branchPadHits = "<<branchPadHits<<endl;
229
230 //
231 //   Loop over tracks
232 //
233
234       Int_t ntracks = (Int_t) treeH->GetEntries();
235
236       for (fTrack = 0; fTrack < ntracks; fTrack++) {
237         if (GetDebug()>2) cerr<<"   fTrack = "<<fTrack<<endl;
238         fHits->Clear();
239         fPadHits->Clear();
240         branchHits->GetEntry(fTrack);
241         branchPadHits->GetEntry(fTrack);
242         
243 //
244 //   Loop over hits
245
246         AliMUONHit* mHit;
247         for(Int_t i = 0; i < fHits->GetEntriesFast(); ++i) {
248           mHit = static_cast<AliMUONHit*>(fHits->At(i));
249           fNch = mHit->Chamber()-1;  // chamber number
250           if (fNch > AliMUONConstants::NCh()-1) {
251             cerr<<"AliMUONDigitizer: ERROR: "
252                 <<"fNch > AliMUONConstants::NCh()-1, fNch, NCh(): "
253                 <<fNch<<", "<< AliMUONConstants::NCh()<<endl;
254             return;
255           }
256           iChamber = &(pMUON->Chamber(fNch));
257 //
258 // Loop over pad hits
259           for (AliMUONPadHit* mPad =
260                  (AliMUONPadHit*)pMUON->FirstPad(mHit,fPadHits);
261                mPad;
262                mPad = (AliMUONPadHit*)pMUON->NextPad(fPadHits))
263           {
264             Int_t cathode  = mPad->Cathode();      // cathode number
265             Int_t ipx      = mPad->PadX();         // pad number on X
266             Int_t ipy      = mPad->PadY();         // pad number on Y
267             Int_t iqpad    = Int_t(mPad->QPad());  // charge per pad
268             if (cathode != (icat+1)) continue;
269
270             fMask = fManager->GetMask(inputFile);
271             fDigits[0] = ipx;
272             fDigits[1] = ipy;
273             fDigits[2] = icat;
274             fDigits[3] = iqpad;
275             if (inputFile == 0) {
276               fDigits[4] = iqpad;
277             } else {
278               fDigits[4] = 0;
279             }
280             if (mHit->Particle() == kMuonPlus ||
281                 mHit->Particle() == kMuonMinus) {
282               fDigits[5] = (mPad->HitNumber()) + fMask; 
283             } else fDigits[5] = -1;
284
285             // build the list of fired pads and update the info, 
286             // fDigits is input for Update(mPad)
287
288             if (!Exists(mPad)) {
289               CreateNew(mPad);
290             } else {
291               Update(mPad);
292             } //  end if Exists(mPad)
293           } //end loop over clusters
294         } // hit loop
295       } // track loop
296     } // end file loop
297     if (GetDebug()>2) cerr<<"END OF FILE LOOP"<<endl;
298
299     Int_t tracks[kMAXTRACKS];
300     Int_t charges[kMAXTRACKS];
301     Int_t nentries = fTDList->GetEntriesFast();
302         
303     for (Int_t nent = 0; nent < nentries; nent++) {
304       AliMUONTransientDigit *address = (AliMUONTransientDigit*)fTDList->At(nent);
305       if (address == 0) continue; 
306       Int_t ich = address->Chamber();
307       Int_t   q = address->Signal(); 
308       iChamber = &(pMUON->Chamber(ich));
309 //
310 //  Digit Response (noise, threshold, saturation, ...)
311       AliMUONResponse * response = iChamber->ResponseModel();
312       q = response->DigitResponse(q);
313             
314       if (!q) continue;
315             
316       fDigits[0] = address->PadX();
317       fDigits[1] = address->PadY();
318       fDigits[2] = address->Cathode();
319       fDigits[3] = q;
320       fDigits[4] = address->Physics();
321       fDigits[5] = address->Hit();
322             
323       Int_t nptracks = address->GetNTracks();
324
325       if (nptracks > kMAXTRACKS) {
326         if (GetDebug() >0) {
327           cerr<<"AliMUONDigitizer: nptracks > 10 "<<nptracks;
328           cerr<<"reset to max value "<<kMAXTRACKS<<endl;
329         }
330         nptracks = kMAXTRACKS;
331       }
332       if (nptracks > 2 && GetDebug() >2) {
333         cerr<<"AliMUONDigitizer: nptracks > 2 "<<nptracks<<endl;
334         printf("cat,ich,ix,iy,q %d %d %d %d %d \n",icat,ich,fDigits[0],fDigits[1],q);
335       }
336       for (Int_t tr = 0; tr < nptracks; tr++) {
337         tracks[tr]   = address->GetTrack(tr);
338         charges[tr]  = address->GetCharge(tr);
339       }      //end loop over list of tracks for one pad
340       // Sort list of tracks according to charge
341       if (nptracks > 1) {
342         SortTracks(tracks,charges,nptracks);
343       }
344       if (nptracks < kMAXTRACKS ) {
345         for (Int_t i = nptracks; i < kMAXTRACKS; i++) {
346           tracks[i]  = 0;
347           charges[i] = 0;
348         }
349       }
350             
351       // fill digits
352       pMUON->AddDigits(ich,tracks,charges,fDigits);
353     }
354
355     fManager->GetTreeD()->Fill();
356
357     pMUON->ResetDigits();  //
358     fTDList->Clear();
359
360         
361     for(Int_t ii = 0; ii < AliMUONConstants::NCh(); ++ii) {
362       if (fHitMap[ii]) {
363         delete fHitMap[ii];
364         fHitMap[ii] = 0;
365       }
366     }
367   } //end loop over cathodes
368   if (GetDebug()>2) 
369     cerr<<"AliMUONDigitizer::Exec: writing the TreeD: "
370         <<fManager->GetTreeD()->GetName()<<endl;
371   fManager->GetTreeD()->Write(0,TObject::kOverwrite);
372   delete [] fHitMap;
373   delete fTDList;
374     
375   if (fHits)    fHits->Delete();
376   if (fPadHits) fPadHits->Delete();
377 }
378
379
380 //------------------------------------------------------------------------
381 void AliMUONDigitizer::SortTracks(Int_t *tracks,Int_t *charges,Int_t ntr)
382 {
383   //
384   // Sort the list of tracks contributing to a given digit
385   // Only the 3 most significant tracks are acctually sorted
386   //
387   
388   //
389   //  Loop over signals, only 3 times
390   //
391   
392   Int_t qmax;
393   Int_t jmax;
394   Int_t idx[3] = {-2,-2,-2};
395   Int_t jch[3] = {-2,-2,-2};
396   Int_t jtr[3] = {-2,-2,-2};
397   Int_t i,j,imax;
398   
399   if (ntr<3) imax=ntr;
400   else imax=3;
401   for(i=0;i<imax;i++){
402     qmax=0;
403     jmax=0;
404     
405     for(j=0;j<ntr;j++){
406       
407       if((i == 1 && j == idx[i-1]) 
408          ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
409       
410       if(charges[j] > qmax) {
411         qmax = charges[j];
412         jmax=j;
413       }       
414     } 
415     
416     if(qmax > 0) {
417       idx[i]=jmax;
418       jch[i]=charges[jmax]; 
419       jtr[i]=tracks[jmax]; 
420     }
421     
422   } 
423   
424   for(i=0;i<3;i++){
425     if (jtr[i] == -2) {
426       charges[i]=0;
427       tracks[i]=0;
428     } else {
429       charges[i]=jch[i];
430       tracks[i]=jtr[i];
431     }
432   }
433 }