]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliModule.cxx
Replacing Header with Id
[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 #include <TNode.h>
36 #include <TObjArray.h>
37 #include <TSystem.h>
38
39 #include "AliModule.h"
40 #include "AliRun.h"
41 #include "AliMagF.h"
42 #include "AliConfig.h"
43
44 ClassImp(AliModule)
45  
46 //_______________________________________________________________________
47 AliModule::AliModule():
48   fEuclidMaterial(""),
49   fEuclidGeometry(""),
50   fIdtmed(0),
51   fIdmate(0),
52   fLoMedium(0),
53   fHiMedium(0),
54   fActive(0),
55   fHistograms(0),
56   fNodes(0),
57   fDebug(0),
58   fEnable(1)
59 {
60   //
61   // Default constructor for the AliModule class
62   //
63 }
64  
65 //_______________________________________________________________________
66 AliModule::AliModule(const char* name,const char *title):
67   TNamed(name,title),
68   fEuclidMaterial(""),
69   fEuclidGeometry(""),
70   fIdtmed(new TArrayI(100)),
71   fIdmate(new TArrayI(100)),
72   fLoMedium(65536),
73   fHiMedium(0),
74   fActive(0),
75   fHistograms(new TList()),
76   fNodes(new TList()),
77   fDebug(0),
78   fEnable(1)
79 {
80   //
81   // Normal constructor invoked by all Modules.
82   // Create the list for Module specific histograms
83   // Add this Module to the global list of Modules in Run.
84   //
85   // Get the Module numeric ID
86   Int_t id = gAlice->GetModuleID(name);
87   if (id>=0) {
88     // Module already added !
89      Warning("Ctor","Module: %s already present at %d\n",name,id);
90      return;
91   }
92   //
93   // Add this Module to the list of Modules
94   gAlice->Modules()->Add(this);
95   //
96   //
97   SetMarkerColor(3);
98   //
99   // Clear space for tracking media and material indexes
100
101   for(Int_t i=0;i<100;i++) (*fIdmate)[i]=(*fIdtmed)[i]=0;
102
103   AliConfig::Instance()->Add(this);    
104     
105   SetDebug(gAlice->GetDebug());
106 }
107  
108 //_______________________________________________________________________
109 AliModule::AliModule(const AliModule &mod):
110   TNamed(mod),
111   TAttLine(mod),
112   TAttMarker(mod),
113   AliRndm(mod),
114   fEuclidMaterial(""),
115   fEuclidGeometry(""),
116   fIdtmed(0),
117   fIdmate(0),
118   fLoMedium(0),
119   fHiMedium(0),
120   fActive(0),
121   fHistograms(0),
122   fNodes(0),
123   fDebug(0),
124   fEnable(0)
125 {
126   //
127   // Copy constructor
128   //
129   mod.Copy(*this);
130 }
131
132 //_______________________________________________________________________
133 AliModule::~AliModule()
134 {
135   //
136   // Destructor
137   //
138
139   // Delete ROOT geometry
140   if(fNodes) {
141     fNodes->Clear();
142     delete fNodes;
143   }
144   //
145   // Delete TArray objects
146   delete fIdtmed;
147   delete fIdmate;
148 }
149  
150 //_______________________________________________________________________
151 void AliModule::Copy(AliModule & /* mod */) const
152 {
153   //
154   // Copy *this onto mod, not implemented for AliModule
155   //
156   Fatal("Copy","Not implemented!\n");
157 }
158
159 //_______________________________________________________________________
160 void AliModule::Disable()
161 {
162   //
163   // Disable Module on viewer
164   //
165   fActive = kFALSE;
166   TIter next(fNodes);
167   TNode *node;
168   //
169   // Loop through geometry to disable all
170   // nodes for this Module
171   while((node = dynamic_cast<TNode*>(next()))) {
172     node->SetVisibility(-1);
173   }   
174 }
175
176 //_______________________________________________________________________
177 Int_t AliModule::DistancetoPrimitive(Int_t, Int_t) const
178 {
179   //
180   // Return distance from mouse pointer to object
181   // Dummy routine for the moment
182   //
183   return 9999;
184 }
185
186 //_______________________________________________________________________
187 void AliModule::Enable()
188 {
189   //
190   // Enable Module on the viewver
191   //
192   fActive = kTRUE;
193   TIter next(fNodes);
194   TNode *node;
195   //
196   // Loop through geometry to enable all
197   // nodes for this Module
198   while((node = dynamic_cast<TNode*>(next()))) {
199     node->SetVisibility(3);
200   }   
201 }
202
203 //_______________________________________________________________________
204 void AliModule::AliMaterial(Int_t imat, const char* name, Float_t a, 
205                             Float_t z, Float_t dens, Float_t radl,
206                             Float_t absl, Float_t *buf, Int_t nwbuf) const
207 {
208   //
209   // Store the parameters for a material
210   //
211   // imat        the material index will be stored in (*fIdmate)[imat]
212   // name        material name
213   // a           atomic mass
214   // z           atomic number
215   // dens        density
216   // radl        radiation length
217   // absl        absorbtion length
218   // buf         adress of an array user words
219   // nwbuf       number of user words
220   //
221   Int_t kmat;
222   gMC->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
223   (*fIdmate)[imat]=kmat;
224 }
225   
226 //_______________________________________________________________________
227 void AliModule::AliGetMaterial(Int_t imat, char* name, Float_t &a, 
228                                Float_t &z, Float_t &dens, Float_t &radl,
229                                Float_t &absl) const
230 {
231   //
232   // Store the parameters for a material
233   //
234   // imat        the material index will be stored in (*fIdmate)[imat]
235   // name        material name
236   // a           atomic mass
237   // z           atomic number
238   // dens        density
239   // radl        radiation length
240   // absl        absorbtion length
241   // buf         adress of an array user words
242   // nwbuf       number of user words
243   //
244
245   Float_t buf[10];
246   Int_t nwbuf, kmat;
247   kmat=(*fIdmate)[imat];
248   gMC->Gfmate(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
249 }
250   
251
252 //_______________________________________________________________________
253 void AliModule::AliMixture(Int_t imat, const char *name, Float_t *a,
254                            Float_t *z, Float_t dens, Int_t nlmat,
255                            Float_t *wmat) const
256
257   //
258   // Defines mixture or compound imat as composed by 
259   // nlmat materials defined by arrays a, z and wmat
260   // 
261   // If nlmat > 0 wmat contains the proportion by
262   // weights of each basic material in the mixture  
263   // 
264   // If nlmat < 0 wmat contains the number of atoms 
265   // of eack kind in the molecule of the compound
266   // In this case, wmat is changed on output to the relative weigths.
267   //
268   // imat        the material index will be stored in (*fIdmate)[imat]
269   // name        material name
270   // a           array of atomic masses
271   // z           array of atomic numbers
272   // dens        density
273   // nlmat       number of components
274   // wmat        array of concentrations
275   //
276   Int_t kmat;
277   gMC->Mixture(kmat, name, a, z, dens, nlmat, wmat);
278   (*fIdmate)[imat]=kmat;
279
280  
281 //_______________________________________________________________________
282 void AliModule::AliMedium(Int_t numed, const char *name, Int_t nmat,
283                           Int_t isvol, Int_t ifield, Float_t fieldm,
284                           Float_t tmaxfd, Float_t stemax, Float_t deemax,
285                           Float_t epsil, Float_t stmin, Float_t *ubuf,
286                           Int_t nbuf) const
287
288   //
289   // Store the parameters of a tracking medium
290   //
291   // numed       the medium number is stored into (*fIdtmed)[numed]
292   // name        medium name
293   // nmat        the material number is stored into (*fIdmate)[nmat]
294   // isvol       sensitive volume if isvol!=0
295   // ifield      magnetic field flag (see below)
296   // fieldm      maximum magnetic field
297   // tmaxfd      maximum deflection angle due to magnetic field
298   // stemax      maximum step allowed
299   // deemax      maximum fractional energy loss in one step
300   // epsil       tracking precision in cm
301   // stmin       minimum step due to continuous processes
302   //
303   // ifield =  0       no magnetic field
304   //        = -1       user decision in guswim
305   //        =  1       tracking performed with Runge Kutta
306   //        =  2       tracking performed with helix
307   //        =  3       constant magnetic field along z
308   //  
309   Int_t kmed;
310   gMC->Medium(kmed,name, (*fIdmate)[nmat], isvol, ifield, fieldm,
311                          tmaxfd, stemax, deemax, epsil, stmin, ubuf, nbuf); 
312   (*fIdtmed)[numed]=kmed;
313
314  
315 //_______________________________________________________________________
316 void AliModule::AliMatrix(Int_t &nmat, Float_t theta1, Float_t phi1,
317                           Float_t theta2, Float_t phi2, Float_t theta3,
318                           Float_t phi3) const
319 {
320   // 
321   // Define a rotation matrix. Angles are in degrees.
322   //
323   // nmat        on output contains the number assigned to the rotation matrix
324   // theta1      polar angle for axis I
325   // phi1        azimuthal angle for axis I
326   // theta2      polar angle for axis II
327   // phi2        azimuthal angle for axis II
328   // theta3      polar angle for axis III
329   // phi3        azimuthal angle for axis III
330   //
331   gMC->Matrix(nmat, theta1, phi1, theta2, phi2, theta3, phi3); 
332
333
334 //_______________________________________________________________________
335 Float_t AliModule::ZMin() const
336 {
337   return -500;
338 }
339
340 //_______________________________________________________________________
341 Float_t AliModule::ZMax() const
342 {
343   return 500;
344 }
345
346 //_______________________________________________________________________
347 void AliModule::SetEuclidFile(char* material, char* geometry)
348 {
349   //
350   // Sets the name of the Euclid file
351   //
352   fEuclidMaterial=material;
353   if(geometry) {
354     fEuclidGeometry=geometry;
355   } else {
356     char* name = new char[strlen(material)];
357     strcpy(name,material);
358     strcpy(&name[strlen(name)-4],".euc");
359     fEuclidGeometry=name;
360     delete [] name;
361   }
362 }
363  
364 //_______________________________________________________________________
365 void AliModule::ReadEuclid(const char* filnam, char* topvol)
366 {
367   //                                                                     
368   //       read in the geometry of the detector in euclid file format    
369   //                                                                     
370   //        id_det : the detector identification (2=its,...)            
371   //        topvol : return parameter describing the name of the top    
372   //        volume of geometry.                                          
373   //                                                                     
374   //            author : m. maire                                        
375   //                                                                     
376   //     28.07.98
377   //     several changes have been made by miroslav helbich
378   //     subroutine is rewrited to follow the new established way of memory
379   //     booking for tracking medias and rotation matrices.
380   //     all used tracking media have to be defined first, for this you can use
381   //     subroutine  greutmed.
382   //     top volume is searched as only volume not positioned into another 
383   //
384
385   Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
386   Int_t ndvmx, nr, flag;
387   char key[5], card[77], natmed[21];
388   char name[5], mother[5], shape[5], konly[5], volst[7000][5];
389   char *filtmp;
390   Float_t par[100];
391   Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
392   Float_t xo, yo, zo;
393   const Int_t kMaxRot=5000;
394   Int_t idrot[kMaxRot],istop[7000];
395   FILE *lun;
396   //
397   // *** The input filnam name will be with extension '.euc'
398   filtmp=gSystem->ExpandPathName(filnam);
399   lun=fopen(filtmp,"r");
400   delete [] filtmp;
401   if(!lun) {
402     Error("ReadEuclid","Could not open file %s\n",filnam);
403     return;
404   }
405   //* --- definition of rotation matrix 0 ---  
406   TArrayI &idtmed = *fIdtmed;
407   for(i=1; i<kMaxRot; ++i) idrot[i]=-99;
408   idrot[0]=0;
409   nvol=0;
410  L10:
411   for(i=0;i<77;i++) card[i]=0;
412   iret=fscanf(lun,"%77[^\n]",card);
413   if(iret<=0) goto L20;
414   fscanf(lun,"%*c");
415   //*
416   strncpy(key,card,4);
417   key[4]='\0';
418   if (!strcmp(key,"TMED")) {
419     sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
420     if( itmed<0 || itmed>=100 ) {
421       Error("ReadEuclid","TMED illegal medium number %d for %s\n",itmed,natmed);
422       exit(1);
423     }
424     //Pad the string with blanks
425     i=-1;
426     while(natmed[++i]);
427     while(i<20) natmed[i++]=' ';
428     natmed[i]='\0';
429     //
430     if( idtmed[itmed]<=0 ) {
431       Error("ReadEuclid","TMED undefined medium number %d for %s\n",itmed,natmed);
432       exit(1);
433     }
434     gMC->Gckmat(idtmed[itmed],natmed);
435     //*
436   } else if (!strcmp(key,"ROTM")) {
437     sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
438     if( irot<=0 || irot>=kMaxRot ) {
439       Error("ReadEuclid","ROTM rotation matrix number %d illegal\n",irot);
440       exit(1);
441     }
442     AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
443     //*
444   } else if (!strcmp(key,"VOLU")) {
445     sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
446     if (npar>0) {
447       for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
448       fscanf(lun,"%*c");
449     }
450     gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
451     //*     save the defined volumes
452     strcpy(volst[++nvol],name);
453     istop[nvol]=1;
454     //*
455   } else if (!strcmp(key,"DIVN")) {
456     sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
457     gMC->Gsdvn  ( name, mother, ndiv, iaxe );
458     //*
459   } else if (!strcmp(key,"DVN2")) {
460     sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
461     gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
462     //*
463   } else if (!strcmp(key,"DIVT")) {
464     sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
465     gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
466     //*
467   } else if (!strcmp(key,"DVT2")) {
468     sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
469     gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
470     //*
471   } else if (!strcmp(key,"POSI")) {
472     sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
473     if( irot<0 || irot>=kMaxRot ) {
474       Error("ReadEuclid","POSI %s#%d rotation matrix number %d illegal\n",name,nr,irot);
475       exit(1);
476     }
477     if( idrot[irot] == -99) {
478       Error("ReadEuclid","POSI %s#%d undefined matrix number %d\n",name,nr,irot);
479       exit(1);
480     }
481     //*** volume name cannot be the top volume
482     for(i=1;i<=nvol;i++) {
483       if (!strcmp(volst[i],name)) istop[i]=0;
484     }
485     //*
486     gMC->Gspos  ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
487     //*
488   } else if (!strcmp(key,"POSP")) {
489     sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
490     if( irot<0 || irot>=kMaxRot ) {
491       Error("ReadEuclid","POSP %s#%d rotation matrix number %d illegal\n",name,nr,irot);
492       exit(1);
493     }
494     if( idrot[irot] == -99) {
495       Error("ReadEuclid","POSP %s#%d undefined matrix number %d\n",name,nr,irot);
496       exit(1);
497     }
498     if (npar > 0) {
499       for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
500       fscanf(lun,"%*c");
501     }
502     //*** volume name cannot be the top volume
503     for(i=1;i<=nvol;i++) {
504       if (!strcmp(volst[i],name)) istop[i]=0;
505     }
506     //*
507     gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
508   }
509   //*
510   if (strcmp(key,"END")) goto L10;
511   //* find top volume in the geometry
512   flag=0;
513   for(i=1;i<=nvol;i++) {
514     if (istop[i] && flag) {
515       Warning("ReadEuclid"," %s is another possible top volume\n",volst[i]);
516     }
517     if (istop[i] && !flag) {
518       strcpy(topvol,volst[i]);
519       if(fDebug) printf("%s::ReadEuclid: volume %s taken as a top volume\n",ClassName(),topvol);
520       flag=1;
521     }
522   }
523   if (!flag) {
524     Warning("ReadEuclid","top volume not found\n");
525   }
526   fclose (lun);
527   //*
528   //*     commented out only for the not cernlib version
529   if(fDebug) printf("%s::ReadEuclid: file: %s is now read in\n",ClassName(),filnam);
530   //
531   return;
532   //*
533   L20:
534   Error("ReadEuclid","reading error or premature end of file\n");
535 }
536
537 //_______________________________________________________________________
538 void AliModule::ReadEuclidMedia(const char* filnam)
539 {
540   //                                                                     
541   //       read in the materials and tracking media for the detector     
542   //                   in euclid file format                             
543   //                                                                     
544   //       filnam: name of the input file                                
545   //       id_det: id_det is the detector identification (2=its,...)     
546   //                                                                     
547   //            author : miroslav helbich                                
548   //
549   Float_t sxmgmx = gAlice->Field()->Max();
550   Int_t   isxfld = gAlice->Field()->Integ();
551   Int_t end, i, iret, itmed;
552   char key[5], card[130], natmed[21], namate[21];
553   Float_t ubuf[50];
554   char* filtmp;
555   FILE *lun;
556   Int_t imate;
557   Int_t nwbuf, isvol, ifield, nmat;
558   Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
559   //
560   end=strlen(filnam);
561   for(i=0;i<end;i++) if(filnam[i]=='.') {
562     end=i;
563     break;
564   }
565   //
566   // *** The input filnam name will be with extension '.euc'
567   if(fDebug) printf("%s::ReadEuclid: The file name is %s\n",ClassName(),filnam); //Debug
568   filtmp=gSystem->ExpandPathName(filnam);
569   lun=fopen(filtmp,"r");
570   delete [] filtmp;
571   if(!lun) {
572     Warning("ReadEuclidMedia","Could not open file %s\n",filnam);
573     return;
574   }
575   //
576   // Retrieve Mag Field parameters
577   Int_t globField=gAlice->Field()->Integ();
578   Float_t globMaxField=gAlice->Field()->Max();
579   //  TArrayI &idtmed = *fIdtmed;
580   //
581  L10:
582   for(i=0;i<130;i++) card[i]=0;
583   iret=fscanf(lun,"%4s %[^\n]",key,card);
584   if(iret<=0) goto L20;
585   fscanf(lun,"%*c");
586   //*
587   //* read material
588   if (!strcmp(key,"MATE")) {
589     sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
590     if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
591     //Pad the string with blanks
592     i=-1;
593     while(namate[++i]);
594     while(i<20) namate[i++]=' ';
595     namate[i]='\0';
596     //
597     AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
598     //* read tracking medium
599   } else if (!strcmp(key,"TMED")) {
600     sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
601            &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
602            &stemax,&deemax,&epsil,&stmin,&nwbuf);
603     if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
604     if (ifield<0) ifield=isxfld;
605     if (fieldm<0) fieldm=sxmgmx;
606     //Pad the string with blanks
607     i=-1;
608     while(natmed[++i]);
609     while(i<20) natmed[i++]=' ';
610     natmed[i]='\0';
611     //
612     AliMedium(itmed,natmed,nmat,isvol,globField,globMaxField,tmaxfd,
613                    stemax,deemax,epsil,stmin,ubuf,nwbuf);
614     //    (*fImedia)[idtmed[itmed]-1]=id_det;
615     //*
616   }
617   //*
618   if (strcmp(key,"END")) goto L10;
619   fclose (lun);
620   //*
621   //*     commented out only for the not cernlib version
622   if(fDebug) printf("%s::ReadEuclidMedia: file %s is now read in\n",
623         ClassName(),filnam);
624   //*
625   return;
626   //*
627  L20:
628   Warning("ReadEuclidMedia","reading error or premature end of file\n");
629
630  
631