ROOT logo
//-------------------------------------------------------------------------
//           Implementation of Class AliHEPDataParser
//
//  This class is used to save the content of hisograms/graphs in the
//  HEP data format and viceversa. The HEP data format is not strictly
//  defined and there are many variants, the class only support a few
//  of them. More will be added as needed.  The input can be a set of
//  2 TH1, TGraphAsymmErrors or TGraphErrors (one for the stat and one
//  for the syst error). If the second one is a null pointer, only the
//  stat error is printed. The class can also import hepdata ascii
//  file (very preliminary)
// 
//  Author: Michele Floris, CERN
//-------------------------------------------------------------------------


#include "AliHEPDataParser.h"
#include "AliLog.h"
#include "TGraphAsymmErrors.h"
#include "TGraph.h"
#include "TGraphErrors.h"
#include "TH1.h"
#include "TObjArray.h"
#include "TObjString.h"
#include "TMath.h"
#include <fstream>
#include <iostream>
using namespace std;

ClassImp(AliHEPDataParser)

AliHEPDataParser::AliHEPDataParser() : TObject(), fHistStat(0),  fHistSyst(0),  fGraphStat(0),  fGraphSyst(0),  fHEPDataFileLines(0), fValueName(""), fXaxisName(""), fTitle(""), fReaction(""), fEnergy(""), fRapidityRange(""), fPrecision(5)
{
  // default ctor

}

AliHEPDataParser::AliHEPDataParser(TH1 * hStat, TH1 * hSyst): TObject(), fHistStat(0),  fHistSyst(0),  fGraphStat(0),  fGraphSyst(0),  fHEPDataFileLines(0), fValueName("y"), fXaxisName(""), fTitle(""), fReaction(""), fEnergy(""), fRapidityRange(""), fPrecision(5)
{
  //ctor
  fHistStat = hStat;
  fHistSyst = hSyst;
  fHEPDataFileLines = new TObjArray;

}
AliHEPDataParser::AliHEPDataParser(TGraph * grStat, TGraph * grSyst): TObject(), fHistStat(0),  fHistSyst(0),  fGraphStat(0),  fGraphSyst(0),  fHEPDataFileLines(0), fValueName(""), fXaxisName(""), fTitle(""), fReaction(""), fEnergy(""), fRapidityRange(""), fPrecision(5)
{
  // ctor
  fGraphStat = grStat;
  fGraphSyst = grSyst;
  fHEPDataFileLines = new TObjArray;
}

AliHEPDataParser::AliHEPDataParser(const char * hepfileName): TObject(), fHistStat(0),  fHistSyst(0),  fGraphStat(0),  fGraphSyst(0),  fHEPDataFileLines(0), fValueName("y"), fXaxisName(""), fTitle(""), fReaction(""), fEnergy(""), fRapidityRange(""), fPrecision(5)
{
  //ctor
  // Put results in graphs
  fGraphSyst = new TGraphAsymmErrors();
  fGraphStat = new TGraphAsymmErrors();
  ifstream infile;
  infile.open(hepfileName);
  TString line;
  Int_t ipoints = 0;
  while (line.ReadLine(infile)) {
    TObjArray * tokens = line.Tokenize(" \t");
    if(!tokens) continue;
    if(! ((TObjString*) tokens->At(0))->String().Atof()) {
      //The first column is not a number: those are the headers: skip!
      delete tokens;
      continue;
    }
    if(tokens->GetEntries() != 8) {
      // this should never happen!
      delete tokens;
      AliError(Form("Wrong number of columns %d! Assuming [x xlow xhigh y dystat+ dystat- dysyst+ dysyst-]", tokens->GetEntries()));
      cout << line.Data() << endl;
      return;      
      //continue;
    }
    // FIXME: Assumes the format
    // x xlow xhigh y dystat+ dystat- dysyst+ dysyst-
    TObjString * xbin   = (TObjString*) tokens->At(0);
    TObjString * xlow   = (TObjString*) tokens->At(1);
    TObjString * xhigh  = (TObjString*) tokens->At(2);
    TObjString * value  = (TObjString*) tokens->At(3);
    TObjString * statp  = (TObjString*) tokens->At(4);
    TObjString * statm  = (TObjString*) tokens->At(5);
    TObjString * systp  = (TObjString*) tokens->At(6);
    TObjString * systm  = (TObjString*) tokens->At(7);
    statm->String().ReplaceAll("+-","");
    statp->String().ReplaceAll("+-","");
    if(systp) systp->String().ReplaceAll("+-","");
    if(systm) systm->String().ReplaceAll("+-","");
    //    if (!binMin->String().Atof()) {delete tokens; continue;} // skip headers
    Float_t binCenter = xbin->String().Atof();
    Float_t binWidth  =  (xlow->String().Atof() - xhigh->String().Atof())/2;


    fGraphStat->SetPoint(ipoints, binCenter, value->String().Atof());
    fGraphSyst->SetPoint(ipoints, binCenter, value->String().Atof());
    ((TGraphAsymmErrors*)fGraphStat)->SetPointError(ipoints, 
						    binWidth,
						    binWidth,
						    statp->String().Atof(),
						    statm->String().Atof());
    if(systp && systm) ((TGraphAsymmErrors*)fGraphSyst)->SetPointError(ipoints, 
							      binWidth,
							      binWidth,
							      systp->String().Atof(),
							      systm->String().Atof());
    ipoints++;
    delete tokens;
  }
  infile.close();
    

}

