File indexing completed on 2024-09-27 07:02:49
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012 #include "eASTPrimaryGeneratorAction.hh"
0013 #include "eASTPrimGenActionMessenger.hh"
0014
0015 #include "G4Event.hh"
0016 #include "G4ParticleGun.hh"
0017 #include "G4GeneralParticleSource.hh"
0018 #include "G4ParticleTable.hh"
0019 #include "G4ParticleDefinition.hh"
0020 #include "G4SystemOfUnits.hh"
0021 #include "Randomize.hh"
0022
0023 #ifdef eAST_USE_HepMC3
0024 #include "eASTHepMC3Interface.hh"
0025 #endif
0026
0027 eASTPrimaryGeneratorAction::eASTPrimaryGeneratorAction(
0028 G4bool useParticleGun, G4bool useParticleSource,
0029 G4bool useHepMC3Interface)
0030 : G4VUserPrimaryGeneratorAction()
0031 {
0032 if(useParticleGun)
0033 {
0034 fParticleGun = new G4ParticleGun(1);
0035
0036 auto particleTable = G4ParticleTable::GetParticleTable();
0037 auto fPion = particleTable->FindParticle("pi+");
0038 fParticleGun->SetParticleDefinition(fPion);
0039
0040
0041 fParticleGun->SetParticlePosition(G4ThreeVector(0.,0.,0.));
0042 fParticleGun->SetParticleMomentumDirection(G4ThreeVector(1.,0.,0.));
0043 fParticleGun->SetParticleEnergy(1.*GeV);
0044 }
0045
0046 if(useParticleSource)
0047 { fParticleSource = new G4GeneralParticleSource(); }
0048
0049 #ifdef eAST_USE_HepMC3
0050 if(useHepMC3Interface)
0051 { fHepMC3Interface = eASTHepMC3Interface::GetInstance(); }
0052 #endif
0053
0054 messenger = new eASTPrimGenActionMessenger(this);
0055 }
0056
0057 eASTPrimaryGeneratorAction::~eASTPrimaryGeneratorAction()
0058 {
0059 if(fParticleGun!=nullptr) delete fParticleGun;
0060 if(fParticleSource!=nullptr) delete fParticleSource;
0061
0062 #ifdef eAST_USE_HepMC3
0063 if(fHepMC3Interface!=nullptr && eASTHepMC3Interface::GetInstance()!=nullptr)
0064 { delete fHepMC3Interface; }
0065 #endif
0066
0067 delete messenger;
0068 }
0069
0070 void eASTPrimaryGeneratorAction::GeneratePrimaries(G4Event* event)
0071 {
0072 if(fParticleGun) fParticleGun->GeneratePrimaryVertex(event);
0073 if(fParticleSource) fParticleSource->GeneratePrimaryVertex(event);
0074 #ifdef eAST_USE_HepMC3
0075 if(fHepMC3Interface) fHepMC3Interface->GeneratePrimaryVertex(event);
0076 #endif
0077
0078 G4int evId = event->GetEventID();
0079 G4double tVtx = deltaT * (G4double)evId + T0;
0080 auto* pv = event->GetPrimaryVertex();
0081 while(pv!=nullptr)
0082 {
0083 pv->SetT0(pv->GetT0() + tVtx);
0084 pv = pv->GetNext();
0085 }
0086
0087 }
0088