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