StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
starsim.hijing.pHe3.C
1 // macro to instantiate the Geant3 from within
2 // STAR C++ framework and get the starsim prompt
3 // To use it do
4 // root4star starsim.C
5 
6 class St_geant_Maker;
7 St_geant_Maker *geant_maker = 0;
8 
9 class StarGenEvent;
10 StarGenEvent *event = 0;
11 
12 class StarPrimaryMaker;
13 StarPrimaryMaker *_primary = 0;
14 
15 // ----------------------------------------------------------------------------
16 void geometry( TString tag, Bool_t agml=true )
17 {
18  TString cmd = "DETP GEOM "; cmd += tag;
19  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
20  geant_maker -> LoadGeometry(cmd);
21  // if ( agml ) command("gexec $STAR_LIB/libxgeometry.so");
22 }
23 // ----------------------------------------------------------------------------
24 void command( TString cmd )
25 {
26  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
27  geant_maker -> Do( cmd );
28 }
29 // ----------------------------------------------------------------------------
30 void trig( Int_t n=0 )
31 {
32  for ( Int_t i=0; i<n+1; i++ ) {
33  chain->Clear();
34  chain->Make();
35  }
36 }
37 // ----------------------------------------------------------------------------
38 void Hijing()
39 {
40  StarHijing *hijing = new StarHijing("hijing");
41  hijing->SetTitle("Hijing 1.383");
42 
43  // Setup collision frame, energy and beam species
44  hijing->SetFrame("CMS",200.0);
45  hijing->SetBlue("p");
46  hijing->SetYell("He3");
47  hijing->SetImpact(0.0, 30.0); // Impact parameter min/max (fm) 0. 30.
48  HiParnt_t &hiparnt = hijing->hiparnt();
49  {
50  hiparnt.ihpr2(4) = 0; // Jet quenching (1=yes/0=no) 0
51  hiparnt.ihpr2(3) = 0; // Hard scattering (1=yes/0=no)
52  hiparnt.hipr1(10) = 2.0; // pT jet
53  hiparnt.ihpr2(8) = 10; // Max number of jets / nucleon
54  hiparnt.ihpr2(11) = 1; // Set baryon production
55  hiparnt.ihpr2(12) = 1; // Turn on/off decay of particles [1=recommended]
56  hiparnt.ihpr2(18) = 0; // 0=charm production, 1=bottom production
57  hiparnt.hipr1(7) = 5.35; // Set B production ???? Not really used... Really ????
58  }
59 
60  // For more configuration options, see the HIJING manual
61  // http://ntc0.lbl.gov/~xnwang/hijing/doc.html
62 
63  _primary -> AddGenerator(hijing);
64  _primary -> SetCuts( 1.0E-6 , -1., -2.5, +2.5 );
65 
66 }
67 // ----------------------------------------------------------------------------
68 // ----------------------------------------------------------------------------
69 // ----------------------------------------------------------------------------
70 void starsim( Int_t nevents=50, Int_t rngSeed=1234 )
71 {
72 
73  gROOT->ProcessLine(".L bfc.C");
74  {
75  TString simple = "y2012 geant gstar usexgeom agml ";
76  //TString full = "tpcrs TpxRaw y2010a MakeEvent ITTF NoSvtIt NoSsdIt Idst IAna l0 ftpc Sti Tree logger genvtx tpcDB TpcHitMover TpxClu pmd bbcSim tofsim tags emcY2 EEfs evout IdTruth geantout -dstout big fzin MiniMcMk clearmem";
77  // TString full = "y2012 geant gstar tpcrs genvtx tpcDb tpxclu dedx event sdt20120224 ";
78  bfc(0, simple );
79  }
80 
81  gSystem->Load( "libVMC.so");
82 
83  gSystem->Load( "StarGeneratorUtil.so" );
84  gSystem->Load( "StarGeneratorEvent.so" );
85  gSystem->Load( "StarGeneratorBase.so" );
86  gSystem->Load( "libMathMore.so" );
87  gSystem->Load( "libHijing1_383.so");
88  gSystem->Load( "xgeometry.so" );
89 
90  // Setup RNG seed and map all ROOT TRandom here
91  StarRandom::seed( rngSeed );
93 
94  //
95  // Create the primary event generator and insert it
96  // before the geant maker
97  //
98  _primary = new StarPrimaryMaker();
99  {
100  _primary -> SetFileName( "hijing.starsim.root");
101  chain -> AddBefore( "geant", _primary );
102  }
103 
104 
105  //
106  // Setup an event generator
107  //
108  Hijing();
109 
110  //
111  // Initialize primary event generator and all sub makers
112  //
113  _primary -> Init();
114 
115  //
116  // Setup geometry and set starsim to use agusread for input
117  //
118  //geometry("y2012");
119  command("gkine -4 0");
120  command("gfile o hijing.starsim.fzd");
121 
122  //
123  // Trigger on nevents
124  //
125  trig( nevents );
126  command("gprint kine");
127 
128  // command("call agexit"); // Make sure that STARSIM exits properly
129 
130 }
131 // ----------------------------------------------------------------------------
132 
HIJING /HIPARNT/ common block/.
Definition: Hijing.h:79
void SetFrame(const Char_t *frame, const Double_t val)
void SetFileName(const Char_t *name)
Set the filename of the output TTree.
void SetImpact(Double_t bmin, Double_t bmax)
Set the minimum and maximum impact parameters for the collision (HI generators)
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
void SetCuts(Double_t ptmin, Double_t ptmax=-1, Double_t ymin=0, Double_t ymax=-1, Double_t phimin=0, Double_t phimax=-1, Double_t zmin=-999, Double_t zmax=+999)
Set particle cuts.
void AddGenerator(StarGenerator *gener)
void SetBlue(const Char_t *b)
Sets the particle species for the blue beam.
int & ihpr2(int i)
Definition: Hijing.h:94
HiParnt_t & hiparnt()
Returns a reference to the hijing parameters.
Definition: StarHijing.h:58
float & hipr1(int i)
Definition: Hijing.h:90
Interface to the HIJING event generator.
Definition: StarHijing.h:48
virtual Int_t Make()
Definition: StChain.cxx:110
static void seed(UInt_t s)
Definition: StarRandom.cxx:119
Base class for event records.
Definition: StarGenEvent.h:81
Main steering class for event generation.
void SetYell(const Char_t *y)
Sets the particle species for the yellow beam.
static void capture()
Capture gRandom random number generator.
Definition: StarRandom.cxx:57