Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4PrimaryGeneratorAction.h
Go to the documentation of this file.
1#ifndef TG4_PRIMARY_GENERATOR_ACTION_H
2#define TG4_PRIMARY_GENERATOR_ACTION_H
3
4//------------------------------------------------
5// The Geant4 Virtual Monte Carlo package
6// Copyright (C) 2007 - 2014 Ivana Hrivnacova
7// All rights reserved.
8//
9// For the licensing terms see geant4_vmc/LICENSE.
10// Contact: root-vmc@cern.ch
11//-------------------------------------------------
12
17
18#include "TG4Verbose.h"
19
20#include <G4VUserPrimaryGeneratorAction.hh>
21#include <globals.hh>
22
23class TVirtualMCStack;
24class TMCManagerStack;
25class TParticle;
26
29class TG4TrackManager;
30
31class G4Event;
33class G4PrimaryVertex;
34
40
42 public TG4Verbose
43{
44 public:
47
48 // methods
49 virtual void GeneratePrimaries(G4Event* event);
50
51 // set methods
52 void SetSkipUnknownParticles(G4bool value);
53
54 // get methods
55 G4bool GetSkipUnknownParticles() const;
56
57 private:
58 // methods
59
60 G4bool CheckVMCStack(TVirtualMCStack* stack) const;
61 G4bool CheckParticleDefinition(const G4ParticleDefinition* particleDefinition,
62 const TParticle* particle) const;
63 G4double GetProperCharge(const G4ParticleDefinition* particleDefinition,
64 const TParticle* particle) const;
65 G4PrimaryVertex* AddParticleToVertex(G4Event* event, G4PrimaryVertex* vertex,
66 const G4ParticleDefinition* particleDefinition,
67 const G4ThreeVector& position, G4double time, const G4ThreeVector& momentum,
68 G4double energy, const G4ThreeVector& polarization, G4double charge,
69 G4double weight) const;
70 void TransformPrimaries(G4Event* event);
71 void TransformTracks(G4Event* event);
72
73 // data members
81 TVirtualMCStack* fMCStack;
82 TMCManagerStack* fMCManagerStack;
84 G4bool fCached;
87};
88
89// inline functions
90
96
102
103#endif // TG4_PRIMARY_GENERATOR_ACTION_H
Definition of the TG4Verbose class.
Provides mapping between TDatabasePDG and Geant4 particles.
Primary generator action defined via TVirtualMCStack and TVirtualMCApplication.
G4bool fSkipUnknownParticles
Option to skip particles which do not exist in Geant4.
G4bool CheckParticleDefinition(const G4ParticleDefinition *particleDefinition, const TParticle *particle) const
void SetSkipUnknownParticles(G4bool value)
Set the option to skip particles which do not exist in Geant4.
G4bool GetSkipUnknownParticles() const
Return the option to skip particles which do not exist in Geant4.
TG4ParticlesManager * fParticlesManager
Thread-local particles manager.
TG4PrimaryGeneratorMessenger * fMessenger
Messenger.
TVirtualMCStack * fMCStack
Thread-local stacks.
G4bool fCached
Flag whether thread-local variables have been cached.
G4PrimaryVertex * AddParticleToVertex(G4Event *event, G4PrimaryVertex *vertex, const G4ParticleDefinition *particleDefinition, const G4ThreeVector &position, G4double time, const G4ThreeVector &momentum, G4double energy, const G4ThreeVector &polarization, G4double charge, G4double weight) const
virtual void GeneratePrimaries(G4Event *event)
G4bool CheckVMCStack(TVirtualMCStack *stack) const
TG4TrackManager * fTrackManager
Thread-local track manager.
G4double GetProperCharge(const G4ParticleDefinition *particleDefinition, const TParticle *particle) const
Messenger class that defines commands for TG4PrimaryGeneratorAction.
The class for storing G4 tracks in VMC sack.
Base class for defining the verbose level and a common messenger.
Definition TG4Verbose.h:36