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