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