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