Geant4-11
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 // Default constructor
53
55 // Destructor
56
57 G4bool operator==(const G4ParticleChangeForDecay& right) const;
58 G4bool operator!=(const G4ParticleChangeForDecay& right) const;
59 // Equality operators
60
61 // --- the following methods are for updating G4Step -----
62 // Return the pointer to the G4Step after updating the step information
63 // by using final state information of the track given by a physics process
64 // !!! No effect for AlongSteyp
65
66 virtual G4Step* UpdateStepForAtRest(G4Step* Step);
67 virtual G4Step* UpdateStepForPostStep(G4Step* Step);
68
69 virtual void Initialize(const G4Track&);
70 // Initialize all properties by using G4Track information
71
74 // Get/Propose the final global/local time
75 // NOTE: DO NOT INVOKE both methods in a step
76 // Each method affects both local and global time
77
78 G4double GetGlobalTime(G4double timeDelay = 0.0) const;
79 G4double GetLocalTime(G4double timeDelay = 0.0) const;
80 // Convert the time delay to the glocbal/local time.
81 // Can get the final global/local time without argument
82
83 const G4ThreeVector* GetPolarization() const;
85 void ProposePolarization(const G4ThreeVector& finalPoralization);
86 // Get/Propose the final Polarization vector
87
88 // --- Dump and debug methods ---
89
90 virtual void DumpInfo() const;
91
92 virtual G4bool CheckIt(const G4Track&);
93
94 protected:
95
98 // Hidden copy constructor and assignment operator
99
101 // The global time at Initial
103 // The local time at Initial
104
106 // The change of local time of a given particle
107
109 // The changed (final) polarization of a given track
110};
111
112// ----------------------
113// Inline methods
114// ----------------------
115
116inline
118{
120}
121
122inline
124{
125 // Convert the time delay to the global time.
126 return theGlobalTime0 + (theTimeChange - theLocalTime0) + timeDelay;
127}
128
129inline
131{
132 theTimeChange = t;
133}
134
135inline
137{
138 // Convert the time delay to the local time.
139 return theTimeChange + timeDelay;
140}
141
142inline
144{
145 return &thePolarizationChange;
146}
147
148inline
150 const G4ThreeVector& finalPoralization)
151{
152 thePolarizationChange = finalPoralization;
153}
154
155inline
157 G4double Py,
158 G4double Pz)
159{
163}
164
165#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)
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
const G4ThreeVector * GetPolarization() const
G4double GetGlobalTime(G4double timeDelay=0.0) const
G4double GetLocalTime(G4double timeDelay=0.0) const
G4ParticleChangeForDecay & operator=(const G4ParticleChangeForDecay &right)
G4bool operator==(const G4ParticleChangeForDecay &right) const
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual void Initialize(const G4Track &)
virtual G4bool CheckIt(const G4Track &)
G4bool operator!=(const G4ParticleChangeForDecay &right) const
Definition: G4Step.hh:62