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