StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
starsim.hijing.D0.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 
7 
8 class St_geant_Maker;
9 St_geant_Maker *geant_maker = 0;
10 
11 class StarGenEvent;
12 StarGenEvent *event = 0;
13 
14 class StarPrimaryMaker;
15 StarPrimaryMaker *_primary = 0;
16 
17 class StarKinematics;
19 
20 Float_t minPt = +0.5;
21 Float_t maxPt = +5.0;
22 Float_t minEta = -1.5;
23 Float_t maxEta = +1.5;
24 
25 // ----------------------------------------------------------------------------
26 void geometry( TString tag, Bool_t agml=true )
27 {
28  TString cmd = "DETP GEOM "; cmd += tag;
29  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
30  geant_maker -> LoadGeometry(cmd);
31  // if ( agml ) command("gexec $STAR_LIB/libxgeometry.so");
32 }
33 // ----------------------------------------------------------------------------
34 void command( TString cmd )
35 {
36  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
37  geant_maker -> Do( cmd );
38 }
39 // ----------------------------------------------------------------------------
40 void trig( Int_t n=0 )
41 {
42  for ( Int_t i=0; i<n+1; i++ ) {
43  chain->Clear();
44  // Throw one hypertriton flat in -1 to 1 w/ momentum btwn 0.1 and 1 GeV
45  if (kinematics) kinematics->Kine( 1, "D0", minPt, maxPt, minEta, maxEta );
46  chain->Make();
47  // command("gprint kine");
48  }
49 }
50 // ----------------------------------------------------------------------------
51 void HyperTritons()
52 {
53  kinematics = new StarKinematics("D0");
54  _primary -> AddGenerator(kinematics);
55 }
56 // ----------------------------------------------------------------------------
57 void Hijing()
58 {
59  StarHijing *hijing = new StarHijing("hijing");
60  hijing->SetTitle("Hijing 1.383");
61 
62  // Setup collision frame, energy and beam species
63  hijing->SetFrame("CMS",200.0);
64  hijing->SetBlue("Au");
65  hijing->SetYell("Au");
66  hijing->SetImpact(0.0, 30.0); // Impact parameter min/max (fm) 0. 30.
67  HiParnt_t &hiparnt = hijing->hiparnt();
68  {
69  hijing->hiparnt().ihpr2(4) = 0; // Jet quenching (1=yes/0=no) 0
70  hijing->hiparnt().ihpr2(3) = 0; // Hard scattering (1=yes/0=no)
71  hijing->hiparnt().hipr1(10) = 2.0; // pT jet
72  hijing->hiparnt().ihpr2(8) = 10; // Max number of jets / nucleon
73  hijing->hiparnt().ihpr2(11) = 1; // Set baryon production
74  hijing->hiparnt().ihpr2(12) = 1; // Turn on/off decay of particles [1=recommended]
75  hijing->hiparnt().ihpr2(18) = 0; // 0=Charm, 1=Bottom production
76  hijing->hiparnt().hipr1(7) = 5.35; // Set B production ???? Not really used... Really ????
77  };
78 
79 
80  hijing->hiparnt().ihpr2(50) = 12345; // unused entry
81 
82  // For more configuration options, see the HIJING manual
83  // http://ntc0.lbl.gov/~xnwang/hijing/doc.html
84 
85  _primary -> AddGenerator(hijing);
86  _primary -> SetCuts( 1.0E-6 , -1., -2.5, +2.5 );
87 
88 }
89 // ----------------------------------------------------------------------------
90 // ----------------------------------------------------------------------------
91 // ----------------------------------------------------------------------------
92 void starsim( Int_t nevents=1, Int_t rngSeed=4321 )
93 {
94 
95  gROOT->ProcessLine(".L bfc.C");
96  {
97  TString simple = "y2014a geant gstar usexgeom agml ";
98  bfc(0, simple );
99  }
100 
101  gSystem->Load( "libVMC.so");
102 
103  gSystem->Load( "StarGeneratorUtil.so" );
104  gSystem->Load( "StarGeneratorEvent.so" );
105  gSystem->Load( "StarGeneratorBase.so" );
106  gSystem->Load( "libMathMore.so" );
107  gSystem->Load( "libHijing1_383.so");
108  gSystem->Load( "libKinematics.so");
109  gSystem->Load( "xgeometry.so" );
110 
111  // force gstar load/call
112  gSystem->Load( "gstar.so" );
113  command("call gstar");
114 
115  // Setup RNG seed and map all ROOT TRandom here
116  StarRandom::seed( rngSeed );
118 
119  //
120  // Create the primary event generator and insert it
121  // before the geant maker
122  //
123  _primary = new StarPrimaryMaker();
124  {
125  _primary -> SetFileName( "hijing.starsim.root");
126  chain -> AddBefore( "geant", _primary );
127  }
128 
129  //
130  // These should be adjusted to your best vertex estimates
131  //
132  _primary -> SetVertex( 0., 0., 0. );
133  _primary -> SetSigma( 0.3, 0.3, 60.0 );
134 
135 
136 
137 
138  // Setup an event generator
139  //
140  Hijing();
141  //
142  // Setup single hypertritons
143  //
144  HyperTritons();
145 
146 
147  //
148  // Initialize primary event generator and all sub makers
149  //
150  _primary -> Init();
151 
152 
153  return;
154 
155  //
156  // Setup geometry and set starsim to use agusread for input
157  //
158  //geometry("y2014a");
159  command("gkine -4 0");
160  command("gfile o hijing.starsim.fzd");
161 
162  //
163  // Trigger on nevents
164  //
165  trig( nevents );
166 
167  _primary->event()->Print();
168 
169  // command("gprint kine");
170 
171  command("call agexit"); // Make sure that STARSIM exits properly
172 
173 }
174 // ----------------------------------------------------------------------------
175 
HIJING /HIPARNT/ common block/.
Definition: Hijing.h:79
void SetFrame(const Char_t *frame, const Double_t val)
void Print(const Option_t *opts="head") const
Star Simple Kinematics Generator.
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 SetBlue(const Char_t *b)
Sets the particle species for the blue beam.
int & ihpr2(int i)
Definition: Hijing.h:94
Int_t Init()
Initialize generator.
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
StarGenEvent * event()
Return a pointer to the event.
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
Sparse class to hold track kinematics.
void Kine(Int_t ntrack, const Char_t *type="pi+,pi-,K+,K-,proton,antiproton", Double_t ptlow=0.0, Double_t pthigh=500.0, Double_t ylow=-10.0, Double_t yhigh=+10.0, Double_t philow=0.0, Double_t phihigh=TMath::TwoPi())