]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STARLIGHT/starlight/TStarLight/testsl.C
104ff8f9e61ecbb877b878effbe65d3190e42b13
[u/mrichter/AliRoot.git] / STARLIGHT / starlight / TStarLight / testsl.C
1 // -*- C++ -*-
2 // $Id$
3
4 // creates histogram of generated rho0 kinematical data
5 // writes a text file "rho0.txt" which should be identical to ../../test/rho0.txt
6
7 void testsl() {
8   gSystem->Load("libStarLight");
9   gSystem->Load("libAliStarLight.so");
10
11   TStarLight* sl = new TStarLight("starlight generator", "title", "");
12
13   sl->SetParameter("BEAM_1_Z     =    1    #Z of projectile");
14   sl->SetParameter("BEAM_1_A     =    1    #A of projectile");
15   sl->SetParameter("BEAM_2_Z     =   82    #Z of target");
16   sl->SetParameter("BEAM_2_A     =  208    #A of target");
17   sl->SetParameter("BEAM_1_GAMMA = 4264.4  #Gamma of the colliding ions");
18   sl->SetParameter("BEAM_2_GAMMA = 1682.4  #Gamma of the colliding ions");
19   sl->SetParameter("W_MAX        =   15    #Max value of w");
20   sl->SetParameter("W_MIN        =    -1   #Min value of w");
21   sl->SetParameter("W_N_BINS     =   40    #Bins i w");
22   sl->SetParameter("RAP_MAX      =    8.   #max y");
23   sl->SetParameter("RAP_N_BINS   =   80    #Bins i y");
24   sl->SetParameter("CUT_PT       =    0    #Cut in pT? 0 = (no, 1 = yes)");
25   sl->SetParameter("PT_MIN       =    1.0  #Minimum pT in GeV");
26   sl->SetParameter("PT_MAX       =    3.0  #Maximum pT in GeV");
27   sl->SetParameter("CUT_ETA      =    0    #Cut in pseudorapidity? (0 = no, 1 = yes)");
28   sl->SetParameter("ETA_MIN      =  -10    #Minimum pseudorapidity");
29   sl->SetParameter("ETA_MAX      =   10    #Maximum pseudorapidity");
30   sl->SetParameter("PROD_MODE    =    2    #gg or gP switch (1 = 2-photon, 2 = coherent vector meson (narrow), 3 = coherent vector meson (wide), # 4 = incoherent vector meson, 5 = A+A DPMJet single, 6 = A+A DPMJet double, 7 = p+A DPMJet single, 8 = p+A Pythia single )");
31   sl->SetParameter("PROD_PID     =  113    #Channel of interest (not relevant for photonuclear processes)");
32   sl->SetParameter("RND_SEED     = 34533   #Random number seed");
33   sl->SetParameter("BREAKUP_MODE  =   5    #Controls the nuclear breakup");
34   sl->SetParameter("INTERFERENCE  =   0    #Interference (0 = off, 1 = on)");
35   sl->SetParameter("IF_STRENGTH   =   1.   #% of intefernce (0.0 - 0.1)");
36   sl->SetParameter("COHERENT      =   1    #Coherent=1,Incoherent=0");
37   sl->SetParameter("INCO_FACTOR   =   1.   #percentage of incoherence");
38   sl->SetParameter("INT_PT_MAX    =   0.24 #Maximum pt considered, when interference is turned on");
39   sl->SetParameter("INT_PT_N_BINS = 120    #Number of pt bins when interference is turned on");
40
41   sl->InitStarLight();
42   sl->PrintInputs(std::cout);
43   TClonesArray tca("TParticle", 1000);
44
45   TLorentzVector v[2], vSum;
46   TH1* hM  = new TH1D("hM",  "STARLIGHT;M#(){#pi^{+}#pi^{-}}",     100, 0., 2.);
47   TH1* hPt = new TH1D("hPt", "STARLIGHT;P_{T}#(){#pi^{+}#pi^{-}}", 100, 0., 1.);
48   TH1* hY  = new TH1D("hY",  "STARLIGHT;Y#(){#pi^{+}#pi^{-}}",     160,-8., 8.);
49
50   std::ofstream ofs("rho0.txt");
51   TParticle *p;
52   for (Int_t counter(0); counter<100;) {
53     sl->GenerateEvent();    
54     sl->BoostEvent();    
55     sl->ImportParticles(&tca, "ALL");
56     Bool_t isOK = kTRUE;
57     for (Int_t i=0; i<tca.GetEntries() && isOK; ++i) {
58       p = (TParticle*)tca.At(i);
59       p->Momentum(v[i]);
60 //       isOK = TMath::Abs(v[i].Rapidity()) <= 0.9;
61     }
62     tca.Clear();
63     if (!isOK) continue;
64     Printf("counter, %d", counter, isOK);
65     ++counter;
66     vSum = v[0] + v[1];
67     ofs << std::fixed << std::setprecision(4)
68         << vSum.M() << " " << vSum.Perp() << " " << vSum.Rapidity() << " "
69         << v[0].Eta() << " " << v[0].Px() << " " << v[0].Py() << " " << v[0].Pz() << " "
70         << v[1].Eta() << " " << v[1].Px() << " " << v[1].Py() << " " << v[1].Pz()
71         << std::endl;
72     hM->Fill(vSum.M());
73     hPt->Fill(vSum.Perp());
74     hY->Fill(vSum.Rapidity());
75   }
76   TFile::Open("sl.root", "RECREATE");
77   sl->Write();
78   gFile->Write();
79
80   hM->Draw();
81   c1->SaveAs("SL.pdf(");
82   hPt->Draw();
83   c1->SaveAs("SL.pdf");
84   hY->Draw();
85   c1->SaveAs("SL.pdf)");
86 }