StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
starsim.hijing.dalitz.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 // Example of how to add a single pi0 decaying 100% through dalitz decay (e+ e- gamm) to a standard
7 // hijing event.
8 
9 
10 
11 
12 class St_geant_Maker;
13 St_geant_Maker *geant_maker = 0;
14 
15 class StarGenEvent;
16 StarGenEvent *event = 0;
17 
18 class StarPrimaryMaker;
19 StarPrimaryMaker *_primary = 0;
20 
21 class StarKinematics;
23 
24 Float_t minPt = +0.5;
25 Float_t maxPt = +5.0;
26 Float_t minEta = -1.5;
27 Float_t maxEta = +1.5;
28 
29 // ----------------------------------------------------------------------------
30 void geometry( TString tag, Bool_t agml=true )
31 {
32  TString cmd = "DETP GEOM "; cmd += tag;
33  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
34  geant_maker -> LoadGeometry(cmd);
35  // if ( agml ) command("gexec $STAR_LIB/libxgeometry.so");
36 }
37 // ----------------------------------------------------------------------------
38 void command( TString cmd )
39 {
40  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
41  geant_maker -> Do( cmd );
42 }
43 // ----------------------------------------------------------------------------
44 void trig( Int_t n=0 )
45 {
46  for ( Int_t i=0; i<n+1; i++ ) {
47  chain->Clear();
48  // Throw one hypertriton flat in -1 to 1 w/ momentum btwn 0.1 and 1 GeV
49  if (kinematics) kinematics->Kine( 1, "Dalitz", minPt, maxPt, minEta, maxEta );
50  chain->Make();
51  // command("gprint kine");
52  }
53 }
54 // ----------------------------------------------------------------------------
55 void Dalitzdecay()
56 {
57 
58  // The Dalitz particle is defined in starsim in gstar_part.g --
59  //
60  // Particle Dalitz code=10007 TrkTyp=4 mass=0.135 charge=0 tlife=8.4e-17, |
61  // pdg=100111 bratio= { 1.0,} mode= {10203,} |
62  //
63  // The particle database does not know about this particle, so we have to add it. It is
64  // important that we do not overwrite PDG id = 111, the ID of the standard pi0. Otherwise,
65  // we will have all pi0's in the event decaying by dalitz.
66  //
68  data.AddParticle("Dalitz","pi0-->e+e-gamma 100%",0.135,0,0,0,"meson",100111,-100111,10007);
69 
70  kinematics = new StarKinematics("Dalitz Decay");
71  _primary -> AddGenerator(kinematics);
72 }
73 // ----------------------------------------------------------------------------
74 void Hijing()
75 {
76  StarHijing *hijing = new StarHijing("hijing");
77  hijing->SetTitle("Hijing 1.383");
78 
79  // Setup collision frame, energy and beam species
80  hijing->SetFrame("CMS",200.0);
81  hijing->SetBlue("Au");
82  hijing->SetYell("Au");
83  hijing->SetImpact(0.0, 30.0); // Impact parameter min/max (fm) 0. 30.
84  HiParnt_t &hiparnt = hijing->hiparnt();
85  {
86  hijing->hiparnt().ihpr2(4) = 0; // Jet quenching (1=yes/0=no) 0
87  hijing->hiparnt().ihpr2(3) = 0; // Hard scattering (1=yes/0=no)
88  hijing->hiparnt().hipr1(10) = 2.0; // pT jet
89  hijing->hiparnt().ihpr2(8) = 10; // Max number of jets / nucleon
90  hijing->hiparnt().ihpr2(11) = 1; // Set baryon production
91  hijing->hiparnt().ihpr2(12) = 1; // Turn on/off decay of particles [1=recommended]
92  hijing->hiparnt().ihpr2(18) = 0; // 0=Charm, 1=Bottom production
93  hijing->hiparnt().hipr1(7) = 5.35; // Set B production ???? Not really used... Really ????
94  };
95 
96 
97  hijing->hiparnt().ihpr2(50) = 12345; // unused entry
98 
99  // For more configuration options, see the HIJING manual
100  // http://ntc0.lbl.gov/~xnwang/hijing/doc.html
101 
102  _primary -> AddGenerator(hijing);
103  _primary -> SetCuts( 1.0E-6 , -1., -2.5, +2.5 );
104 
105 }
106 // ----------------------------------------------------------------------------
107 // ----------------------------------------------------------------------------
108 // ----------------------------------------------------------------------------
109 void starsim( Int_t nevents=10, Int_t rngSeed=4321 )
110 {
111 
112  gROOT->ProcessLine(".L bfc.C");
113  {
114  TString simple = "y2014a geant gstar usexgeom agml ";
115  bfc(0, simple );
116  }
117 
118  gSystem->Load( "libVMC.so");
119 
120  gSystem->Load( "StarGeneratorUtil.so" );
121  gSystem->Load( "StarGeneratorEvent.so" );
122  gSystem->Load( "StarGeneratorBase.so" );
123  gSystem->Load( "libMathMore.so" );
124  gSystem->Load( "libHijing1_383.so");
125  gSystem->Load( "libKinematics.so");
126  gSystem->Load( "xgeometry.so" );
127 
128  // force gstar load/call
129  gSystem->Load( "gstar.so" );
130  command("call gstar");
131 
132  // Setup RNG seed and map all ROOT TRandom here
133  StarRandom::seed( rngSeed );
135 
136  //
137  // Create the primary event generator and insert it
138  // before the geant maker
139  //
140  _primary = new StarPrimaryMaker();
141  {
142  _primary -> SetFileName( "hijing.starsim.root");
143  chain -> AddBefore( "geant", _primary );
144  }
145 
146  //
147  // These should be adjusted to your best vertex estimates
148  //
149  _primary -> SetVertex( 0., 0., 0. );
150  _primary -> SetSigma( 0.3, 0.3, 60.0 );
151 
152 
153 
154 
155  // Setup an event generator
156  //
157  Hijing();
158  //
159  // Setup single dalitzdecay
160  //
161  Dalitzdecay();
162 
163 
164  //
165  // Initialize primary event generator and all sub makers
166  //
167  _primary -> Init();
168 
169 
170  //
171  // Setup geometry and set starsim to use agusread for input
172  //
173 
174  command("gkine -4 0");
175  command("gfile o hijing.starsim.fzd");
176 
177  //
178  // Trigger on nevents
179  //
180  trig( nevents );
181 
182  _primary->event()->Print();
183 
184  // command("gprint kine");
185 
186  command("call agexit"); // Make sure that STARSIM exits properly
187 
188 }
189 // ----------------------------------------------------------------------------
190 
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
static StarParticleData & instance()
Returns a reference to the single instance of this class.
Interface to PDG information.
void SetBlue(const Char_t *b)
Sets the particle species for the blue beam.
int & ihpr2(int i)
Definition: Hijing.h:94
void AddParticle(const Char_t *name, TParticlePDG *particle)
Add a particle to the database.
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())