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