StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
lambda.C
1 //usr/bin/env root4star -l -b -q $0'('$1', '$2')'; exit $?
2 // macro to instantiate the Geant3 from within
3 // STAR C++ framework and get the starsim prompt
4 // To use it do
5 // root4star starsim.C
6 
7 class St_geant_Maker;
8 St_geant_Maker *geant_maker = 0;
9 
10 class StarGenEvent;
11 StarGenEvent *event = 0;
12 
13 class StarPrimaryMaker;
14 StarPrimaryMaker *_primary = 0;
15 
16 class StarKinematics;
18 
19 
20 TH1F* hMll = 0;
21 float numParticles = 1;
22 
23 // ----------------------------------------------------------------------------
24 void geometry( TString tag, Bool_t agml=true )
25 {
26  TString cmd = "DETP GEOM "; cmd += tag + " field=-5.0";
27  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
28  geant_maker -> LoadGeometry(cmd);
29  // if ( agml ) command("gexec $STAR_LIB/libxgeometry.so");
30 }
31 // ----------------------------------------------------------------------------
32 void command( TString cmd )
33 {
34  if ( !geant_maker ) geant_maker = (St_geant_Maker *)chain->GetMaker("geant");
35  geant_maker -> Do( cmd );
36 }
37 // ----------------------------------------------------------------------------
38 void trig( Int_t n=1 )
39 {
40 
41 
42  for ( Int_t i=0; i<n; i++ ) {
43 
44  // Clear the chain from the previous event
45  chain->Clear();
46 
47  //(Momentum, Energy units are Gev/C, GeV)
48  Double_t masses[2] = { 0.13957, 0.938} ;
49 
50  TGenPhaseSpace genEvent;
51  TLorentzVector W;
52  // W.SetPtEtaPhiM( 0.0, 100.0, 0, 3.096 );
53  W.SetXYZM( 0, 0, 5, 1.11568 );
54  genEvent.SetDecay(W, 2, masses);
55 
56  TLorentzVector lv;
57  for ( int j = 0; j < numParticles; j++ ){
58  Double_t weight = genEvent.Generate();
59  TLorentzVector *pPion = genEvent.GetDecay(0);
60  TLorentzVector *pProton = genEvent.GetDecay(1);
61  lv = *pPion + *pProton;
62 
63  StarGenParticle *pion;
64  pion = kinematics->AddParticle( "pi-" );
65 
66  pion->SetPx(pPion->Px());
67  pion->SetPy(pPion->Py());
68  pion->SetPz(pPion->Pz());
69  pion->SetMass( masses[0] );
70 
71  StarGenParticle *proton;
72  proton = kinematics->AddParticle( "p" );
73 
74 
75  proton->SetPx(pProton->Px());
76  proton->SetPy(pProton->Py());
77  proton->SetPz(pProton->Pz());
78  proton->SetMass( masses[1] );
79 
80  hMll->Fill( lv.M() );
81 
82  cout << "pion eta = " << pPion->Eta() << endl;
83  cout << "proton eta = " << pProton->Eta() << endl;
84  }
85 
86 
87  // kinematics->Kine( numParticles, nameParticle.Data(), 10.2, 12.0, 2.5, 4.00 );
88 
89  // Generate the event
90  chain->Make();
91 
92  // TTable* hits = chain->GetDataSet("bfc/.make/geant/.data/g2t_stg_hit");
93  // if ( hits ) {
94  // double nhits = hits->GetNRows();
95  // hNumHits->Fill( double(i), nhits / 4.0 / numParticles );
96  // std::cout << "N hits = " << nhits << std::endl;
97  // }
98 
99  // Print the event
100  // command("gprint hits stgh");
101 
102  }
103 }
104 // ----------------------------------------------------------------------------
105 // ----------------------------------------------------------------------------
106 // ----------------------------------------------------------------------------
107 void Kinematics()
108 {
109 
110  // gSystem->Load( "libStarGeneratorPoolPythia6_4_23.so" );
111  gSystem->Load( "libKinematics.so");
112  kinematics = new StarKinematics();
113 
114  _primary->AddGenerator(kinematics);
115 }
116 // ----------------------------------------------------------------------------
117 // ----------------------------------------------------------------------------
118 // ----------------------------------------------------------------------------
119 void lambda( Int_t nevents=100, Int_t rngSeed=12352342 )
120 {
121  hMll = new TH1F("hMll",";Mll;counts [10MeV]", 200, 2.0, 4.0 );
122  cout << "Generating: " << nevents << " events with seed: " << rngSeed << endl;
123  cout << "Simulating J/psi->e+e-" << endl;
124 
125  gSystem->Load( "libStarRoot.so" );
126  gROOT->SetMacroPath(".:/star-sw/StRoot/macros/:./StRoot/macros:./StRoot/macros/graphics:./StRoot/macros/analysis:./StRoot/macros/test:./StRoot/macros/examples:./StRoot/macros/html:./StRoot/macros/qa:./StRoot/macros/calib:./StRoot/macros/mudst:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/graphics:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/analysis:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/test:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/examples:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/html:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/qa:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/calib:/afs/rhic.bnl.gov/star/packages/DEV/StRoot/macros/mudst:/afs/rhic.bnl.gov/star/ROOT/36/5.34.38/.sl73_x8664_gcc485/rootdeb/macros:/afs/rhic.bnl.gov/star/ROOT/36/5.34.38/.sl73_x8664_gcc485/rootdeb/tutorials");
127 
128  gROOT->ProcessLine(".L bfc.C");
129  {
130  TString simple = "sdt20211016 y2023 geant gstar usexgeom agml ";
131  bfc(0, simple );
132  }
133 
134  gSystem->Load( "libVMC.so");
135 
136  gSystem->Load( "StarGeneratorUtil.so" );
137  gSystem->Load( "StarGeneratorEvent.so" );
138  gSystem->Load( "StarGeneratorBase.so" );
139 
140  gSystem->Load( "libMathMore.so" );
141  gSystem->Load( "xgeometry.so" );
142 
143  // Setup RNG seed and map all ROOT TRandom here
144  StarRandom::seed( rngSeed );
146 
147  //
148  // Create the primary event generator and insert it
149  // before the geant maker
150  //
151  // StarPrimaryMaker *
152  _primary = new StarPrimaryMaker();
153  {
154  _primary -> SetFileName( "lambda_fwd_gun.root");
155  chain -> AddBefore( "geant", _primary );
156  }
157 
158  Kinematics();
159 
160  //
161  // Initialize primary event generator and all sub makers
162  //
163  _primary -> Init();
164  _primary->SetSigma( 0.1, 0.1, 0.1 ); // 1mm x 1mm x 1mm smearing at the vertex
165  _primary->SetVertex(0.0, 0.0, 0.0 );
166 
167  //
168  // Setup geometry and set starsim to use agusread for input
169  //
170  //geometry("y2012");
171  command("gkine -4 0");
172  command("gfile o lambda_fwd_gun.fzd");
173 
174 
175  hNumHits = new TH1F("hNumEvents","Nhits/plane/incident track vs event number",nevents + 1, -0.5, (float)( nevents ) + 0.5 );
176 
177  //
178  // Trigger on nevents
179  //
180  trig( nevents );
181 
182  TFile * f = new TFile( "lambda_gen.root", "RECREATE" );
183  f->cd();
184  hMll->Write();
185  f->Write();
186 
187  command("call agexit"); // Make sure that STARSIM exits properly
188 
189 }
190 // ----------------------------------------------------------------------------
191 
StarGenParticle * AddParticle()
void SetPx(Float_t px)
Set the x-component of the momentum.
void SetSigma(Double_t sx, Double_t sy, Double_t sz, Double_t rho=0)
Star Simple Kinematics Generator.
Yet another particle class.
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
void AddGenerator(StarGenerator *gener)
void SetPz(Float_t pz)
Set the z-component of the momentum.
Int_t Init()
Initialize generator.
void SetMass(Float_t mass)
Set the mass.
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 SetPy(Float_t py)
Set the y-component of the momentum.
static void capture()
Capture gRandom random number generator.
Definition: StarRandom.cxx:57
Sparse class to hold track kinematics.
void SetVertex(Double_t x, Double_t y, Double_t z)
Set the x, y and z vertex position.