]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliModule.cxx
321db45087926bc9a86c4db0ac555e15edecfcff
[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 #include <TString.h>
45
46 #include "AliLog.h"
47 #include "AliConfig.h"
48 #include "AliLoader.h"
49 #include "AliMagF.h"
50 #include "AliModule.h"
51 #include "AliRun.h"
52 #include "AliTrackReference.h"
53 #include "AliMC.h"
54 #include "AliRawDataHeader.h"
55
56 #include "AliDAQ.h"
57
58 ClassImp(AliModule)
59  
60 //_______________________________________________________________________
61 AliModule::AliModule():
62   fEuclidMaterial(""),
63   fEuclidGeometry(""),
64   fIdtmed(0),
65   fIdmate(0),
66   fLoMedium(0),
67   fHiMedium(0),
68   fActive(0),
69   fHistograms(0),
70   fNodes(0),
71   fEnable(1),
72   fTrackReferences(0),
73   fMaxIterTrackRef(0),
74   fCurrentIterTrackRef(0),
75   fRunLoader(0)
76 {
77   //
78   // Default constructor for the AliModule class
79   //
80 }
81  
82 //_______________________________________________________________________
83 AliModule::AliModule(const char* name,const char *title):
84   TNamed(name,title),
85   fEuclidMaterial(""),
86   fEuclidGeometry(""),
87   fIdtmed(new TArrayI(100)),
88   fIdmate(new TArrayI(100)),
89   fLoMedium(65536),
90   fHiMedium(0),
91   fActive(0),
92   fHistograms(new TList()),
93   fNodes(new TList()),
94   fEnable(1),
95   fTrackReferences(new TClonesArray("AliTrackReference", 100)),
96   fMaxIterTrackRef(0),
97   fCurrentIterTrackRef(0),
98   fRunLoader(0)
99 {
100   //
101   // Normal constructor invoked by all Modules.
102   // Create the list for Module specific histograms
103   // Add this Module to the global list of Modules in Run.
104   //
105   // Get the Module numeric ID
106   Int_t id = gAlice->GetModuleID(name);
107   if (id>=0) {
108     // Module already added !
109      AliWarning(Form("Module: %s already present at %d",name,id));
110      return;
111   }
112   //
113   // Add this Module to the list of Modules
114
115   gAlice->AddModule(this);
116
117   SetMarkerColor(3);
118   //
119   // Clear space for tracking media and material indexes
120
121   for(Int_t i=0;i<100;i++) (*fIdmate)[i]=(*fIdtmed)[i]=0;
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   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   //Build the string uniquename as "DET_materialname"
250   TString uniquename = GetName();
251   uniquename.Append("_");
252   uniquename.Append(name);
253   //if geometry loaded from file only fill fIdmate, else create material too
254   if(gAlice->IsRootGeometry()){
255     TGeoMaterial *mat = gGeoManager->GetMaterial(uniquename.Data());
256     kmat = mat->GetUniqueID();
257     (*fIdmate)[imat]=kmat;
258   }else{
259     gMC->Material(kmat, uniquename.Data(), a, z, dens, radl, absl, buf, nwbuf);
260     (*fIdmate)[imat]=kmat;
261   }
262 }
263   
264 //_______________________________________________________________________
265 void AliModule::AliGetMaterial(Int_t imat, char* name, Float_t &a, 
266                                Float_t &z, Float_t &dens, Float_t &radl,
267                                Float_t &absl) const
268 {
269   //
270   // Store the parameters for a material
271   //
272   // imat        the material index will be stored in (*fIdmate)[imat]
273   // name        material name
274   // a           atomic mass
275   // z           atomic number
276   // dens        density
277   // radl        radiation length
278   // absl        absorbtion length
279   // buf         adress of an array user words
280   // nwbuf       number of user words
281   //
282
283   Float_t buf[10];
284   Int_t nwbuf, kmat;
285   kmat=(*fIdmate)[imat];
286   gMC->Gfmate(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
287 }
288   
289
290 //_______________________________________________________________________
291 void AliModule::AliMixture(Int_t imat, const char *name, Float_t *a,
292                            Float_t *z, Float_t dens, Int_t nlmat,
293                            Float_t *wmat) const
294
295   //
296   // Defines mixture or compound imat as composed by 
297   // nlmat materials defined by arrays a, z and wmat
298   // 
299   // If nlmat > 0 wmat contains the proportion by
300   // weights of each basic material in the mixture  
301   // 
302   // If nlmat < 0 wmat contains the number of atoms 
303   // of eack kind in the molecule of the compound
304   // In this case, wmat is changed on output to the relative weigths.
305   //
306   // imat        the material index will be stored in (*fIdmate)[imat]
307   // name        material name
308   // a           array of atomic masses
309   // z           array of atomic numbers
310   // dens        density
311   // nlmat       number of components
312   // wmat        array of concentrations
313   //
314   Int_t kmat;
315   //Build the string uniquename as "DET_mixturename"
316   TString uniquename = GetName();
317   uniquename.Append("_");
318   uniquename.Append(name);
319   //if geometry loaded from file only fill fIdmate, else create mixture too
320   if(gAlice->IsRootGeometry()){
321     TGeoMaterial *mat = gGeoManager->GetMaterial(uniquename.Data());
322     kmat = mat->GetUniqueID();
323     (*fIdmate)[imat]=kmat;
324   }else{
325     gMC->Mixture(kmat, uniquename.Data(), a, z, dens, nlmat, wmat);
326     (*fIdmate)[imat]=kmat;
327   }
328
329  
330 //_______________________________________________________________________
331 void AliModule::AliMedium(Int_t numed, const char *name, Int_t nmat,
332                           Int_t isvol, Int_t ifield, Float_t fieldm,
333                           Float_t tmaxfd, Float_t stemax, Float_t deemax,
334                           Float_t epsil, Float_t stmin, Float_t *ubuf,
335                           Int_t nbuf) const
336
337   //
338   // Store the parameters of a tracking medium
339   //
340   // numed       the medium number is stored into (*fIdtmed)[numed]
341   // name        medium name
342   // nmat        the material number is stored into (*fIdmate)[nmat]
343   // isvol       sensitive volume if isvol!=0
344   // ifield      magnetic field flag (see below)
345   // fieldm      maximum magnetic field
346   // tmaxfd      maximum deflection angle due to magnetic field
347   // stemax      maximum step allowed
348   // deemax      maximum fractional energy loss in one step
349   // epsil       tracking precision in cm
350   // stmin       minimum step due to continuous processes
351   //
352   // ifield =  0       no magnetic field
353   //        = -1       user decision in guswim
354   //        =  1       tracking performed with Runge Kutta
355   //        =  2       tracking performed with helix
356   //        =  3       constant magnetic field along z
357   //  
358   Int_t kmed;
359   //Build the string uniquename as "DET_mediumname"
360   TString uniquename = GetName();
361   uniquename.Append("_");
362   uniquename.Append(name);
363   //if geometry loaded from file only fill fIdtmed, else create medium too
364   if(gAlice->IsRootGeometry()){
365     TGeoMedium *med = gGeoManager->GetMedium(uniquename.Data());
366     kmed = med->GetId();
367     (*fIdtmed)[numed]=kmed;
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::AddAlignableVolumes() const
692 {
693   // 
694   AliWarning(Form(" %s still has to implement the AddAlignableVolumes method!",GetName()));
695 }
696
697 //_______________________________________________________________________
698 void AliModule::RemapTrackReferencesIDs(Int_t *map)
699 {
700   // 
701   // Remapping track reference
702   // Called at finish primary
703   //
704   if (!fTrackReferences) return;
705   Int_t nEntries = fTrackReferences->GetEntries();
706   
707   for (Int_t i=0;i<nEntries;i++){
708     AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
709     if (ref) {
710       Int_t newID = map[ref->GetTrack()];
711       if (newID>=0) ref->SetTrack(newID);
712       else {
713         ref->SetBit(kNotDeleted,kFALSE);
714         fTrackReferences->RemoveAt(i);  
715       }      
716     }
717   }
718   fTrackReferences->Compress();
719
720 }
721
722
723 //_______________________________________________________________________
724 AliTrackReference* AliModule::FirstTrackReference(Int_t track)
725 {
726   //
727   // Initialise the hit iterator
728   // Return the address of the first hit for track
729   // If track>=0 the track is read from disk
730   // while if track<0 the first hit of the current
731   // track is returned
732   // 
733   if(track>=0) 
734    {
735      if (fRunLoader == 0x0)
736        AliFatal("AliRunLoader not initialized. Can not proceed");
737      fRunLoader->GetAliRun()->GetMCApp()->ResetTrackReferences();
738      fRunLoader->TreeTR()->GetEvent(track);
739    }
740   //
741   fMaxIterTrackRef     = fTrackReferences->GetEntriesFast();
742   fCurrentIterTrackRef = 0;
743   if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
744   else            return 0;
745 }
746
747 //_______________________________________________________________________
748 AliTrackReference* AliModule::NextTrackReference()
749 {
750   //
751   // Return the next hit for the current track
752   //
753   if(fMaxIterTrackRef) {
754     if(++fCurrentIterTrackRef<fMaxIterTrackRef) 
755       return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
756     else        
757       return 0;
758   } else {
759     AliWarning("Iterator called without calling FistTrackReference before");
760     return 0;
761   }
762 }
763
764
765 //_______________________________________________________________________
766 void AliModule::ResetTrackReferences()
767 {
768   //
769   // Reset number of hits and the hits array
770   //
771   fMaxIterTrackRef   = 0;
772   if (fTrackReferences)   fTrackReferences->Clear();
773 }
774  
775 //_____________________________________________________________________________
776
777 AliLoader*  AliModule::MakeLoader(const char* /*topfoldername*/)
778 {
779   return 0x0;
780 }
781  
782 //PH Merged with v3-09-08 |
783 //                        V
784 //_____________________________________________________________________________
785
786 void AliModule::SetTreeAddress()
787 {
788   //
789   // Set branch address for track reference Tree
790   //
791
792   TBranch *branch;
793
794   // Branch address for track reference tree
795   TTree *treeTR = TreeTR();
796
797   if (treeTR && fTrackReferences) {
798      branch = treeTR->GetBranch(GetName());
799     if (branch) 
800      {
801        AliDebug(3, Form("(%s) Setting for TrackRefs",GetName()));
802        branch->SetAddress(&fTrackReferences);
803      }
804     else
805      { 
806      //can be called before MakeBranch and than does not make sense to issue the warning
807        AliDebug(1, Form("(%s) Failed for Track References. Can not find branch in tree.",
808                  GetName()));
809      }
810   }
811 }
812
813 //_____________________________________________________________________________
814 AliTrackReference*  AliModule::AddTrackReference(Int_t label){
815   //
816   // add a trackrefernce to the list
817   if (!fTrackReferences) {
818     AliError("Container trackrefernce not active");
819     return 0;
820   }
821   Int_t nref = fTrackReferences->GetEntriesFast();
822   TClonesArray &lref = *fTrackReferences;
823   return new(lref[nref]) AliTrackReference(label);
824 }
825
826
827 //_____________________________________________________________________________
828 void AliModule::MakeBranchTR(Option_t */*option*/)
829
830     //
831     // Makes branch in treeTR
832     //  
833   AliDebug(2,Form("Making Track Refs. Branch for %s",GetName()));
834   TTree * tree = TreeTR();
835   if (fTrackReferences && tree) 
836    {
837       TBranch *branch = tree->GetBranch(GetName());
838      if (branch) 
839        {  
840          AliDebug(2,Form("Branch %s is already in tree.",GetName()));
841          return;
842        }
843   
844      branch = tree->Branch(GetName(),&fTrackReferences);
845    }
846   else
847     {
848       AliDebug(2,Form("FAILED for %s: tree=%#x fTrackReferences=%#x",
849              GetName(),tree,fTrackReferences));
850     }
851 }
852
853 //_____________________________________________________________________________
854 TTree* AliModule::TreeTR()
855 {
856   //
857   // Return TR tree pointer
858   //
859   if ( fRunLoader == 0x0)
860    {
861      AliError("Can not get the run loader");
862      return 0x0;
863    }
864
865   TTree* tree = fRunLoader->TreeTR();
866   return tree;
867 }
868
869
870 //_____________________________________________________________________________
871 void AliModule::Digits2Raw()
872 {
873 // This is a dummy version that just copies the digits file contents
874 // to a raw data file.
875
876   AliWarning(Form("Dummy version called for %s", GetName()));
877
878   Int_t nDDLs = AliDAQ::NumberOfDdls(GetName());
879
880   if (!GetLoader()) return;
881   fstream digitsFile(GetLoader()->GetDigitsFileName(), ios::in);
882   if (!digitsFile) return;
883
884   digitsFile.seekg(0, ios::end);
885   UInt_t size = digitsFile.tellg();
886   UInt_t ddlSize = 4 * (size / (4*nDDLs));
887   Char_t* buffer = new Char_t[ddlSize+1];
888
889   for (Int_t iDDL = 0; iDDL < nDDLs; iDDL++) {
890     char fileName[20];
891     strcpy(fileName,AliDAQ::DdlFileName(GetName(),iDDL));
892     fstream rawFile(fileName, ios::out);
893     if (!rawFile) return;
894
895     AliRawDataHeader header;
896     header.fSize = ddlSize + sizeof(header);
897     rawFile.write((char*) &header, sizeof(header));
898
899     digitsFile.read(buffer, ddlSize);
900     rawFile.write(buffer, ddlSize);
901     rawFile.close();
902   }
903
904   digitsFile.close();
905   delete[] buffer;
906 }