Version: SMASH-3.1
ROOT format

SMASH ROOT output is a fast and disk-space efficient, but not human-readable output. It is a custom format making information about the SMASH calculation accessible with ROOT, mostly mirroring the information of the OSCAR format formats. This output is distinct from the standarized HepMC Output that is also available in ROOT format and more widely adopted.

SMASH ROOT output files can be viewed using ROOT's TBrowser. One can also access them using ROOT functions. The full memory structure of the ROOT files can be found here: http://root.cern.ch/root/html/TFile.html. We only describe the logical structure of the SMASH ROOT output. Knowing the logical structure is enough to read and write ROOT files, be able to view them in TBrowser, or write a ROOT macro to analyze them.

Producing ROOT output requires ROOT installed (see http://root.cern.ch).

Depending on configuration (see Output) SMASH can produces up to two ROOT files per run: Particles.root and Collisions.root. These files contain a TTree called particles and a TTree called collisions. The particles tree contains information about the parameters of the run (such as the number of testparticles and event number), information relating to individual particles (such as their position or charge), and information about bulk observables in the system (kinetic energy, mean field energy, and total energy). The collisions tree contains information about each collision, such as number of incoming and outgoing particles. It also has the full information about the incoming and outgoing particles of each collision.

In case that the ROOT format is used for dilepton output (see Dileptons), the ROOT file is called Dileptons.root and only contains a collisions tree with all the dilepton decays.

Every physical quantity corresponds to a separate TBranch. One entry in the particles TTree is:

ev tcounter npart test_p modus_l current_t impact_b empty_event
pdgcode[npart] charge[npart] t[npart] x[npart] y[npart] z[npart] p0[npart]
px[npart] py[npart] pz[npart] E_kinetic_tot E_fields_tot E_tot

The maximal number of particles in one entry is limited to 500000. This is done to limit the buffer size needed for ROOT output. If the number of particles in one block exceeds 500000, then they are written in separate blocks with the same tcounter and ev. The fields have the following meaning:

  • ev is event number
  • tcounter is number of output block in a given event in terms of OSCAR
  • npart is number of particles in the block
  • test_p is number of testparticles per particle
  • modus_l is modus length
  • current_t is time associated with the output block, in fm
  • impact_b is the impact parameter of the event
  • empty_event indicates whether the projectile did not interact with the target
  • pdgcode is PDG id array
  • charge is the electric charge array
  • p0, px, py, pz are 4-momenta arrays
  • t, x, y, z are position arrays
  • E_kinetic_tot is total kinetic energy in the system
  • E_fields_tot is total mean field energy * test_p
  • E_total is the sum of E_kinetic_tot and E_fields_tot

In case of extended output (see input_output_content_specific_) more fields are added. Their description is the same that in case of OSCAR format, see extended_output_format_.

The entries in the collisions tree are organized in the same way, but a few additional fields are present:

  • nin and nout are added to characterize number of incoming and outgoing particles in the reaction, with nin + nout = npart.
  • weight is an action weight, whose meaning depends on the type of action: For collisions it is the total cross section, for decays it is the total decay width and for dilepton decays it is the shining weight.

Currently writing initial and final configuration to collisions tree is not supported.

See also Collision output in box modus.

Here is an example of a basic ROOT macro to read the ROOT output of SMASH:

// file name: basic_macro.C
#include <TFile.h>
#include <TTree.h>
int rootfile_basic_example() {
// open SMASH output file to be read in
TFile *input_file = TFile::Open("../build/data/0/Particles.root");
if (input_file->IsOpen()) {
printf("Successfully opened file %s\n", input_file->GetName());
} else {
printf("Error at opening file %s\n", input_file->GetName());
}
// Get a tree from file
TTree *tree = static_cast<TTree*>(input_file->Get("particles"));
// Get number of entries in a tree
Int_t nentries = tree->GetEntries();
printf("Number of entries in a tree is %d\n", nentries);
// This draws p_T distribution at initialization
// tree->Draw("sqrt(px*px + py*py)","tcounter==0");
// This draws 3D momentum space distribution at initialization
tree->Draw("px:py:pz","tcounter==0");
return 0;
}

To execute this macro:

root -l
.L basic_macro.C+
rootfile_basic_example()

To generate a ROOT macro for a more involved analysis:

root -l Particles.root
particles->MakeClass("analysis")

This creates analysis.h and analysis.C, the latter of which provides an example of a basic loop over entries. The user can modify the Loop() function and build a more complicated analysis. The following is an example of a macro using an array of histograms to obtain the radial momentum distribution for each of the entries in the Particles.root:

#define analysis_cxx
#include "analysis.h"
#include <TH2.h>
#include <TStyle.h>
#include <TCanvas.h>
#include <iostream>
void analysis::Loop()
{
if (fChain == 0) return;
Long64_t n_entries = fChain->GetEntriesFast();
// an array of histograms
TH1D *h_p_avg[n_entries];
// Each histogram needs to be declared with a unique name and title
for (int i = 0; i < n_entries; i++){
char h_p_avg_name[256];
char h_p_avg_title[256];
sprintf(h_p_avg_name, "h_p_avg_entry_%d", i);
sprintf(h_p_avg_title, "momentum distribution at entry %d", i);
h_p_avg[i] = new TH1D(h_p_avg_name, h_p_avg_title, 50, 0, 1.0);
h_p_avg[i]->Sumw2();
h_p_avg[i]->SetStats(1);
}
Long64_t nb = 0;
// A loop over all entries
for (Long64_t j_entry = 0; j_entry < n_entries; j_entry++){
//Load the TTree data for that entry
Long64_t i_entry = LoadTree(j_entry);
if (i_entry < 0){
std::cout << "Failed to load the TTree at j_entry = "
<< j_entry << std::endl;
break;
}
nb = fChain->GetEntry(j_entry);
// A loop over all particles in the entry
for (int i = 0; i < npart; i++) {
const double p = sqrt( px[i] * px[i] + py[i] * py[i] + pz[i] * pz[i]
);
// filling the j_entry-th histogram
h_p_avg[j_entry]->Fill(p, 1.0/(p*p));
}
}
// drawing the histogram corresponding to the j_entry = 0 entry
h_p_avg[0]->Draw();
}
/endcode
To run this analysis:
\code
root -l
.L analysis.C+
analysis t
t.Loop()
constexpr int p
Proton.

To quickly view a ROOT file from the command line:

root -l Particles.root // attaches the .root file
.ls // lists objects contained in the .root file; here: particles
particles->Draw("p0", "tcounter == 0")

Viewing a ROOT file can be also done in a TBrowser:

root -l
new TBrowser

For more examples of extracting info from .root file see root.cern.ch