91 #define FILENAME "pythia.root"
92 #define TREENAME "tree"
93 #define BRANCHNAME "particles"
94 #define HISTNAME "ptSpectra"
102 gSystem->Load(
"libpythia6_dummy");
107 int makeEventSample(Int_t nEvents)
113 TPythia6* pythia =
new TPythia6;
116 pythia->Initialize(
"cms",
"p",
"p", 200);
119 TFile* file = TFile::Open(FILENAME,
"RECREATE");
120 if (!file || !file->IsOpen()) {
121 Error(
"makeEventSample",
"Couldn;t open file %s", FILENAME);
126 TTree* tree =
new TTree(TREENAME,
"Pythia 6 tree");
130 TClonesArray* particles = (TClonesArray*)pythia->GetListOfParticles();
131 tree->Branch(BRANCHNAME, &particles);
134 for (Int_t i = 0; i < nEvents; i++) {
137 cout <<
"Event # " << i << endl;
140 pythia->GenerateEvent();
154 TH1D* hist =
new TH1D(HISTNAME,
"p_{#perp} spectrum for #pi^{+}",
156 hist->SetXTitle(
"p_{#perp}");
157 hist->SetYTitle(
"dN/dp_{#perp}");
159 sprintf(expression,
"sqrt(pow(%s.fPx,2)+pow(%s.fPy,2))>>%s",
160 BRANCHNAME, BRANCHNAME, HISTNAME);
162 sprintf(selection,
"%s.fKF==%d", BRANCHNAME, PDGNUMBER);
163 tree->Draw(expression,selection);
167 hist->Scale(3 / 100. / hist->Integral());
168 hist->Fit(
"expo",
"QO+",
"", .25, 1.75);
169 TF1* func = hist->GetFunction(
"expo");
170 func->SetParNames(
"A",
"-1/T");
179 int showEventSample()
185 TFile* file = TFile::Open(FILENAME,
"READ");
186 if (!file || !file->IsOpen()) {
187 Error(
"showEventSample",
"Couldn;t open file %s", FILENAME);
192 TTree* tree = (TTree*)file->Get(TREENAME);
194 Error(
"showEventSample",
"couldn't get TTree %s", TREENAME);
202 TH1D* hist = (TH1D*)file->Get(HISTNAME);
204 Error(
"showEventSample",
"couldn't get TH1D %s", HISTNAME);
209 gStyle->SetOptStat(1);
210 TCanvas* canvas =
new TCanvas(
"canvas",
"canvas");
213 TF1* func = hist->GetFunction(
"expo");
216 sprintf(expression,
"T #approx %5.1f", -1000 / func->GetParameter(1));
217 TLatex* latex =
new TLatex(1.5, 1e-4, expression);
218 latex->SetTextSize(.1);
219 latex->SetTextColor(4);
225 void pythiaExample(Int_t n=1000) {
231 int main(
int argc,
char** argv)
233 TApplication app(
"app", &argc, argv);
237 n = strtol(argv[1], NULL, 0);
241 retVal = makeEventSample(n);
243 retVal = showEventSample();