Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleChangeForDecay.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// G4ParticleChangeForDecay
27//
28// Class description:
29//
30// Concrete class for ParticleChange which has functionality for G4Decay.
31//
32// This class contains the results after invocation of the decay process.
33// This includes secondary particles generated by the interaction.
34
35// Author: Hisaya Kurashige, 23 March 1998
36// --------------------------------------------------------------------
37#ifndef G4ParticleChangeForDecay_hh
38#define G4ParticleChangeForDecay_hh 1
39
40#include "globals.hh"
41#include "G4ios.hh"
42#include "G4ThreeVector.hh"
43#include "G4VParticleChange.hh"
44
46
48{
49 public:
50
52
53 ~G4ParticleChangeForDecay() override = default;
54
57
58 // --- the following methods are for updating G4Step -----
59 // Return the pointer to the G4Step after updating the step information
60 // by using final state information of the track given by a physics process
61 // !!! No effect for AlongSteyp
62
63 G4Step* UpdateStepForAtRest(G4Step* Step) final;
65
66 void Initialize(const G4Track&) final;
67 // Initialize all properties by using G4Track information
68
71 // Get/Propose the final global/local time
72 // NOTE: DO NOT INVOKE both methods in a step
73 // Each method affects both local and global time
74
75 G4double GetGlobalTime(G4double timeDelay = 0.0) const;
76 G4double GetLocalTime(G4double timeDelay = 0.0) const;
77 // Convert the time delay to the glocbal/local time.
78 // Can get the final global/local time without argument
79
80 const G4ThreeVector* GetPolarization() const;
82 void ProposePolarization(const G4ThreeVector& finalPoralization);
83 // Get/Propose the final Polarization vector
84
85 // --- Dump and debug methods ---
86
87 void DumpInfo() const final;
88
89 G4bool CheckIt(const G4Track&) final;
90
91 private:
92
93 G4double theGlobalTime0 = 0.0;
94 // The global time at Initial
95 G4double theLocalTime0 = 0.0;
96 // The local time at Initial
97
98 G4double theTimeChange = 0.0;
99 // The change of local time of a given particle
100
101 G4ThreeVector thePolarizationChange;
102 // The changed (final) polarization of a given track
103};
104
105// ----------------------
106// Inline methods
107// ----------------------
108
109inline
111{
112 theTimeChange = (t - theGlobalTime0) + theLocalTime0;
113}
114
115inline
117{
118 // Convert the time delay to the global time.
119 return theGlobalTime0 + (theTimeChange - theLocalTime0) + timeDelay;
120}
121
122inline
124{
125 theTimeChange = t;
126}
127
128inline
130{
131 // Convert the time delay to the local time.
132 return theTimeChange + timeDelay;
133}
134
135inline
137{
138 return &thePolarizationChange;
139}
140
141inline
143 const G4ThreeVector& finalPoralization)
144{
145 thePolarizationChange = finalPoralization;
146}
147
148inline
150 G4double Py,
151 G4double Pz)
152{
153 thePolarizationChange.setX(Px);
154 thePolarizationChange.setY(Py);
155 thePolarizationChange.setZ(Pz);
156}
157
158#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
void setY(double)
void setZ(double)
void setX(double)
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
const G4ThreeVector * GetPolarization() const
void Initialize(const G4Track &) final
~G4ParticleChangeForDecay() override=default
G4Step * UpdateStepForPostStep(G4Step *Step) final
G4double GetGlobalTime(G4double timeDelay=0.0) const
G4double GetLocalTime(G4double timeDelay=0.0) const
G4bool CheckIt(const G4Track &) final
G4ParticleChangeForDecay & operator=(const G4ParticleChangeForDecay &right)=delete
G4ParticleChangeForDecay(const G4ParticleChangeForDecay &right)=delete
G4Step * UpdateStepForAtRest(G4Step *Step) final
Definition: G4Step.hh:62
#define inline
Definition: internal.h:104