1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
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. //
25 // This class is also responsible for building the geometry of the //
30 <img src="picts/AliModuleClass.gif">
34 ///////////////////////////////////////////////////////////////////////////////
37 #include <TObjArray.h>
38 #include <TClonesArray.h>
41 #include <TDirectory.h>
42 #include <TVirtualMC.h>
44 #include "AliConfig.h"
45 #include "AliLoader.h"
47 #include "AliModule.h"
49 #include "AliTrackReference.h"
54 //_______________________________________________________________________
55 AliModule::AliModule():
69 fCurrentIterTrackRef(0)
72 // Default constructor for the AliModule class
76 //_______________________________________________________________________
77 AliModule::AliModule(const char* name,const char *title):
81 fIdtmed(new TArrayI(100)),
82 fIdmate(new TArrayI(100)),
86 fHistograms(new TList()),
90 fTrackReferences(new TClonesArray("AliTrackReference", 100)),
92 fCurrentIterTrackRef(0)
95 // Normal constructor invoked by all Modules.
96 // Create the list for Module specific histograms
97 // Add this Module to the global list of Modules in Run.
99 // Get the Module numeric ID
100 Int_t id = gAlice->GetModuleID(name);
102 // Module already added !
103 Warning("Ctor","Module: %s already present at %d\n",name,id);
107 // Add this Module to the list of Modules
109 gAlice->AddModule(this);
113 // Clear space for tracking media and material indexes
115 for(Int_t i=0;i<100;i++) (*fIdmate)[i]=(*fIdtmed)[i]=0;
118 SetDebug(gAlice->GetDebug());
121 //_______________________________________________________________________
122 AliModule::AliModule(const AliModule &mod):
140 fCurrentIterTrackRef(0)
148 //_______________________________________________________________________
149 AliModule::~AliModule()
155 // Remove this Module from the list of Modules
157 TObjArray * modules = gAlice->Modules();
158 if (modules) modules->Remove(this);
160 // Delete ROOT geometry
168 fHistograms->Clear();
172 // Delete track references
173 if (fTrackReferences) {
174 fTrackReferences->Delete();
175 delete fTrackReferences;
176 fTrackReferences = 0;
178 // Delete TArray objects
184 //_______________________________________________________________________
185 void AliModule::Copy(AliModule & /* mod */) const
188 // Copy *this onto mod, not implemented for AliModule
190 Fatal("Copy","Not implemented!\n");
193 //_______________________________________________________________________
194 void AliModule::Disable()
197 // Disable Module on viewer
203 // Loop through geometry to disable all
204 // nodes for this Module
205 while((node = dynamic_cast<TNode*>(next()))) {
206 node->SetVisibility(-1);
210 //_______________________________________________________________________
211 Int_t AliModule::DistancetoPrimitive(Int_t, Int_t) const
214 // Return distance from mouse pointer to object
215 // Dummy routine for the moment
220 //_______________________________________________________________________
221 void AliModule::Enable()
224 // Enable Module on the viewver
230 // Loop through geometry to enable all
231 // nodes for this Module
232 while((node = dynamic_cast<TNode*>(next()))) {
233 node->SetVisibility(3);
237 //_______________________________________________________________________
238 void AliModule::AliMaterial(Int_t imat, const char* name, Float_t a,
239 Float_t z, Float_t dens, Float_t radl,
240 Float_t absl, Float_t *buf, Int_t nwbuf) const
243 // Store the parameters for a material
245 // imat the material index will be stored in (*fIdmate)[imat]
246 // name material name
250 // radl radiation length
251 // absl absorbtion length
252 // buf adress of an array user words
253 // nwbuf number of user words
256 gMC->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
257 (*fIdmate)[imat]=kmat;
260 //_______________________________________________________________________
261 void AliModule::AliGetMaterial(Int_t imat, char* name, Float_t &a,
262 Float_t &z, Float_t &dens, Float_t &radl,
266 // Store the parameters for a material
268 // imat the material index will be stored in (*fIdmate)[imat]
269 // name material name
273 // radl radiation length
274 // absl absorbtion length
275 // buf adress of an array user words
276 // nwbuf number of user words
281 kmat=(*fIdmate)[imat];
282 gMC->Gfmate(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
286 //_______________________________________________________________________
287 void AliModule::AliMixture(Int_t imat, const char *name, Float_t *a,
288 Float_t *z, Float_t dens, Int_t nlmat,
292 // Defines mixture or compound imat as composed by
293 // nlmat materials defined by arrays a, z and wmat
295 // If nlmat > 0 wmat contains the proportion by
296 // weights of each basic material in the mixture
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.
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
307 // nlmat number of components
308 // wmat array of concentrations
311 gMC->Mixture(kmat, name, a, z, dens, nlmat, wmat);
312 (*fIdmate)[imat]=kmat;
315 //_______________________________________________________________________
316 void AliModule::AliMedium(Int_t numed, const char *name, Int_t nmat,
317 Int_t isvol, Int_t ifield, Float_t fieldm,
318 Float_t tmaxfd, Float_t stemax, Float_t deemax,
319 Float_t epsil, Float_t stmin, Float_t *ubuf,
323 // Store the parameters of a tracking medium
325 // numed the medium number is stored into (*fIdtmed)[numed]
327 // nmat the material number is stored into (*fIdmate)[nmat]
328 // isvol sensitive volume if isvol!=0
329 // ifield magnetic field flag (see below)
330 // fieldm maximum magnetic field
331 // tmaxfd maximum deflection angle due to magnetic field
332 // stemax maximum step allowed
333 // deemax maximum fractional energy loss in one step
334 // epsil tracking precision in cm
335 // stmin minimum step due to continuous processes
337 // ifield = 0 no magnetic field
338 // = -1 user decision in guswim
339 // = 1 tracking performed with Runge Kutta
340 // = 2 tracking performed with helix
341 // = 3 constant magnetic field along z
344 gMC->Medium(kmed,name, (*fIdmate)[nmat], isvol, ifield, fieldm,
345 tmaxfd, stemax, deemax, epsil, stmin, ubuf, nbuf);
346 (*fIdtmed)[numed]=kmed;
349 //_______________________________________________________________________
350 void AliModule::AliMatrix(Int_t &nmat, Float_t theta1, Float_t phi1,
351 Float_t theta2, Float_t phi2, Float_t theta3,
355 // Define a rotation matrix. Angles are in degrees.
357 // nmat on output contains the number assigned to the rotation matrix
358 // theta1 polar angle for axis I
359 // phi1 azimuthal angle for axis I
360 // theta2 polar angle for axis II
361 // phi2 azimuthal angle for axis II
362 // theta3 polar angle for axis III
363 // phi3 azimuthal angle for axis III
365 gMC->Matrix(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
368 //_______________________________________________________________________
369 Float_t AliModule::ZMin() const
374 //_______________________________________________________________________
375 Float_t AliModule::ZMax() const
380 //_______________________________________________________________________
381 void AliModule::SetEuclidFile(char* material, char* geometry)
384 // Sets the name of the Euclid file
386 fEuclidMaterial=material;
388 fEuclidGeometry=geometry;
390 char* name = new char[strlen(material)];
391 strcpy(name,material);
392 strcpy(&name[strlen(name)-4],".euc");
393 fEuclidGeometry=name;
398 //_______________________________________________________________________
399 void AliModule::ReadEuclid(const char* filnam, char* topvol)
402 // read in the geometry of the detector in euclid file format
404 // id_det : the detector identification (2=its,...)
405 // topvol : return parameter describing the name of the top
406 // volume of geometry.
411 // several changes have been made by miroslav helbich
412 // subroutine is rewrited to follow the new established way of memory
413 // booking for tracking medias and rotation matrices.
414 // all used tracking media have to be defined first, for this you can use
415 // subroutine greutmed.
416 // top volume is searched as only volume not positioned into another
419 Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
420 Int_t ndvmx, nr, flag;
421 char key[5], card[77], natmed[21];
422 char name[5], mother[5], shape[5], konly[5], volst[7000][5];
425 Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
427 const Int_t kMaxRot=5000;
428 Int_t idrot[kMaxRot],istop[7000];
431 // *** The input filnam name will be with extension '.euc'
432 filtmp=gSystem->ExpandPathName(filnam);
433 lun=fopen(filtmp,"r");
436 Error("ReadEuclid","Could not open file %s\n",filnam);
439 //* --- definition of rotation matrix 0 ---
440 TArrayI &idtmed = *fIdtmed;
441 for(i=1; i<kMaxRot; ++i) idrot[i]=-99;
445 for(i=0;i<77;i++) card[i]=0;
446 iret=fscanf(lun,"%77[^\n]",card);
447 if(iret<=0) goto L20;
452 if (!strcmp(key,"TMED")) {
453 sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
454 if( itmed<0 || itmed>=100 ) {
455 Error("ReadEuclid","TMED illegal medium number %d for %s\n",itmed,natmed);
458 //Pad the string with blanks
461 while(i<20) natmed[i++]=' ';
464 if( idtmed[itmed]<=0 ) {
465 Error("ReadEuclid","TMED undefined medium number %d for %s\n",itmed,natmed);
468 gMC->Gckmat(idtmed[itmed],natmed);
470 } else if (!strcmp(key,"ROTM")) {
471 sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
472 if( irot<=0 || irot>=kMaxRot ) {
473 Error("ReadEuclid","ROTM rotation matrix number %d illegal\n",irot);
476 AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
478 } else if (!strcmp(key,"VOLU")) {
479 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
481 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
484 gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
485 //* save the defined volumes
486 strcpy(volst[++nvol],name);
489 } else if (!strcmp(key,"DIVN")) {
490 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
491 gMC->Gsdvn ( name, mother, ndiv, iaxe );
493 } else if (!strcmp(key,"DVN2")) {
494 sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
495 gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
497 } else if (!strcmp(key,"DIVT")) {
498 sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
499 gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
501 } else if (!strcmp(key,"DVT2")) {
502 sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
503 gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
505 } else if (!strcmp(key,"POSI")) {
506 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
507 if( irot<0 || irot>=kMaxRot ) {
508 Error("ReadEuclid","POSI %s#%d rotation matrix number %d illegal\n",name,nr,irot);
511 if( idrot[irot] == -99) {
512 Error("ReadEuclid","POSI %s#%d undefined matrix number %d\n",name,nr,irot);
515 //*** volume name cannot be the top volume
516 for(i=1;i<=nvol;i++) {
517 if (!strcmp(volst[i],name)) istop[i]=0;
520 gMC->Gspos ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
522 } else if (!strcmp(key,"POSP")) {
523 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
524 if( irot<0 || irot>=kMaxRot ) {
525 Error("ReadEuclid","POSP %s#%d rotation matrix number %d illegal\n",name,nr,irot);
528 if( idrot[irot] == -99) {
529 Error("ReadEuclid","POSP %s#%d undefined matrix number %d\n",name,nr,irot);
533 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
536 //*** volume name cannot be the top volume
537 for(i=1;i<=nvol;i++) {
538 if (!strcmp(volst[i],name)) istop[i]=0;
541 gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
544 if (strcmp(key,"END")) goto L10;
545 //* find top volume in the geometry
547 for(i=1;i<=nvol;i++) {
548 if (istop[i] && flag) {
549 Warning("ReadEuclid"," %s is another possible top volume\n",volst[i]);
551 if (istop[i] && !flag) {
552 strcpy(topvol,volst[i]);
553 if(fDebug) printf("%s::ReadEuclid: volume %s taken as a top volume\n",ClassName(),topvol);
558 Warning("ReadEuclid","top volume not found\n");
562 //* commented out only for the not cernlib version
563 if(fDebug) printf("%s::ReadEuclid: file: %s is now read in\n",ClassName(),filnam);
568 Error("ReadEuclid","reading error or premature end of file\n");
571 //_______________________________________________________________________
572 void AliModule::ReadEuclidMedia(const char* filnam)
575 // read in the materials and tracking media for the detector
576 // in euclid file format
578 // filnam: name of the input file
579 // id_det: id_det is the detector identification (2=its,...)
581 // author : miroslav helbich
583 Float_t sxmgmx = gAlice->Field()->Max();
584 Int_t isxfld = gAlice->Field()->Integ();
585 Int_t end, i, iret, itmed;
586 char key[5], card[130], natmed[21], namate[21];
591 Int_t nwbuf, isvol, ifield, nmat;
592 Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
595 for(i=0;i<end;i++) if(filnam[i]=='.') {
600 // *** The input filnam name will be with extension '.euc'
601 if(fDebug) printf("%s::ReadEuclid: The file name is %s\n",ClassName(),filnam); //Debug
602 filtmp=gSystem->ExpandPathName(filnam);
603 lun=fopen(filtmp,"r");
606 Warning("ReadEuclidMedia","Could not open file %s\n",filnam);
610 // Retrieve Mag Field parameters
611 Int_t globField=gAlice->Field()->Integ();
612 Float_t globMaxField=gAlice->Field()->Max();
613 // TArrayI &idtmed = *fIdtmed;
616 for(i=0;i<130;i++) card[i]=0;
617 iret=fscanf(lun,"%4s %[^\n]",key,card);
618 if(iret<=0) goto L20;
622 if (!strcmp(key,"MATE")) {
623 sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
624 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
625 //Pad the string with blanks
628 while(i<20) namate[i++]=' ';
631 AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
632 //* read tracking medium
633 } else if (!strcmp(key,"TMED")) {
634 sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
635 &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
636 &stemax,&deemax,&epsil,&stmin,&nwbuf);
637 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
638 if (ifield<0) ifield=isxfld;
639 if (fieldm<0) fieldm=sxmgmx;
640 //Pad the string with blanks
643 while(i<20) natmed[i++]=' ';
646 AliMedium(itmed,natmed,nmat,isvol,globField,globMaxField,tmaxfd,
647 stemax,deemax,epsil,stmin,ubuf,nwbuf);
648 // (*fImedia)[idtmed[itmed]-1]=id_det;
652 if (strcmp(key,"END")) goto L10;
655 //* commented out only for the not cernlib version
656 if(fDebug) printf("%s::ReadEuclidMedia: file %s is now read in\n",
662 Warning("ReadEuclidMedia","reading error or premature end of file\n");
665 //_______________________________________________________________________
666 void AliModule::RemapTrackReferencesIDs(Int_t *map)
669 // Remapping track reference
670 // Called at finish primary
672 if (!fTrackReferences) return;
673 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
674 AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
676 Int_t newID = map[ref->GetTrack()];
677 if (newID>=0) ref->SetTrack(newID);
680 ref->SetBit(kNotDeleted,kFALSE);
681 fTrackReferences->RemoveAt(i);
685 fTrackReferences->Compress();
690 //_______________________________________________________________________
691 AliTrackReference* AliModule::FirstTrackReference(Int_t track)
694 // Initialise the hit iterator
695 // Return the address of the first hit for track
696 // If track>=0 the track is read from disk
697 // while if track<0 the first hit of the current
702 AliRunLoader* rl = AliRunLoader::GetRunLoader();
704 rl->GetAliRun()->GetMCApp()->ResetTrackReferences();
705 rl->TreeTR()->GetEvent(track);
707 Fatal("FirstTrackReference","AliRunLoader not initialized. Can not proceed");
710 fMaxIterTrackRef = fTrackReferences->GetEntriesFast();
711 fCurrentIterTrackRef = 0;
712 if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
716 //_______________________________________________________________________
717 AliTrackReference* AliModule::NextTrackReference()
720 // Return the next hit for the current track
722 if(fMaxIterTrackRef) {
723 if(++fCurrentIterTrackRef<fMaxIterTrackRef)
724 return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
728 printf("* AliModule::NextTrackReference * TrackReference Iterator called without calling FistTrackReference before\n");
734 //_______________________________________________________________________
735 void AliModule::ResetTrackReferences()
738 // Reset number of hits and the hits array
740 fMaxIterTrackRef = 0;
741 if (fTrackReferences) fTrackReferences->Clear();
744 //_____________________________________________________________________________
746 AliLoader* AliModule::MakeLoader(const char* /*topfoldername*/)
751 //PH Merged with v3-09-08 |
753 //_____________________________________________________________________________
755 void AliModule::SetTreeAddress()
758 // Set branch address for track reference Tree
763 // Branch address for track reference tree
764 TTree *treeTR = TreeTR();
766 if (treeTR && fTrackReferences) {
767 branch = treeTR->GetBranch(GetName());
771 Info("SetTreeAddress","(%s) Setting for TrackRefs",GetName());
772 branch->SetAddress(&fTrackReferences);
776 //can be called before MakeBranch and than does not make sense to issue the warning
778 Warning("SetTreeAddress",
779 "(%s) Failed for Track References. Can not find branch in tree.",
785 //_____________________________________________________________________________
786 void AliModule::AddTrackReference(Int_t label){
788 // add a trackrefernce to the list
789 if (!fTrackReferences) {
790 cerr<<"Container trackrefernce not active\n";
793 Int_t nref = fTrackReferences->GetEntriesFast();
794 TClonesArray &lref = *fTrackReferences;
795 new(lref[nref]) AliTrackReference(label);
799 //_____________________________________________________________________________
800 void AliModule::MakeBranchTR(Option_t */*option*/)
803 // Makes branch in treeTR
805 if(GetDebug()) Info("MakeBranchTR","Making Track Refs. Branch for %s",GetName());
806 TTree * tree = TreeTR();
807 if (fTrackReferences && tree)
809 TBranch *branch = tree->GetBranch(GetName());
812 if(GetDebug()) Info("MakeBranch","Branch %s is already in tree.",GetName());
816 branch = tree->Branch(GetName(),&fTrackReferences);
821 Info("MakeBranchTR","FAILED for %s: tree=%#x fTrackReferences=%#x",
822 GetName(),tree,fTrackReferences);
826 //_____________________________________________________________________________
827 TTree* AliModule::TreeTR()
829 AliRunLoader* rl = AliRunLoader::GetRunLoader();
833 Error("TreeTR","Can not get the run loader");
837 TTree* tree = rl->TreeTR();