Geant4-11
G4ParticleChangeForDecay.cc
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 class implementation
27//
28// Author: Hisaya Kurashige, 23 March 1998
29// --------------------------------------------------------------------
30
32#include "G4SystemOfUnits.hh"
33#include "G4Track.hh"
34#include "G4Step.hh"
35#include "G4TrackFastVector.hh"
36#include "G4DynamicParticle.hh"
38
39// --------------------------------------------------------------------
42{
43}
44
45// --------------------------------------------------------------------
47{
48}
49
50// --------------------------------------------------------------------
53 : G4VParticleChange(right)
54{
59}
60
61// --------------------------------------------------------------------
64{
65 if(this != &right)
66 {
68 {
69 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
70 {
71 if((*theListOfSecondaries)[index])
72 delete(*theListOfSecondaries)[index];
73 }
74 }
78 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
79 {
80 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
81 theListOfSecondaries->SetElement(index, newTrack);
82 }
83
88
93 }
94 return *this;
95}
96
97// --------------------------------------------------------------------
99 const G4ParticleChangeForDecay& right) const
100{
101 return ((G4VParticleChange*) this == (G4VParticleChange*) &right);
102}
103
104// --------------------------------------------------------------------
106 const G4ParticleChangeForDecay& right) const
107{
108 return ((G4VParticleChange*) this != (G4VParticleChange*) &right);
109}
110
111// --------------------------------------------------------------------
113{
114 // use base class's method at first
116
117 const G4DynamicParticle* pParticle = track.GetDynamicParticle();
118
119 // set TimeChange equal to local time of the parent track
120 theTimeChange = track.GetLocalTime();
121
122 // set initial Local/Global time of the parent track
123 theLocalTime0 = track.GetLocalTime();
125
126 // set the Polarization equal to those of the parent track
128}
129
130// --------------------------------------------------------------------
132{
133 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
134
136 {
137 pPostStepPoint->SetWeight(theParentWeight);
138 }
139 // update polarization
140 pPostStepPoint->SetPolarization(thePolarizationChange);
141
142 // update the G4Step specific attributes
143 return UpdateStepInfo(pStep);
144}
145
146// --------------------------------------------------------------------
148{
149 // A physics process always calculates the final state of the particle
150
151 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
152
153 // update polarization
154 pPostStepPoint->SetPolarization(thePolarizationChange);
155
156 // update time
157 pPostStepPoint->SetGlobalTime(GetGlobalTime());
158 pPostStepPoint->SetLocalTime(theTimeChange);
159 pPostStepPoint->AddProperTime(theTimeChange - theLocalTime0);
160
161#ifdef G4VERBOSE
162 G4Track* aTrack = pStep->GetTrack();
163 if(debugFlag) { CheckIt(*aTrack); }
164#endif
165
167 pPostStepPoint->SetWeight(theParentWeight);
168
169 // Update the G4Step specific attributes
170 return UpdateStepInfo(pStep);
171}
172
173// --------------------------------------------------------------------
175{
176 // Show header
178
179 G4int oldprc = G4cout.precision(3);
180 G4cout << " proposed local Time (ns) : " << std::setw(20)
181 << theTimeChange / ns << G4endl;
182 G4cout << " initial local Time (ns) : " << std::setw(20)
183 << theLocalTime0 / ns << G4endl;
184 G4cout << " initial global Time (ns) : " << std::setw(20)
185 << theGlobalTime0 / ns << G4endl;
186 G4cout.precision(oldprc);
187}
188
189// --------------------------------------------------------------------
191{
192 G4bool exitWithError = false;
193
194 G4double accuracy;
195
196 // local time should not go back
197 G4bool itsOK = true;
198 accuracy = -1.0 * (theTimeChange - theLocalTime0) / ns;
199 if(accuracy > accuracyForWarning)
200 {
201 itsOK = false;
202 exitWithError = (accuracy > accuracyForException);
203#ifdef G4VERBOSE
204 G4cout << " G4ParticleChangeForDecay::CheckIt : ";
205 G4cout << "the local time goes back !!"
206 << " Difference: " << accuracy << "[ns] " << G4endl;
207 G4cout << "initial local time " << theLocalTime0 / ns << "[ns] "
208 << "initial global time " << theGlobalTime0 / ns << "[ns] "
209 << G4endl;
211 << " E=" << aTrack.GetKineticEnergy() / MeV
212 << " pos=" << aTrack.GetPosition().x() / m << ", "
213 << aTrack.GetPosition().y() / m << ", "
214 << aTrack.GetPosition().z() / m << G4endl;
215#endif
216 }
217
218 // dump out information of this particle change
219#ifdef G4VERBOSE
220 if(!itsOK) { DumpInfo(); }
221#endif
222
223 // exit with error
224 if(exitWithError)
225 {
226 G4Exception("G4ParticleChangeForDecay::CheckIt()", "TRACK005",
227 EventMustBeAborted, "time was illegal");
228 }
229
230 // correction
231 if(!itsOK)
232 {
233 theTimeChange = aTrack.GetLocalTime();
234 }
235
236 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
237 return itsOK;
238}
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double MeV
Definition: G4SIunits.hh:200
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
const G4ThreeVector & GetPolarization() const
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
G4double GetGlobalTime(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
const G4String & GetParticleName() const
void SetLocalTime(const G4double aValue)
void SetWeight(G4double aValue)
void AddProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetPolarization(const G4ThreeVector &aValue)
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4double GetLocalTime() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
static const G4double accuracyForException
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
static const G4double accuracyForWarning
virtual void Initialize(const G4Track &)
virtual void DumpInfo() const
G4Step * UpdateStepInfo(G4Step *Step)
#define ns
Definition: xmlparse.cc:614