]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUON.cxx
Use TLorentzVector for position and momentum
[u/mrichter/AliRoot.git] / MUON / AliMUON.cxx
1 ////////////////////////////////////////////////
2 //  Manager and hits classes for set:MUON     //
3 ////////////////////////////////////////////////
4
5 #include <TTUBE.h>
6 #include <TNode.h> 
7 #include <TRandom.h> 
8 #include <TObject.h>
9 #include <TVector.h>
10 #include <TObjArray.h>
11
12 #include "AliMUON.h"
13 #include "AliRun.h"
14 #include "AliMC.h"
15 #include "iostream.h"
16 #include "AliCallf77.h" 
17
18 // Static variables for the pad-hit iterator routines
19 static Int_t sMaxIterPad=0;
20 static Int_t sCurIterPad=0;
21  
22 ClassImp(AliMUON)
23  
24 //___________________________________________
25 AliMUON::AliMUON()
26 {
27    fIshunt     = 0;
28    fHits       = 0;
29    fClusters   = 0;
30    fNclusters  = 0;
31    fDchambers  = 0;
32    fRecClusters= 0;
33    fNdch       = 0;
34 }
35  
36 //___________________________________________
37 AliMUON::AliMUON(const char *name, const char *title)
38        : AliDetector(name,title)
39 {
40 //Begin_Html
41 /*
42 <img src="picts/alimuon.gif">
43 */
44 //End_Html
45  
46    fHits     = new TClonesArray("AliMUONhit",1000  );
47    fClusters = new TClonesArray("AliMUONcluster",10000);
48    fNclusters  =  0;
49    fIshunt     =  0;
50
51    fNdch      = new Int_t[11];
52
53    fDchambers = new TObjArray(11);
54
55    Int_t i;
56    
57    for (i=0; i<11 ;i++) {
58        (*fDchambers)[i] = new TClonesArray("AliMUONdigit",10000); 
59        fNdch[i]=0;
60    }
61
62    fRecClusters=new TObjArray(20);
63    for (i=0; i<20;i++)
64      (*fRecClusters)[i] = new TObjArray(1000);
65
66 //   
67 // Transport angular cut
68    fAccCut=0;
69    fAccMin=2;
70    fAccMax=9;
71
72    SetMarkerColor(kRed);
73 }
74  
75 //___________________________________________
76 AliMUON::~AliMUON()
77 {
78   fIshunt  = 0;
79   delete fHits;
80   delete fClusters;
81 //  for (Int_t i=0;i<10;i++) {
82 //      delete (*fDchambers)[i];
83 //      fNdch[i]=0;
84 //  }
85 //  delete fDchambers;
86   for (Int_t i=0;i<20;i++) 
87       delete (*fRecClusters)[i];
88   delete fRecClusters;
89
90 }
91  
92 //___________________________________________
93 void AliMUON::AddHit(Int_t track, Int_t *vol, Float_t *hits)
94 {
95   TClonesArray &lhits = *fHits;
96   new(lhits[fNhits++]) AliMUONhit(fIshunt,track,vol,hits);
97 }
98 //___________________________________________
99 void AliMUON::AddCluster(Int_t *clhits)
100 {
101    TClonesArray &lclusters = *fClusters;
102    new(lclusters[fNclusters++]) AliMUONcluster(clhits);
103 }
104 //_____________________________________________________________________________
105 void AliMUON::AddDigits(Int_t id, Int_t *tracks, Int_t *charges, Int_t *digits)
106 {
107     //
108     // Add a MUON digit to the list
109     //
110
111     TClonesArray &ldigits = *((TClonesArray*)(*fDchambers)[id]);
112     new(ldigits[fNdch[id]++]) AliMUONdigit(tracks,charges,digits);
113 }
114
115
116 //_____________________________________________________________________________
117 void AliMUON::AddRecCluster(Int_t iCh, Int_t iCat, AliMUONRecCluster* Cluster)
118 {
119     //
120     // Add a MUON reconstructed cluster to the list
121     //
122     TObjArray* ClusterList = RecClusters(iCh,iCat);
123     ClusterList->Add(Cluster);
124 }
125
126 //___________________________________________
127 void AliMUON::BuildGeometry()
128 {
129   TNode *Node, *Top;
130   const int kColorMUON = kBlue;
131   //
132   Top=gAlice->GetGeometry()->GetNode("alice");
133
134   // MUON
135   const float cz[10] = { 511, 519, 686, 694, 971, 979, 1245, 1253, 1445, 1453};
136   const float acc_min = TMath::Tan( 2*.0174532925199432958);
137   const float acc_max = TMath::Tan(9*.0174532925199432958);
138   float rmin, rmax;
139   
140   // Chamber 1
141   rmin = (cz[0]+0.25)*acc_min;
142   rmax = cz[0]*acc_max;
143   new TTUBE("S_MUON1","MUON chamber 1","void",rmin,rmax,0.25);
144   Top->cd();
145   Node = new TNode("MUON1","MUON chamber 1","S_MUON1",0,0,cz[0],"");
146   Node->SetLineColor(kColorMUON);
147   fNodes->Add(Node);
148   
149   // Chamber 2
150   rmin = (cz[1]+0.25)*acc_min;
151   rmax = cz[1]*acc_max;
152   new TTUBE("S_MUON2","MUON chamber 2","void",rmin,rmax,0.25);
153   Top->cd();
154   Node = new TNode("MUON2","MUON chamber 2","S_MUON2",0,0,cz[1],"");
155   fNodes->Add(Node);
156   Node->SetLineColor(kColorMUON);
157   
158   // Chamber 3
159   rmin = (cz[2]+0.25)*acc_min;
160   rmax =  cz[2]*acc_max;
161   new TTUBE("S_MUON3","MUON chamber 3","void",rmin,rmax,0.25);
162   Top->cd();
163   Node = new TNode("MUON3","MUON chamber 3","S_MUON3",0,0,cz[2],"");
164   Node->SetLineColor(kColorMUON);
165   fNodes->Add(Node);
166   
167   // Chamber 4
168   rmin = (cz[3]+0.25)*acc_min;
169   rmax = cz[3]*acc_max;
170   new TTUBE("S_MUON4","MUON chamber 4","void",rmin,rmax,0.25);
171   Top->cd();
172   Node = new TNode("MUON4","MUON chamber 4","S_MUON4",0,0,cz[3],"");
173   Node->SetLineColor(kColorMUON);
174   fNodes->Add(Node);
175
176   // Chamber 5
177   rmin = 30;
178   rmax = cz[4]*acc_max;
179   new TTUBE("S_MUON5","MUON chamber 5","void",rmin,rmax,0.25);
180   Top->cd();
181   Node = new TNode("MUON5","MUON chamber 5","S_MUON5",0,0,cz[4],"");
182   Node->SetLineColor(kColorMUON);
183   fNodes->Add(Node);
184   
185   // Chamber 6
186   rmin = 30;
187   rmax = cz[5]*acc_max;
188   new TTUBE("S_MUON6","MUON chamber 6","void",rmin,rmax,0.25);
189   Top->cd();
190   Node = new TNode("MUON6","MUON chamber 6","S_MUON6",0,0,cz[5],"");
191   fNodes->Add(Node);
192   Node->SetLineColor(kColorMUON);
193   
194   // Chamber 7
195   rmin = 30;
196   rmax = cz[6]*acc_max;
197   new TTUBE("S_MUON7","MUON chamber 7","void",rmin,rmax,0.25);
198   Top->cd();
199   Node = new TNode("MUON7","MUON chamber 7","S_MUON7",0,0,cz[6],"");
200   Node->SetLineColor(kColorMUON);
201   fNodes->Add(Node);
202
203   // Chamber 8
204   rmin = 30;
205   rmax = cz[7]*acc_max;
206   new TTUBE("S_MUON8","MUON chamber 8","void",rmin,rmax,0.25);
207   Top->cd();
208   Node = new TNode("MUON8","MUON chamber 8","S_MUON8",0,0,cz[7],"");
209   Node->SetLineColor(kColorMUON);
210   fNodes->Add(Node);
211
212   // Chamber 9
213   rmin = 30;
214   rmax = cz[8]*acc_max;
215   new TTUBE("S_MUON9","MUON chamber 9","void",rmin,rmax,0.25);
216   Top->cd();
217   Node = new TNode("MUON9","MUON chamber 9","S_MUON9",0,0,cz[8],"");
218   Node->SetLineColor(kColorMUON);
219   fNodes->Add(Node);
220
221   // Chamber 10
222   rmin = 30;
223   rmax = cz[9]*acc_max;
224   new TTUBE("S_MUON10","MUON chamber 10","void",rmin,rmax,0.25);
225   Top->cd();
226   Node = new TNode("MUON10","MUON chamber 10","S_MUON10",0,0,cz[9],"");
227   Node->SetLineColor(kColorMUON);
228   fNodes->Add(Node);
229
230 }
231
232 //___________________________________________
233 Int_t AliMUON::DistancetoPrimitive(Int_t , Int_t )
234 {
235    return 9999;
236 }
237
238 //___________________________________________
239 void AliMUON::MakeBranch(Option_t* option)
240 {
241   // Create Tree branches for the MUON.
242   
243   const Int_t buffersize = 4000;
244   char branchname[20];
245   sprintf(branchname,"%sCluster",GetName());
246
247   AliDetector::MakeBranch(option);
248
249   if (fClusters   && gAlice->TreeH()) {
250     gAlice->TreeH()->Branch(branchname,&fClusters, buffersize);
251     printf("Making Branch %s for clusters\n",branchname);
252   }
253
254 // one branch for digits per chamber
255   Int_t i;
256   
257   for (i=0; i<10 ;i++) {
258       sprintf(branchname,"%sDigits%d",GetName(),i+1);
259       
260       if (fDchambers   && gAlice->TreeD()) {
261           gAlice->TreeD()->Branch(branchname,&((*fDchambers)[i]), buffersize);
262           printf("Making Branch %s for digits in chamber %d\n",branchname,i+1);
263       } 
264   }
265 // one branch for rec clusters
266   for (i=0; i<20; i++) {
267       sprintf(branchname,"%sRecClus%d",GetName(),i+1);
268       if (fRecClusters   && gAlice->TreeD()) {
269           gAlice->TreeR()
270               ->Branch(branchname,"TObjArray", 
271                        &((*fRecClusters)[i]), buffersize,0);
272           printf("Making Branch %s for clusters in chamber %d\n",
273                  branchname,i+1);
274       }
275   }
276 }
277
278 //___________________________________________
279 void AliMUON::SetTreeAddress()
280 {
281   // Set branch address for the Hits and Digits Tree.
282   char branchname[20];
283   AliDetector::SetTreeAddress();
284
285   TBranch *branch;
286   TTree *treeH = gAlice->TreeH();
287   TTree *treeD = gAlice->TreeD();
288
289   if (treeH) {
290     if (fClusters) {
291       branch = treeH->GetBranch("MUONCluster");
292       if (branch) branch->SetAddress(&fClusters);
293     }
294   }
295
296   if (treeD) {
297       for (int i=0; i<10; i++) {
298           sprintf(branchname,"%sDigits%d",GetName(),i+1);
299           if (fDchambers) {
300               branch = treeD->GetBranch(branchname);
301               if (branch) branch->SetAddress(&((*fDchambers)[i]));
302           }
303       }
304   }
305 }
306 //___________________________________________
307 void AliMUON::ResetHits()
308 {
309   // Reset number of clusters and the cluster array for this detector
310   AliDetector::ResetHits();
311   fNclusters = 0;
312   if (fClusters) fClusters->Clear();
313 }
314
315 //____________________________________________
316 void AliMUON::ResetDigits()
317 {
318     //
319     // Reset number of digits and the digits array for this detector
320     //
321     for ( int i=0;i<10;i++ ) {
322         if ((*fDchambers)[i])   (*fDchambers)[i]->Clear();
323         if (fNdch)  fNdch[i]=0;
324     }
325 }
326 //____________________________________________
327 void AliMUON::ResetRecClusters()
328 {
329     //
330     // Reset the rec clusters
331     //
332     for ( int i=0;i<20;i++ ) {
333         if ((*fRecClusters)[i])   (*fRecClusters)[i]->Clear();
334     }
335 }
336 //___________________________________________
337
338 void AliMUON::SetPADSIZ(Int_t id, Int_t isec, Float_t p1, Float_t p2)
339 {
340     Int_t i=2*(id-1);
341     ((AliMUONchamber*) (*fChambers)[i])  ->SetPADSIZ(isec,p1,p2);
342     ((AliMUONchamber*) (*fChambers)[i+1])->SetPADSIZ(isec,p1,p2);
343 }
344
345 //___________________________________________
346 void AliMUON::SetMUCHSP(Int_t id, Float_t p1)
347 {
348     Int_t i=2*(id-1);
349     ((AliMUONchamber*) (*fChambers)[i])->SetMUCHSP(p1);
350     ((AliMUONchamber*) (*fChambers)[i+1])->SetMUCHSP(p1);
351 }
352
353 //___________________________________________
354 void AliMUON::SetMUSIGM(Int_t id, Float_t p1, Float_t p2)
355 {
356     Int_t i=2*(id-1);
357     ((AliMUONchamber*) (*fChambers)[i])->SetMUSIGM(p1,p2);
358     ((AliMUONchamber*) (*fChambers)[i+1])->SetMUSIGM(p1,p2);
359 }
360
361 //___________________________________________
362 void AliMUON::SetRSIGM(Int_t id, Float_t p1)
363 {
364     Int_t i=2*(id-1);
365     ((AliMUONchamber*) (*fChambers)[i])->SetRSIGM(p1);
366     ((AliMUONchamber*) (*fChambers)[i+1])->SetRSIGM(p1);
367 }
368
369 //___________________________________________
370 void AliMUON::SetMAXADC(Int_t id, Float_t p1)
371 {
372     Int_t i=2*(id-1);
373     ((AliMUONchamber*) (*fChambers)[i])->SetMAXADC(p1);
374     ((AliMUONchamber*) (*fChambers)[i+1])->SetMAXADC(p1);
375 }
376
377 //___________________________________________
378 void AliMUON::SetSMAXAR(Float_t p1)
379 {
380      fMaxStepGas=p1;
381 }
382
383 //___________________________________________
384 void AliMUON::SetSMAXAL(Float_t p1)
385 {
386     fMaxStepAlu=p1;
387 }
388
389 //___________________________________________
390 void AliMUON::SetDMAXAR(Float_t p1)
391 {
392     fMaxDestepGas=p1;
393 }
394
395 //___________________________________________
396 void AliMUON::SetDMAXAL(Float_t p1)
397 {
398     fMaxDestepAlu=p1;
399 }
400 //___________________________________________
401 void AliMUON::SetMUONACC(Bool_t acc, Float_t angmin, Float_t angmax)
402 {
403    fAccCut=acc;
404    fAccMin=angmin;
405    fAccMax=angmax;
406 }
407 //___________________________________________
408 void   AliMUON::SetSegmentationModel(Int_t id, Int_t isec, AliMUONsegmentation *segmentation)
409 {
410     ((AliMUONchamber*) (*fChambers)[id])->SegmentationModel(isec, segmentation);
411
412 }
413 //___________________________________________
414 void   AliMUON::SetResponseModel(Int_t id, AliMUONresponse *response)
415 {
416     ((AliMUONchamber*) (*fChambers)[id])->ResponseModel(response);
417 }
418
419 void   AliMUON::SetNsec(Int_t id, Int_t nsec)
420 {
421     ((AliMUONchamber*) (*fChambers)[id])->SetNsec(nsec);
422 }
423
424
425 //___________________________________________
426
427 void AliMUON::StepManager()
428 {
429     printf("Dummy version of muon step -- it should never happen!!\n");
430     const Float_t kRaddeg = 180/TMath::Pi();
431     Int_t nsec, ipart;
432     Float_t x[4], p[4];
433     Float_t pt, th0, th1;
434     char proc[5];
435     if(fAccCut) {
436         if((nsec=gMC->NSecondaries())>0) {
437             gMC->ProdProcess(proc);
438             if((gMC->TrackPid()==113 || gMC->TrackPid()==114) && !strcmp(proc,"DCAY")) {
439                 //
440                 // Check angular acceptance
441                 //* --- and have muons from resonance decays in the wanted window --- 
442                 if(nsec != 2) {
443                     printf(" AliMUON::StepManager: Strange resonance Decay into %d particles\n",nsec);
444                     gMC->StopEvent();
445                 } else {
446                     gMC->GetSecondary(0,ipart,x,p);
447                     pt  = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
448                     th0 = TMath::ATan2(pt,p[2])*kRaddeg;
449                     gMC->GetSecondary(1,ipart,x,p);
450                     pt  = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
451                     th1 = TMath::ATan2(pt,p[2])*kRaddeg;
452                     if(!(fAccMin < th0 && th0 < fAccMax) ||
453                        !(fAccMin < th1 && th1 < fAccMax)) 
454                         gMC->StopEvent();
455                 }
456             }
457         }
458     }
459 }
460 void AliMUON::ReconstructClusters()
461 {
462 //
463 // Initialize the necessary correspondance table
464 //
465     static const Int_t kMaxNpadx = 600;
466     static const Int_t kMaxNpady = 600;
467     Int_t elem[kMaxNpadx*2][kMaxNpady*2];
468 //
469 // Loop on chambers and on cathode planes
470 //
471     for (Int_t ich=0;ich<10;ich++)
472         for (Int_t icat=0;icat<2;icat++) {
473             //
474             // Get ready the current chamber stuff
475             //
476             AliMUONchamber*  iChamber= &(this->Chamber(ich));
477             AliMUONsegmentation*  segmentation = 
478                 iChamber->GetSegmentationModel(icat+1);
479             if (!segmentation) 
480                 continue;
481             TClonesArray *MUONdigits  = this->DigitsAddress(ich);
482             if (MUONdigits == 0) 
483                 continue;
484             cout << "Npx " << segmentation->Npx() 
485                  << " Npy " << segmentation->Npy() << endl;
486             //      
487             // Ready the digits
488             //  
489             gAlice->ResetDigits();
490             gAlice->TreeD()->GetEvent(icat+1); // spurious +1 ...
491             Int_t ndigits = MUONdigits->GetEntriesFast();
492             if (ndigits == 0) 
493                 continue;
494             printf("Found %d digits for cathode %d in chamber %d \n",
495                    ndigits,icat,ich+1);
496             AliMUONdigit  *mdig;
497             AliMUONRecCluster *Cluster;
498             //
499             // Build the correspondance table
500             //
501             memset(elem,0,sizeof(Int_t)*kMaxNpadx*kMaxNpady*4);
502             Int_t digit;
503             for (digit=0; digit<ndigits; digit++) 
504             {
505                 mdig    = (AliMUONdigit*)MUONdigits->UncheckedAt(digit);
506                 elem[kMaxNpadx+mdig->fPadX][kMaxNpady+mdig->fPadY] = digit+1;
507                 // because default is 0
508             }
509             //
510             // Declare some useful variables
511             //
512             Int_t Xlist[10];
513             Int_t Ylist[10];
514             Int_t Nlist;
515             Int_t nclust=0;
516             //
517             // loop over digits
518             //
519             for (digit=0;digit<ndigits;digit++) {
520                 mdig    = (AliMUONdigit*)MUONdigits->UncheckedAt(digit);
521                 //
522                 // if digit still available, start clustering
523                 //
524                 if (elem[kMaxNpadx+mdig->fPadX][kMaxNpady+mdig->fPadY]) {
525                     Cluster = new AliMUONRecCluster(digit, ich, icat);
526                     elem[kMaxNpadx+mdig->fPadX][kMaxNpady+mdig->fPadY]=0;
527                     //
528                     // loop over the current list of digits 
529                     // and look for neighbours
530                     //
531                     for(Int_t clusDigit=Cluster->FirstDigitIndex();
532                         clusDigit!=Cluster->InvalidDigitIndex();
533                         clusDigit=Cluster->NextDigitIndex()) {
534                         AliMUONdigit* pDigit=(AliMUONdigit*)MUONdigits
535                             ->UncheckedAt(clusDigit);
536                         segmentation->Neighbours(pDigit->fPadX,pDigit->fPadY, 
537                                                  &Nlist, Xlist, Ylist);
538                         for (Int_t Ilist=0;Ilist<Nlist;Ilist++) {
539                             if (elem[kMaxNpadx+Xlist[Ilist]][kMaxNpady
540                                                             +Ylist[Ilist]]) {
541                                 //
542                                 // Add the digit at the end and reset the table
543                                 //
544                                 Cluster->AddDigit(elem[kMaxNpadx+Xlist[Ilist]][kMaxNpady+Ylist[Ilist]]-1);
545                                 elem[kMaxNpadx+Xlist[Ilist]][kMaxNpady
546                                                             +Ylist[Ilist]] =0;
547                             } // if elem
548                         } // for Ilist
549                     } // for pDigit
550                     //
551                     // Store the cluster (good time to do Cluster polishing)
552                     //
553                     segmentation->FitXY(Cluster,MUONdigits);
554                     nclust ++;
555                     AddRecCluster(ich,icat,Cluster);
556                 }
557             }
558             printf("===> %d Clusters\n",nclust); 
559         } // for icat
560 }
561
562  
563 //______________________________________________________________________________
564 void AliMUON::Streamer(TBuffer &R__b)
565 {
566    // Stream an object of class AliMUON.
567       AliMUONchamber       *iChamber;
568       AliMUONsegmentation  *segmentation;
569       AliMUONresponse      *response;
570       TClonesArray         *digitsaddress;
571       
572    if (R__b.IsReading()) {
573       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
574       AliDetector::Streamer(R__b);
575       R__b >> fNclusters;
576       R__b >> fClusters; // diff
577       R__b >> fDchambers;
578       R__b.ReadArray(fNdch);
579       //
580       R__b >> fAccCut;
581       R__b >> fAccMin;
582       R__b >> fAccMax; 
583       //   
584       R__b >> fChambers;
585 // Stream chamber related information
586       for (Int_t i =0; i<10; i++) {
587           iChamber=(AliMUONchamber*) (*fChambers)[i];
588           iChamber->Streamer(R__b);
589           if (iChamber->Nsec()==1) {
590               segmentation=iChamber->GetSegmentationModel(1);
591               segmentation->Streamer(R__b);
592           } else {
593               segmentation=iChamber->GetSegmentationModel(1);
594               segmentation->Streamer(R__b);
595               segmentation=iChamber->GetSegmentationModel(2);
596               segmentation->Streamer(R__b);
597           }
598           response=iChamber->GetResponseModel();
599           response->Streamer(R__b);       
600           digitsaddress=(TClonesArray*) (*fDchambers)[i];
601           digitsaddress->Streamer(R__b);
602       }
603       
604    } else {
605       R__b.WriteVersion(AliMUON::IsA());
606       AliDetector::Streamer(R__b);
607       R__b << fNclusters;
608       R__b << fClusters; // diff
609       R__b << fDchambers;
610       R__b.WriteArray(fNdch, 10);
611       //
612       R__b << fAccCut;
613       R__b << fAccMin;
614       R__b << fAccMax; 
615       //   
616       R__b << fChambers;
617 //  Stream chamber related information
618       for (Int_t i =0; i<10; i++) {
619           iChamber=(AliMUONchamber*) (*fChambers)[i];
620           iChamber->Streamer(R__b);
621           if (iChamber->Nsec()==1) {
622               segmentation=iChamber->GetSegmentationModel(1);
623               segmentation->Streamer(R__b);
624           } else {
625               segmentation=iChamber->GetSegmentationModel(1);
626               segmentation->Streamer(R__b);
627               segmentation=iChamber->GetSegmentationModel(2);
628               segmentation->Streamer(R__b);
629           }
630           response=iChamber->GetResponseModel();
631           response->Streamer(R__b);
632          
633           digitsaddress=(TClonesArray*) (*fDchambers)[i];
634           digitsaddress->Streamer(R__b);
635       }
636    }
637 }
638 AliMUONcluster* AliMUON::FirstPad(AliMUONhit*  hit) 
639 {
640 //
641     // Initialise the pad iterator
642     // Return the address of the first padhit for hit
643     TClonesArray *theClusters = Clusters();
644     Int_t nclust = theClusters->GetEntriesFast();
645     if (nclust && hit->fPHlast > 0) {
646         sMaxIterPad=hit->fPHlast;
647         sCurIterPad=hit->fPHfirst;
648         return (AliMUONcluster*) fClusters->UncheckedAt(sCurIterPad-1);
649     } else {
650         return 0;
651     }
652 }
653
654 AliMUONcluster* AliMUON::NextPad() 
655 {
656     sCurIterPad++;
657     if (sCurIterPad <= sMaxIterPad) {
658         return (AliMUONcluster*) fClusters->UncheckedAt(sCurIterPad-1);
659     } else {
660         return 0;
661     }
662 }
663
664 ClassImp(AliMUONcluster)
665  
666 //___________________________________________
667 AliMUONcluster::AliMUONcluster(Int_t *clhits)
668 {
669    fHitNumber=clhits[0];
670    fCathode=clhits[1];
671    fQ=clhits[2];
672    fPadX=clhits[3];
673    fPadY=clhits[4];
674    fQpad=clhits[5];
675    fRSec=clhits[6];
676 }
677 ClassImp(AliMUONdigit)
678 //_____________________________________________________________________________
679 AliMUONdigit::AliMUONdigit(Int_t *digits)
680 {
681   //
682   // Creates a MUON digit object to be updated
683   //
684   fPadX        = digits[0];
685   fPadY        = digits[1];
686   fSignal      = digits[2];
687
688 }
689 //_____________________________________________________________________________
690 AliMUONdigit::AliMUONdigit(Int_t *tracks, Int_t *charges, Int_t *digits)
691 {
692     //
693     // Creates a MUON digit object
694     //
695     fPadX        = digits[0];
696     fPadY        = digits[1];
697     fSignal      = digits[2];
698     for(Int_t i=0; i<10; i++) {
699         fTcharges[i]  = charges[i];
700         fTracks[i]    = tracks[i];
701     }
702 }
703
704 ClassImp(AliMUONlist)
705  
706 //____________________________________________________________________________
707     AliMUONlist::AliMUONlist(Int_t rpad, Int_t *digits): 
708         AliMUONdigit(digits)
709 {
710     //
711     // Creates a MUON digit list object
712     //
713
714     fRpad = rpad;
715     fTrackList   = new TObjArray;
716  
717 }
718 //_____________________________________________________________________________
719
720
721 ClassImp(AliMUONhit)
722  
723 //___________________________________________
724     AliMUONhit::AliMUONhit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
725         AliHit(shunt, track)
726 {
727     fChamber=vol[0];
728     fParticle=(Int_t) hits[0];
729     fX=hits[1];
730     fY=hits[2];
731     fZ=hits[3];
732     fTheta=hits[4];
733     fPhi=hits[5];
734     fTlength=hits[6];
735     fEloss=hits[7];
736     fPHfirst=(Int_t) hits[8];
737     fPHlast=(Int_t) hits[9];
738 }
739 ClassImp(AliMUONreccluster) 
740
741 ClassImp(AliMUONRecCluster)
742
743 //_____________________________________________________________________
744 AliMUONRecCluster::AliMUONRecCluster()
745 {
746     fDigits=0;
747     fNdigit=-1;
748 }
749
750 AliMUONRecCluster::
751 AliMUONRecCluster(Int_t FirstDigit,Int_t Ichamber, Int_t Icathod)
752 {
753     fX = 0.;
754     fY = 0.;
755     fDigits = new TArrayI(10);
756     fNdigit=0;
757     AddDigit(FirstDigit);
758     fChamber=Ichamber;
759     fCathod=Icathod;
760 }
761
762 void AliMUONRecCluster::AddDigit(Int_t Digit)
763 {
764     if (fNdigit==fDigits->GetSize()) {
765         //enlarge the list by hand!
766         Int_t *array= new Int_t[fNdigit*2];
767         for(Int_t i=0;i<fNdigit;i++)
768             array[i] = fDigits->At(i);
769         fDigits->Adopt(fNdigit*2,array);
770     }
771     fDigits->AddAt(Digit,fNdigit);
772     fNdigit++;
773 }
774
775
776 AliMUONRecCluster::~AliMUONRecCluster()
777 {
778     if (fDigits)
779         delete fDigits;
780 }
781
782 Int_t AliMUONRecCluster::FirstDigitIndex()
783 {
784     fCurrentDigit=0;
785     return fDigits->At(fCurrentDigit);
786 }
787
788 Int_t AliMUONRecCluster::NextDigitIndex()
789 {
790     fCurrentDigit++;
791     if (fCurrentDigit<fNdigit)
792         return fDigits->At(fCurrentDigit);
793     else 
794         return InvalidDigitIndex();
795 }
796
797 Int_t AliMUONRecCluster::NDigits()
798 {
799     return fNdigit;
800 }
801
802 void AliMUONRecCluster::Finish()
803 {
804     // In order to reconstruct coordinates, one has to
805     // get back to the digits which is not trivial here,
806     // because we don't know where digits are stored!
807     // Center Of Gravity, or other method should be
808     // a property of AliMUON class!
809 }
810
811
812
813
814
815
816
817
818
819
820