Geant4  10.02.p03
B1PrimaryGeneratorAction Class Reference

#include <B1PrimaryGeneratorAction.hh>

Inheritance diagram for B1PrimaryGeneratorAction:
Collaboration diagram for B1PrimaryGeneratorAction:

Public Member Functions

 B1PrimaryGeneratorAction ()
 
virtual ~B1PrimaryGeneratorAction ()
 
virtual void GeneratePrimaries (G4Event *)
 
const G4ParticleGunGetParticleGun () const
 
 B1PrimaryGeneratorAction ()
 
virtual ~B1PrimaryGeneratorAction ()
 
virtual void GeneratePrimaries (G4Event *)
 
const G4ParticleGunGetParticleGun () const
 
 B1PrimaryGeneratorAction ()
 
virtual ~B1PrimaryGeneratorAction ()
 
virtual void GeneratePrimaries (G4Event *)
 
const G4ParticleGunGetParticleGun () const
 
- Public Member Functions inherited from G4VUserPrimaryGeneratorAction
 G4VUserPrimaryGeneratorAction ()
 
virtual ~G4VUserPrimaryGeneratorAction ()
 

Private Attributes

G4ParticleGunfParticleGun
 
G4BoxfEnvelopeBox
 

Detailed Description

The primary generator action class with particle gun.

The default kinematic is a 6 MeV gamma, randomly distribued in front of the phantom across 80% of the (X,Y) phantom size.

Definition at line 47 of file basic/B1/include/B1PrimaryGeneratorAction.hh.

Constructor & Destructor Documentation

◆ B1PrimaryGeneratorAction() [1/3]

B1PrimaryGeneratorAction::B1PrimaryGeneratorAction ( )

Definition at line 45 of file basic/B1/src/B1PrimaryGeneratorAction.cc.

47  fParticleGun(0),
48  fEnvelopeBox(0)
49 {
50  G4int n_particle = 1;
51  fParticleGun = new G4ParticleGun(n_particle);
52 
53  // default particle kinematic
55  G4String particleName;
56  G4ParticleDefinition* particle
57  = particleTable->FindParticle(particleName="gamma");
61 }
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
CLHEP::Hep3Vector G4ThreeVector
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
int G4int
Definition: G4Types.hh:78
void SetParticleEnergy(G4double aKineticEnergy)
static G4ParticleTable * GetParticleTable()
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
Here is the call graph for this function:

◆ ~B1PrimaryGeneratorAction() [1/3]

B1PrimaryGeneratorAction::~B1PrimaryGeneratorAction ( )
virtual

Definition at line 65 of file basic/B1/src/B1PrimaryGeneratorAction.cc.

66 {
67  delete fParticleGun;
68 }

◆ B1PrimaryGeneratorAction() [2/3]

B1PrimaryGeneratorAction::B1PrimaryGeneratorAction ( )

◆ ~B1PrimaryGeneratorAction() [2/3]

virtual B1PrimaryGeneratorAction::~B1PrimaryGeneratorAction ( )
virtual

◆ B1PrimaryGeneratorAction() [3/3]

B1PrimaryGeneratorAction::B1PrimaryGeneratorAction ( )

◆ ~B1PrimaryGeneratorAction() [3/3]

virtual B1PrimaryGeneratorAction::~B1PrimaryGeneratorAction ( )
virtual

Member Function Documentation

◆ GeneratePrimaries() [1/3]

void B1PrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Implements G4VUserPrimaryGeneratorAction.

Definition at line 72 of file basic/B1/src/B1PrimaryGeneratorAction.cc.

73 {
74  //this function is called at the begining of ecah event
75  //
76 
77  // In order to avoid dependence of PrimaryGeneratorAction
78  // on DetectorConstruction class we get Envelope volume
79  // from G4LogicalVolumeStore.
80 
81  G4double envSizeXY = 0;
82  G4double envSizeZ = 0;
83 
84  if (!fEnvelopeBox)
85  {
86  G4LogicalVolume* envLV
88  if ( envLV ) fEnvelopeBox = dynamic_cast<G4Box*>(envLV->GetSolid());
89  }
90 
91  if ( fEnvelopeBox ) {
92  envSizeXY = fEnvelopeBox->GetXHalfLength()*2.;
93  envSizeZ = fEnvelopeBox->GetZHalfLength()*2.;
94  }
95  else {
97  msg << "Envelope volume of box shape not found.\n";
98  msg << "Perhaps you have changed geometry.\n";
99  msg << "The gun will be place at the center.";
100  G4Exception("B1PrimaryGeneratorAction::GeneratePrimaries()",
101  "MyCode0002",JustWarning,msg);
102  }
103 
104  G4double size = 0.8;
105  G4double x0 = size * envSizeXY * (G4UniformRand()-0.5);
106  G4double y0 = size * envSizeXY * (G4UniformRand()-0.5);
107  G4double z0 = -0.5 * envSizeZ;
108 
110 
112 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:64
virtual void GeneratePrimaryVertex(G4Event *evt)
G4double GetXHalfLength() const
void SetParticlePosition(G4ThreeVector aPosition)
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetZHalfLength() const
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4LogicalVolumeStore * GetInstance()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76
G4VSolid * GetSolid() const
Here is the call graph for this function:

◆ GeneratePrimaries() [2/3]

virtual void B1PrimaryGeneratorAction::GeneratePrimaries ( G4Event )
virtual

◆ GeneratePrimaries() [3/3]

virtual void B1PrimaryGeneratorAction::GeneratePrimaries ( G4Event )
virtual

◆ GetParticleGun() [1/3]

const G4ParticleGun* B1PrimaryGeneratorAction::GetParticleGun ( void  ) const
inline

◆ GetParticleGun() [2/3]

const G4ParticleGun* B1PrimaryGeneratorAction::GetParticleGun ( void  ) const
inline

Definition at line 57 of file basic/B1/include/B1PrimaryGeneratorAction.hh.

Here is the caller graph for this function:

◆ GetParticleGun() [3/3]

const G4ParticleGun* B1PrimaryGeneratorAction::GetParticleGun ( void  ) const
inline

Member Data Documentation

◆ fEnvelopeBox

G4Box * B1PrimaryGeneratorAction::fEnvelopeBox
private

Definition at line 61 of file basic/B1/include/B1PrimaryGeneratorAction.hh.

◆ fParticleGun

G4ParticleGun * B1PrimaryGeneratorAction::fParticleGun
private

Definition at line 60 of file basic/B1/include/B1PrimaryGeneratorAction.hh.


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