]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliModule.cxx
Deleting the list of histograms in the destructor
[u/mrichter/AliRoot.git] / STEER / AliModule.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 // Base class for ALICE modules. Both sensitive modules (Modules) and      //
21 // non-sensitive ones are described by this base class. This class           //
22 // supports the hit and digit trees produced by the simulation and also      //
23 // the objects produced by the reconstruction.                               //
24 //                                                                           //
25 // This class is also responsible for building the geometry of the           //
26 // Modules.                                                                //
27 //                                                                           //
28 //Begin_Html
29 /*
30 <img src="picts/AliModuleClass.gif">
31 */
32 //End_Html
33 //                                                                           //
34 ///////////////////////////////////////////////////////////////////////////////
35 #include <TNode.h>
36 #include <TObjArray.h>
37 #include <TClonesArray.h>
38 #include <TTree.h>
39 #include <TSystem.h>
40 #include <TDirectory.h>
41
42 #include "AliModule.h"
43 #include "AliRun.h"
44 #include "AliMagF.h"
45 #include "AliConfig.h"
46 #include "AliTrackReference.h"
47
48 ClassImp(AliModule)
49  
50 //_______________________________________________________________________
51 AliModule::AliModule():
52   fEuclidMaterial(""),
53   fEuclidGeometry(""),
54   fIdtmed(0),
55   fIdmate(0),
56   fLoMedium(0),
57   fHiMedium(0),
58   fActive(0),
59   fHistograms(0),
60   fNodes(0),
61   fDebug(0),
62   fEnable(1),
63   fTrackReferences(0),
64   fMaxIterTrackRef(0),
65   fCurrentIterTrackRef(0)
66 {
67   //
68   // Default constructor for the AliModule class
69   //
70 }
71  
72 //_______________________________________________________________________
73 AliModule::AliModule(const char* name,const char *title):
74   TNamed(name,title),
75   fEuclidMaterial(""),
76   fEuclidGeometry(""),
77   fIdtmed(new TArrayI(100)),
78   fIdmate(new TArrayI(100)),
79   fLoMedium(65536),
80   fHiMedium(0),
81   fActive(0),
82   fHistograms(new TList()),
83   fNodes(new TList()),
84   fDebug(0),
85   fEnable(1),
86   fTrackReferences(new TClonesArray("AliTrackReference", 100)),
87   fMaxIterTrackRef(0),
88   fCurrentIterTrackRef(0)
89 {
90   //
91   // Normal constructor invoked by all Modules.
92   // Create the list for Module specific histograms
93   // Add this Module to the global list of Modules in Run.
94   //
95   // Get the Module numeric ID
96   Int_t id = gAlice->GetModuleID(name);
97   if (id>=0) {
98     // Module already added !
99      Warning("Ctor","Module: %s already present at %d\n",name,id);
100      return;
101   }
102   //
103   // Add this Module to the list of Modules
104   gAlice->Modules()->Add(this);
105   //
106   //
107   SetMarkerColor(3);
108   //
109   // Clear space for tracking media and material indexes
110
111   for(Int_t i=0;i<100;i++) (*fIdmate)[i]=(*fIdtmed)[i]=0;
112
113   AliConfig::Instance()->Add(this);    
114     
115   SetDebug(gAlice->GetDebug());
116 }
117  
118 //_______________________________________________________________________
119 AliModule::AliModule(const AliModule &mod):
120   TNamed(mod),
121   TAttLine(mod),
122   TAttMarker(mod),
123   AliRndm(mod),
124   fEuclidMaterial(""),
125   fEuclidGeometry(""),
126   fIdtmed(0),
127   fIdmate(0),
128   fLoMedium(0),
129   fHiMedium(0),
130   fActive(0),
131   fHistograms(0),
132   fNodes(0),
133   fDebug(0),
134   fEnable(0),
135   fTrackReferences(0),
136   fMaxIterTrackRef(0),
137   fCurrentIterTrackRef(0)
138 {
139   //
140   // Copy constructor
141   //
142   mod.Copy(*this);
143 }
144
145 //_______________________________________________________________________
146 AliModule::~AliModule()
147 {
148   //
149   // Destructor
150   //
151
152   // Remove this Module from the list of Modules
153   if (gAlice) {
154     TObjArray * modules = gAlice->Modules();
155     if (modules) modules->Remove(this);
156   }
157   // Delete ROOT geometry
158   if(fNodes) {
159     fNodes->Clear();
160     delete fNodes;
161   }
162   // Delete histograms
163   if(fHistograms) {
164     fHistograms->Clear();
165     delete fHistograms;
166   }
167   //
168   // Delete TArray objects
169   delete fIdtmed;
170   delete fIdmate;
171 }
172  
173 //_______________________________________________________________________
174 void AliModule::Copy(AliModule & /* mod */) const
175 {
176   //
177   // Copy *this onto mod, not implemented for AliModule
178   //
179   Fatal("Copy","Not implemented!\n");
180 }
181
182 //_______________________________________________________________________
183 void AliModule::Disable()
184 {
185   //
186   // Disable Module on viewer
187   //
188   fActive = kFALSE;
189   TIter next(fNodes);
190   TNode *node;
191   //
192   // Loop through geometry to disable all
193   // nodes for this Module
194   while((node = dynamic_cast<TNode*>(next()))) {
195     node->SetVisibility(-1);
196   }   
197 }
198
199 //_______________________________________________________________________
200 Int_t AliModule::DistancetoPrimitive(Int_t, Int_t) const
201 {
202   //
203   // Return distance from mouse pointer to object
204   // Dummy routine for the moment
205   //
206   return 9999;
207 }
208
209 //_______________________________________________________________________
210 void AliModule::Enable()
211 {
212   //
213   // Enable Module on the viewver
214   //
215   fActive = kTRUE;
216   TIter next(fNodes);
217   TNode *node;
218   //
219   // Loop through geometry to enable all
220   // nodes for this Module
221   while((node = dynamic_cast<TNode*>(next()))) {
222     node->SetVisibility(3);
223   }   
224 }
225
226 //_______________________________________________________________________
227 void AliModule::AliMaterial(Int_t imat, const char* name, Float_t a, 
228                             Float_t z, Float_t dens, Float_t radl,
229                             Float_t absl, Float_t *buf, Int_t nwbuf) const
230 {
231   //
232   // Store the parameters for a material
233   //
234   // imat        the material index will be stored in (*fIdmate)[imat]
235   // name        material name
236   // a           atomic mass
237   // z           atomic number
238   // dens        density
239   // radl        radiation length
240   // absl        absorbtion length
241   // buf         adress of an array user words
242   // nwbuf       number of user words
243   //
244   Int_t kmat;
245   gMC->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
246   (*fIdmate)[imat]=kmat;
247 }
248   
249 //_______________________________________________________________________
250 void AliModule::AliGetMaterial(Int_t imat, char* name, Float_t &a, 
251                                Float_t &z, Float_t &dens, Float_t &radl,
252                                Float_t &absl) const
253 {
254   //
255   // Store the parameters for a material
256   //
257   // imat        the material index will be stored in (*fIdmate)[imat]
258   // name        material name
259   // a           atomic mass
260   // z           atomic number
261   // dens        density
262   // radl        radiation length
263   // absl        absorbtion length
264   // buf         adress of an array user words
265   // nwbuf       number of user words
266   //
267
268   Float_t buf[10];
269   Int_t nwbuf, kmat;
270   kmat=(*fIdmate)[imat];
271   gMC->Gfmate(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
272 }
273   
274
275 //_______________________________________________________________________
276 void AliModule::AliMixture(Int_t imat, const char *name, Float_t *a,
277                            Float_t *z, Float_t dens, Int_t nlmat,
278                            Float_t *wmat) const
279
280   //
281   // Defines mixture or compound imat as composed by 
282   // nlmat materials defined by arrays a, z and wmat
283   // 
284   // If nlmat > 0 wmat contains the proportion by
285   // weights of each basic material in the mixture  
286   // 
287   // If nlmat < 0 wmat contains the number of atoms 
288   // of eack kind in the molecule of the compound
289   // In this case, wmat is changed on output to the relative weigths.
290   //
291   // imat        the material index will be stored in (*fIdmate)[imat]
292   // name        material name
293   // a           array of atomic masses
294   // z           array of atomic numbers
295   // dens        density
296   // nlmat       number of components
297   // wmat        array of concentrations
298   //
299   Int_t kmat;
300   gMC->Mixture(kmat, name, a, z, dens, nlmat, wmat);
301   (*fIdmate)[imat]=kmat;
302
303  
304 //_______________________________________________________________________
305 void AliModule::AliMedium(Int_t numed, const char *name, Int_t nmat,
306                           Int_t isvol, Int_t ifield, Float_t fieldm,
307                           Float_t tmaxfd, Float_t stemax, Float_t deemax,
308                           Float_t epsil, Float_t stmin, Float_t *ubuf,
309                           Int_t nbuf) const
310
311   //
312   // Store the parameters of a tracking medium
313   //
314   // numed       the medium number is stored into (*fIdtmed)[numed]
315   // name        medium name
316   // nmat        the material number is stored into (*fIdmate)[nmat]
317   // isvol       sensitive volume if isvol!=0
318   // ifield      magnetic field flag (see below)
319   // fieldm      maximum magnetic field
320   // tmaxfd      maximum deflection angle due to magnetic field
321   // stemax      maximum step allowed
322   // deemax      maximum fractional energy loss in one step
323   // epsil       tracking precision in cm
324   // stmin       minimum step due to continuous processes
325   //
326   // ifield =  0       no magnetic field
327   //        = -1       user decision in guswim
328   //        =  1       tracking performed with Runge Kutta
329   //        =  2       tracking performed with helix
330   //        =  3       constant magnetic field along z
331   //  
332   Int_t kmed;
333   gMC->Medium(kmed,name, (*fIdmate)[nmat], isvol, ifield, fieldm,
334                          tmaxfd, stemax, deemax, epsil, stmin, ubuf, nbuf); 
335   (*fIdtmed)[numed]=kmed;
336
337  
338 //_______________________________________________________________________
339 void AliModule::AliMatrix(Int_t &nmat, Float_t theta1, Float_t phi1,
340                           Float_t theta2, Float_t phi2, Float_t theta3,
341                           Float_t phi3) const
342 {
343   // 
344   // Define a rotation matrix. Angles are in degrees.
345   //
346   // nmat        on output contains the number assigned to the rotation matrix
347   // theta1      polar angle for axis I
348   // phi1        azimuthal angle for axis I
349   // theta2      polar angle for axis II
350   // phi2        azimuthal angle for axis II
351   // theta3      polar angle for axis III
352   // phi3        azimuthal angle for axis III
353   //
354   gMC->Matrix(nmat, theta1, phi1, theta2, phi2, theta3, phi3); 
355
356
357 //_______________________________________________________________________
358 Float_t AliModule::ZMin() const
359 {
360   return -500;
361 }
362
363 //_______________________________________________________________________
364 Float_t AliModule::ZMax() const
365 {
366   return 500;
367 }
368
369 //_______________________________________________________________________
370 void AliModule::SetEuclidFile(char* material, char* geometry)
371 {
372   //
373   // Sets the name of the Euclid file
374   //
375   fEuclidMaterial=material;
376   if(geometry) {
377     fEuclidGeometry=geometry;
378   } else {
379     char* name = new char[strlen(material)];
380     strcpy(name,material);
381     strcpy(&name[strlen(name)-4],".euc");
382     fEuclidGeometry=name;
383     delete [] name;
384   }
385 }
386  
387 //_______________________________________________________________________
388 void AliModule::ReadEuclid(const char* filnam, char* topvol)
389 {
390   //                                                                     
391   //       read in the geometry of the detector in euclid file format    
392   //                                                                     
393   //        id_det : the detector identification (2=its,...)            
394   //        topvol : return parameter describing the name of the top    
395   //        volume of geometry.                                          
396   //                                                                     
397   //            author : m. maire                                        
398   //                                                                     
399   //     28.07.98
400   //     several changes have been made by miroslav helbich
401   //     subroutine is rewrited to follow the new established way of memory
402   //     booking for tracking medias and rotation matrices.
403   //     all used tracking media have to be defined first, for this you can use
404   //     subroutine  greutmed.
405   //     top volume is searched as only volume not positioned into another 
406   //
407
408   Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
409   Int_t ndvmx, nr, flag;
410   char key[5], card[77], natmed[21];
411   char name[5], mother[5], shape[5], konly[5], volst[7000][5];
412   char *filtmp;
413   Float_t par[100];
414   Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
415   Float_t xo, yo, zo;
416   const Int_t kMaxRot=5000;
417   Int_t idrot[kMaxRot],istop[7000];
418   FILE *lun;
419   //
420   // *** The input filnam name will be with extension '.euc'
421   filtmp=gSystem->ExpandPathName(filnam);
422   lun=fopen(filtmp,"r");
423   delete [] filtmp;
424   if(!lun) {
425     Error("ReadEuclid","Could not open file %s\n",filnam);
426     return;
427   }
428   //* --- definition of rotation matrix 0 ---  
429   TArrayI &idtmed = *fIdtmed;
430   for(i=1; i<kMaxRot; ++i) idrot[i]=-99;
431   idrot[0]=0;
432   nvol=0;
433  L10:
434   for(i=0;i<77;i++) card[i]=0;
435   iret=fscanf(lun,"%77[^\n]",card);
436   if(iret<=0) goto L20;
437   fscanf(lun,"%*c");
438   //*
439   strncpy(key,card,4);
440   key[4]='\0';
441   if (!strcmp(key,"TMED")) {
442     sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
443     if( itmed<0 || itmed>=100 ) {
444       Error("ReadEuclid","TMED illegal medium number %d for %s\n",itmed,natmed);
445       exit(1);
446     }
447     //Pad the string with blanks
448     i=-1;
449     while(natmed[++i]);
450     while(i<20) natmed[i++]=' ';
451     natmed[i]='\0';
452     //
453     if( idtmed[itmed]<=0 ) {
454       Error("ReadEuclid","TMED undefined medium number %d for %s\n",itmed,natmed);
455       exit(1);
456     }
457     gMC->Gckmat(idtmed[itmed],natmed);
458     //*
459   } else if (!strcmp(key,"ROTM")) {
460     sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
461     if( irot<=0 || irot>=kMaxRot ) {
462       Error("ReadEuclid","ROTM rotation matrix number %d illegal\n",irot);
463       exit(1);
464     }
465     AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
466     //*
467   } else if (!strcmp(key,"VOLU")) {
468     sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
469     if (npar>0) {
470       for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
471       fscanf(lun,"%*c");
472     }
473     gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
474     //*     save the defined volumes
475     strcpy(volst[++nvol],name);
476     istop[nvol]=1;
477     //*
478   } else if (!strcmp(key,"DIVN")) {
479     sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
480     gMC->Gsdvn  ( name, mother, ndiv, iaxe );
481     //*
482   } else if (!strcmp(key,"DVN2")) {
483     sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
484     gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
485     //*
486   } else if (!strcmp(key,"DIVT")) {
487     sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
488     gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
489     //*
490   } else if (!strcmp(key,"DVT2")) {
491     sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
492     gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
493     //*
494   } else if (!strcmp(key,"POSI")) {
495     sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
496     if( irot<0 || irot>=kMaxRot ) {
497       Error("ReadEuclid","POSI %s#%d rotation matrix number %d illegal\n",name,nr,irot);
498       exit(1);
499     }
500     if( idrot[irot] == -99) {
501       Error("ReadEuclid","POSI %s#%d undefined matrix number %d\n",name,nr,irot);
502       exit(1);
503     }
504     //*** volume name cannot be the top volume
505     for(i=1;i<=nvol;i++) {
506       if (!strcmp(volst[i],name)) istop[i]=0;
507     }
508     //*
509     gMC->Gspos  ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
510     //*
511   } else if (!strcmp(key,"POSP")) {
512     sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
513     if( irot<0 || irot>=kMaxRot ) {
514       Error("ReadEuclid","POSP %s#%d rotation matrix number %d illegal\n",name,nr,irot);
515       exit(1);
516     }
517     if( idrot[irot] == -99) {
518       Error("ReadEuclid","POSP %s#%d undefined matrix number %d\n",name,nr,irot);
519       exit(1);
520     }
521     if (npar > 0) {
522       for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
523       fscanf(lun,"%*c");
524     }
525     //*** volume name cannot be the top volume
526     for(i=1;i<=nvol;i++) {
527       if (!strcmp(volst[i],name)) istop[i]=0;
528     }
529     //*
530     gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
531   }
532   //*
533   if (strcmp(key,"END")) goto L10;
534   //* find top volume in the geometry
535   flag=0;
536   for(i=1;i<=nvol;i++) {
537     if (istop[i] && flag) {
538       Warning("ReadEuclid"," %s is another possible top volume\n",volst[i]);
539     }
540     if (istop[i] && !flag) {
541       strcpy(topvol,volst[i]);
542       if(fDebug) printf("%s::ReadEuclid: volume %s taken as a top volume\n",ClassName(),topvol);
543       flag=1;
544     }
545   }
546   if (!flag) {
547     Warning("ReadEuclid","top volume not found\n");
548   }
549   fclose (lun);
550   //*
551   //*     commented out only for the not cernlib version
552   if(fDebug) printf("%s::ReadEuclid: file: %s is now read in\n",ClassName(),filnam);
553   //
554   return;
555   //*
556   L20:
557   Error("ReadEuclid","reading error or premature end of file\n");
558 }
559
560 //_______________________________________________________________________
561 void AliModule::ReadEuclidMedia(const char* filnam)
562 {
563   //                                                                     
564   //       read in the materials and tracking media for the detector     
565   //                   in euclid file format                             
566   //                                                                     
567   //       filnam: name of the input file                                
568   //       id_det: id_det is the detector identification (2=its,...)     
569   //                                                                     
570   //            author : miroslav helbich                                
571   //
572   Float_t sxmgmx = gAlice->Field()->Max();
573   Int_t   isxfld = gAlice->Field()->Integ();
574   Int_t end, i, iret, itmed;
575   char key[5], card[130], natmed[21], namate[21];
576   Float_t ubuf[50];
577   char* filtmp;
578   FILE *lun;
579   Int_t imate;
580   Int_t nwbuf, isvol, ifield, nmat;
581   Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
582   //
583   end=strlen(filnam);
584   for(i=0;i<end;i++) if(filnam[i]=='.') {
585     end=i;
586     break;
587   }
588   //
589   // *** The input filnam name will be with extension '.euc'
590   if(fDebug) printf("%s::ReadEuclid: The file name is %s\n",ClassName(),filnam); //Debug
591   filtmp=gSystem->ExpandPathName(filnam);
592   lun=fopen(filtmp,"r");
593   delete [] filtmp;
594   if(!lun) {
595     Warning("ReadEuclidMedia","Could not open file %s\n",filnam);
596     return;
597   }
598   //
599   // Retrieve Mag Field parameters
600   Int_t globField=gAlice->Field()->Integ();
601   Float_t globMaxField=gAlice->Field()->Max();
602   //  TArrayI &idtmed = *fIdtmed;
603   //
604  L10:
605   for(i=0;i<130;i++) card[i]=0;
606   iret=fscanf(lun,"%4s %[^\n]",key,card);
607   if(iret<=0) goto L20;
608   fscanf(lun,"%*c");
609   //*
610   //* read material
611   if (!strcmp(key,"MATE")) {
612     sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
613     if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
614     //Pad the string with blanks
615     i=-1;
616     while(namate[++i]);
617     while(i<20) namate[i++]=' ';
618     namate[i]='\0';
619     //
620     AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
621     //* read tracking medium
622   } else if (!strcmp(key,"TMED")) {
623     sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
624            &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
625            &stemax,&deemax,&epsil,&stmin,&nwbuf);
626     if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
627     if (ifield<0) ifield=isxfld;
628     if (fieldm<0) fieldm=sxmgmx;
629     //Pad the string with blanks
630     i=-1;
631     while(natmed[++i]);
632     while(i<20) natmed[i++]=' ';
633     natmed[i]='\0';
634     //
635     AliMedium(itmed,natmed,nmat,isvol,globField,globMaxField,tmaxfd,
636                    stemax,deemax,epsil,stmin,ubuf,nwbuf);
637     //    (*fImedia)[idtmed[itmed]-1]=id_det;
638     //*
639   }
640   //*
641   if (strcmp(key,"END")) goto L10;
642   fclose (lun);
643   //*
644   //*     commented out only for the not cernlib version
645   if(fDebug) printf("%s::ReadEuclidMedia: file %s is now read in\n",
646         ClassName(),filnam);
647   //*
648   return;
649   //*
650  L20:
651   Warning("ReadEuclidMedia","reading error or premature end of file\n");
652
653
654 //_______________________________________________________________________
655 void AliModule::RemapTrackReferencesIDs(Int_t *map)
656 {
657   // 
658   // Remapping track reference
659   // Called at finish primary
660   //
661   if (!fTrackReferences) return;
662   for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
663     AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
664     if (ref) {
665       Int_t newID = map[ref->GetTrack()];
666       if (newID>=0) ref->SetTrack(newID);
667       else ref->SetTrack(-1);
668       
669     }
670   }
671 }
672
673
674 //_______________________________________________________________________
675 AliTrackReference* AliModule::FirstTrackReference(Int_t track)
676 {
677   //
678   // Initialise the hit iterator
679   // Return the address of the first hit for track
680   // If track>=0 the track is read from disk
681   // while if track<0 the first hit of the current
682   // track is returned
683   // 
684   if(track>=0) {
685     gAlice->ResetTrackReferences();
686     gAlice->TreeTR()->GetEvent(track);
687   }
688   //
689   fMaxIterTrackRef     = fTrackReferences->GetEntriesFast();
690   fCurrentIterTrackRef = 0;
691   if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
692   else            return 0;
693 }
694
695 //_______________________________________________________________________
696 AliTrackReference* AliModule::NextTrackReference()
697 {
698   //
699   // Return the next hit for the current track
700   //
701   if(fMaxIterTrackRef) {
702     if(++fCurrentIterTrackRef<fMaxIterTrackRef) 
703       return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
704     else        
705       return 0;
706   } else {
707     printf("* AliDetector::NextTrackReference * TrackReference  Iterator called without calling FistTrackReference before\n");
708     return 0;
709   }
710 }
711
712
713 //_______________________________________________________________________
714 void AliModule::ResetTrackReferences()
715 {
716   //
717   // Reset number of hits and the hits array
718   //
719   fMaxIterTrackRef   = 0;
720   if (fTrackReferences)   fTrackReferences->Clear();
721 }
722  
723  
724
725 void AliModule::SetTreeAddress()
726 {
727   //
728   // Set branch address for the Hits and Digits Trees
729   //
730     TBranch *branch;
731     char branchname[20];
732     sprintf(branchname,"%s",GetName());
733     // Branch address for track reference tree
734     TTree *treeTR = gAlice->TreeTR();
735     if (treeTR && fTrackReferences) {
736         branch = treeTR->GetBranch(branchname);
737         if (branch) branch->SetAddress(&fTrackReferences);
738     }
739 }
740
741 void  AliModule::AddTrackReference(Int_t label){
742   //
743   // add a trackrefernce to the list
744   if (!fTrackReferences) {
745     cerr<<"Container trackrefernce not active\n";
746     return;
747   }
748   Int_t nref = fTrackReferences->GetEntriesFast();
749   TClonesArray &lref = *fTrackReferences;
750   new(lref[nref]) AliTrackReference(label);
751 }
752
753
754 void AliModule::MakeBranchTR(Option_t *option, const char *file)
755
756     //
757     // Makes branch in treeTR
758     //  
759     char name[10];
760     sprintf(name,"%s",GetName());
761
762     if (GetDebug()>1)
763         printf("* MakeBranch * Making Branch %s \n",name);
764     
765     TDirectory *cwd = gDirectory;
766     TBranch *branch = 0;
767     TTree* tree = gAlice->TreeTR();
768     if (tree) {
769         branch = tree->Branch(name, &fTrackReferences, 1600);
770         if (file) {
771             char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
772             sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
773             branch->SetFile(outFile);
774             TIter next( branch->GetListOfBranches());
775             while ((branch=dynamic_cast<TBranch*>(next()))) {
776                 branch->SetFile(outFile);
777             } 
778             delete outFile;
779             
780             cwd->cd();
781             
782             if (GetDebug()>1)
783                 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
784         }
785     }
786 }