Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
Par02PrimaryGeneratorAction Class Reference

Construction of a primary generation action. More...

#include <Doxymodules_parameterisations.h>

Inheritance diagram for Par02PrimaryGeneratorAction:
G4VUserPrimaryGeneratorAction

Public Member Functions

 Par02PrimaryGeneratorAction ()
 
 ~Par02PrimaryGeneratorAction ()
 
virtual void GeneratePrimaries (G4Event *anEvent)
 
G4ParticleGunGetParticleGun ()
 

Private Attributes

G4ParticleGunfParticleGun
 

Detailed Description

Construction of a primary generation action.

For simplicity, we use here the particle gun, but in the original application for FCC (developed by Anna Zaborowska), the Monte Carlo event generator Pythia8 is used as generator and it is interfaced to Geant4 via HepMC.

Definition at line 46 of file Doxymodules_parameterisations.h.

Constructor & Destructor Documentation

◆ Par02PrimaryGeneratorAction()

Par02PrimaryGeneratorAction::Par02PrimaryGeneratorAction ( )

Definition at line 41 of file Par02PrimaryGeneratorAction.cc.

41 {
42 G4int n_particle = 1;
43 fParticleGun = new G4ParticleGun( n_particle );
44
45 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
46 G4String particleName;
47 G4ParticleDefinition* particle =
48 particleTable->FindParticle( particleName = "geantino" );
49 fParticleGun->SetParticleDefinition( particle );
50
51 fParticleGun->SetParticleMomentumDirection( G4ThreeVector( 0.0, 1.0, 0.0 ) );
52 fParticleGun->SetParticleEnergy( 100.0*GeV );
53 fParticleGun->SetParticlePosition( G4ThreeVector( 0.0, 0.0, 0.0 ) );
54}

◆ ~Par02PrimaryGeneratorAction()

Par02PrimaryGeneratorAction::~Par02PrimaryGeneratorAction ( )

Definition at line 58 of file Par02PrimaryGeneratorAction.cc.

58 {
59 delete fParticleGun;
60}

Member Function Documentation

◆ GeneratePrimaries()

void Par02PrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Definition at line 64 of file Par02PrimaryGeneratorAction.cc.

64 {
65 fParticleGun->GeneratePrimaryVertex( anEvent );
66
67 // Loop over the vertices, and then over primary particles,
68 // and for each primary particle create an info object, in
69 // which to store "Monte Carlo true" information.
70 // This approach could appear unnecessarily heavy in the present case
71 // of a trivial particle gun generator, but it is useful in the more
72 // realistic case of a Monte Carlo event generator like Pythia8.
73 G4int count_particles = 0;
74 for ( G4int ivtx = 0; ivtx < anEvent->GetNumberOfPrimaryVertex(); ivtx++ ) {
75 for ( G4int ipp = 0; ipp < anEvent->GetPrimaryVertex( ivtx )->GetNumberOfParticle();
76 ipp++ ) {
77 G4PrimaryParticle* primary_particle =
78 anEvent->GetPrimaryVertex( ivtx )->GetPrimary( ipp );
79 if ( primary_particle ) {
80 primary_particle->SetUserInformation( new Par02PrimaryParticleInformation(
81 count_particles, primary_particle->GetPDGcode(),
82 primary_particle->GetMomentum() ) );
83 count_particles++;
84 }
85 }
86 }
87}

◆ GetParticleGun()

G4ParticleGun * Par02PrimaryGeneratorAction::GetParticleGun ( )

Definition at line 91 of file Par02PrimaryGeneratorAction.cc.

91 {
92 return fParticleGun;
93}

Member Data Documentation

◆ fParticleGun

G4ParticleGun* Par02PrimaryGeneratorAction::fParticleGun
private

Definition at line 54 of file Par02PrimaryGeneratorAction.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration