StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarPrimaryMaker.h
1 #ifndef __StarPrimaryMaker_h__
2 #define __StarPrimaryMaker_h__
3 
4 #include "StMaker.h"
5 #include "TLorentzVector.h"
6 #include "TVector3.h"
7 #include "TClonesArray.h"
8 #include "TDatabasePDG.h"
9 #include "StarParticleStack.h"
10 #include "StMessMgr.h"
11 
12 #include "StarGenerator/UTIL/StarParticleData.h"
13 
14 #include "TEventList.h"
15 
16 #ifndef __CINT__
17 #include <functional>
18 #include <map>
19 #endif
20 
21 class StarGenerator;
22 class StarGenEvent;
23 class StarFilterMaker;
24 
25 class TTree;
26 class TFile;
27 
103 class StarPrimaryMaker : public StMaker
104 {
105  public:
107  virtual ~StarPrimaryMaker();
108 
109  Int_t InitRun( Int_t runnumber );
110  Int_t Init();
111  Int_t Make();
112  void Clear( const Option_t *opts="" );
113  Int_t Finish();
114 
116  void SetFileName( const Char_t *name ){ mFileName = name; }
117 
120  static TParticlePDG *pdg( Int_t pdgid );
121 
124  void AddGenerator( StarGenerator *gener );
125 
128  void AddFilter( StarFilterMaker *filt );
129 
131  void SetVertex( Double_t x, Double_t y, Double_t z ){ mVx=x; mVy=y; mVz=z; }
132 
138  void SetSigma( Double_t sx, Double_t sy, Double_t sz, Double_t rho=0 ){ mSx=sx; mSy=sy; mSz=sz; mRho=rho; }
139 
143  void SetSlope( Double_t dxdz, Double_t dydz ){ mVdxdz=dxdz; mVdydz=dydz; }
144 
150  void SetBeamline( Int_t beamline=1 ){ assert(beamline>=0); SetAttr("beamline",beamline); }
151 
153  void SetCuts( Double_t ptmin, Double_t ptmax=-1,
154  Double_t ymin=0, Double_t ymax=-1,
155  Double_t phimin=0, Double_t phimax=-1,
156  Double_t zmin=-999, Double_t zmax=+999 );
157 
159  void SetPtRange( Double_t ptmin, Double_t ptmax=-1 ){ LOG_INFO << " -- StarPrimaryMaker -- ptmin=" << (mPtMin = ptmin) << " ptmax=" << (mPtMax = ptmax) << endm; }
161  void SetEtaRange( Double_t etamin, Double_t etamax ){ LOG_INFO << " -- StarPrimaryMaker -- ymin=" << (mRapidityMin = etamin) << " ymax= " << (mRapidityMax = etamax) << endm; }
163  void SetPhiRange( Double_t phimin, Double_t phimax ){ LOG_INFO << " -- StarPrimaryMaker -- phimn=" << (mPhiMin = phimin) << " phimx=" << (mPhiMax = phimax) << endm; }
164 
166  void SetZvertexRange( Double_t zmin, Double_t zmax ){ LOG_INFO << " -- StarPrimaryMaker -- zmin=" << (mZMin = zmin) << " zmax=" << (mZMax = zmax) << endm; }
167 
169  StarGenEvent *event() { return mPrimaryEvent; }
170 
171  virtual const char *GetCVS() const
172  {static const char cvs[]="Tag $Name: $ $Id: StarPrimaryMaker.h,v 1.8 2015/06/15 13:23:00 jwebb Exp $ built " __DATE__ " " __TIME__ ; return cvs;}
173 
174  private:
175  protected:
176 
177 
178  Int_t PreGenerate();
179  Int_t Generate();
180  Int_t PostGenerate();
181  Int_t Finalize();
182  void BuildTables();
183 
186 
188  TTree *mTree;
189  TFile *mFile;
190  TString mTreeName;
191  TString mFileName;
192  StarParticleStack *mStack;
193  StarGenEvent *mPrimaryEvent;
194 
195  // Vertex position, width and x,y correlation
196  Double_t mVx, mVy, mVz, mSx, mSy, mSz, mRho;
197 
198  // Vertex slope for beamline constraint
199  Double_t mVdxdz, mVdydz;
200 
201  // Handle boost onto the beamline
202  Bool_t mDoBeamline;
203 
214  void RotateBeamline( Double_t &px, Double_t &py, Double_t &pz, Double_t &E, Double_t &M, Double_t &vx, Double_t &vy, Double_t &vz, Double_t &vt );
215 
216  Double_t mPtMin, mPtMax, mRapidityMin, mRapidityMax, mPhiMin, mPhiMax, mZMin, mZMax;
217 
219  Bool_t Simulate( StarGenParticle *p );
220 
221  Int_t mRunNumber;
222 
223  TLorentzVector mPrimaryVertex;
224 
225  StarFilterMaker *mFilter;
226  TEventList *mAccepted; //*< event list containing accepted events
227 
228 #ifndef __CINT__
229  std::function< TLorentzVector() > mVertexFunction;
230  std::map< std::string, std::function< TLorentzVector() > > mVertexFunctionMap;
231  std::function< TLorentzVector() > GetVertexFunction(const char* name);
232 #endif
233 
234  TLorentzVector vertexGaussXYZ(); // gaussian XYZ
235  TLorentzVector vertexFlatZ(); // uniform distribution along z (point in XY)
236  TLorentzVector vertexFlatXYZ(); // uniform distribution in rectangle xyz
237  TLorentzVector vertexFlatRZ(); // uniform distribution in cylinder rz
238  TLorentzVector vertexFlatABZ(); // uniform distribution in ellipical cylinder xyz
239 
240  ClassDef(StarPrimaryMaker,1);
241 
242 };
243 
244 #endif
static TParticlePDG * pdg(Int_t pdgid)
TTree * mTree
The output tree.
Bool_t Simulate(StarGenParticle *p)
Tests to see whether the particle passes all appropriate cuts to be passed to the simulator...
void SetSigma(Double_t sx, Double_t sy, Double_t sz, Double_t rho=0)
void SetFileName(const Char_t *name)
Set the filename of the output TTree.
void SetSlope(Double_t dxdz, Double_t dydz)
Yet another particle class.
void SetPhiRange(Double_t phimin, Double_t phimax)
Set phi range. Particles falling outside this range will be dropped from simulation.
ABC for defining event generator interfaces.
Definition: StarGenerator.h:34
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 SetBeamline(Int_t beamline=1)
Int_t mNumParticles
Total number of particles.
Main filter class. Goes anywhere in the chain, filters StarGenEvent objects.
void AddGenerator(StarGenerator *gener)
Implementation of the VMC particle stack for use in STAR.
Base class for event records.
Definition: StarGenEvent.h:81
StarGenEvent * event()
Return a pointer to the event.
Main steering class for event generation.
void SetZvertexRange(Double_t zmin, Double_t zmax)
Set z-vertex range. Primary vertices outside these bounds will be rejected.
void SetPtRange(Double_t ptmin, Double_t ptmax=-1)
Set PT range. Particles falling outside this range will be dropped from simulation.
void RotateBeamline(Double_t &px, Double_t &py, Double_t &pz, Double_t &E, Double_t &M, Double_t &vx, Double_t &vy, Double_t &vz, Double_t &vt)
void SetEtaRange(Double_t etamin, Double_t etamax)
Set rapidity range. Particles falling outside this range will be dropped from simulation.
void AddFilter(StarFilterMaker *filt)
void SetVertex(Double_t x, Double_t y, Double_t z)
Set the x, y and z vertex position.