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