StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarGenEvent.cxx
1 #include "StarGenEvent.h"
2 ClassImp(StarGenEvent);
3 #include "StarGenParticle.h"
4 
5 #include "TDatabasePDG.h"
6 #include "TParticlePDG.h"
7 #include <assert.h>
8 #include "TMath.h"
9 #include <iostream>
10 
11 using namespace std;
12 
13 // ----------------------------------------------------------------------------
14 StarGenEvent::StarGenEvent(const Char_t *name, const Char_t *title ) : TObject(),
15  mName(name),
16  mTitle(title),
17  mParticles(0),
18  mGeneratorId(0),
19  mProcessId(0),
20  mOffset(0),
21  mEventNumber(0),
22  mRunNumber(0),
23  mDaqRunNumber(0),
24  mDaqFileNumber(0),
25  mBlueId(0),
26  mYellId(0),
27  mCmsEnergy(0),
28  mFilterResult( StarGenEvent::kAccept ),// default accept
29  mNumParticles(0)
30 {
31  for ( Int_t i=0;i<3;i++ ) mNumRejected[i]=0;
32  mParticles = new TClonesArray( "StarGenParticle", 1000 );
33 }
34 // ----------------------------------------------------------------------------
36 {
37  if ( mParticles ) delete mParticles; mParticles=0;
38 }
39 // ----------------------------------------------------------------------------
40 void StarGenEvent::Clear( Option_t *opts )
41 {
42  TString Opts = opts;
43  if ( Opts.Contains("part") ) {
44  mParticles->Clear();
45  mNumParticles = 0;
46  }
47  if ( Opts.Contains("data") ) {
48  mWeights.clear();
49  mProcessId = 0;
50  mFilterResult = StarGenEvent::kAccept; // Default is always accept
51  mCmsEnergy = 0;
52  for ( Int_t i=0;i<3;i++ ) mNumRejected[i]=0;
53  }
54 }
55 // ----------------------------------------------------------------------------
57 {
58  TClonesArray &particles = *mParticles;
59  StarGenParticle *particle = new( particles[ mNumParticles ] ) StarGenParticle();
60  particle->SetIndex( mNumParticles++ );
61  return particle;
62 }
63 StarGenParticle *StarGenEvent::AddParticle( Int_t status, Int_t pdg, Int_t m1, Int_t m2, Int_t d1, Int_t d2,
64  Double_t px, Double_t py, Double_t pz, Double_t E, Double_t M,
65  Double_t vx, Double_t vy, Double_t vz, Double_t vt )
66 {
67 
68  StarGenParticle *particle = AddParticle();
69  particle->SetStatus(status);
70  particle->SetId(pdg);
71  particle->SetFirstMother(m1);
72  particle->SetLastMother(m2);
73  particle->SetFirstDaughter(d1);
74  particle->SetLastDaughter(d2);
75  particle->SetPx(px);
76  particle->SetPy(py);
77  particle->SetPz(pz);
78  particle->SetEnergy(E);
79  particle->SetMass(M);
80  particle->SetVx(vx);
81  particle->SetVy(vy);
82  particle->SetVz(vz);
83  particle->SetTof(vt);
84 
85  return particle;
86 }
87 // ----------------------------------------------------------------------------
89 {
90  // Double_t P[]={p->GetPx(),p->GetPy(),p->GetPz()};
91  // Double_t V[]={p->GetVx(),p->GetVy(),p->GetVz()};
92  // Double_t tof=p->GetTof();
93  return AddParticle( p->GetStatus(), p->GetId(), p->GetFirstMother(), p->GetLastMother(), p->GetFirstDaughter(), p->GetLastDaughter(),
94  p->GetPx(), p->GetPy(), p->GetPz(), p->GetEnergy(), p->GetMass(),
95  p->GetVx(), p->GetVy(), p->GetVz(), p->GetTof() );
96 }
97 // ----------------------------------------------------------------------------
98 void StarGenEvent::Print( const Option_t *opts )const
99 {
100  TString Opts = opts;
101 
102  if ( Opts.Contains("head") )
103  {
104  cout << "----------------------------------------------------------------------------- "
105  << mName.Data() << endl;
106  cout << "Run number: " << mRunNumber << endl;
107  cout << "Event number: " << mEventNumber << endl;
108  cout << "Generator: " << mGeneratorId << endl;
109  cout << "Offset: " << mOffset << endl;
110  cout << Form("Filter: 0x%04x",int(mFilterResult)) << mOffset << endl;
111  cout << "----------------------------------------------------------------------------- "
112  << endl;
113  }
114 
115  for ( Int_t i=0;i<mParticles->GetEntriesFast();i++ )
116  {
117  StarGenParticle *part = (StarGenParticle *)(*mParticles)[i];
118  if ( !Opts.Contains("simu")||part->Simulate() )
119  part->Print();
120  }
121 
122 }
123 // ----------------------------------------------------------------------------
Float_t GetVz()
Get the z-component of the start vertex.
~StarGenEvent()
Destructor.
void SetPx(Float_t px)
Set the x-component of the momentum.
StarGenEvent(const Char_t *name="event", const Char_t *title="")
Constructor.
Int_t mEventNumber
Event number.
Definition: StarGenEvent.h:233
void SetVy(Float_t vy)
Set the y-component of the start vertex.
Int_t mProcessId
Event generator process ID.
Definition: StarGenEvent.h:230
void Print(const Option_t *opts="head") const
UInt_t mFilterResult
Result of filter.
Definition: StarGenEvent.h:244
void SetLastMother(Int_t last)
Set the last mother particle in the array of particles.
Yet another particle class.
Int_t GetStatus()
Get the status code of the particle according to the HEPEVT standard.
Float_t GetTof()
Get the tof.
Float_t GetPz()
Get the z-component of the momentum.
Float_t GetEnergy()
Get the energy.
void SetVz(Float_t vz)
Set the z-component of the start vertex.
void SetVx(Float_t vx)
Set the x-component of the start vertex.
void SetIndex(Int_t i)
Set the line number in the event record.
virtual void Clear(const Option_t *opts="part,data")
Clear the event.
void SetId(Int_t id)
Int_t GetId()
Get the id code of the particle according to the PDG standard.
std::vector< Double_t > mWeights
User weights.
Definition: StarGenEvent.h:246
void SetStatus(Int_t status)
Set the status code of the particle according to the HEPEVT standard.
void SetPz(Float_t pz)
Set the z-component of the momentum.
Int_t mOffset
Event generator offset.
Definition: StarGenEvent.h:231
Double_t mCmsEnergy
aka sqrt(s)
Definition: StarGenEvent.h:241
void SetMass(Float_t mass)
Set the mass.
void SetTof(Float_t tof)
Set the tof.
Float_t GetPx()
Get the x-component of the momentum.
Int_t mGeneratorId
Generator Id.
Definition: StarGenEvent.h:229
Int_t GetLastMother()
Get the last mother particle.
Float_t GetVy()
Get the y-component of the start vertex.
Float_t GetPy()
Get the y-component of the momentum.
Int_t GetFirstMother()
Get the first mother particle.
Base class for event records.
Definition: StarGenEvent.h:81
void SetLastDaughter(Int_t last)
Set the last daughter particle in the array of particles.
Float_t GetVx()
Get the x-component of the start vertex.
Int_t GetLastDaughter()
Get the last daughter particle.
void SetFirstDaughter(Int_t first)
Set the first daughter particle in the array of particles.
void SetPy(Float_t py)
Set the y-component of the momentum.
void SetEnergy(Float_t energy)
Set the energy.
Int_t mNumRejected[3]
0=total, 1=EG, 2=filter
Definition: StarGenEvent.h:243
void Print(const Option_t *opts="") const
Print the particle.
TClonesArray * mParticles
Array of particles.
Definition: StarGenEvent.h:227
void SetFirstMother(Int_t first)
Set the first mother particle in the array of particles.
Int_t mNumParticles
Number of particles in the record.
Definition: StarGenEvent.h:248
Int_t mRunNumber
Monte Carlo run number.
Definition: StarGenEvent.h:234
Float_t GetMass()
Get the mass.
Int_t GetFirstDaughter()
Get the first daughter particle.
StarGenParticle * AddParticle()