///////////////////////////////////////////////////////////////////////////////
// //
// Control class for Alice C++ //
// Only one single instance of this class exists. //
// The object is created in main program aliroot //
// and is pointed by the global gAlice. //
// //
// -Supports the list of all Alice Detectors (fModules). //
// -Supports the list of particles (fParticles). //
// -Supports the Trees. //
// -Supports the geometry. //
// -Supports the event display. //
//Begin_Html
/*
*/
//End_Html
//Begin_Html
/*
*/
//End_Html
// //
///////////////////////////////////////////////////////////////////////////////
#include
#include
#include
#include
#include
#include
#include
#include "GParticle.h"
#include "AliRun.h"
#include "AliModule.h"
#include "AliDisplay.h"
#include "AliCallf77.h"
#include
#include
#include
AliRun *gAlice;
static AliHeader *header;
#ifndef WIN32
# define rxgtrak rxgtrak_
# define rxstrak rxstrak_
# define rxkeep rxkeep_
# define rxouth rxouth_
# define sxpart sxpart_
#else
# define rxgtrak RXGTRAK
# define rxstrak RXSTRAK
# define rxkeep RXKEEP
# define rxouth RXOUTH
# define sxpart SXPART
#endif
static TArrayF sEventEnergy;
static TArrayF sSummEnergy;
static TArrayF sSum2Energy;
extern "C" void type_of_call sxpart();
ClassImp(AliRun)
//_____________________________________________________________________________
AliRun::AliRun()
{
//
// Default constructor for AliRun
//
header=&fHeader;
fRun = 0;
fEvent = 0;
fCurrent = -1;
fModules = 0;
fGenerator = 0;
fTreeD = 0;
fTreeK = 0;
fTreeH = 0;
fTreeE = 0;
fTreeR = 0;
fParticles = 0;
fGeometry = 0;
fDisplay = 0;
fField = 0;
fMC = 0;
fNdets = 0;
fImedia = 0;
fTrRmax = 1.e10;
fTrZmax = 1.e10;
fIdtmed = 0;
fInitDone = kFALSE;
fLego = 0;
}
//_____________________________________________________________________________
AliRun::AliRun(const char *name, const char *title)
: TNamed(name,title)
{
//
// Constructor for the main processor.
// Creates the geometry
// Creates the list of Detectors.
// Creates the list of particles.
//
Int_t i;
gAlice = this;
fTreeD = 0;
fTreeK = 0;
fTreeH = 0;
fTreeE = 0;
fTreeR = 0;
fTrRmax = 1.e10;
fTrZmax = 1.e10;
fGenerator = 0;
fInitDone = kFALSE;
fLego = 0;
fField = 0;
gROOT->GetListOfBrowsables()->Add(this,name);
//
// create the support list for the various Detectors
fModules = new TObjArray(77);
//
// Create the TNode geometry for the event display
BuildSimpleGeometry();
fNtrack=0;
fHgwmk=0;
fCurrent=-1;
header=&fHeader;
fRun = 0;
fEvent = 0;
//
// Create the particle stack
fParticles = new TClonesArray("GParticle",100);
fDisplay = 0;
//
// Create default mag field
SetField();
//
fMC = AliMC::GetMC();
//
//---------------Load detector names
fNdets=21;
strcpy(fDnames[0],"BODY");
strcpy(fDnames[1],"NULL");
strcpy(fDnames[2],"ITS");
strcpy(fDnames[3],"MAG");
strcpy(fDnames[4],"TPC");
strcpy(fDnames[5],"TOF");
strcpy(fDnames[6],"PMD");
strcpy(fDnames[7],"PHOS");
strcpy(fDnames[8],"ZDC");
strcpy(fDnames[9],"FMD");
strcpy(fDnames[10],"RICH");
strcpy(fDnames[11],"MUON");
strcpy(fDnames[12],"FRAME");
strcpy(fDnames[13],"TRD");
strcpy(fDnames[14],"NULL");
strcpy(fDnames[15],"CASTOR");
strcpy(fDnames[16],"ABSO");
strcpy(fDnames[17],"SHIL");
strcpy(fDnames[18],"DIPO");
strcpy(fDnames[19],"HALL");
strcpy(fDnames[20],"PIPE");
//
// Prepare the tracking medium lists
fImedia = new TArrayI(1000);
for(i=0;i<1000;i++) (*fImedia)[i]=-99;
fIdtmed = new Int_t[fNdets*100];
for(i=0;iDelete();
delete fModules;
}
if (fParticles) {
fParticles->Delete();
delete fParticles;
}
}
//_____________________________________________________________________________
void AliRun::AddHit(Int_t id, Int_t track, Int_t *vol, Float_t *hits) const
{
//
// Add a hit to detector id
//
TObjArray &dets = *fModules;
if(dets[id]) ((AliModule*) dets[id])->AddHit(track,vol,hits);
}
//_____________________________________________________________________________
void AliRun::AddDigit(Int_t id, Int_t *tracks, Int_t *digits) const
{
//
// Add digit to detector id
//
TObjArray &dets = *fModules;
if(dets[id]) ((AliModule*) dets[id])->AddDigit(tracks,digits);
}
//_____________________________________________________________________________
void AliRun::Browse(TBrowser *b)
{
//
// Called when the item "Run" is clicked on the left pane
// of the Root browser.
// It displays the Root Trees and all detectors.
//
if (fTreeK) b->Add(fTreeK,fTreeK->GetName());
if (fTreeH) b->Add(fTreeH,fTreeH->GetName());
if (fTreeD) b->Add(fTreeD,fTreeD->GetName());
if (fTreeE) b->Add(fTreeE,fTreeE->GetName());
if (fTreeR) b->Add(fTreeR,fTreeR->GetName());
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
b->Add(detector,detector->GetName());
}
}
//_____________________________________________________________________________
void AliRun::Build()
{
//
// Initialize Alice geometry
// Dummy routine
//
}
//_____________________________________________________________________________
void AliRun::BuildSimpleGeometry()
{
//
// Create a simple TNode geometry used by Root display engine
//
// Initialise geometry
//
fGeometry = new TGeometry("AliceGeom","Galice Geometry for Hits");
new TMaterial("void","Vacuum",0,0,0); //Everything is void
TBRIK *brik = new TBRIK("S_alice","alice volume","void",2000,2000,3000);
brik->SetVisibility(0);
new TNode("alice","alice","S_alice");
}
//_____________________________________________________________________________
void AliRun::CleanDetectors()
{
//
// Clean Detectors at the end of event
//
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->FinishEvent();
}
}
//_____________________________________________________________________________
void AliRun::CleanParents()
{
//
// Clean Particles stack.
// Set parent/child relations
//
TClonesArray &particles = *(gAlice->Particles());
GParticle *part;
int i;
for(i=0; iTestBit(Children_Bit)) {
part->SetFirstChild(-1);
part->SetLastChild(-1);
}
}
}
//_____________________________________________________________________________
Int_t AliRun::DistancetoPrimitive(Int_t, Int_t)
{
//
// Return the distance from the mouse to the AliRun object
// Dummy routine
//
return 9999;
}
//_____________________________________________________________________________
void AliRun::DumpPart (Int_t i)
{
//
// Dumps particle i in the stack
//
TClonesArray &particles = *fParticles;
((GParticle*) particles[i])->Dump();
}
//_____________________________________________________________________________
void AliRun::DumpPStack ()
{
//
// Dumps the particle stack
//
TClonesArray &particles = *fParticles;
printf(
"\n\n=======================================================================\n");
for (Int_t i=0;i %d ",i); ((GParticle*) particles[i])->Dump();
printf("--------------------------------------------------------------\n");
}
printf(
"\n=======================================================================\n\n");
}
//_____________________________________________________________________________
void AliRun::SetField(Int_t type, Int_t version, Float_t scale,
Float_t maxField, char* filename)
{
//
// Set magnetic field parameters
// type Magnetic field transport flag 0=no field, 2=helix, 3=Runge Kutta
// version Magnetic field map version (only 1 active now)
// scale Scale factor for the magnetic field
// maxField Maximum value for the magnetic field
//
// --- Sanity check on mag field flags
if(type<0 || type > 2) {
printf(" Invalid magnetic field flag: %5d; Helix tracking chosen instead\n"
,type);
type=2;
}
if(fField) delete fField;
if(version==1) {
fField = new AliMagFC("Map1"," ",type,version,scale,maxField);
} else if(version<=3) {
fField = new AliMagFCM("Map2-3",filename,type,version,scale,maxField);
fField->ReadField();
} else {
printf("Invalid map %d\n",version);
}
}
//_____________________________________________________________________________
void AliRun::FillTree()
{
//
// Fills all AliRun TTrees
//
if (fTreeK) fTreeK->Fill();
if (fTreeH) fTreeH->Fill();
if (fTreeD) fTreeD->Fill();
if (fTreeR) fTreeR->Fill();
}
//_____________________________________________________________________________
void AliRun::FinishPrimary()
{
//
// Called at the end of each primary track
//
// This primary is finished, purify stack
gAlice->PurifyKine();
// Write out hits if any
if (gAlice->TreeH()) {
gAlice->TreeH()->Fill();
}
// Reset Hits info
gAlice->ResetHits();
}
//_____________________________________________________________________________
void AliRun::FinishEvent()
{
//
// Called at the end of the event.
//
//Update the energy deposit tables
Int_t i;
for(i=0;iFill();
}
// Write out the digits
if (fTreeD) {
fTreeD->Fill();
ResetDigits();
}
// Write out reconstructed clusters
if (fTreeR) {
fTreeR->Fill();
}
// Write out the event Header information
if (fTreeE) fTreeE->Fill();
// Reset stack info
ResetStack();
// Write Tree headers
Int_t ievent = fHeader.GetEvent();
char hname[30];
sprintf(hname,"TreeK%d",ievent);
if (fTreeK) fTreeK->Write(hname);
sprintf(hname,"TreeH%d",ievent);
if (fTreeH) fTreeH->Write(hname);
sprintf(hname,"TreeD%d",ievent);
if (fTreeD) fTreeD->Write(hname);
sprintf(hname,"TreeR%d",ievent);
if (fTreeR) fTreeR->Write(hname);
}
//_____________________________________________________________________________
void AliRun::FinishRun()
{
//
// Called at the end of the run.
//
// Clean detector information
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->FinishRun();
}
//Output energy summary tables
EnergySummary();
// file is retrieved from whatever tree
TFile *File = 0;
if (fTreeK) File = fTreeK->GetCurrentFile();
if ((!File) && (fTreeH)) File = fTreeH->GetCurrentFile();
if ((!File) && (fTreeD)) File = fTreeD->GetCurrentFile();
if ((!File) && (fTreeE)) File = fTreeE->GetCurrentFile();
if( NULL==File ) {
Error("FinishRun","There isn't root file!");
exit(1);
}
File->cd();
fTreeE->Write();
// Clean tree information
delete fTreeK; fTreeK = 0;
delete fTreeH; fTreeH = 0;
delete fTreeD; fTreeD = 0;
delete fTreeR; fTreeR = 0;
delete fTreeE; fTreeE = 0;
// Write AliRun info and all detectors parameters
Write();
// Close output file
File->Write();
File->Close();
}
//_____________________________________________________________________________
void AliRun::FlagTrack(Int_t track)
{
//
// Flags a track and all its family tree to be kept
//
int curr;
GParticle *particle;
curr=track;
while(1) {
particle=(GParticle*)fParticles->UncheckedAt(curr);
// If the particle is flagged the three from here upward is saved already
if(particle->TestBit(Keep_Bit)) return;
// Save this particle
particle->SetBit(Keep_Bit);
// Move to father if any
if((curr=particle->GetParent())==-1) return;
}
}
//_____________________________________________________________________________
void AliRun::EnergySummary()
{
//
// Print summary of deposited energy
//
AliMC* pMC = AliMC::GetMC();
Int_t ndep=0;
Float_t edtot=0;
Float_t ed, ed2;
Int_t kn, i, left, j, id;
const Float_t zero=0;
Int_t ievent=fHeader.GetEvent()+1;
//
// Energy loss information
if(ievent) {
printf("***************** Energy Loss Information per event (GEV) *****************\n");
for(kn=1;kn0) {
sEventEnergy[ndep]=kn;
if(ievent>1) {
ed=ed/ievent;
ed2=sSum2Energy[kn];
ed2=ed2/ievent;
ed2=100*TMath::Sqrt(TMath::Max(ed2-ed*ed,zero))/ed;
} else
ed2=99;
sSummEnergy[ndep]=ed;
sSum2Energy[ndep]=TMath::Min((Float_t) 99.,TMath::Max(ed2,zero));
edtot+=ed;
ndep++;
}
}
for(kn=0;kn<(ndep-1)/3+1;kn++) {
left=ndep-kn*3;
for(i=0;i<(3VolName(id),sSummEnergy[j],sSum2Energy[j]);
}
printf("\n");
}
//
// Relative energy loss in different detectors
printf("******************** Relative Energy Loss per event ********************\n");
printf("Total energy loss per event %10.3f GeV\n",edtot);
for(kn=0;kn<(ndep-1)/5+1;kn++) {
left=ndep-kn*5;
for(i=0;i<(5VolName(id),100*sSummEnergy[j]/edtot);
}
printf("\n");
}
for(kn=0;kn<75;kn++) printf("*");
printf("\n");
}
//
// Reset the TArray's
sEventEnergy.Set(0);
sSummEnergy.Set(0);
sSum2Energy.Set(0);
}
//_____________________________________________________________________________
AliModule *AliRun::GetModule(const char *name)
{
//
// Return pointer to detector from name
//
return (AliModule*)fModules->FindObject(name);
}
//_____________________________________________________________________________
AliDetector *AliRun::GetDetector(const char *name)
{
//
// Return pointer to detector from name
//
return (AliDetector*)fModules->FindObject(name);
}
//_____________________________________________________________________________
Int_t AliRun::GetModuleID(const char *name)
{
//
// Return galice internal detector identifier from name
//
Int_t i;
for(i=0;iGet(treeName);
if (fTreeK) fTreeK->SetBranchAddress("Particles", &fParticles);
else printf("ERROR: cannot find Kine Tree for event:%d\n",event);
// Get Hits Tree header from file
sprintf(treeName,"TreeH%d",event);
fTreeH = (TTree*)gDirectory->Get(treeName);
if (!fTreeH) {
printf("ERROR: cannot find Hits Tree for event:%d\n",event);
return 0;
}
// Get Digits Tree header from file
sprintf(treeName,"TreeD%d",event);
fTreeD = (TTree*)gDirectory->Get(treeName);
if (!fTreeD) {
printf("WARNING: cannot find Digits Tree for event:%d\n",event);
}
// Get Reconstruct Tree header from file
sprintf(treeName,"TreeR%d",event);
fTreeR = (TTree*)gDirectory->Get(treeName);
if (!fTreeR) {
// printf("WARNING: cannot find Reconstructed Tree for event:%d\n",event);
}
// Set Trees branch addresses
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->SetTreeAddress();
}
if (fTreeK) fTreeK->GetEvent(0);
fNtrack = Int_t (fParticles->GetEntries());
return fNtrack;
}
//_____________________________________________________________________________
TGeometry *AliRun::GetGeometry()
{
//
// Import Alice geometry from current file
// Return pointer to geometry object
//
if (!fGeometry) fGeometry = (TGeometry*)gDirectory->Get("AliceGeom");
//
// Unlink and relink nodes in detectors
// This is bad and there must be a better way...
//
TList *tnodes=fGeometry->GetListOfNodes();
TNode *alice=(TNode*)tnodes->At(0);
TList *gnodes=alice->GetListOfNodes();
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->SetTreeAddress();
TList *dnodes=detector->Nodes();
Int_t j;
TNode *node, *node1;
for ( j=0; jGetSize(); j++) {
node = (TNode*) dnodes->At(j);
node1 = (TNode*) gnodes->FindObject(node->GetName());
dnodes->Remove(node);
dnodes->AddAt(node1,j);
}
}
return fGeometry;
}
//_____________________________________________________________________________
void AliRun::GetNextTrack(Int_t &mtrack, Int_t &ipart, Float_t *pmom,
Float_t &e, Float_t *vpos, Float_t *polar,
Float_t &tof)
{
//
// Return next track from stack of particles
//
fCurrent=-1;
GParticle *track;
for(Int_t i=fNtrack-1; i>=0; i--) {
track=(GParticle*) fParticles->UncheckedAt(i);
if(!track->TestBit(Done_Bit)) {
//
// The track has not yet been processed
fCurrent=i;
ipart=track->GetKF();
pmom[0]=track->GetPx();
pmom[1]=track->GetPy();
pmom[2]=track->GetPz();
e =track->GetEnergy();
vpos[0]=track->GetVx();
vpos[1]=track->GetVy();
vpos[2]=track->GetVz();
polar[0]=track->GetPolx();
polar[1]=track->GetPoly();
polar[2]=track->GetPolz();
tof=track->GetTime();
track->SetBit(Done_Bit);
break;
}
}
mtrack=fCurrent;
//
// stop and start timer when we start a primary track
Int_t nprimaries = fHeader.GetNprimary();
if (fCurrent >= nprimaries) return;
if (fCurrent < nprimaries-1) {
fTimer.Stop();
track=(GParticle*) fParticles->UncheckedAt(fCurrent+1);
track->SetProcessTime(fTimer.CpuTime());
}
fTimer.Start();
}
//_____________________________________________________________________________
Int_t AliRun::GetPrimary(Int_t track)
{
//
// return number of primary that has generated track
//
int current, parent;
GParticle *part;
//
parent=track;
while (1) {
current=parent;
part = (GParticle *)fParticles->UncheckedAt(current);
parent=part->GetParent();
if(parent<0) return current;
}
}
//_____________________________________________________________________________
void AliRun::Init(const char *setup)
{
//
// Initialize the Alice setup
//
gROOT->LoadMacro(setup);
gInterpreter->ProcessLine("Config();");
AliMC* pMC = AliMC::GetMC();
pMC->Gpart(); //Create standard Geant particles
sxpart(); //Define additional particles
TObject *objfirst, *objlast;
//
//=================Create Materials, geometry, histograms, etc
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->SetTreeAddress();
objlast = gDirectory->GetList()->Last();
// Initialise detector materials, geometry, histograms,etc
detector->CreateMaterials();
detector->CreateGeometry();
detector->BuildGeometry();
detector->Init();
// Add Detector histograms in Detector list of histograms
if (objlast) objfirst = gDirectory->GetList()->After(objlast);
else objfirst = gDirectory->GetList()->First();
while (objfirst) {
detector->Histograms()->Add(objfirst);
objfirst = gDirectory->GetList()->After(objfirst);
}
}
SetTransPar(); //Read the cuts for all materials
MediaTable(); //Build the special IMEDIA table
//Close the geometry structure
pMC->Ggclos();
//Initialise geometry deposition table
sEventEnergy.Set(pMC->NofVolumes()+1);
sSummEnergy.Set(pMC->NofVolumes()+1);
sSum2Energy.Set(pMC->NofVolumes()+1);
//Create the color table
pMC->SetColors();
//Compute cross-sections
pMC->Gphysi();
//Write Geometry object to current file.
fGeometry->Write();
fInitDone = kTRUE;
}
//_____________________________________________________________________________
void AliRun::MediaTable()
{
//
// Built media table to get from the media number to
// the detector id
//
Int_t kz, ibeg, nz, idt, lz, i, k, ind;
TObjArray &dets = *gAlice->Detectors();
AliModule *det;
//
// For all detectors
for (kz=0;kzLoMedium() = det->LoMedium() < idt ? det->LoMedium() : idt;
det->HiMedium() = det->HiMedium() > idt ? det->HiMedium() : idt;
}
}
if(det->LoMedium() > det->HiMedium()) {
det->LoMedium() = 0;
det->HiMedium() = 0;
} else {
if(det->HiMedium() > fImedia->GetSize()) {
Error("MediaTable","Increase fImedia");
return;
}
// Tag all materials in rage as belonging to detector kz
for(lz=det->LoMedium(); lz<= det->HiMedium(); lz++) {
(*fImedia)[lz]=kz;
}
}
}
}
//
// Print summary table
printf(" Traking media ranges:\n");
for(i=0;i<(fNdets-1)/6+1;i++) {
for(k=0;k< (6 %3d;",det->GetName(),det->LoMedium(),
det->HiMedium());
else
printf(" %6s: %3d -> %3d;","NULL",0,0);
}
printf("\n");
}
}
//____________________________________________________________________________
void AliRun::SetGenerator(AliGenerator *generator)
{
//
// Load the event generator
//
if(!fGenerator) fGenerator = generator;
}
//____________________________________________________________________________
void AliRun::SetTransPar(char* filename)
{
//
// Read filename to set the transport parameters
//
AliMC* pMC = AliMC::GetMC();
const Int_t ncuts=10;
const Int_t nflags=11;
const Int_t npars=ncuts+nflags;
const char pars[npars][7] = {"CUTGAM" ,"CUTELE","CUTNEU","CUTHAD","CUTMUO",
"BCUTE","BCUTM","DCUTE","DCUTM","PPCUTM","ANNI",
"BREM","COMP","DCAY","DRAY","HADR","LOSS",
"MULS","PAIR","PHOT","RAYL"};
char line[256];
char* filtmp;
Float_t cut[ncuts];
Int_t flag[nflags];
Int_t i, itmed, iret, ktmed, kz;
FILE *lun;
//
// See whether the file is there
filtmp=gSystem->ExpandPathName(filename);
lun=fopen(filtmp,"r");
delete [] filtmp;
if(!lun) {
printf(" * AliRun::SetTransPar * file %s does not exist!\n",filename);
return;
}
//
printf(" "); for(i=0;i<60;i++) printf("*"); printf("\n");
printf(" *%59s\n","*");
printf(" * Please check carefully what you are doing!%10s\n","*");
printf(" *%59s\n","*");
//
while(1) {
// Initialise cuts and flags
for(i=0;i=0) {
printf(" * %-6s set to %10.3E for tracking medium code %4d *\n",pars[kz],cut[kz],itmed);
pMC->Gstpar(ktmed,pars[kz],cut[kz]);
}
}
// Set transport mechanisms
for(kz=0;kz=0) {
printf(" * %-6s set to %10d for tracking medium code %4d *\n",pars[ncuts+kz],flag[kz],itmed);
pMC->Gstpar(ktmed,pars[ncuts+kz],Float_t(flag[kz]));
}
}
} else {
printf(" * Invalid tracking medium code %d *\n",itmed);
continue;
}
}
}
//_____________________________________________________________________________
void AliRun::MakeTree(Option_t *option)
{
//
// Create the ROOT trees
// Loop on all detectors to create the Root branch (if any)
//
//
// Analyse options
char *K = strstr(option,"K");
char *H = strstr(option,"H");
char *E = strstr(option,"E");
char *D = strstr(option,"D");
char *R = strstr(option,"R");
//
if (K && !fTreeK) fTreeK = new TTree("TK","Kinematics");
if (H && !fTreeH) fTreeH = new TTree("TH","Hits");
if (D && !fTreeD) fTreeD = new TTree("TD","Digits");
if (E && !fTreeE) fTreeE = new TTree("TE","Header");
if (R && !fTreeR) fTreeR = new TTree("TR","Reconstruction");
if (fTreeH) fTreeH->SetAutoSave(1000000000); //no autosave
//
// Create a branch for hits/digits for each detector
// Each branch is a TClonesArray. Each data member of the Hits classes
// will be in turn a subbranch of the detector master branch
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
if (H || D || R) detector->MakeBranch(option);
}
// Create a branch for particles
if (fTreeK && K) fTreeK->Branch("Particles",&fParticles,4000);
// Create a branch for Header
if (fTreeE && E) fTreeE->Branch("Header","AliHeader",&header,4000);
}
//_____________________________________________________________________________
Int_t AliRun::PurifyKine(Int_t lastSavedTrack, Int_t nofTracks)
{
//
// PurifyKine with external parameters
//
fHgwmk = lastSavedTrack;
fNtrack = nofTracks;
PurifyKine();
return fHgwmk;
}
//_____________________________________________________________________________
void AliRun::PurifyKine()
{
//
// Compress kinematic tree keeping only flagged particles
// and renaming the particle id's in all the hits
//
TClonesArray &particles = *fParticles;
int nkeep=fHgwmk+1, parent, i;
GParticle *part, *partnew, *father;
AliHit *OneHit;
int *map = new int[particles.GetEntries()];
// Save in Header total number of tracks before compression
fHeader.SetNtrack(fHeader.GetNtrack()+fNtrack-fHgwmk);
// Preset map, to be removed later
for(i=0; iTestBit(Keep_Bit)) {
// This particle has to be kept
map[i]=nkeep;
if(i!=nkeep) {
// Old and new are different, have to copy
partnew = (GParticle *)particles.UncheckedAt(nkeep);
*partnew = *part;
} else partnew = part;
// as the parent is always *before*, it must be already
// in place. This is what we are checking anyway!
if((parent=partnew->GetParent())>fHgwmk) {
if(map[parent]==-99) printf("map[%d] = -99!\n",parent);
partnew->SetParent(map[parent]);
}
nkeep++;
}
}
fNtrack=nkeep;
// Fix children information
for (i=fHgwmk+1; iGetParent();
father = (GParticle *)particles.UncheckedAt(parent);
if(father->TestBit(Children_Bit)) {
if(iGetFirstChild()) father->SetFirstChild(i);
if(i>father->GetLastChild()) father->SetLastChild(i);
} else {
// Iitialise children info for first pass
father->SetFirstChild(i);
father->SetLastChild(i);
father->SetBit(Children_Bit);
}
}
// Now loop on all detectors and reset the hits
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
if (!detector->Hits()) continue;
TClonesArray &vHits=*(detector->Hits());
if(vHits.GetEntries() != detector->GetNhits())
printf("vHits.GetEntries()!=detector->GetNhits(): %d != %d\n",
vHits.GetEntries(),detector->GetNhits());
for (i=0; iGetNhits(); i++) {
OneHit = (AliHit *)vHits.UncheckedAt(i);
OneHit->SetTrack(map[OneHit->GetTrack()]);
}
}
fHgwmk=nkeep-1;
particles.SetLast(fHgwmk);
delete [] map;
}
//_____________________________________________________________________________
void AliRun::Reset(Int_t run, Int_t idevent)
{
//
// Reset all Detectors & kinematics & trees
//
ResetStack();
ResetHits();
ResetDigits();
// Initialise event header
fHeader.Reset(run,idevent);
if(fTreeK) fTreeK->Reset();
if(fTreeH) fTreeH->Reset();
if(fTreeD) fTreeD->Reset();
}
//_____________________________________________________________________________
void AliRun::ResetDigits()
{
//
// Reset all Detectors digits
//
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->ResetDigits();
}
}
//_____________________________________________________________________________
void AliRun::ResetHits()
{
//
// Reset all Detectors hits
//
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->ResetHits();
}
}
//_____________________________________________________________________________
void AliRun::ResetPoints()
{
//
// Reset all Detectors points
//
TIter next(fModules);
AliModule *detector;
while((detector = (AliModule*)next())) {
detector->ResetPoints();
}
}
//_____________________________________________________________________________
void AliRun::Run(Int_t nevent, const char *setup)
{
//
// Main function to be called to process a galice run
// example
// Root > gAlice.Run();
// a positive number of events will cause the finish routine
// to be called
//
Int_t i, todo;
// check if initialisation has been done
if (!fInitDone) Init(setup);
AliMC* pMC = AliMC::GetMC();
// Create the Root Tree with one branch per detector
if(!fEvent) {
gAlice->MakeTree("KHDER");
}
todo = TMath::Abs(nevent);
for (i=0; iReset(fRun, fEvent);
pMC->Gtrigi();
pMC->Gtrigc();
pMC->Gtrig();
gAlice->FinishEvent();
fEvent++;
}
// End of this run, close files
if(nevent>0) gAlice->FinishRun();
}
//_____________________________________________________________________________
void AliRun::RunLego(const char *setup,Int_t ntheta,Float_t themin,
Float_t themax,Int_t nphi,Float_t phimin,Float_t phimax,
Float_t rmin,Float_t rmax,Float_t zmax)
{
//
// Generates lego plots of:
// - radiation length map phi vs theta
// - radiation length map phi vs eta
// - interaction length map
// - g/cm2 length map
//
// ntheta bins in theta, eta
// themin minimum angle in theta (degrees)
// themax maximum angle in theta (degrees)
// nphi bins in phi
// phimin minimum angle in phi (degrees)
// phimax maximum angle in phi (degrees)
// rmin minimum radius
// rmax maximum radius
//
//
// The number of events generated = ntheta*nphi
// run input parameters in macro setup (default="Config.C")
//
// Use macro "lego.C" to visualize the 3 lego plots in spherical coordinates
//Begin_Html
/*
*/
//End_Html
//Begin_Html
/*
*/
//End_Html
//Begin_Html
/*
*/
//End_Html
//
// check if initialisation has been done
if (!fInitDone) Init(setup);
fLego = new AliLego("lego","lego");
fLego->Init(ntheta,themin,themax,nphi,phimin,phimax,rmin,rmax,zmax);
fLego->Run();
// Create only the Root event Tree
gAlice->MakeTree("E");
// End of this run, close files
gAlice->FinishRun();
}
//_____________________________________________________________________________
void AliRun::SetCurrentTrack(Int_t track)
{
//
// Set current track number
//
fCurrent = track;
}
//_____________________________________________________________________________
void AliRun::SetTrack(Int_t done, Int_t parent, Int_t ipart, Float_t *pmom,
Float_t *vpos, Float_t *polar, Float_t tof,
const char *mecha, Int_t &ntr, Float_t weight)
{
//
// Load a track on the stack
//
// done 0 if the track has to be transported
// 1 if not
// parent identifier of the parent track. -1 for a primary
// ipart particle code
// pmom momentum GeV/c
// vpos position
// polar polarisation
// tof time of flight in seconds
// mecha production mechanism
// ntr on output the number of the track stored
//
TClonesArray &particles = *fParticles;
GParticle *particle;
Float_t mass;
char pname[21];
const Int_t firstchild=-1;
const Int_t lastchild=-1;
const Int_t KS=0;
const Float_t tlife=0;
AliMC::GetMC()->GetParticle(ipart,pname,mass);
Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
pmom[1]*pmom[1]+pmom[2]*pmom[2]);
//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",
//pname,mass,e,fNtrack,ipart,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],KS,mecha);
particle=new(particles[fNtrack]) GParticle(KS,ipart,parent,firstchild,
lastchild,pmom[0],pmom[1],pmom[2],
e,mass,vpos[0],vpos[1],vpos[2],
polar[0],polar[1],polar[2],tof,
tlife,mecha,weight);
if(!done) particle->SetBit(Done_Bit);
if(parent>=0) {
particle=(GParticle*) fParticles->UncheckedAt(parent);
particle->SetLastChild(fNtrack);
if(particle->GetFirstChild()<0) particle->SetFirstChild(fNtrack);
} else {
//
// This is a primary track. Set high water mark for this event
fHgwmk=fNtrack;
//
// Set also number if primary tracks
fHeader.SetNprimary(fHgwmk+1);
fHeader.SetNtrack(fHgwmk+1);
}
ntr = fNtrack++;
}
//_____________________________________________________________________________
void AliRun::KeepTrack(const Int_t track)
{
//
// flags a track to be kept
//
TClonesArray &particles = *fParticles;
((GParticle*)particles[track])->SetBit(Keep_Bit);
}
//_____________________________________________________________________________
void AliRun::StepManager(Int_t id) const
{
//
// Called at every step during transport
//
AliMC* pMC = AliMC::GetMC();
Int_t copy;
//
// --- If lego option, do it and leave
if (fLego) {
fLego->StepManager();
return;
}
//Update energy deposition tables
sEventEnergy[pMC->CurrentVol(0,copy)]+=pMC->Edep();
//Call the appropriate stepping routine;
AliModule *det = (AliModule*)fModules->At(id);
if(det) det->StepManager();
}
//_____________________________________________________________________________
void AliRun::ReadEuclid(const char* filnam, Int_t id_det, const char* topvol)
{
//
// read in the geometry of the detector in euclid file format
//
// id_det : the detector identification (2=its,...)
// topvol : return parameter describing the name of the top
// volume of geometry.
//
// author : m. maire
//
// 28.07.98
// several changes have been made by miroslav helbich
// subroutine is rewrited to follow the new established way of memory
// booking for tracking medias and rotation matrices.
// all used tracking media have to be defined first, for this you can use
// subroutine greutmed.
// top volume is searched as only volume not positioned into another
//
AliMC* pMC = AliMC::GetMC();
Int_t i, nvol, iret, itmed, irot, numed, npar, ndiv, iaxe;
Int_t ndvmx, nr, flag;
char key[5], card[77], natmed[21];
char name[5], mother[5], shape[5], konly[5], volst[7000][5];
char *filtmp;
Float_t par[50];
Float_t teta1, phi1, teta2, phi2, teta3, phi3, orig, step;
Float_t xo, yo, zo;
Int_t idrot[5000],istop[7000];
FILE *lun;
AliModule *det;
//
TObjArray &dets = *fModules;
if(!dets[id_det]) {
printf(" *** GREUTMED *** Detector %d not defined\n",id_det);
return;
} else {
det = (AliModule*) dets[id_det];
}
//
// *** The input filnam name will be with extension '.euc'
filtmp=gSystem->ExpandPathName(filnam);
lun=fopen(filtmp,"r");
delete [] filtmp;
if(!lun) {
printf(" *** GREUCL *** Could not open file %s\n",filnam);
return;
}
//* --- definition of rotation matrix 0 ---
idrot[0]=0;
nvol=0;
L10:
for(i=0;i<77;i++) card[i]=0;
iret=fscanf(lun,"%77[^\n]",card);
if(iret<=0) goto L20;
fscanf(lun,"%*c");
//*
strncpy(key,card,4);
key[4]='\0';
if (!strcmp(key,"TMED")) {
sscanf(&card[5],"%d '%[^']'",&itmed,natmed);
//Pad the string with blanks
i=-1;
while(natmed[++i]);
while(i<20) natmed[i++]=' ';
natmed[i]='\0';
//
pMC->Gckmat(fIdtmed[itmed+id_det*100-1],natmed);
//*
} else if (!strcmp(key,"ROTM")) {
sscanf(&card[4],"%d %f %f %f %f %f %f",&irot,&teta1,&phi1,&teta2,&phi2,&teta3,&phi3);
det->AliMatrix(idrot[irot],teta1,phi1,teta2,phi2,teta3,phi3);
//*
} else if (!strcmp(key,"VOLU")) {
sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, shape, &numed, &npar);
if (npar>0) {
for(i=0;iGsvolu( name, shape, fIdtmed[numed+id_det*100-1], par, npar);
//* save the defined volumes
strcpy(volst[++nvol],name);
istop[nvol]=1;
//*
} else if (!strcmp(key,"DIVN")) {
sscanf(&card[5],"'%[^']' '%[^']' %d %d", name, mother, &ndiv, &iaxe);
pMC->Gsdvn ( name, mother, ndiv, iaxe );
//*
} else if (!strcmp(key,"DVN2")) {
sscanf(&card[5],"'%[^']' '%[^']' %d %d %f %d",name, mother, &ndiv, &iaxe, &orig, &numed);
pMC->Gsdvn2( name, mother, ndiv, iaxe, orig,fIdtmed[numed+id_det*100-1]);
//*
} else if (!strcmp(key,"DIVT")) {
sscanf(&card[5],"'%[^']' '%[^']' %f %d %d %d", name, mother, &step, &iaxe, &numed, &ndvmx);
pMC->Gsdvt ( name, mother, step, iaxe, fIdtmed[numed+id_det*100-1], ndvmx);
//*
} else if (!strcmp(key,"DVT2")) {
sscanf(&card[5],"'%[^']' '%[^']' %f %d %f %d %d", name, mother, &step, &iaxe, &orig, &numed, &ndvmx);
pMC->Gsdvt2 ( name, mother, step, iaxe, orig, fIdtmed[numed+id_det*100-1], ndvmx );
//*
} else if (!strcmp(key,"POSI")) {
sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']'", name, &nr, mother, &xo, &yo, &zo, &irot, konly);
//*** volume name cannot be the top volume
for(i=1;i<=nvol;i++) {
if (!strcmp(volst[i],name)) istop[i]=0;
}
//*
pMC->Gspos ( name, nr, mother, xo, yo, zo, idrot[irot], konly );
//*
} else if (!strcmp(key,"POSP")) {
sscanf(&card[5],"'%[^']' %d '%[^']' %f %f %f %d '%[^']' %d", name, &nr, mother, &xo, &yo, &zo, &irot, konly, &npar);
if (npar > 0) {
for(i=0;iGsposp ( name, nr, mother, xo,yo,zo, idrot[irot], konly, par, npar);
}
//*
if (strcmp(key,"END")) goto L10;
//* find top volume in the geometry
flag=0;
for(i=1;i<=nvol;i++) {
if (istop[i] && flag) {
printf(" *** GREUCL *** warning: %s is another possible top volume\n",volst[i]);
}
if (istop[i] && !flag) {
topvol=volst[i];
printf(" *** GREUCL *** volume %s taken as a top volume\n",topvol);
flag=1;
}
}
if (!flag) {
printf("*** GREUCL *** warning: top volume not found\n");
}
fclose (lun);
//*
//* commented out only for the not cernlib version
printf(" *** GREUCL *** file: %s is now read in\n",filnam);
//
return;
//*
L20:
printf(" *** GREUCL *** reading error or premature end of file\n");
}
//_____________________________________________________________________________
void AliRun::ReadEuclidMedia(const char* filnam, Int_t id_det)
{
//
// read in the materials and tracking media for the detector
// in euclid file format
//
// filnam: name of the input file
// id_det: id_det is the detector identification (2=its,...)
//
// author : miroslav helbich
//
Float_t sxmgmx = gAlice->Field()->Max();
Int_t isxfld = gAlice->Field()->Integ();
Int_t end, i, iret, itmed;
char key[5], card[130], natmed[21], namate[21];
Float_t ubuf[50];
char* filtmp;
FILE *lun;
Int_t imate;
Int_t nwbuf, isvol, ifield, nmat;
Float_t a, z, dens, radl, absl, fieldm, tmaxfd, stemax, deemax, epsil, stmin;
AliModule* det;
//
TObjArray &dets = *fModules;
if(!dets[id_det]) {
printf(" *** GREUTMED *** Detector %d not defined\n",id_det);
return;
} else {
det = (AliModule*) dets[id_det];
}
end=strlen(filnam);
for(i=0;iExpandPathName(filnam);
lun=fopen(filtmp,"r");
delete [] filtmp;
if(!lun) {
printf(" *** GREUTMED *** Could not open file %s\n",filnam);
return;
}
//
// Retrieve Mag Field parameters
Int_t ISXFLD=gAlice->Field()->Integ();
Float_t SXMGMX=gAlice->Field()->Max();
//
L10:
for(i=0;i<130;i++) card[i]=0;
iret=fscanf(lun,"%4s %[^\n]",key,card);
if(iret<=0) goto L20;
fscanf(lun,"%*c");
//*
//* read material
if (!strcmp(key,"MATE")) {
sscanf(card,"%d '%[^']' %f %f %f %f %f %d",&imate,namate,&a,&z,&dens,&radl,&absl,&nwbuf);
if (nwbuf>0) for(i=0;iAliMaterial(imate,namate,a,z,dens,radl,absl,ubuf,nwbuf);
//* read tracking medium
} else if (!strcmp(key,"TMED")) {
sscanf(card,"%d '%[^']' %d %d %d %f %f %f %f %f %f %d",
&itmed,natmed,&nmat,&isvol,&ifield,&fieldm,&tmaxfd,
&stemax,&deemax,&epsil,&stmin,&nwbuf);
if (nwbuf>0) for(i=0;iAliMedium(itmed+id_det*100,natmed,nmat,isvol,ISXFLD,SXMGMX,tmaxfd,
stemax,deemax,epsil,stmin,ubuf,nwbuf);
(*fImedia)[fIdtmed[itmed+id_det*100-1]-1]=id_det;
//*
}
//*
if (strcmp(key,"END")) goto L10;
fclose (lun);
//*
//* commented out only for the not cernlib version
printf(" *** GREUTMED *** file: %s is now read in\n",filnam);
//*
return;
//*
L20:
printf(" *** GREUTMED *** reading error or premature end of file\n");
}
//_____________________________________________________________________________
void AliRun::Streamer(TBuffer &R__b)
{
//
// Stream an object of class AliRun.
//
if (R__b.IsReading()) {
Version_t R__v = R__b.ReadVersion(); if (R__v) { }
TNamed::Streamer(R__b);
if (!gAlice) gAlice = this;
gROOT->GetListOfBrowsables()->Add(this,"Run");
R__b >> fNtrack;
R__b >> fHgwmk;
R__b >> fDebug;
fHeader.Streamer(R__b);
R__b >> fModules;
R__b >> fParticles;
R__b >> fField;
// R__b >> fMC;
R__b >> fNdets;
R__b >> fTrRmax;
R__b >> fTrZmax;
R__b >> fGenerator;
} else {
R__b.WriteVersion(AliRun::IsA());
TNamed::Streamer(R__b);
R__b << fNtrack;
R__b << fHgwmk;
R__b << fDebug;
fHeader.Streamer(R__b);
R__b << fModules;
R__b << fParticles;
R__b << fField;
// R__b << fMC;
R__b << fNdets;
R__b << fTrRmax;
R__b << fTrZmax;
R__b << fGenerator;
}
}
//_____________________________________________________________________________
//
// Interfaces to Fortran
//
//_____________________________________________________________________________
extern "C" void type_of_call rxgtrak (Int_t &mtrack, Int_t &ipart, Float_t *pmom,
Float_t &e, Float_t *vpos, Float_t &tof)
{
//
// Fetches next track from the ROOT stack for transport. Called by the
// modified version of GTREVE.
//
// Track number in the ROOT stack. If MTRACK=0 no
// mtrack more tracks are left in the stack to be
// transported.
// ipart Particle code in the GEANT conventions.
// pmom[3] Particle momentum in GeV/c
// e Particle energy in GeV
// vpos[3] Particle position
// tof Particle time of flight in seconds
//
Float_t polar[3];
gAlice->GetNextTrack(mtrack, ipart, pmom, e, vpos, polar, tof);
mtrack++;
}
//_____________________________________________________________________________
extern "C" void type_of_call
#ifndef WIN32
rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
Float_t *vpos, Float_t &tof, const char* cmech, Int_t &ntr, const int cmlen)
#else
rxstrak (Int_t &keep, Int_t &parent, Int_t &ipart, Float_t *pmom,
Float_t *vpos, Float_t &tof, const char* cmech, const int cmlen,
Int_t &ntr)
#endif
{
//
// Fetches next track from the ROOT stack for transport. Called by GUKINE
// and GUSTEP.
//
// Status of the track. If keep=0 the track is put
// keep on the ROOT stack but it is not fetched for
// transport.
// parent Parent track. If parent=0 the track is a primary.
// In GUSTEP the routine is normally called to store
// secondaries generated by the current track whose
// ROOT stack number is MTRACK (common SCKINE.
// ipart Particle code in the GEANT conventions.
// pmom[3] Particle momentum in GeV/c
// vpos[3] Particle position
// tof Particle time of flight in seconds
//
// cmech (CHARACTER*10) Particle origin. This field is user
// defined and it is not used inside the GALICE code.
// ntr Number assigned to the particle in the ROOT stack.
//
char mecha[11];
Float_t polar[3]={0.,0.,0.};
for(int i=0; i<10 && iSetTrack(keep, parent-1, ipart, pmom, vpos, polar, tof, mecha, ntr);
ntr++;
}
//_____________________________________________________________________________
extern "C" void type_of_call rxkeep(const Int_t &n)
{
if( NULL==gAlice ) exit(1);
if( n<=0 || n>gAlice->Particles()->GetEntries() )
{
printf(" Bad index n=%d must be 0Particles()->GetEntries());
exit(1);
}
((GParticle*)(gAlice->Particles()->UncheckedAt(n-1)))->SetBit(Keep_Bit);
}
//_____________________________________________________________________________
extern "C" void type_of_call rxouth ()
{
//
// Called by Gtreve at the end of each primary track
//
gAlice->FinishPrimary();
}