Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PrimaryTransformer.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29
30#ifndef G4PromaryTransformer_h
31#define G4PromaryTransformer_h 1
32
33#include "G4TrackVector.hh"
34#include "G4ParticleTable.hh"
35#include "globals.hh"
36class G4Event;
37class G4PrimaryVertex;
38#include "G4PrimaryParticle.hh"
39#include "G4DynamicParticle.hh"
40
41// class description:
42//
43// This class is exclusively used by G4EventManager for the conversion
44// from G4PrimaryVertex/G4PrimaryParticle to G4DynamicParticle/G4Track.
45//
46
48{
49 public:
51 virtual ~G4PrimaryTransformer();
52
53 G4TrackVector* GimmePrimaries(G4Event* anEvent, G4int trackIDCounter=0);
54 void CheckUnknown();
55
56 protected:
61
64
65 public:
66 inline void SetVerboseLevel(G4int vl)
67 { verboseLevel = vl; };
68
69 public: //with description
71 // By invoking this Set method, the user can alter the treatment of unknown
72 // particle. The ideal place to invoke this method is in the BeginOfRunAction.
73
75 { return unknownParticleDefined; }
76
77 protected:
78 void GenerateTracks(G4PrimaryVertex* primaryVertex);
79 void GenerateSingleTrack(G4PrimaryParticle* primaryParticle,
82 G4DynamicParticle* motherDP);
84
85 protected: //with description
86 // Following two virtual methods are provided to customize the use of PrimaryTransformer
87 // for particle types exotic to Geant4.
88
90 // Return appropriate G4ParticleDefinition w.r.t. the primary particle.
91 // If NULL is returned, the particle will not be treated as a track, but its daughters
92 // will be examined in case it has "pre-assigned decay products".
93
95 // Return true if a primary particle should be converted into a track.
96 // By default, all particles of non-shortlived and shortlived with valid decay
97 // tables are converted.
98};
99
100#endif
101
102
std::vector< G4Track * > G4TrackVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4bool CheckDynamicParticle(G4DynamicParticle *DP)
void SetUnknnownParticleDefined(G4bool vl)
void SetVerboseLevel(G4int vl)
virtual G4ParticleDefinition * GetDefinition(G4PrimaryParticle *pp)
G4ParticleTable * particleTable
G4ParticleDefinition * unknown
G4bool GetUnknownParticleDefined() const
virtual G4bool IsGoodForTrack(G4ParticleDefinition *pd)
void SetDecayProducts(G4PrimaryParticle *mother, G4DynamicParticle *motherDP)
G4TrackVector * GimmePrimaries(G4Event *anEvent, G4int trackIDCounter=0)
void GenerateSingleTrack(G4PrimaryParticle *primaryParticle, G4double x0, G4double y0, G4double z0, G4double t0, G4double wv)
void GenerateTracks(G4PrimaryVertex *primaryVertex)