AliHEPDataParser::~AliHEPDataParser(){
  // dtor
  if(fHistStat) delete fHistStat;
  if(fHistSyst) delete fHistSyst;
  if(fGraphStat) delete fGraphStat;
  if(fGraphSyst) delete fGraphSyst;
  if(fHEPDataFileLines) delete fHEPDataFileLines;
}
  
void AliHEPDataParser::SaveHEPDataFile(const char * hepfileName, Bool_t trueUseGraphFalesUseHisto) {

  // Fills fHEPDataFileLines and saves its content to a file
  if(!fHEPDataFileLines)   fHEPDataFileLines = new TObjArray;
  // Write headers if relevant
  if(fTitle.Length())         fHEPDataFileLines->Add(new TObjString(fTitle));
  if(fReaction.Length())      fHEPDataFileLines->Add(new TObjString(fReaction));
  if(fEnergy.Length())        fHEPDataFileLines->Add(new TObjString(fEnergy));
  if(fRapidityRange.Length()) fHEPDataFileLines->Add(new TObjString(fRapidityRange));
  if(!fValueName.Length() || !fXaxisName.Length()) AliFatal("At least x and y titles should be given!");
  fHEPDataFileLines->Add(new TObjString(Form("x: %s", fXaxisName.Data())));
  fHEPDataFileLines->Add(new TObjString(Form("y: %s", fValueName.Data())));


  if(trueUseGraphFalesUseHisto) {
    AliWarning("Graph saving not thoroughly tested!!");
    if(!fGraphStat) {
      AliError("Graph not set");
      return;
    }
    Bool_t asym = kFALSE; // check if this has asymmetric errors
    if       (!strcmp(fGraphStat->ClassName(), "TGraphErrors")) asym = kFALSE;
    else if  (!strcmp(fGraphStat->ClassName(), "TGraphAsymmErrors")) asym = kTRUE;
    else     {AliError("Unsupported graph type"); return;}
    Int_t npoint = fGraphStat->GetN();
    if(asym) AliInfo("Assymmetric errors");
    for(Int_t ipoint = 0; ipoint < npoint; ipoint++){            
      if(ipoint == 0) {
	if(fGraphSyst) {
	  fHEPDataFileLines->Add(new TObjString("x\txlow\txhigh\t+stat\t-stat\t+syst\t-syst"));
	}
	else {
	  fHEPDataFileLines->Add(new TObjString("x\txlow\txhigh\t+stat\t-stat"));
	}
      }
      // Skip empty bins
      if(!fGraphStat->GetY()[ipoint]) continue;
      TObjString * line = new TObjString;    
      if(fGraphSyst) {
	if (asym)
	  line->String().Form("%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s",
			      //	line->String().Form("%10s %10s %10s %10s %10s %10s %10s %10s", 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint],                                                                   10).Data(), 			    
			      GetFixWidthCol(RoundToSignificantFigures(fGraphStat->GetX()[ipoint]-((TGraphAsymmErrors*)fGraphStat)->GetEXlow()[ipoint],  fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(fGraphStat->GetX()[ipoint]+((TGraphAsymmErrors*)fGraphStat)->GetEXhigh()[ipoint], fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(fGraphStat->GetY()[ipoint],                            fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphAsymmErrors*)fGraphStat)->GetEYhigh()[ipoint], fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphAsymmErrors*)fGraphStat)->GetEYlow()[ipoint] , fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphAsymmErrors*)fGraphSyst)->GetEYhigh()[ipoint], fPrecision), 10).Data(),
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphAsymmErrors*)fGraphSyst)->GetEYlow()[ipoint] , fPrecision), 10).Data());
	else 
	  line->String().Form("%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s",
			      //	line->String().Form("%10s %10s %10s %10s %10s %10s %10s %10s", 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint],                                                                   10).Data(), 			    
			      GetFixWidthCol(fGraphStat->GetX()[ipoint]-fGraphStat->GetEX()[ipoint],                                    10).Data(), 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint]+fGraphStat->GetEX()[ipoint],                                10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(fGraphStat->GetY()[ipoint],                            fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphErrors*)fGraphStat)->GetEY()[ipoint], fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphErrors*)fGraphStat)->GetEY()[ipoint], fPrecision), 10).Data(),
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphErrors*)fGraphSyst)->GetEY()[ipoint], fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphErrors*)fGraphSyst)->GetEY()[ipoint], fPrecision), 10).Data());
	  // line->String().Form("%f %f +-%f +-%f", 
	  // 		      fGraphStat->GetX()[ipoint], RoundToSignificantFigures(fGraphStat->GetY()[ipoint],fPrecision),
	  // 		      RoundToSignificantFigures(((TGraphErrors*)fGraphStat)->GetEY()[ipoint],fPrecision), 
	  // 		      RoundToSignificantFigures(((TGraphErrors*)fGraphSyst)->GetEY()[ipoint],fPrecision));

	fHEPDataFileLines->Add(line);
      } else {
	if (asym)
	  line->String().Form("%s\t%s\t%s\t%s\t%s\t%s",
			      //	line->String().Form("%10s %10s %10s %10s %10s %10s %10s %10s", 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint],                                                                   10).Data(), 			    
			      GetFixWidthCol(fGraphStat->GetX()[ipoint]-fGraphStat->GetEXlow()[ipoint],                                    10).Data(), 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint]+fGraphStat->GetEXhigh()[ipoint],                                   10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(fGraphStat->GetY()[ipoint],                            fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphAsymmErrors*)fGraphStat)->GetEYhigh()[ipoint], fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphAsymmErrors*)fGraphStat)->GetEYlow()[ipoint] , fPrecision), 10).Data()); 
	else 
	  line->String().Form("%s\t%s\t%s\t%s\t%s\t%s",
			      //	line->String().Form("%10s %10s %10s %10s %10s %10s %10s %10s", 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint],                                                                   10).Data(), 			    
			      GetFixWidthCol(fGraphStat->GetX()[ipoint]-fGraphStat->GetEX()[ipoint],                                    10).Data(), 
			      GetFixWidthCol(fGraphStat->GetX()[ipoint]+fGraphStat->GetEX()[ipoint],                                10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(fGraphStat->GetY()[ipoint],                            fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphErrors*)fGraphStat)->GetEY()[ipoint], fPrecision), 10).Data(), 
			      GetFixWidthCol(RoundToSignificantFigures(((TGraphErrors*)fGraphStat)->GetEY()[ipoint], fPrecision), 10).Data());

	fHEPDataFileLines->Add(line);
      }
    }    
  }
  else {
    if(!fHistStat) {
      AliError("Hist not set");
      return;
    }
    Int_t nbin = fHistStat->GetNbinsX();
    
    for(Int_t ibin = 1; ibin <= nbin; ibin++){
      if(ibin == 1) {
	if(fHistSyst)
	  fHEPDataFileLines->Add(new TObjString("x\t\txlow\t\txhigh\t\ty\t\tdystat+\t\tdystat-\t\tdysyst+\t\tdysyst-")); 
	else 
	  fHEPDataFileLines->Add(new TObjString("x\t\txlow\t\txhigh\t\ty\t\tdystat+\t\tdystat-"));       	
      }
      // Skip empty bins
      if(!fHistStat->GetBinContent(ibin)) continue;
      TObjString * line = new TObjString;      
      if(fHistSyst) {
	
	line->String().Form("%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s",
			    //	line->String().Form("%10s %10s %10s %10s %10s %10s %10s %10s", 
			    GetFixWidthCol(fHistStat->GetBinLowEdge(ibin)+fHistStat->GetBinWidth(ibin)/2,         12).Data(), 			    
			    GetFixWidthCol(fHistStat->GetBinLowEdge(ibin),                                        12).Data(), 
			    GetFixWidthCol(fHistStat->GetBinLowEdge(ibin)+fHistStat->GetBinWidth(ibin),           12).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistStat->GetBinContent(ibin),fPrecision),  12).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistStat->GetBinError(ibin),  fPrecision),  12).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistStat->GetBinError(ibin),  fPrecision),  12).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistSyst->GetBinError(ibin),  fPrecision),  12).Data(),
			    GetFixWidthCol(RoundToSignificantFigures(fHistSyst->GetBinError(ibin),  fPrecision),  12).Data());

	fHEPDataFileLines->Add(line);
      } else {
	line->String().Form("%s\t%s\t%s\t%s\t%s\t%s", 
			    GetFixWidthCol(fHistStat->GetBinLowEdge(ibin)+fHistStat->GetBinWidth(ibin)/2,         10).Data(), 			    
			    GetFixWidthCol(fHistStat->GetBinLowEdge(ibin),                                        10).Data(), 
			    GetFixWidthCol(fHistStat->GetBinLowEdge(ibin)+fHistStat->GetBinWidth(ibin),           10).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistStat->GetBinContent(ibin),fPrecision),  10).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistStat->GetBinError(ibin),  fPrecision),  10).Data(), 
			    GetFixWidthCol(RoundToSignificantFigures(fHistStat->GetBinError(ibin),  fPrecision),  10).Data()); 
	fHEPDataFileLines->Add(line);
      }
      //      delete line;      
    }           
  }

  TIterator * lineIter = fHEPDataFileLines->MakeIterator();
  TObjString * obj = 0;
  ofstream outfile;
  outfile.open (hepfileName);
  cout << "Saving HEP File " << hepfileName << endl;
  
  while ((obj = (TObjString*) lineIter->Next())) {    
    cout << obj->String().Data() << endl;    
    outfile << obj->String().Data() << endl;
  }
  outfile.close();
}


Double_t AliHEPDataParser::RoundToSignificantFigures(double num, int n) {
  // Rounds num to n significant digits.
  // Recipe from :http://stackoverflow.com/questions/202302/rounding-to-an-arbitrary-number-of-significant-digits
  // Basically the log is used to determine the number of leading 0s, than convert to an integer by multipliing by the expo, 
  // round the integer and shift back.
  if(num == 0) {
    return 0;
  }

  Double_t d = TMath::Ceil(TMath::Log10(num < 0 ? -num: num));
  Int_t power = n - (int) d;

  Double_t magnitude = TMath::Power(10, power);
  Long_t shifted = TMath::Nint(num*magnitude);
  return shifted/magnitude;

}

TString AliHEPDataParser::GetFixWidthCol(Double_t number, Int_t width) {

  // Formats a column to fixed width
  TString col;
  char format[100];
  snprintf(format,100,"%%%d#g", fPrecision);
  col.Form(format, number);
  if(col.Length()>width) AliError("larger than width, cannot align!");

  if(col.Contains("e"))
    while (col.Length() < width) col.Append(" ");
  else
    while (col.Length() < width) col.Append("0");
  
  return col;
}
 AliHEPDataParser.cxx:1
 AliHEPDataParser.cxx:2
 AliHEPDataParser.cxx:3
 AliHEPDataParser.cxx:4
 AliHEPDataParser.cxx:5
 AliHEPDataParser.cxx:6
 AliHEPDataParser.cxx:7
 AliHEPDataParser.cxx:8
 AliHEPDataParser.cxx:9
 AliHEPDataParser.cxx:10
 AliHEPDataParser.cxx:11
 AliHEPDataParser.cxx:12
 AliHEPDataParser.cxx:13
 AliHEPDataParser.cxx:14
 AliHEPDataParser.cxx:15
 AliHEPDataParser.cxx:16
 AliHEPDataParser.cxx:17
 AliHEPDataParser.cxx:18
 AliHEPDataParser.cxx:19
 AliHEPDataParser.cxx:20
 AliHEPDataParser.cxx:21
 AliHEPDataParser.cxx:22
 AliHEPDataParser.cxx:23
 AliHEPDataParser.cxx:24
 AliHEPDataParser.cxx:25
 AliHEPDataParser.cxx:26
 AliHEPDataParser.cxx:27
 AliHEPDataParser.cxx:28
 AliHEPDataParser.cxx:29
 AliHEPDataParser.cxx:30
 AliHEPDataParser.cxx:31
 AliHEPDataParser.cxx:32
 AliHEPDataParser.cxx:33
 AliHEPDataParser.cxx:34
 AliHEPDataParser.cxx:35
 AliHEPDataParser.cxx:36
 AliHEPDataParser.cxx:37
 AliHEPDataParser.cxx:38
 AliHEPDataParser.cxx:39
 AliHEPDataParser.cxx:40
 AliHEPDataParser.cxx:41
 AliHEPDataParser.cxx:42
 AliHEPDataParser.cxx:43
 AliHEPDataParser.cxx:44
 AliHEPDataParser.cxx:45
 AliHEPDataParser.cxx:46
 AliHEPDataParser.cxx:47
 AliHEPDataParser.cxx:48
 AliHEPDataParser.cxx:49
 AliHEPDataParser.cxx:50
 AliHEPDataParser.cxx:51
 AliHEPDataParser.cxx:52
 AliHEPDataParser.cxx:53
 AliHEPDataParser.cxx:54
 AliHEPDataParser.cxx:55
 AliHEPDataParser.cxx:56
 AliHEPDataParser.cxx:57
 AliHEPDataParser.cxx:58
 AliHEPDataParser.cxx:59
 AliHEPDataParser.cxx:60
 AliHEPDataParser.cxx:61
 AliHEPDataParser.cxx:62
 AliHEPDataParser.cxx:63
 AliHEPDataParser.cxx:64
 AliHEPDataParser.cxx:65
 AliHEPDataParser.cxx:66
 AliHEPDataParser.cxx:67
 AliHEPDataParser.cxx:68
 AliHEPDataParser.cxx:69
 AliHEPDataParser.cxx:70
 AliHEPDataParser.cxx:71
 AliHEPDataParser.cxx:72
 AliHEPDataParser.cxx:73
 AliHEPDataParser.cxx:74
 AliHEPDataParser.cxx:75
 AliHEPDataParser.cxx:76
 AliHEPDataParser.cxx:77
 AliHEPDataParser.cxx:78
 AliHEPDataParser.cxx:79
 AliHEPDataParser.cxx:80
 AliHEPDataParser.cxx:81
 AliHEPDataParser.cxx:82
 AliHEPDataParser.cxx:83
 AliHEPDataParser.cxx:84
 AliHEPDataParser.cxx:85
 AliHEPDataParser.cxx:86
 AliHEPDataParser.cxx:87
 AliHEPDataParser.cxx:88
 AliHEPDataParser.cxx:89
 AliHEPDataParser.cxx:90
 AliHEPDataParser.cxx:91
 AliHEPDataParser.cxx:92
 AliHEPDataParser.cxx:93
 AliHEPDataParser.cxx:94
 AliHEPDataParser.cxx:95
 AliHEPDataParser.cxx:96
 AliHEPDataParser.cxx:97
 AliHEPDataParser.cxx:98
 AliHEPDataParser.cxx:99
 AliHEPDataParser.cxx:100
 AliHEPDataParser.cxx:101
 AliHEPDataParser.cxx:102
 AliHEPDataParser.cxx:103
 AliHEPDataParser.cxx:104
 AliHEPDataParser.cxx:105
 AliHEPDataParser.cxx:106
 AliHEPDataParser.cxx:107
 AliHEPDataParser.cxx:108
 AliHEPDataParser.cxx:109
 AliHEPDataParser.cxx:110
 AliHEPDataParser.cxx:111
 AliHEPDataParser.cxx:112
 AliHEPDataParser.cxx:113
 AliHEPDataParser.cxx:114
 AliHEPDataParser.cxx:115
 AliHEPDataParser.cxx:116
 AliHEPDataParser.cxx:117
 AliHEPDataParser.cxx:118
 AliHEPDataParser.cxx:119
 AliHEPDataParser.cxx:120
 AliHEPDataParser.cxx:121
 AliHEPDataParser.cxx:122
 AliHEPDataParser.cxx:123
 AliHEPDataParser.cxx:124
 AliHEPDataParser.cxx:125
 AliHEPDataParser.cxx:126
 AliHEPDataParser.cxx:127
 AliHEPDataParser.cxx:128
 AliHEPDataParser.cxx:129
 AliHEPDataParser.cxx:130
 AliHEPDataParser.cxx:131
 AliHEPDataParser.cxx:132
 AliHEPDataParser.cxx:133
 AliHEPDataParser.cxx:134
 AliHEPDataParser.cxx:135
 AliHEPDataParser.cxx:136
 AliHEPDataParser.cxx:137
 AliHEPDataParser.cxx:138
 AliHEPDataParser.cxx:139
 AliHEPDataParser.cxx:140
 AliHEPDataParser.cxx:141
 AliHEPDataParser.cxx:142
 AliHEPDataParser.cxx:143
 AliHEPDataParser.cxx:144
 AliHEPDataParser.cxx:145
 AliHEPDataParser.cxx:146
 AliHEPDataParser.cxx:147
 AliHEPDataParser.cxx:148
 AliHEPDataParser.cxx:149
 AliHEPDataParser.cxx:150
 AliHEPDataParser.cxx:151
 AliHEPDataParser.cxx:152
 AliHEPDataParser.cxx:153
 AliHEPDataParser.cxx:154
 AliHEPDataParser.cxx:155
 AliHEPDataParser.cxx:156
 AliHEPDataParser.cxx:157
 AliHEPDataParser.cxx:158
 AliHEPDataParser.cxx:159
 AliHEPDataParser.cxx:160
 AliHEPDataParser.cxx:161
 AliHEPDataParser.cxx:162
 AliHEPDataParser.cxx:163
 AliHEPDataParser.cxx:164
 AliHEPDataParser.cxx:165
 AliHEPDataParser.cxx:166
 AliHEPDataParser.cxx:167
 AliHEPDataParser.cxx:168
 AliHEPDataParser.cxx:169
 AliHEPDataParser.cxx:170
 AliHEPDataParser.cxx:171
 AliHEPDataParser.cxx:172
 AliHEPDataParser.cxx:173
 AliHEPDataParser.cxx:174
 AliHEPDataParser.cxx:175
 AliHEPDataParser.cxx:176
 AliHEPDataParser.cxx:177
 AliHEPDataParser.cxx:178
 AliHEPDataParser.cxx:179
 AliHEPDataParser.cxx:180
 AliHEPDataParser.cxx:181
 AliHEPDataParser.cxx:182
 AliHEPDataParser.cxx:183
 AliHEPDataParser.cxx:184
 AliHEPDataParser.cxx:185
 AliHEPDataParser.cxx:186
 AliHEPDataParser.cxx:187
 AliHEPDataParser.cxx:188
 AliHEPDataParser.cxx:189
 AliHEPDataParser.cxx:190
 AliHEPDataParser.cxx:191
 AliHEPDataParser.cxx:192
 AliHEPDataParser.cxx:193
 AliHEPDataParser.cxx:194
 AliHEPDataParser.cxx:195
 AliHEPDataParser.cxx:196
 AliHEPDataParser.cxx:197
 AliHEPDataParser.cxx:198
 AliHEPDataParser.cxx:199
 AliHEPDataParser.cxx:200
 AliHEPDataParser.cxx:201
 AliHEPDataParser.cxx:202
 AliHEPDataParser.cxx:203
 AliHEPDataParser.cxx:204
 AliHEPDataParser.cxx:205
 AliHEPDataParser.cxx:206
 AliHEPDataParser.cxx:207
 AliHEPDataParser.cxx:208
 AliHEPDataParser.cxx:209
 AliHEPDataParser.cxx:210
 AliHEPDataParser.cxx:211
 AliHEPDataParser.cxx:212
 AliHEPDataParser.cxx:213
 AliHEPDataParser.cxx:214
 AliHEPDataParser.cxx:215
 AliHEPDataParser.cxx:216
 AliHEPDataParser.cxx:217
 AliHEPDataParser.cxx:218
 AliHEPDataParser.cxx:219
 AliHEPDataParser.cxx:220
 AliHEPDataParser.cxx:221
 AliHEPDataParser.cxx:222
 AliHEPDataParser.cxx:223
 AliHEPDataParser.cxx:224
 AliHEPDataParser.cxx:225
 AliHEPDataParser.cxx:226
 AliHEPDataParser.cxx:227
 AliHEPDataParser.cxx:228
 AliHEPDataParser.cxx:229
 AliHEPDataParser.cxx:230
 AliHEPDataParser.cxx:231
 AliHEPDataParser.cxx:232
 AliHEPDataParser.cxx:233
 AliHEPDataParser.cxx:234
 AliHEPDataParser.cxx:235
 AliHEPDataParser.cxx:236
 AliHEPDataParser.cxx:237
 AliHEPDataParser.cxx:238
 AliHEPDataParser.cxx:239
 AliHEPDataParser.cxx:240
 AliHEPDataParser.cxx:241
 AliHEPDataParser.cxx:242
 AliHEPDataParser.cxx:243
 AliHEPDataParser.cxx:244
 AliHEPDataParser.cxx:245
 AliHEPDataParser.cxx:246
 AliHEPDataParser.cxx:247
 AliHEPDataParser.cxx:248
 AliHEPDataParser.cxx:249
 AliHEPDataParser.cxx:250
 AliHEPDataParser.cxx:251
 AliHEPDataParser.cxx:252
 AliHEPDataParser.cxx:253
 AliHEPDataParser.cxx:254
 AliHEPDataParser.cxx:255
 AliHEPDataParser.cxx:256
 AliHEPDataParser.cxx:257
 AliHEPDataParser.cxx:258
 AliHEPDataParser.cxx:259
 AliHEPDataParser.cxx:260
 AliHEPDataParser.cxx:261
 AliHEPDataParser.cxx:262
 AliHEPDataParser.cxx:263
 AliHEPDataParser.cxx:264
 AliHEPDataParser.cxx:265
 AliHEPDataParser.cxx:266
 AliHEPDataParser.cxx:267
 AliHEPDataParser.cxx:268
 AliHEPDataParser.cxx:269
 AliHEPDataParser.cxx:270
 AliHEPDataParser.cxx:271
 AliHEPDataParser.cxx:272
 AliHEPDataParser.cxx:273
 AliHEPDataParser.cxx:274
 AliHEPDataParser.cxx:275
 AliHEPDataParser.cxx:276
 AliHEPDataParser.cxx:277
 AliHEPDataParser.cxx:278
 AliHEPDataParser.cxx:279
 AliHEPDataParser.cxx:280
 AliHEPDataParser.cxx:281
 AliHEPDataParser.cxx:282
 AliHEPDataParser.cxx:283
 AliHEPDataParser.cxx:284
 AliHEPDataParser.cxx:285
 AliHEPDataParser.cxx:286
 AliHEPDataParser.cxx:287
 AliHEPDataParser.cxx:288
 AliHEPDataParser.cxx:289
 AliHEPDataParser.cxx:290
 AliHEPDataParser.cxx:291
 AliHEPDataParser.cxx:292
 AliHEPDataParser.cxx:293
 AliHEPDataParser.cxx:294
 AliHEPDataParser.cxx:295
 AliHEPDataParser.cxx:296
 AliHEPDataParser.cxx:297
 AliHEPDataParser.cxx:298
 AliHEPDataParser.cxx:299
 AliHEPDataParser.cxx:300
 AliHEPDataParser.cxx:301
 AliHEPDataParser.cxx:302
 AliHEPDataParser.cxx:303
 AliHEPDataParser.cxx:304
 AliHEPDataParser.cxx:305
 AliHEPDataParser.cxx:306
 AliHEPDataParser.cxx:307
 AliHEPDataParser.cxx:308
 AliHEPDataParser.cxx:309
 AliHEPDataParser.cxx:310
 AliHEPDataParser.cxx:311