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 Revision 1.21 1999/11/25 10:40:08 fca
19 Fixing daughters information also in primary tracks
21 Revision 1.20 1999/10/04 18:08:49 fca
22 Adding protection against inconsistent Euclid files
24 Revision 1.19 1999/09/29 07:50:40 fca
25 Introduction of the Copyright and cvs Log
29 ///////////////////////////////////////////////////////////////////////////////
31 // Control class for Alice C++ //
32 // Only one single instance of this class exists. //
33 // The object is created in main program aliroot //
34 // and is pointed by the global gAlice. //
36 // -Supports the list of all Alice Detectors (fModules). //
37 // -Supports the list of particles (fParticles). //
38 // -Supports the Trees. //
39 // -Supports the geometry. //
40 // -Supports the event display. //
43 <img src="picts/AliRunClass.gif">
48 <img src="picts/alirun.gif">
52 ///////////////////////////////////////////////////////////////////////////////
60 #include <TObjectTable.h>
62 #include "TParticle.h"
64 #include "AliDisplay.h"
66 #include "AliCallf77.h"
74 static AliHeader *header;
78 # define rxgtrak rxgtrak_
79 # define rxstrak rxstrak_
80 # define rxkeep rxkeep_
81 # define rxouth rxouth_
84 # define rxgtrak RXGTRAK
85 # define rxstrak RXSTRAK
86 # define rxkeep RXKEEP
87 # define rxouth RXOUTH
90 static TArrayF sEventEnergy;
91 static TArrayF sSummEnergy;
92 static TArrayF sSum2Energy;
96 //_____________________________________________________________________________
100 // Default constructor for AliRun
124 fPDGDB = 0; //Particle factory object!
127 //_____________________________________________________________________________
128 AliRun::AliRun(const char *name, const char *title)
132 // Constructor for the main processor.
133 // Creates the geometry
134 // Creates the list of Detectors.
135 // Creates the list of particles.
152 gROOT->GetListOfBrowsables()->Add(this,name);
154 // create the support list for the various Detectors
155 fModules = new TObjArray(77);
157 // Create the TNode geometry for the event display
159 BuildSimpleGeometry();
169 // Create the particle stack
170 fParticles = new TClonesArray("TParticle",100);
174 // Create default mag field
179 // Prepare the tracking medium lists
180 fImedia = new TArrayI(1000);
181 for(i=0;i<1000;i++) (*fImedia)[i]=-99;
184 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
187 //_____________________________________________________________________________
191 // Defaullt AliRun destructor
210 fParticles->Delete();
215 //_____________________________________________________________________________
216 void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
219 // Add a hit to detector id
221 TObjArray &dets = *fModules;
222 if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
225 //_____________________________________________________________________________
226 void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
229 // Add digit to detector id
231 TObjArray &dets = *fModules;
232 if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
235 //_____________________________________________________________________________
236 void AliRun::Browse(TBrowser *b)
239 // Called when the item "Run" is clicked on the left pane
240 // of the Root browser.
241 // It displays the Root Trees and all detectors.
243 if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
244 if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
245 if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
246 if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
247 if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
249 TIter next(fModules);
251 while((detector = (AliModule*)next())) {
252 b->Add(detector,detector->GetName());
256 //_____________________________________________________________________________
260 // Initialize Alice geometry
265 //_____________________________________________________________________________
266 void AliRun::BuildSimpleGeometry()
269 // Create a simple TNode geometry used by Root display engine
271 // Initialise geometry
273 fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
274 new TMaterial("void","Vacuum",0,0,0); //Everything is void
275 TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
276 brik->SetVisibility(0);
277 new TNode("alice","alice","S_alice");
280 //_____________________________________________________________________________
281 void AliRun::CleanDetectors()
284 // Clean Detectors at the end of event
286 TIter next(fModules);
288 while((detector = (AliModule*)next())) {
289 detector->FinishEvent();
293 //_____________________________________________________________________________
294 void AliRun::CleanParents()
297 // Clean Particles stack.
298 // Set parent/daughter relations
300 TClonesArray &particles = *(gAlice->Particles());
303 for(i=0; i<fNtrack; i++) {
304 part = (TParticle *)particles.UncheckedAt(i);
305 if(!part->TestBit(Daughters_Bit)) {
306 part->SetFirstDaughter(-1);
307 part->SetLastDaughter(-1);
312 //_____________________________________________________________________________
313 Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
316 // Return the distance from the mouse to the AliRun object
322 //_____________________________________________________________________________
323 void AliRun::DumpPart (Int_t i)
326 // Dumps particle i in the stack
328 TClonesArray &particles = *fParticles;
329 ((TParticle*) particles[i])->Print();
332 //_____________________________________________________________________________
333 void AliRun::DumpPStack ()
336 // Dumps the particle stack
338 TClonesArray &particles = *fParticles;
340 "\n\n=======================================================================\n");
341 for (Int_t i=0;i<fNtrack;i++)
343 printf("-> %d ",i); ((TParticle*) particles[i])->Print();
344 printf("--------------------------------------------------------------\n");
347 "\n=======================================================================\n\n");
350 //_____________________________________________________________________________
351 void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
352 Float_t maxField, char* filename)
355 // Set magnetic field parameters
356 // type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
357 // version Magnetic field map version (only 1 active now)
358 // scale Scale factor for the magnetic field
359 // maxField Maximum value for the magnetic field
362 // --- Sanity check on mag field flags
363 if(type<0 || type > 2) {
365 "Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
369 if(fField) delete fField;
371 fField = new AliMagFC("Map1"," ",type,version,scale,maxField);
372 } else if(version<=3) {
373 fField = new AliMagFCM("Map2-3",filename,type,version,scale,maxField);
376 Warning("SetField","Invalid map %d\n",version);
380 //_____________________________________________________________________________
381 void AliRun::FillTree()
384 // Fills all AliRun TTrees
386 if (fTreeK) fTreeK->Fill();
387 if (fTreeH) fTreeH->Fill();
388 if (fTreeD) fTreeD->Fill();
389 if (fTreeR) fTreeR->Fill();
392 //_____________________________________________________________________________
393 void AliRun::FinishPrimary()
396 // Called at the end of each primary track
399 // static Int_t count=0;
400 // const Int_t times=10;
401 // This primary is finished, purify stack
402 gAlice->PurifyKine();
404 // Write out hits if any
405 if (gAlice->TreeH()) {
406 gAlice->TreeH()->Fill();
413 // if(++count%times==1) gObjectTable->Print();
416 //_____________________________________________________________________________
417 void AliRun::FinishEvent()
420 // Called at the end of the event.
422 printf("\n\n\n\nEntering FinishEvent \n\n\n");
423 //Update the energy deposit tables
425 for(i=0;i<sEventEnergy.GetSize();i++) {
426 sSummEnergy[i]+=sEventEnergy[i];
427 sSum2Energy[i]+=sEventEnergy[i]*sEventEnergy[i];
429 sEventEnergy.Reset();
431 // Clean detector information
434 // Write out the kinematics
440 // Write out the digits
446 // Write out reconstructed clusters
451 // Write out the event Header information
452 if (fTreeE) fTreeE->Fill();
454 Int_t np=fParticles->GetEntriesFast();
455 printf("Checking %d\n",np);
456 for (Int_t inpart=0; inpart<np;++inpart) {
458 TParticle *particle = (TParticle*)fParticles->UncheckedAt(inpart);
459 Int_t icode = particle->GetPdgCode();
460 Double_t chg = particle->GetPDG()->Charge();
462 if (icode==111 && chg) {
463 printf("%s charge %f %p %s\n",
464 particle->GetName(),chg,particle->GetPDG(),
465 particle->GetPDG()->GetName());
472 // Write Tree headers
473 // Int_t ievent = fHeader.GetEvent();
475 // sprintf(hname,"TreeK%d",ievent);
476 if (fTreeK) fTreeK->Write();
477 // sprintf(hname,"TreeH%d",ievent);
478 if (fTreeH) fTreeH->Write();
479 // sprintf(hname,"TreeD%d",ievent);
480 if (fTreeD) fTreeD->Write();
481 // sprintf(hname,"TreeR%d",ievent);
482 if (fTreeR) fTreeR->Write();
485 //_____________________________________________________________________________
486 void AliRun::FinishRun()
489 // Called at the end of the run.
492 // Clean detector information
493 TIter next(fModules);
495 while((detector = (AliModule*)next())) {
496 detector->FinishRun();
499 //Output energy summary tables
502 // file is retrieved from whatever tree
504 if (fTreeK) File = fTreeK->GetCurrentFile();
505 if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
506 if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
507 if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
509 Error("FinishRun","There isn't root file!");
515 // Clean tree information
516 delete fTreeK; fTreeK = 0;
517 delete fTreeH; fTreeH = 0;
518 delete fTreeD; fTreeD = 0;
519 delete fTreeR; fTreeR = 0;
520 delete fTreeE; fTreeE = 0;
522 // Write AliRun info and all detectors parameters
530 //_____________________________________________________________________________
531 void AliRun::FlagTrack(Int_t track)
534 // Flags a track and all its family tree to be kept
541 particle=(TParticle*)fParticles->UncheckedAt(curr);
543 // If the particle is flagged the three from here upward is saved already
544 if(particle->TestBit(Keep_Bit)) return;
546 // Save this particle
547 particle->SetBit(Keep_Bit);
549 // Move to father if any
550 if((curr=particle->GetFirstMother())==-1) return;
554 //_____________________________________________________________________________
555 void AliRun::EnergySummary()
558 // Print summary of deposited energy
564 Int_t kn, i, left, j, id;
565 const Float_t zero=0;
566 Int_t ievent=fHeader.GetEvent()+1;
568 // Energy loss information
570 printf("***************** Energy Loss Information per event (GEV) *****************\n");
571 for(kn=1;kn<sEventEnergy.GetSize();kn++) {
574 sEventEnergy[ndep]=kn;
579 ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
582 sSummEnergy[ndep]=ed;
583 sSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
588 for(kn=0;kn<(ndep-1)/3+1;kn++) {
590 for(i=0;i<(3<left?3:left);i++) {
592 id=Int_t (sEventEnergy[j]+0.1);
593 printf(" %s %10.3f +- %10.3f%%;",gMC->VolName(id),sSummEnergy[j],sSum2Energy[j]);
598 // Relative energy loss in different detectors
599 printf("******************** Relative Energy Loss per event ********************\n");
600 printf("Total energy loss per event %10.3f GeV\n",edtot);
601 for(kn=0;kn<(ndep-1)/5+1;kn++) {
603 for(i=0;i<(5<left?5:left);i++) {
605 id=Int_t (sEventEnergy[j]+0.1);
606 printf(" %s %10.3f%%;",gMC->VolName(id),100*sSummEnergy[j]/edtot);
610 for(kn=0;kn<75;kn++) printf("*");
614 // Reset the TArray's
620 //_____________________________________________________________________________
621 AliModule *AliRun::GetModule(const char *name)
624 // Return pointer to detector from name
626 return (AliModule*)fModules->FindObject(name);
629 //_____________________________________________________________________________
630 AliDetector *AliRun::GetDetector(const char *name)
633 // Return pointer to detector from name
635 return (AliDetector*)fModules->FindObject(name);
638 //_____________________________________________________________________________
639 Int_t AliRun::GetModuleID(const char *name)
642 // Return galice internal detector identifier from name
645 TObject *mod=fModules->FindObject(name);
646 if(mod) i=fModules->IndexOf(mod);
650 //_____________________________________________________________________________
651 Int_t AliRun::GetEvent(Int_t event)
654 // Connect the Trees Kinematics and Hits for event # event
655 // Set branch addresses
658 // Reset existing structures
663 // Delete Trees already connected
664 if (fTreeK) delete fTreeK;
665 if (fTreeH) delete fTreeH;
666 if (fTreeD) delete fTreeD;
667 if (fTreeR) delete fTreeR;
669 // Get header from file
670 if(fTreeE) fTreeE->GetEntry(event);
671 else Error("GetEvent","Cannot file Header Tree\n");
673 // Get Kine Tree from file
675 sprintf(treeName,"TreeK%d",event);
676 fTreeK = (TTree*)gDirectory->Get(treeName);
677 if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
678 else Error("GetEvent","cannot find Kine Tree for event:%d\n",event);
680 // Get Hits Tree header from file
681 sprintf(treeName,"TreeH%d",event);
682 fTreeH = (TTree*)gDirectory->Get(treeName);
684 Error("GetEvent","cannot find Hits Tree for event:%d\n",event);
687 // Get Digits Tree header from file
688 sprintf(treeName,"TreeD%d",event);
689 fTreeD = (TTree*)gDirectory->Get(treeName);
691 Warning("GetEvent","cannot find Digits Tree for event:%d\n",event);
695 // Get Reconstruct Tree header from file
696 sprintf(treeName,"TreeR%d",event);
697 fTreeR = (TTree*)gDirectory->Get(treeName);
699 // printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
702 // Set Trees branch addresses
703 TIter next(fModules);
705 while((detector = (AliModule*)next())) {
706 detector->SetTreeAddress();
709 if (fTreeK) fTreeK->GetEvent(0);
710 fNtrack = Int_t (fParticles->GetEntries());
714 //_____________________________________________________________________________
715 TGeometry *AliRun::GetGeometry()
718 // Import Alice geometry from current file
719 // Return pointer to geometry object
721 if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
723 // Unlink and relink nodes in detectors
724 // This is bad and there must be a better way...
727 TIter next(fModules);
729 while((detector = (AliModule*)next())) {
730 detector->SetTreeAddress();
731 TList *dnodes=detector->Nodes();
734 for ( j=0; j<dnodes->GetSize(); j++) {
735 node = (TNode*) dnodes->At(j);
736 node1 = fGeometry->GetNode(node->GetName());
737 dnodes->Remove(node);
738 dnodes->AddAt(node1,j);
744 //_____________________________________________________________________________
745 void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
746 Float_t &e, Float_t *vpos, Float_t *polar,
750 // Return next track from stack of particles
755 for(Int_t i=fNtrack-1; i>=0; i--) {
756 track=(TParticle*) fParticles->UncheckedAt(i);
757 if(!track->TestBit(Done_Bit)) {
759 // The track has not yet been processed
761 ipart=track->GetPdgCode();
769 track->GetPolarisation(pol);
774 track->SetBit(Done_Bit);
780 // stop and start timer when we start a primary track
781 Int_t nprimaries = fHeader.GetNprimary();
782 if (fCurrent >= nprimaries) return;
783 if (fCurrent < nprimaries-1) {
785 track=(TParticle*) fParticles->UncheckedAt(fCurrent+1);
786 // track->SetProcessTime(fTimer.CpuTime());
791 //_____________________________________________________________________________
792 Int_t AliRun::GetPrimary(Int_t track)
795 // return number of primary that has generated track
803 part = (TParticle *)fParticles->UncheckedAt(current);
804 parent=part->GetFirstMother();
805 if(parent<0) return current;
809 //_____________________________________________________________________________
810 void AliRun::Init(const char *setup)
813 // Initialize the Alice setup
816 gROOT->LoadMacro(setup);
817 gInterpreter->ProcessLine("Config();");
819 gMC->DefineParticles(); //Create standard MC particles
821 TObject *objfirst, *objlast;
823 fNdets = fModules->GetLast()+1;
826 //=================Create Materials, geometry, histograms, etc
827 TIter next(fModules);
829 while((detector = (AliModule*)next())) {
830 detector->SetTreeAddress();
831 objlast = gDirectory->GetList()->Last();
833 // Initialise detector materials, geometry, histograms,etc
834 detector->CreateMaterials();
835 detector->CreateGeometry();
836 detector->BuildGeometry();
839 // Add Detector histograms in Detector list of histograms
840 if (objlast) objfirst = gDirectory->GetList()->After(objlast);
841 else objfirst = gDirectory->GetList()->First();
843 detector->Histograms()->Add(objfirst);
844 objfirst = gDirectory->GetList()->After(objfirst);
847 SetTransPar(); //Read the cuts for all materials
849 MediaTable(); //Build the special IMEDIA table
851 //Close the geometry structure
854 //Initialise geometry deposition table
855 sEventEnergy.Set(gMC->NofVolumes()+1);
856 sSummEnergy.Set(gMC->NofVolumes()+1);
857 sSum2Energy.Set(gMC->NofVolumes()+1);
859 //Create the color table
862 //Compute cross-sections
865 //Write Geometry object to current file.
871 //_____________________________________________________________________________
872 void AliRun::MediaTable()
875 // Built media table to get from the media number to
878 Int_t kz, nz, idt, lz, i, k, ind;
880 TObjArray &dets = *gAlice->Detectors();
884 for (kz=0;kz<fNdets;kz++) {
885 // If detector is defined
886 if((det=(AliModule*) dets[kz])) {
887 TArrayI &idtmed = *(det->GetIdtmed());
888 for(nz=0;nz<100;nz++) {
889 // Find max and min material number
890 if((idt=idtmed[nz])) {
891 det->LoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
892 det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
895 if(det->LoMedium() > det->HiMedium()) {
899 if(det->HiMedium() > fImedia->GetSize()) {
900 Error("MediaTable","Increase fImedia from %d to %d",
901 fImedia->GetSize(),det->HiMedium());
904 // Tag all materials in rage as belonging to detector kz
905 for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
912 // Print summary table
913 printf(" Traking media ranges:\n");
914 for(i=0;i<(fNdets-1)/6+1;i++) {
915 for(k=0;k< (6<fNdets-i*6?6:fNdets-i*6);k++) {
917 det=(AliModule*)dets[ind];
919 printf(" %6s: %3d -> %3d;",det->GetName(),det->LoMedium(),
922 printf(" %6s: %3d -> %3d;","NULL",0,0);
928 //____________________________________________________________________________
929 void AliRun::SetGenerator(AliGenerator *generator)
932 // Load the event generator
934 if(!fGenerator) fGenerator = generator;
937 //____________________________________________________________________________
938 void AliRun::SetTransPar(char* filename)
941 // Read filename to set the transport parameters
945 const Int_t ncuts=10;
946 const Int_t nflags=11;
947 const Int_t npars=ncuts+nflags;
948 const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
949 "BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
950 "BREM","COMP","DCAY","DRAY","HADR","LOSS",
951 "MULS","PAIR","PHOT","RAYL"};
957 Int_t i, itmed, iret, ktmed, kz;
960 // See whether the file is there
961 filtmp=gSystem->ExpandPathName(filename);
962 lun=fopen(filtmp,"r");
965 Warning("SetTransPar","File %s does not exist!\n",filename);
969 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
970 printf(" *%59s\n","*");
971 printf(" * Please check carefully what you are doing!%10s\n","*");
972 printf(" *%59s\n","*");
975 // Initialise cuts and flags
976 for(i=0;i<ncuts;i++) cut[i]=-99;
977 for(i=0;i<nflags;i++) flag[i]=-99;
979 for(i=0;i<256;i++) line[i]='\0';
980 // Read up to the end of line excluded
981 iret=fscanf(lun,"%[^\n]",line);
985 printf(" *%59s\n","*");
986 printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
989 // Read the end of line
992 if(line[0]=='*') continue;
994 iret=sscanf(line,"%s %d %f %f %f %f %f %f %f %f %f %f %d %d %d %d %d %d %d %d %d %d %d",
995 detName,&itmed,&cut[0],&cut[1],&cut[2],&cut[3],&cut[4],&cut[5],&cut[6],&cut[7],&cut[8],
996 &cut[9],&flag[0],&flag[1],&flag[2],&flag[3],&flag[4],&flag[5],&flag[6],&flag[7],
997 &flag[8],&flag[9],&flag[10]);
1001 Warning("SetTransPar","Error reading file %s\n",filename);
1004 // Check that the module exist
1005 AliModule *mod = GetModule(detName);
1007 // Get the array of media numbers
1008 TArrayI &idtmed = *mod->GetIdtmed();
1009 // Check that the tracking medium code is valid
1010 if(0<=itmed && itmed < 100) {
1011 ktmed=idtmed[itmed];
1013 Warning("SetTransPar","Invalid tracking medium code %d for %s\n",itmed,mod->GetName());
1016 // Set energy thresholds
1017 for(kz=0;kz<ncuts;kz++) {
1019 printf(" * %-6s set to %10.3E for tracking medium code %4d for %s\n",
1020 pars[kz],cut[kz],itmed,mod->GetName());
1021 gMC->Gstpar(ktmed,pars[kz],cut[kz]);
1024 // Set transport mechanisms
1025 for(kz=0;kz<nflags;kz++) {
1027 printf(" * %-6s set to %10d for tracking medium code %4d for %s\n",
1028 pars[ncuts+kz],flag[kz],itmed,mod->GetName());
1029 gMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
1033 Warning("SetTransPar","Invalid medium code %d *\n",itmed);
1037 Warning("SetTransPar","Module %s not present\n",detName);
1043 //_____________________________________________________________________________
1044 void AliRun::MakeTree(Option_t *option)
1047 // Create the ROOT trees
1048 // Loop on all detectors to create the Root branch (if any)
1053 char *K = strstr(option,"K");
1054 char *H = strstr(option,"H");
1055 char *E = strstr(option,"E");
1056 char *D = strstr(option,"D");
1057 char *R = strstr(option,"R");
1059 if (K && !fTreeK) fTreeK = new TTree("TreeK0","Kinematics");
1060 if (H && !fTreeH) fTreeH = new TTree("TreeH0","Hits");
1061 if (D && !fTreeD) fTreeD = new TTree("TreeD0","Digits");
1062 if (E && !fTreeE) fTreeE = new TTree("TE","Header");
1063 if (R && !fTreeR) fTreeR = new TTree("TreeR0","Reconstruction");
1064 if (fTreeH) fTreeH->SetAutoSave(1000000000); //no autosave
1066 // Create a branch for hits/digits for each detector
1067 // Each branch is a TClonesArray. Each data member of the Hits classes
1068 // will be in turn a subbranch of the detector master branch
1069 TIter next(fModules);
1070 AliModule *detector;
1071 while((detector = (AliModule*)next())) {
1072 if (H || D || R) detector->MakeBranch(option);
1074 // Create a branch for particles
1075 if (fTreeK && K) fTreeK->Branch("Particles",&fParticles,4000);
1077 // Create a branch for Header
1078 if (fTreeE && E) fTreeE->Branch("Header","AliHeader",&header,4000);
1081 //_____________________________________________________________________________
1082 Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
1085 // PurifyKine with external parameters
1087 fHgwmk = lastSavedTrack;
1088 fNtrack = nofTracks;
1093 //_____________________________________________________________________________
1094 void AliRun::PurifyKine()
1097 // Compress kinematic tree keeping only flagged particles
1098 // and renaming the particle id's in all the hits
1100 TClonesArray &particles = *fParticles;
1101 int nkeep=fHgwmk+1, parent, i;
1102 TParticle *part, *partnew, *father;
1104 int *map = new int[particles.GetEntries()];
1106 // Save in Header total number of tracks before compression
1107 fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
1109 // Preset map, to be removed later
1110 for(i=0; i<fNtrack; i++) {
1111 if(i<=fHgwmk) map[i]=i ; else map[i] = -99 ;}
1112 // Second pass, build map between old and new numbering
1113 for(i=fHgwmk+1; i<fNtrack; i++) {
1114 part = (TParticle *)particles.UncheckedAt(i);
1115 if(part->TestBit(Keep_Bit)) {
1117 // This particle has to be kept
1121 // Old and new are different, have to copy
1122 partnew = (TParticle *)particles.UncheckedAt(nkeep);
1124 } else partnew = part;
1126 // as the parent is always *before*, it must be already
1127 // in place. This is what we are checking anyway!
1128 if((parent=partnew->GetFirstMother())>fHgwmk) {
1129 if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
1130 partnew->SetFirstMother(map[parent]);
1137 // Fix daughters information
1138 for (i=0; i<fNtrack; i++) {
1139 part = (TParticle *)particles.UncheckedAt(i);
1140 parent = part->GetFirstMother();
1142 father = (TParticle *)particles.UncheckedAt(parent);
1143 if(father->TestBit(Daughters_Bit)) {
1145 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
1146 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
1148 // Iitialise daughters info for first pass
1149 father->SetFirstDaughter(i);
1150 father->SetLastDaughter(i);
1151 father->SetBit(Daughters_Bit);
1156 // Now loop on all detectors and reset the hits
1157 TIter next(fModules);
1158 AliModule *detector;
1159 while((detector = (AliModule*)next())) {
1160 if (!detector->Hits()) continue;
1161 TClonesArray &vHits=*(detector->Hits());
1162 if(vHits.GetEntries() != detector->GetNhits())
1163 printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
1164 vHits.GetEntries(),detector->GetNhits());
1165 for (i=0; i<detector->GetNhits(); i++) {
1166 OneHit = (AliHit *)vHits.UncheckedAt(i);
1167 OneHit->SetTrack(map[OneHit->GetTrack()]);
1172 particles.SetLast(fHgwmk);
1176 //_____________________________________________________________________________
1177 void AliRun::Reset(Int_t run, Int_t idevent)
1180 // Reset all Detectors & kinematics & trees
1188 // Initialise event header
1189 fHeader.Reset(run,idevent);
1193 sprintf(hname,"TreeK%d",idevent);
1194 fTreeK->SetName(hname);
1198 sprintf(hname,"TreeH%d",idevent);
1199 fTreeH->SetName(hname);
1203 sprintf(hname,"TreeD%d",idevent);
1204 fTreeD->SetName(hname);
1208 sprintf(hname,"TreeR%d",idevent);
1209 fTreeR->SetName(hname);
1213 //_____________________________________________________________________________
1214 void AliRun::ResetDigits()
1217 // Reset all Detectors digits
1219 TIter next(fModules);
1220 AliModule *detector;
1221 while((detector = (AliModule*)next())) {
1222 detector->ResetDigits();
1226 //_____________________________________________________________________________
1227 void AliRun::ResetHits()
1230 // Reset all Detectors hits
1232 TIter next(fModules);
1233 AliModule *detector;
1234 while((detector = (AliModule*)next())) {
1235 detector->ResetHits();
1239 //_____________________________________________________________________________
1240 void AliRun::ResetPoints()
1243 // Reset all Detectors points
1245 TIter next(fModules);
1246 AliModule *detector;
1247 while((detector = (AliModule*)next())) {
1248 detector->ResetPoints();
1252 //_____________________________________________________________________________
1253 void AliRun::Run(Int_t nevent, const char *setup)
1256 // Main function to be called to process a galice run
1258 // Root > gAlice.Run();
1259 // a positive number of events will cause the finish routine
1264 // check if initialisation has been done
1265 if (!fInitDone) Init(setup);
1267 // Create the Root Tree with one branch per detector
1269 gAlice->MakeTree("KHDER");
1272 todo = TMath::Abs(nevent);
1273 for (i=0; i<todo; i++) {
1274 // Process one run (one run = one event)
1275 gAlice->Reset(fRun, fEvent);
1279 gAlice->FinishEvent();
1283 // End of this run, close files
1284 if(nevent>0) gAlice->FinishRun();
1287 //_____________________________________________________________________________
1288 void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
1289 Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
1290 Float_t rmin,Float_t rmax,Float_t zmax)
1293 // Generates lego plots of:
1294 // - radiation length map phi vs theta
1295 // - radiation length map phi vs eta
1296 // - interaction length map
1297 // - g/cm2 length map
1299 // ntheta bins in theta, eta
1300 // themin minimum angle in theta (degrees)
1301 // themax maximum angle in theta (degrees)
1303 // phimin minimum angle in phi (degrees)
1304 // phimax maximum angle in phi (degrees)
1305 // rmin minimum radius
1306 // rmax maximum radius
1309 // The number of events generated = ntheta*nphi
1310 // run input parameters in macro setup (default="Config.C")
1312 // Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
1315 <img src="picts/AliRunLego1.gif">
1320 <img src="picts/AliRunLego2.gif">
1325 <img src="picts/AliRunLego3.gif">
1330 // check if initialisation has been done
1331 if (!fInitDone) Init(setup);
1333 fLego = new AliLego("lego","lego");
1334 fLego->Init(ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
1337 // Create only the Root event Tree
1338 gAlice->MakeTree("E");
1340 // End of this run, close files
1341 gAlice->FinishRun();
1344 //_____________________________________________________________________________
1345 void AliRun::SetCurrentTrack(Int_t track)
1348 // Set current track number
1353 //_____________________________________________________________________________
1354 void AliRun::SetTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
1355 Float_t *vpos, Float_t *polar, Float_t tof,
1356 const char *mecha, Int_t &ntr, Float_t weight)
1359 // Load a track on the stack
1361 // done 0 if the track has to be transported
1363 // parent identifier of the parent track. -1 for a primary
1364 // pdg particle code
1365 // pmom momentum GeV/c
1367 // polar polarisation
1368 // tof time of flight in seconds
1369 // mecha production mechanism
1370 // ntr on output the number of the track stored
1372 TClonesArray &particles = *fParticles;
1373 TParticle *particle;
1375 const Int_t firstdaughter=-1;
1376 const Int_t lastdaughter=-1;
1378 // const Float_t tlife=0;
1381 // Here we get the static mass
1382 // For MC is ok, but a more sophisticated method could be necessary
1383 // if the calculated mass is required
1384 // also, this method is potentially dangerous if the mass
1385 // used in the MC is not the same of the PDG database
1387 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
1388 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
1389 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
1391 //printf("Loading particle %s mass %f ene %f No %d ip %d pos %f %f %f mom %f %f %f KS %d m %s\n",
1392 //pname,mass,e,fNtrack,pdg,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
1394 particle=new(particles[fNtrack]) TParticle(pdg,KS,parent,-1,firstdaughter,
1395 lastdaughter,pmom[0],pmom[1],pmom[2],
1396 e,vpos[0],vpos[1],vpos[2],tof);
1397 // polar[0],polar[1],polar[2],tof,
1399 ((TParticle*)particles[fNtrack])->SetPolarisation(TVector3(polar[0],polar[1],polar[2]));
1400 ((TParticle*)particles[fNtrack])->SetWeight(weight);
1401 if(!done) particle->SetBit(Done_Bit);
1404 particle=(TParticle*) fParticles->UncheckedAt(parent);
1405 particle->SetLastDaughter(fNtrack);
1406 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
1409 // This is a primary track. Set high water mark for this event
1412 // Set also number if primary tracks
1413 fHeader.SetNprimary(fHgwmk+1);
1414 fHeader.SetNtrack(fHgwmk+1);
1419 //_____________________________________________________________________________
1420 void AliRun::KeepTrack(const Int_t track)
1423 // flags a track to be kept
1425 TClonesArray &particles = *fParticles;
1426 ((TParticle*)particles[track])->SetBit(Keep_Bit);
1429 //_____________________________________________________________________________
1430 void AliRun::StepManager(Int_t id) const
1433 // Called at every step during transport
1438 // --- If lego option, do it and leave
1440 fLego->StepManager();
1443 //Update energy deposition tables
1444 sEventEnergy[gMC->CurrentVolID(copy)]+=gMC->Edep();
1446 //Call the appropriate stepping routine;
1447 AliModule *det = (AliModule*)fModules->At(id);
1448 if(det) det->StepManager();
1451 //_____________________________________________________________________________
1452 void AliRun::ReadEuclid(const char* filnam, const AliModule *det, char* topvol)
1455 // read in the geometry of the detector in euclid file format
1457 // id_det : the detector identification (2=its,...)
1458 // topvol : return parameter describing the name of the top
1459 // volume of geometry.
1461 // author : m. maire
1464 // several changes have been made by miroslav helbich
1465 // subroutine is rewrited to follow the new established way of memory
1466 // booking for tracking medias and rotation matrices.
1467 // all used tracking media have to be defined first, for this you can use
1468 // subroutine greutmed.
1469 // top volume is searched as only volume not positioned into another
1472 Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
1473 Int_t ndvmx, nr, flag;
1474 char key[5], card[77], natmed[21];
1475 char name[5], mother[5], shape[5], konly[5], volst[7000][5];
1478 Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
1480 const Int_t maxrot=5000;
1481 Int_t idrot[maxrot],istop[7000];
1484 // *** The input filnam name will be with extension '.euc'
1485 filtmp=gSystem->ExpandPathName(filnam);
1486 lun=fopen(filtmp,"r");
1489 Error("ReadEuclid","Could not open file %s\n",filnam);
1492 //* --- definition of rotation matrix 0 ---
1493 TArrayI &idtmed = *(det->GetIdtmed());
1494 for(i=1; i<maxrot; ++i) idrot[i]=-99;
1498 for(i=0;i<77;i++) card[i]=0;
1499 iret=fscanf(lun,"%77[^\n]",card);
1500 if(iret<=0) goto L20;
1503 strncpy(key,card,4);
1505 if (!strcmp(key,"TMED")) {
1506 sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
1507 if( itmed<0 || itmed>=100 ) {
1508 Error("ReadEuclid","TMED illegal medium number %d for %s\n",itmed,natmed);
1511 //Pad the string with blanks
1514 while(i<20) natmed[i++]=' ';
1517 if( idtmed[itmed]<=0 ) {
1518 Error("ReadEuclid","TMED undefined medium number %d for %s\n",itmed,natmed);
1521 gMC->Gckmat(idtmed[itmed],natmed);
1523 } else if (!strcmp(key,"ROTM")) {
1524 sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
1525 if( irot<=0 || irot>=maxrot ) {
1526 Error("ReadEuclid","ROTM rotation matrix number %d illegal\n",irot);
1529 det->AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
1531 } else if (!strcmp(key,"VOLU")) {
1532 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
1534 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
1537 gMC->Gsvolu( name, shape, idtmed[numed], par, npar);
1538 //* save the defined volumes
1539 strcpy(volst[++nvol],name);
1542 } else if (!strcmp(key,"DIVN")) {
1543 sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
1544 gMC->Gsdvn ( name, mother, ndiv, iaxe );
1546 } else if (!strcmp(key,"DVN2")) {
1547 sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
1548 gMC->Gsdvn2( name, mother, ndiv, iaxe, orig,idtmed[numed]);
1550 } else if (!strcmp(key,"DIVT")) {
1551 sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
1552 gMC->Gsdvt ( name, mother, step, iaxe, idtmed[numed], ndvmx);
1554 } else if (!strcmp(key,"DVT2")) {
1555 sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
1556 gMC->Gsdvt2 ( name, mother, step, iaxe, orig, idtmed[numed], ndvmx );
1558 } else if (!strcmp(key,"POSI")) {
1559 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
1560 if( irot<0 || irot>=maxrot ) {
1561 Error("ReadEuclid","POSI %s#%d rotation matrix number %d illegal\n",name,nr,irot);
1564 if( idrot[irot] == -99) {
1565 Error("ReadEuclid","POSI %s#%d undefined matrix number %d\n",name,nr,irot);
1568 //*** volume name cannot be the top volume
1569 for(i=1;i<=nvol;i++) {
1570 if (!strcmp(volst[i],name)) istop[i]=0;
1573 gMC->Gspos ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
1575 } else if (!strcmp(key,"POSP")) {
1576 sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
1577 if( irot<0 || irot>=maxrot ) {
1578 Error("ReadEuclid","POSP %s#%d rotation matrix number %d illegal\n",name,nr,irot);
1581 if( idrot[irot] == -99) {
1582 Error("ReadEuclid","POSP %s#%d undefined matrix number %d\n",name,nr,irot);
1586 for(i=0;i<npar;i++) fscanf(lun,"%f",&par[i]);
1589 //*** volume name cannot be the top volume
1590 for(i=1;i<=nvol;i++) {
1591 if (!strcmp(volst[i],name)) istop[i]=0;
1594 gMC->Gsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
1597 if (strcmp(key,"END")) goto L10;
1598 //* find top volume in the geometry
1600 for(i=1;i<=nvol;i++) {
1601 if (istop[i] && flag) {
1602 Warning("ReadEuclid"," %s is another possible top volume\n",volst[i]);
1604 if (istop[i] && !flag) {
1605 strcpy(topvol,volst[i]);
1606 printf(" *** GREUCL *** volume %s taken as a top volume\n",topvol);
1611 Warning("ReadEuclid","top volume not found\n");
1615 //* commented out only for the not cernlib version
1616 printf(" *** GREUCL *** file: %s is now read in\n",filnam);
1621 Error("ReadEuclid","reading error or premature end of file\n");
1624 //_____________________________________________________________________________
1625 void AliRun::ReadEuclidMedia(const char* filnam, const AliModule *det)
1628 // read in the materials and tracking media for the detector
1629 // in euclid file format
1631 // filnam: name of the input file
1632 // id_det: id_det is the detector identification (2=its,...)
1634 // author : miroslav helbich
1636 Float_t sxmgmx = gAlice->Field()->Max();
1637 Int_t isxfld = gAlice->Field()->Integ();
1638 Int_t end, i, iret, itmed;
1639 char key[5], card[130], natmed[21], namate[21];
1644 Int_t nwbuf, isvol, ifield, nmat;
1645 Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
1648 for(i=0;i<end;i++) if(filnam[i]=='.') {
1653 // *** The input filnam name will be with extension '.euc'
1654 printf("The file name is %s\n",filnam); //Debug
1655 filtmp=gSystem->ExpandPathName(filnam);
1656 lun=fopen(filtmp,"r");
1659 Warning("ReadEuclidMedia","Could not open file %s\n",filnam);
1663 // Retrieve Mag Field parameters
1664 Int_t ISXFLD=gAlice->Field()->Integ();
1665 Float_t SXMGMX=gAlice->Field()->Max();
1666 // TArrayI &idtmed = *(det->GetIdtmed());
1669 for(i=0;i<130;i++) card[i]=0;
1670 iret=fscanf(lun,"%4s %[^\n]",key,card);
1671 if(iret<=0) goto L20;
1675 if (!strcmp(key,"MATE")) {
1676 sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
1677 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
1678 //Pad the string with blanks
1681 while(i<20) namate[i++]=' ';
1684 det->AliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
1685 //* read tracking medium
1686 } else if (!strcmp(key,"TMED")) {
1687 sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
1688 &itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
1689 &stemax,&deemax,&epsil,&stmin,&nwbuf);
1690 if (nwbuf>0) for(i=0;i<nwbuf;i++) fscanf(lun,"%f",&ubuf[i]);
1691 if (ifield<0) ifield=isxfld;
1692 if (fieldm<0) fieldm=sxmgmx;
1693 //Pad the string with blanks
1696 while(i<20) natmed[i++]=' ';
1699 det->AliMedium(itmed,natmed,nmat,isvol,ISXFLD,SXMGMX,tmaxfd,
1700 stemax,deemax,epsil,stmin,ubuf,nwbuf);
1701 // (*fImedia)[idtmed[itmed]-1]=id_det;
1705 if (strcmp(key,"END")) goto L10;
1708 //* commented out only for the not cernlib version
1709 Warning("ReadEuclidMedia","file: %s is now read in\n",filnam);
1714 Warning("ReadEuclidMedia","reading error or premature end of file\n");
1717 //_____________________________________________________________________________
1718 void AliRun::Streamer(TBuffer &R__b)
1721 // Stream an object of class AliRun.
1723 if (R__b.IsReading()) {
1724 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1725 TNamed::Streamer(R__b);
1726 if (!gAlice) gAlice = this;
1727 gROOT->GetListOfBrowsables()->Add(this,"Run");
1728 fTreeE = (TTree*)gDirectory->Get("TE");
1729 if (fTreeE) fTreeE->SetBranchAddress("Header", &header);
1730 else Error("Streamer","cannot find Header Tree\n");
1734 fHeader.Streamer(R__b);
1744 R__b >> fPDGDB; //Particle factory object!
1745 fTreeE->GetEntry(0);
1747 fHeader.SetEvent(0);
1748 fPDGDB = TDatabasePDG::Instance(); //Particle factory object!
1751 R__b.WriteVersion(AliRun::IsA());
1752 TNamed::Streamer(R__b);
1756 fHeader.Streamer(R__b);
1765 R__b << fPDGDB; //Particle factory object!
1770 //_____________________________________________________________________________
1772 // Interfaces to Fortran
1774 //_____________________________________________________________________________
1776 extern "C" void type_of_call rxgtrak (Int_t &mtrack, Int_t &ipart, Float_t *pmom,
1777 Float_t &e, Float_t *vpos, Float_t *polar,
1781 // Fetches next track from the ROOT stack for transport. Called by the
1782 // modified version of GTREVE.
1784 // Track number in the ROOT stack. If MTRACK=0 no
1785 // mtrack more tracks are left in the stack to be
1787 // ipart Particle code in the GEANT conventions.
1788 // pmom[3] Particle momentum in GeV/c
1789 // e Particle energy in GeV
1790 // vpos[3] Particle position
1791 // tof Particle time of flight in seconds
1794 gAlice->GetNextTrack(mtrack, pdg, pmom, e, vpos, polar, tof);
1795 ipart = gMC->IdFromPDG(pdg);
1799 //_____________________________________________________________________________
1800 extern "C" void type_of_call
1802 rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
1803 Float_t *vpos, Float_t &tof, const char* cmech, Int_t &ntr, const int cmlen)
1805 rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
1806 Float_t *vpos, Float_t &tof, const char* cmech, const int cmlen,
1811 // Fetches next track from the ROOT stack for transport. Called by GUKINE
1814 // Status of the track. If keep=0 the track is put
1815 // keep on the ROOT stack but it is not fetched for
1817 // parent Parent track. If parent=0 the track is a primary.
1818 // In GUSTEP the routine is normally called to store
1819 // secondaries generated by the current track whose
1820 // ROOT stack number is MTRACK (common SCKINE.
1821 // ipart Particle code in the GEANT conventions.
1822 // pmom[3] Particle momentum in GeV/c
1823 // vpos[3] Particle position
1824 // tof Particle time of flight in seconds
1826 // cmech (CHARACTER*10) Particle origin. This field is user
1827 // defined and it is not used inside the GALICE code.
1828 // ntr Number assigned to the particle in the ROOT stack.
1831 Float_t polar[3]={0.,0.,0.};
1832 for(int i=0; i<10 && i<cmlen; i++) mecha[i]=cmech[i];
1834 Int_t pdg=gMC->PDGFromId(ipart);
1835 gAlice->SetTrack(keep, parent-1, pdg, pmom, vpos, polar, tof, mecha, ntr);
1839 //_____________________________________________________________________________
1840 extern "C" void type_of_call rxkeep(const Int_t &n)
1842 if( NULL==gAlice ) exit(1);
1844 if( n<=0 || n>gAlice->Particles()->GetEntries() )
1846 printf(" Bad index n=%d must be 0<n<=%d\n",
1847 n,gAlice->Particles()->GetEntries());
1851 ((TParticle*)(gAlice->Particles()->UncheckedAt(n-1)))->SetBit(Keep_Bit);
1854 //_____________________________________________________________________________
1855 extern "C" void type_of_call rxouth ()
1858 // Called by Gtreve at the end of each primary track
1860 gAlice->FinishPrimary();