Geant4-11
G4ParticleChangeForGamma.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// G4ParticleChangeForGamma
27//
28// Class description:
29//
30// Concrete class for ParticleChange for gamma processes.
31
32// Author: Hisaya Kurashige, 23 March 1998
33// Revision: Vladimir Ivantchenko, 15 April 2005
34// --------------------------------------------------------------------
35#ifndef G4ParticleChangeForGamma_hh
36#define G4ParticleChangeForGamma_hh 1
37
38#include "globals.hh"
39#include "G4ios.hh"
40#include "G4VParticleChange.hh"
41
43
45{
46 public:
47
49 // Default constructor
50
52 // Destructor
53
54 // --- the following methods are for updating G4Step -----
55
58 // A physics process gives the final state of the particle
59 // based on information of G4Track
60
61 inline void InitializeForPostStep(const G4Track&);
62 // Initialize all properties by using G4Track information
63
64 void AddSecondary(G4DynamicParticle* aParticle);
65 // Add next secondary
66
69 // Get/Set the final kinetic energy of the current particle
70
71 inline const G4ThreeVector& GetProposedMomentumDirection() const;
72 inline void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz);
73 inline void ProposeMomentumDirection(const G4ThreeVector& Pfinal);
74 // Get/Propose the MomentumDirection vector: it is the final momentum
75 // direction
76
77 inline const G4ThreeVector& GetProposedPolarization() const;
78 inline void ProposePolarization(const G4ThreeVector& dir);
79 inline void ProposePolarization(G4double Px, G4double Py, G4double Pz);
80
81 inline const G4Track* GetCurrentTrack() const;
82
83 virtual void DumpInfo() const;
84
85 virtual G4bool CheckIt(const G4Track&);
86
87 protected:
88
91 // Hidden copy constructor and assignment operator
92
93 private:
94
95 const G4Track* currentTrack = nullptr;
96 // The pointer to G4Track
97
99 // The final kinetic energy of the current particle
100
102 // The final momentum direction of the current particle
103
105 // The final polarization of the current particle
106};
107
108// ----------------------
109// Inline methods
110// ----------------------
111
112inline
114{
115 return proposedKinEnergy;
116}
117
118inline
120{
122}
123
124inline
127{
129}
130
131inline
134{
136}
137
138inline
140 G4double Py,
141 G4double Pz)
142{
146}
147
148inline
150{
151 return currentTrack;
152}
153
154inline
156{
158}
159
160inline
162{
164}
165
166inline
168 G4double Py,
169 G4double Pz)
170{
174}
175
176inline
178{
183 theParentWeight = track.GetWeight();
188 currentTrack = &track;
189}
190
191#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
void setY(double)
void setZ(double)
void setX(double)
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void InitializeForPostStep(const G4Track &)
void ProposePolarization(const G4ThreeVector &dir)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4Step * UpdateStepForAtRest(G4Step *pStep)
const G4ThreeVector & GetProposedMomentumDirection() const
G4ParticleChangeForGamma & operator=(const G4ParticleChangeForGamma &right)
virtual G4bool CheckIt(const G4Track &)
void AddSecondary(G4DynamicParticle *aParticle)
G4Step * UpdateStepForPostStep(G4Step *Step)
const G4ThreeVector & GetProposedPolarization() const
Definition: G4Step.hh:62
G4TrackStatus GetTrackStatus() const
G4double GetWeight() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit
void InitializeSecondaries(const G4Track &)
G4double energy(const ThreeVector &p, const G4double m)