]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUON.cxx
Adapting FillESD to the AliRecostruction class (Y.Schutz)
[u/mrichter/AliRoot.git] / MUON / AliMUON.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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
19 ///////////////////////////////////////////////
20 //  Manager and hits classes for set:MUON     //
21 ////////////////////////////////////////////////
22
23 #include "Riostream.h"
24
25 #include <AliPDG.h>
26 #include <TBRIK.h>
27 #include <TCanvas.h>
28 #include <TDirectory.h>
29 #include <TFile.h>
30 #include <TGeometry.h>
31 #include <TMinuit.h>
32 #include <TNode.h> 
33 #include <TNtuple.h>
34 #include <TObjArray.h>
35 #include <TObject.h>
36 #include <TObjectTable.h>
37 #include <TPad.h>
38 #include <TParticle.h>
39 #include <TROOT.h>
40 #include <TRandom.h> 
41 #include <TRotMatrix.h>
42 #include <TTUBE.h>
43 #include <TTUBE.h>
44 #include <TTree.h> 
45 #include <TVector.h>
46 #include <TVirtualMC.h>
47
48 #include "AliConst.h" 
49 #include "AliHeader.h"
50 #include "AliHitMap.h"
51 #include "AliLoader.h"
52 #include "AliESD.h"
53 #include "AliESDMuonTrack.h"
54 #include "AliMUONLoader.h"
55 #include "AliMUON.h"
56 #include "AliMUONTriggerTrack.h"
57 #include "AliMUONEventReconstructor.h"
58 #include "AliMUONTrack.h"
59 #include "AliMUONTrackParam.h"
60 #include "AliMUONChamberTrigger.h"
61 #include "AliMUONClusterFinderVS.h"
62 #include "AliMUONClusterInput.h"
63 #include "AliMUONConstants.h"
64 #include "AliMUONDigit.h"
65 #include "AliMUONGlobalTrigger.h"
66 #include "AliMUONHit.h"
67 #include "AliMUONHitMapA1.h"
68 #include "AliMUONLocalTrigger.h"
69 #include "AliMUONMerger.h"      
70 #include "AliMUONPadHit.h"
71 #include "AliMUONRawCluster.h"
72 #include "AliMUONTransientDigit.h"
73 #include "AliMUONTriggerCircuit.h"
74 #include "AliMUONTriggerDecision.h"
75 #include "AliMUONVGeometryBuilder.h"    
76 #include "AliRun.h"     
77 #include "AliMUONDigitizerv1.h"
78
79
80 // Defaults parameters for Z positions of chambers
81 // taken from values for "stations" in AliMUON::AliMUON
82 //     const Float_t zch[7]={528, 690., 975., 1249., 1449., 1610, 1710.};
83 // and from array "dstation" in AliMUONv1::CreateGeometry
84 //          Float_t dstation[5]={20., 20., 20, 20., 20.};
85 //     for tracking chambers,
86 //          according to (Z1 = zch - dstation) and  (Z2 = zch + dstation)
87 //          for the first and second chambers in the station, respectively,
88 // and from "DTPLANES" in AliMUONv1::CreateGeometry
89 //           const Float_t DTPLANES = 15.;
90 //     for trigger chambers,
91 //          according to (Z1 = zch) and  (Z2 = zch + DTPLANES)
92 //          for the first and second chambers in the station, respectively
93
94 ClassImp(AliMUON)
95 //__________________________________________________________________
96 AliMUON::AliMUON()
97 {
98 // Default Constructor
99 //
100     fNCh             = 0;
101     fNTrackingCh     = 0;
102     fIshunt          = 0;
103     fChambers        = 0;
104     fGeometryBuilders = 0; 
105     fTriggerCircuits = 0;
106     fAccMin          = 0.;
107     fAccMax          = 0.;   
108     fAccCut          = kFALSE;
109     fMerger          = 0;
110     fFileName        = 0;
111     fMUONData        = 0;
112     fSplitLevel      = 0;
113 }
114 //__________________________________________________________________
115 AliMUON::AliMUON(const char *name, const char *title)
116   : AliDetector(name,title)
117 {
118 //Begin_Html
119 /*
120 <img src="gif/alimuon.gif">
121 */
122 //End_Html
123   fMUONData  = 0x0;
124   fSplitLevel= 0;
125   fIshunt     =  0;
126
127   fNCh             = AliMUONConstants::NCh(); 
128   fNTrackingCh     = AliMUONConstants::NTrackingCh();
129
130   SetMarkerColor(kRed);//
131 //
132 // Creating List of Chambers
133     Int_t ch;
134     fChambers = new TObjArray(AliMUONConstants::NCh());
135     fGeometryBuilders = new TObjArray(AliMUONConstants::NCh());
136
137     // Loop over stations
138     for (Int_t st = 0; st < AliMUONConstants::NCh() / 2; st++) {
139       // Loop over 2 chambers in the station
140       for (Int_t stCH = 0; stCH < 2; stCH++) {
141         //
142         //    
143         //    Default Parameters for Muon Tracking Stations
144         ch = 2 * st + stCH;
145         if (ch < AliMUONConstants::NTrackingCh()) {
146           fChambers->AddAt(new AliMUONChamber(ch),ch);
147         } else {
148           fChambers->AddAt(new AliMUONChamberTrigger(ch),ch);
149         }
150         AliMUONChamber* chamber = (AliMUONChamber*) fChambers->At(ch);
151         //chamber->SetGid(0);
152         // Default values for Z of chambers
153         chamber->SetZ(AliMUONConstants::DefaultChamberZ(ch));
154         //
155         chamber->InitGeo(AliMUONConstants::DefaultChamberZ(ch));
156         //          Set chamber inner and outer radius to default
157         chamber->SetRInner(AliMUONConstants::Dmin(st)/2);
158         chamber->SetROuter(AliMUONConstants::Dmax(st)/2);
159         //
160       } // Chamber stCH (0, 1) in 
161     }     // Station st (0...)
162     
163     // Negatives values are ignored by geant3 CONS200 in the calculation of the tracking parameters
164     fMaxStepGas=0.1; 
165     fMaxStepAlu=0.1;  
166     fMaxDestepGas=-1;
167     fMaxDestepAlu=-1;
168     
169     fMaxIterPad   = 0;
170     fCurIterPad   = 0;
171     
172     fAccMin          = 0.;
173     fAccMax          = 0.;   
174     fAccCut          = kFALSE;
175     
176     // cp new design of AliMUONTriggerDecision
177     fTriggerCircuits = new TObjArray(AliMUONConstants::NTriggerCircuit());
178     for (Int_t circ=0; circ<AliMUONConstants::NTriggerCircuit(); circ++) {
179       fTriggerCircuits->AddAt(new AliMUONTriggerCircuit(),circ);          
180     }
181     fMerger = 0;
182 }
183 //____________________________________________________________________
184 AliMUON::AliMUON(const AliMUON& rMUON):AliDetector(rMUON)
185 {
186 // Dummy copy constructor
187     ;
188     
189 }
190 //____________________________________________________________________
191 AliMUON::~AliMUON()
192 {
193 // Destructor
194   if(fDebug) printf("%s: Calling AliMUON destructor !!!\n",ClassName());
195   fIshunt  = 0;
196   if (fMerger) delete fMerger;
197
198   if (fGeometryBuilders){
199     fGeometryBuilders->Delete();
200     delete fGeometryBuilders;
201   } 
202 }
203 //_____________________________________________________________________________
204 void AliMUON::AddGeometryBuilder(AliMUONVGeometryBuilder* geomBuilder)
205 {
206 // Adds the geometry builder to the list
207 // ---
208
209   fGeometryBuilders->Add(geomBuilder);
210 }
211 //____________________________________________________________________
212 void AliMUON::BuildGeometry()
213 {
214 // Geometry for event display
215   for (Int_t i=0; i<7; i++) {
216     for (Int_t j=0; j<2; j++) {
217       Int_t id=2*i+j+1;
218       this->Chamber(id-1).SegmentationModel(1)->Draw("eventdisplay");
219     }
220   }
221 }
222 //___________________________________________________________________
223 Int_t AliMUON::DistancetoPrimitive(Int_t , Int_t )
224 {
225   return 9999;
226 }
227 //__________________________________________________________________
228 void  AliMUON::SetTreeAddress()
229 {
230   GetMUONData()->SetLoader(fLoader); 
231   GetMUONData()->SetTreeAddress("H,D,RC");
232   fHits = GetMUONData()->Hits(); // Added by Ivana to use the methods FisrtHit, NextHit of AliDetector
233 }
234
235 //____________________________________________________________________
236 void AliMUON::SetPadSize(Int_t id, Int_t isec, Float_t p1, Float_t p2)
237 {
238 // Set the pad size for chamber id and cathode isec
239     Int_t i=2*(id-1);
240     ((AliMUONChamber*) fChambers->At(i))  ->SetPadSize(isec,p1,p2);
241     ((AliMUONChamber*) fChambers->At(i+1))->SetPadSize(isec,p1,p2);
242 }
243
244 //___________________________________________
245 void AliMUON::SetChambersZ(const Float_t *Z)
246 {
247   // Set Z values for all chambers (tracking and trigger)
248   // from the array pointed to by "Z"
249     for (Int_t ch = 0; ch < AliMUONConstants::NCh(); ch++)
250         ((AliMUONChamber*) fChambers->At(ch))->SetZ(Z[ch]);
251     return;
252 }
253 //_________________________________________________________________
254 void AliMUON::SetChambersZToDefault()
255 {
256   // Set Z values for all chambers (tracking and trigger)
257   // to default values
258   SetChambersZ(AliMUONConstants::DefaultChamberZ());
259   return;
260 }
261 //_________________________________________________________________
262 void AliMUON::SetChargeSlope(Int_t id, Float_t p1)
263 {
264 // Set the inverse charge slope for chamber id
265     Int_t i=2*(id-1);    //PH    ((AliMUONChamber*) (*fChambers)[i])->SetSigmaIntegration(p1);
266     //PH    ((AliMUONChamber*) (*fChambers)[i+1])->SetSigmaIntegration(p1);
267     ((AliMUONChamber*) fChambers->At(i))->SetChargeSlope(p1);
268     ((AliMUONChamber*) fChambers->At(i+1))->SetChargeSlope(p1);
269 }
270 //__________________________________________________________________
271 void AliMUON::SetChargeSpread(Int_t id, Float_t p1, Float_t p2)
272 {
273 // Set sigma of charge spread for chamber id
274     Int_t i=2*(id-1);
275     ((AliMUONChamber*) fChambers->At(i))->SetChargeSpread(p1,p2);
276     ((AliMUONChamber*) fChambers->At(i+1))->SetChargeSpread(p1,p2);
277 }
278 //___________________________________________________________________
279 void AliMUON::SetSigmaIntegration(Int_t id, Float_t p1)
280 {
281 // Set integration limits for charge spread
282     Int_t i=2*(id-1);
283     ((AliMUONChamber*) fChambers->At(i))->SetSigmaIntegration(p1);
284     ((AliMUONChamber*) fChambers->At(i+1))->SetSigmaIntegration(p1);
285 }
286
287 //__________________________________________________________________
288 void AliMUON::SetMaxAdc(Int_t id, Int_t p1)
289 {
290 // Set maximum number for ADCcounts (saturation)
291     Int_t i=2*(id-1);
292     ((AliMUONChamber*) fChambers->At(i))->SetMaxAdc(p1);
293     ((AliMUONChamber*) fChambers->At(i+1))->SetMaxAdc(p1);
294 }
295 //__________________________________________________________________
296 void AliMUON::SetMaxStepGas(Float_t p1)
297 {
298 // Set stepsize in gas
299   fMaxStepGas=p1;
300 }
301 //__________________________________________________________________
302 void AliMUON::SetMaxStepAlu(Float_t p1)
303 {
304 // Set step size in Alu
305     fMaxStepAlu=p1;
306 }
307 //__________________________________________________________________
308 void AliMUON::SetMaxDestepGas(Float_t p1)
309 {
310 // Set maximum step size in Gas
311     fMaxDestepGas=p1;
312 }
313 //__________________________________________________________________
314 void AliMUON::SetMaxDestepAlu(Float_t p1)
315 {
316 // Set maximum step size in Alu
317   fMaxDestepAlu=p1;
318 }
319 //___________________________________________________________________
320 void AliMUON::SetAcceptance(Bool_t acc, Float_t angmin, Float_t angmax)
321 {
322 // Set acceptance cuts 
323   fAccCut=acc;
324   fAccMin=angmin*TMath::Pi()/180;
325   fAccMax=angmax*TMath::Pi()/180;
326   Int_t ch;
327   if (acc) {
328     for (Int_t st = 0; st < AliMUONConstants::NCh() / 2; st++) {
329       // Loop over 2 chambers in the station
330       for (Int_t stCH = 0; stCH < 2; stCH++) {
331         ch = 2 * st + stCH;
332         //         Set chamber inner and outer radius according to acceptance cuts
333         Chamber(ch).SetRInner(TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)*TMath::Tan(fAccMin)));
334         Chamber(ch).SetROuter(TMath::Abs(AliMUONConstants::DefaultChamberZ(ch)*TMath::Tan(fAccMax)));
335       } // chamber loop
336     } // station loop
337   }
338 }
339
340 //____________________________________________________________________
341 Float_t  AliMUON::GetMaxStepGas() const
342 {
343 // Return stepsize in gas
344   
345   return fMaxStepGas;
346 }  
347
348 //____________________________________________________________________
349 Float_t  AliMUON::GetMaxStepAlu() const
350 {
351 // Return step size in Alu
352   
353   return fMaxStepAlu;
354 }
355   
356 //____________________________________________________________________
357 Float_t  AliMUON::GetMaxDestepGas() const
358 {
359 // Return maximum step size in Gas
360   
361   return fMaxDestepGas;
362 }
363   
364 //____________________________________________________________________
365 Float_t  AliMUON::GetMaxDestepAlu() const
366 {
367 // Return maximum step size in Gas
368   
369   return fMaxDestepAlu;
370 }
371
372 //____________________________________________________________________
373 void   AliMUON::SetSegmentationModel(Int_t id, Int_t isec, AliSegmentation *segmentation)
374 {
375 // Set the segmentation for chamber id cathode isec
376     ((AliMUONChamber*) fChambers->At(id))->SetSegmentationModel(isec, segmentation);
377
378 }
379 //____________________________________________________________________
380 void   AliMUON::SetResponseModel(Int_t id, AliMUONResponse *response)
381 {
382 // Set the response for chamber id
383     ((AliMUONChamber*) fChambers->At(id))->SetResponseModel(response);
384 }
385 //____________________________________________________________________
386 void   AliMUON::SetReconstructionModel(Int_t id, AliMUONClusterFinderVS *reconst)
387 {
388 // Set ClusterFinder for chamber id
389     ((AliMUONChamber*) fChambers->At(id))->SetReconstructionModel(reconst);
390 }
391 //____________________________________________________________________
392 void   AliMUON::SetNsec(Int_t id, Int_t nsec)
393 {
394 // Set number of segmented cathods for chamber id
395     ((AliMUONChamber*) fChambers->At(id))->SetNsec(nsec);
396 }
397 //____________________________________________________________________
398 AliDigitizer* AliMUON::CreateDigitizer(AliRunDigitizer* manager) const
399 {
400   return new AliMUONDigitizerv1(manager);
401 }
402 //_____________________________________________________________________
403 void AliMUON::SDigits2Digits()
404 {
405
406 // write TreeD here 
407
408     if (!fMerger) {
409       if (gAlice->GetDebug()>0) {
410         cerr<<"AliMUON::SDigits2Digits: create default AliMUONMerger "<<endl;
411         cerr<<" no merging, just digitization of 1 event will be done"<<endl;
412       }
413       fMerger = new AliMUONMerger();
414     }
415     fMerger->Init();
416     fMerger->Digitise();
417     char hname[30];
418     //    sprintf(hname,"TreeD%d",fLoader->GetHeader()->GetEvent());
419     fLoader->TreeD()->Write(hname,TObject::kOverwrite);
420     fLoader->TreeD()->Reset();
421 }
422
423 //_______________________________________________________________________
424 AliLoader* AliMUON::MakeLoader(const char* topfoldername)
425
426 //builds standard getter (AliLoader type)
427 //if detector wants to use castomized getter, it must overload this method
428
429  if (GetDebug())
430    Info("MakeLoader",
431         "Creating standard getter for detector %s. Top folder is %s.",
432          GetName(),topfoldername);
433      
434  fLoader   = new AliLoader(GetName(),topfoldername);
435  fMUONData = new AliMUONData(fLoader,GetName(),GetName()); 
436  fMUONData->SetSplitLevel(fSplitLevel);
437  return fLoader;
438 }
439
440 //_______________________________________________________________________
441 void AliMUON::Trigger(Int_t nev){
442 // call the Trigger Algorithm and fill TreeR
443
444   Int_t singlePlus[3]  = {0,0,0}; 
445   Int_t singleMinus[3] = {0,0,0}; 
446   Int_t singleUndef[3] = {0,0,0};
447   Int_t pairUnlike[3]  = {0,0,0}; 
448   Int_t pairLike[3]    = {0,0,0};
449   
450   ResetTrigger();
451   AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
452   decision->Trigger();   
453   decision->GetGlobalTrigger(singlePlus, singleMinus, singleUndef,
454                              pairUnlike, pairLike);
455   
456   // add a local trigger in the list 
457   GetMUONData()->AddGlobalTrigger(singlePlus, singleMinus, singleUndef, pairUnlike, pairLike);
458   Int_t i;
459   
460   for (Int_t icirc=0; icirc<AliMUONConstants::NTriggerCircuit(); icirc++) { 
461     if(decision->GetITrigger(icirc)==1) {
462       Int_t localtr[7]={0,0,0,0,0,0,0};      
463       Int_t loLpt[2]={0,0}; Int_t loHpt[2]={0,0}; Int_t loApt[2]={0,0};
464       decision->GetLutOutput(icirc, loLpt, loHpt, loApt);
465       localtr[0] = icirc;
466       localtr[1] = decision->GetStripX11(icirc);
467       localtr[2] = decision->GetDev(icirc);
468       localtr[3] = decision->GetStripY11(icirc);
469       for (i=0; i<2; i++) {    // convert the Lut output in 1 digit 
470         localtr[4] = localtr[4]+Int_t(loLpt[i]*TMath::Power(2,i));
471         localtr[5] = localtr[5]+Int_t(loHpt[i]*TMath::Power(2,i));
472         localtr[6] = localtr[6]+Int_t(loApt[i]*TMath::Power(2,i));
473       }
474       GetMUONData()->AddLocalTrigger(localtr);  // add a local trigger in the list
475     }
476   }
477   
478   delete decision;
479
480   //  fLoader->TreeR()->Fill();
481   GetMUONData()->Fill("GLT"); //Filling Global and Local Trigger GLT
482   //  char hname[30];
483   //  sprintf(hname,"TreeR%d",nev);
484   //  fLoader->TreeR()->Write(hname,TObject::kOverwrite);
485     //  fLoader->TreeR()->Reset();
486   fLoader->WriteRecPoints("OVERWRITE");
487   
488   printf("\n End of trigger for event %d", nev);
489 }
490
491 //____________________________________________________________________
492 void AliMUON::Digits2Reco()
493 {
494   FindClusters();
495   Int_t nev = gAlice->GetHeader()->GetEvent();
496   GetMUONData()->Fill("RC"); //Filling Reconstructed Cluster
497   fLoader->WriteRecPoints("OVERWRITE");
498   GetMUONData()->ResetRawClusters();        
499   Info("Digits2Reco","End of cluster finding for event %d", nev);
500 }
501 //____________________________________________________________________
502 void AliMUON::FindClusters()
503 {
504 //
505 //  Perform cluster finding
506 //
507     TClonesArray *dig1, *dig2;
508     Int_t ndig, k;
509     dig1 = new TClonesArray("AliMUONDigit",1000);
510     dig2 = new TClonesArray("AliMUONDigit",1000);
511     AliMUONDigit *digit;
512 // Loop on chambers and on cathode planes
513 //
514     ResetRawClusters();        
515     TClonesArray * muonDigits;
516
517     for (Int_t ich = 0; ich < 10; ich++) {
518       //PH      AliMUONChamber* iChamber = (AliMUONChamber*) (*fChambers)[ich];
519         AliMUONChamber* iChamber = (AliMUONChamber*) fChambers->At(ich);
520         AliMUONClusterFinderVS* rec = iChamber->ReconstructionModel();
521     
522         ResetDigits();
523         GetMUONData()->GetCathode(0);
524         //TClonesArray *
525         muonDigits = GetMUONData()->Digits(ich); 
526         ndig=muonDigits->GetEntriesFast();
527         printf("\n 1 Found %d digits in %p chamber %d", ndig, muonDigits,ich);
528         TClonesArray &lhits1 = *dig1;
529         Int_t n = 0;
530         for (k = 0; k < ndig; k++) {
531             digit = (AliMUONDigit*) muonDigits->UncheckedAt(k);
532             if (rec->TestTrack(digit->Track(0)))
533                 new(lhits1[n++]) AliMUONDigit(*digit);
534         }
535         GetMUONData()->ResetDigits();
536         GetMUONData()->GetCathode(1);
537         muonDigits =  GetMUONData()->Digits(ich);  
538         ndig=muonDigits->GetEntriesFast();
539         printf("\n 2 Found %d digits in %p %d", ndig, muonDigits, ich);
540         TClonesArray &lhits2 = *dig2;
541         n=0;
542         
543         for (k=0; k<ndig; k++) {
544             digit= (AliMUONDigit*) muonDigits->UncheckedAt(k);
545             if (rec->TestTrack(digit->Track(0)))
546             new(lhits2[n++]) AliMUONDigit(*digit);
547         }
548
549         if (rec) {       
550             AliMUONClusterInput::Instance()->SetDigits(ich, dig1, dig2);
551             rec->FindRawClusters();
552         }
553         dig1->Delete();
554         dig2->Delete();
555     } // for ich
556     delete dig1;
557     delete dig2;
558 }
559 //______________________________________________________________________
560 #ifdef never
561 void AliMUON::Streamer(TBuffer &R__b)_
562 {
563    // Stream an object of class AliMUON.
564       AliMUONChamber        *iChamber;
565       AliMUONTriggerCircuit *iTriggerCircuit;
566       AliSegmentation       *segmentation;
567       AliMUONResponse       *response;
568       TClonesArray          *digitsaddress;
569       TClonesArray          *rawcladdress;
570       Int_t i;
571       if (R__b.IsReading()) {
572           Version_t R__v = R__b.ReadVersion(); if (R__v) { }
573           AliDetector::Streamer(R__b);
574           R__b >> fNPadHits;
575           R__b >> fPadHits; // diff
576           R__b >> fNLocalTrigger;       
577           R__b >> fLocalTrigger;       
578           R__b >> fNGlobalTrigger;       
579           R__b >> fGlobalTrigger;   
580           R__b >> fDchambers;
581           R__b >> fRawClusters;
582           R__b.ReadArray(fNdch);
583           R__b.ReadArray(fNrawch);
584           R__b >> fAccCut;
585           R__b >> fAccMin;
586           R__b >> fAccMax; 
587           R__b >> fChambers;
588           R__b >> fTriggerCircuits;
589           for (i =0; i<AliMUONConstants::NTriggerCircuit(); i++) {
590               iTriggerCircuit=(AliMUONTriggerCircuit*) (*fTriggerCircuits)[i];
591               iTriggerCircuit->Streamer(R__b);
592           }
593 // Stream chamber related information
594           for (i =0; i<AliMUONConstants::NCh(); i++) {
595               iChamber=(AliMUONChamber*) (*fChambers)[i];
596               iChamber->Streamer(R__b);
597               if (iChamber->Nsec()==1) {
598                   segmentation=iChamber->SegmentationModel(1);
599                   if (segmentation)
600                       segmentation->Streamer(R__b);
601               } else {
602                   segmentation=iChamber->SegmentationModel(1);
603                   if (segmentation)
604                       segmentation->Streamer(R__b);
605                   if (segmentation)
606                       segmentation=iChamber->SegmentationModel(2);
607                   segmentation->Streamer(R__b);
608               }
609               response=iChamber->ResponseModel();
610               if (response)
611                   response->Streamer(R__b);       
612               digitsaddress=(TClonesArray*) (*fDchambers)[i];
613               digitsaddress->Streamer(R__b);
614               if (i < AliMUONConstants::NTrackingCh()) {
615                   rawcladdress=(TClonesArray*) (*fRawClusters)[i];
616                   rawcladdress->Streamer(R__b);
617               }
618           }
619           
620       } else {
621           R__b.WriteVersion(AliMUON::IsA());
622           AliDetector::Streamer(R__b);
623           R__b << fNPadHits;
624           R__b << fPadHits; // diff
625           R__b << fNLocalTrigger;       
626           R__b << fLocalTrigger;       
627           R__b << fNGlobalTrigger;       
628           R__b << fGlobalTrigger; 
629           R__b << fDchambers;
630           R__b << fRawClusters;
631           R__b.WriteArray(fNdch, AliMUONConstants::NCh());
632           R__b.WriteArray(fNrawch, AliMUONConstants::NTrackingCh());
633           
634           R__b << fAccCut;
635           R__b << fAccMin;
636           R__b << fAccMax; 
637           
638           R__b << fChambers;
639           R__b << fTriggerCircuits;
640           for (i =0; i<AliMUONConstants::NTriggerCircuit(); i++) {
641               iTriggerCircuit=(AliMUONTriggerCircuit*) (*fTriggerCircuits)[i];
642               iTriggerCircuit->Streamer(R__b);
643           }
644           for (i =0; i<AliMUONConstants::NCh(); i++) {
645               iChamber=(AliMUONChamber*) (*fChambers)[i];
646               iChamber->Streamer(R__b);
647               if (iChamber->Nsec()==1) {
648                   segmentation=iChamber->SegmentationModel(1);
649                   if (segmentation)
650                       segmentation->Streamer(R__b);
651               } else {
652                   segmentation=iChamber->SegmentationModel(1);
653                   if (segmentation)
654                       segmentation->Streamer(R__b);
655                   segmentation=iChamber->SegmentationModel(2);
656                   if (segmentation)
657                       segmentation->Streamer(R__b);
658               }
659               response=iChamber->ResponseModel();
660               if (response)
661                   response->Streamer(R__b);
662               digitsaddress=(TClonesArray*) (*fDchambers)[i];
663               digitsaddress->Streamer(R__b);
664               if (i < AliMUONConstants::NTrackingCh()) {
665                   rawcladdress=(TClonesArray*) (*fRawClusters)[i];
666                   rawcladdress->Streamer(R__b);
667               }
668           }
669       }
670 }
671 #endif
672 //_______________________________________________________________________
673 AliMUONPadHit* AliMUON::FirstPad(AliMUONHit*  hit, TClonesArray *clusters) 
674 {
675 // to be removed
676     // Initialise the pad iterator
677     // Return the address of the first padhit for hit
678     TClonesArray *theClusters = clusters;
679     Int_t nclust = theClusters->GetEntriesFast();
680     if (nclust && hit->PHlast() > 0) {
681         AliMUON::fMaxIterPad=hit->PHlast();
682         AliMUON::fCurIterPad=hit->PHfirst();
683         return (AliMUONPadHit*) clusters->UncheckedAt(AliMUON::fCurIterPad-1);
684     } else {
685         return 0;
686     }
687 }
688 //_______________________________________________________________________
689 AliMUONPadHit* AliMUON::NextPad(TClonesArray *clusters) 
690 {
691   // To be removed
692 // Get next pad (in iterator) 
693 //
694     AliMUON::fCurIterPad++;
695     if (AliMUON::fCurIterPad <= AliMUON::fMaxIterPad) {
696         return (AliMUONPadHit*) clusters->UncheckedAt(AliMUON::fCurIterPad-1);
697     } else {
698         return 0;
699     }
700 }
701 //_______________________________________________________________________
702
703 AliMUONRawCluster *AliMUON::RawCluster(Int_t ichamber, Int_t icathod, Int_t icluster)
704 {
705 //
706 //  Return rawcluster (icluster) for chamber ichamber and cathode icathod
707 //  Obsolete ??
708     TClonesArray *muonRawCluster  = GetMUONData()->RawClusters(ichamber);
709     ResetRawClusters();
710     TTree *treeR = fLoader->TreeR();
711     Int_t nent=(Int_t)treeR->GetEntries();
712     treeR->GetEvent(nent-2+icathod-1);
713     //treeR->GetEvent(icathod);
714     //Int_t nrawcl = (Int_t)muonRawCluster->GetEntriesFast();
715
716     AliMUONRawCluster * mRaw = (AliMUONRawCluster*)muonRawCluster->UncheckedAt(icluster);
717     //printf("RawCluster _ nent nrawcl icluster mRaw %d %d %d%p\n",nent,nrawcl,icluster,mRaw);
718     
719     return  mRaw;
720 }
721 //________________________________________________________________________
722 void   AliMUON::SetMerger(AliMUONMerger* merger)
723 {
724 // Set pointer to merger 
725     fMerger = merger;
726 }
727 //________________________________________________________________________
728 AliMUONMerger*  AliMUON::Merger()
729 {
730 // Return pointer to merger
731     return fMerger;
732 }
733 //________________________________________________________________________
734 AliMUON& AliMUON::operator = (const AliMUON& /*rhs*/)
735 {
736 // copy operator
737 // dummy version
738     return *this;
739 }
740 //________________________________________________________________________
741 void AliMUON::Reconstruct() const
742 {
743 // reconstruct clusters
744
745   AliLoader* loader = GetLoader();
746   AliRunLoader* runLoader = loader->GetRunLoader();
747   Int_t nEvents = runLoader->GetNumberOfEvents();
748
749   for (Int_t i = 0; i < 10; i++) {
750     AliMUONClusterFinderVS *RecModel = new AliMUONClusterFinderVS();
751     RecModel->SetGhostChi2Cut(10);
752     const_cast<AliMUON*>(this)->SetReconstructionModel(i,RecModel);
753   } 
754
755   AliMUONEventReconstructor* reco = new AliMUONEventReconstructor();
756
757   loader->LoadRecPoints("RECREATE");
758   loader->LoadTracks("RECREATE");
759   loader->LoadDigits("READ");
760
761
762   //   Loop over events              
763   for(Int_t ievent = 0; ievent < nEvents; ievent++) {
764     printf("event %d\n",ievent);
765     runLoader->GetEvent(ievent);
766
767     //---------------------------- digit2Reco & Trigger ---------------------
768     if (!loader->TreeR()) loader->MakeRecPointsContainer();
769      
770     // tracking branch
771     fMUONData->MakeBranch("RC");
772     fMUONData->SetTreeAddress("D,RC");
773     const_cast<AliMUON*>(this)->Digits2Reco(); 
774
775     // trigger branch
776     fMUONData->MakeBranch("GLT");
777     fMUONData->SetTreeAddress("D,GLT");
778     const_cast<AliMUON*>(this)->Trigger(ievent); 
779
780     //---------------------------- Track & TriggerTrack ---------------------
781     if (!loader->TreeT()) loader->MakeTracksContainer();
782
783     fMUONData->MakeBranch("RT"); //track
784     fMUONData->SetTreeAddress("RT");
785     reco->EventReconstruct();
786     for(Int_t i=0; i<reco->GetNRecTracks(); i++) {
787       AliMUONTrack * track = (AliMUONTrack*) reco->GetRecTracksPtr()->At(i);
788       fMUONData->AddRecTrack(*track);
789     }
790     fMUONData->Fill("RT");
791
792     fMUONData->MakeBranch("RL"); //trigger track
793     fMUONData->SetTreeAddress("RL");
794     reco->EventReconstructTrigger();
795     for(Int_t i=0; i<reco->GetNRecTriggerTracks(); i++) {
796       AliMUONTriggerTrack * triggertrack = (AliMUONTriggerTrack*) reco->GetRecTriggerTracksPtr()->At(i);
797       fMUONData->AddRecTriggerTrack(*triggertrack);
798     }
799     fMUONData->Fill("RL");
800
801     loader->WriteTracks("OVERWRITE");  
802   
803     //--------------------------- Resetting branches -----------------------
804     fMUONData->ResetDigits();
805     fMUONData->ResetRawClusters();
806     fMUONData->ResetTrigger();
807     fMUONData->ResetRecTracks();
808     fMUONData->ResetRecTriggerTracks();
809   }
810   loader->UnloadDigits();
811   loader->UnloadRecPoints();
812   loader->UnloadTracks();
813 }
814 //________________________________________________________________________
815 void AliMUON::FillESD(AliESD* event) const
816 {
817
818   TClonesArray * recTracksArray;
819   
820   //YS AliLoader* loader = GetLoader();
821   AliRunLoader* runLoader = fLoader->GetRunLoader(); //YS loader->GetRunLoader();
822   fLoader->LoadTracks("READ"); //YS
823
824
825   // declaration  
826   Int_t ievent;
827   Int_t ntrackhits;
828   Double_t fitfmin;
829   Int_t nrectracks;
830
831   Double_t bendingSlope, nonBendingSlope, fInverseBendingMomentum;
832   Double_t fXRec, fYRec, fZRec;
833
834   Float_t  x11, y11, thetaX,thetaY ;
835
836   //YS Int_t nEvents = runLoader->GetNumberOfEvents();
837
838   // setting pointer for tracks, triggertracks& trackparam at vertex
839   AliMUONTrack * rectrack;
840   AliMUONTriggerTrack * rectriggertrack;
841   AliMUONTrackParam *trackParam;
842   
843   ievent = runLoader->GetEventNumber() ; //YS 
844   //YS for (ievent = 0; ievent < nEvents; ievent++) {
845     runLoader->GetEvent(ievent);
846
847     // setting ESD MUON class
848     AliESDMuonTrack* ESDTrack = new  AliESDMuonTrack() ;
849
850     // -------------------- tracks-------------
851     fMUONData->SetTreeAddress("RT");
852     fMUONData->GetRecTracks();
853     recTracksArray = fMUONData->RecTracks();
854         
855     nrectracks = (Int_t) recTracksArray->GetEntriesFast(); //
856  
857     // printf(">>> Event %d Number of Recconstructed tracks %d \n",ievent, nrectracks);
858    
859     // read track infos
860     for (Int_t irectracks = 0; irectracks <  nrectracks;  irectracks++) {
861
862       rectrack = (AliMUONTrack*) recTracksArray->At(irectracks);
863
864       trackParam = rectrack->GetTrackParamAtVertex();
865
866       bendingSlope            = trackParam->GetBendingSlope();
867       nonBendingSlope         = trackParam->GetNonBendingSlope();
868       fInverseBendingMomentum = trackParam->GetInverseBendingMomentum();
869       fXRec  = trackParam->GetNonBendingCoor();
870       fYRec  = trackParam->GetBendingCoor();
871       fZRec  = trackParam->GetZ();
872
873       ntrackhits = rectrack->GetNTrackHits();
874       fitfmin = rectrack->GetFitFMin();
875
876       // setting data member of ESD MUON
877       ESDTrack->SetInverseBendingMomentum(fInverseBendingMomentum);
878       ESDTrack->SetThetaX(TMath::ATan(nonBendingSlope));
879       ESDTrack->SetThetaY(TMath::ATan(bendingSlope));
880       ESDTrack->SetZ(fZRec);
881       ESDTrack->SetBendingCoor(fYRec);
882       ESDTrack->SetNonBendingCoor(fXRec);
883       ESDTrack->SetChi2(fitfmin);
884       ESDTrack->SetNHit(ntrackhits);
885     }
886
887  //    // -------------------- trigger tracks-------------
888     fMUONData->SetTreeAddress("RL");
889     fMUONData->GetRecTriggerTracks();
890     recTracksArray = fMUONData->RecTriggerTracks();
891         
892     Int_t ntrectracks = (Int_t) recTracksArray->GetEntriesFast(); //YS
893  
894     //  printf(">>> Event %d Number of Recconstructed tracks %d \n",ievent, nrectracks);
895    
896     // read trigger track infos
897     for (Int_t irectracks = 0; irectracks <  ntrectracks;  irectracks++) {
898
899       rectriggertrack = (AliMUONTriggerTrack*) recTracksArray->At(irectracks);
900     
901       x11 = rectriggertrack->GetY11();
902       y11 = rectriggertrack->GetY11();
903       thetaX = rectriggertrack->GetThetax();
904       thetaY = rectriggertrack->GetThetay();
905
906       // setting data member of ESD MUON trigger
907       ESDTrack->SetThetaX11(thetaX);
908       ESDTrack->SetThetaY11(thetaY);
909       ESDTrack->SetX11(x11);
910       ESDTrack->SetY11(y11);
911     }
912
913     // storing ESD MUON Track into ESD Event & reset muondata
914     if (ntrectracks+ntrectracks != 0) //YS 
915       event->AddMuonTrack(ESDTrack);
916     fMUONData->ResetRecTracks();
917     fMUONData->ResetRecTriggerTracks();
918
919     //YS } // end loop on event  
920     fLoader->UnloadTracks(); //YS
921 }
922