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.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 "G4DynamicParticle.hh"
37
38// --------------------------------------------------------------------
40{}
41
42// --------------------------------------------------------------------
44{
45 // use base class's method at first
47
48 // set TimeChange equal to local time of the parent track
49 theLocalTime0 = theTimeChange = track.GetLocalTime();
50
51 // set initial Local/Global time of the parent track
52 theGlobalTime0 = track.GetGlobalTime();
53
54 // set the Polarization equal to those of the parent track
55 thePolarizationChange = track.GetDynamicParticle()->GetPolarization();
56}
57
58// --------------------------------------------------------------------
60{
61 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
62
64 {
65 pPostStepPoint->SetWeight(theParentWeight);
66 }
67 // update polarization
68 pPostStepPoint->SetPolarization(thePolarizationChange);
69
70 // update the G4Step specific attributes
71 return UpdateStepInfo(pStep);
72}
73
74// --------------------------------------------------------------------
76{
77 // A physics process always calculates the final state of the particle
78
79 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
80
81 // update polarization
82 pPostStepPoint->SetPolarization(thePolarizationChange);
83
84 // update time
85 pPostStepPoint->SetGlobalTime(GetGlobalTime());
86 pPostStepPoint->SetLocalTime(theTimeChange);
87 pPostStepPoint->AddProperTime(theTimeChange - theLocalTime0);
88
89#ifdef G4VERBOSE
91#endif
92
94 pPostStepPoint->SetWeight(theParentWeight);
95
96 // Update the G4Step specific attributes
97 return UpdateStepInfo(pStep);
98}
99
100// --------------------------------------------------------------------
102{
103 // Show header
105
106 G4long oldprc = G4cout.precision(8);
107 G4cout << " -----------------------------------------------" << G4endl;
108 G4cout << " G4ParticleChangeForDecay proposes: " << G4endl;
109 G4cout << " Proposed local Time (ns): " << std::setw(20)
110 << theTimeChange / ns << G4endl;
111 G4cout << " Initial local Time (ns) : " << std::setw(20)
112 << theLocalTime0 / ns << G4endl;
113 G4cout << " Initial global Time (ns): " << std::setw(20)
114 << theGlobalTime0 / ns << G4endl;
115 G4cout << " Current global Time (ns): " << std::setw(20)
117 G4cout.precision(oldprc);
118}
119
120// --------------------------------------------------------------------
122{
123 // local time should not go back
124 G4bool itsOK = true;
125 if(theTimeChange < theLocalTime0)
126 {
127 itsOK = false;
128 ++nError;
129#ifdef G4VERBOSE
130 if(nError < maxError)
131 {
132 G4cout << " G4ParticleChangeForDecay::CheckIt : ";
133 G4cout << "the local time goes back !!"
134 << " Difference: " << (theTimeChange - theLocalTime0) / ns
135 << "[ns] " << G4endl;
136 G4cout << "initial local time " << theLocalTime0 / ns << "[ns] "
137 << "initial global time " << theGlobalTime0 / ns << "[ns] "
138 << G4endl;
139 }
140#endif
141 theTimeChange = theLocalTime0;
142 }
143
144 if(!itsOK)
145 {
146 if(nError < maxError)
147 {
148#ifdef G4VERBOSE
149 // dump out information of this particle change
150 DumpInfo();
151#endif
152 G4Exception("G4ParticleChangeForDecay::CheckIt()", "TRACK005",
153 JustWarning, "time is illegal");
154 }
155 }
156
157 return (itsOK) && G4VParticleChange::CheckIt(aTrack);
158}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
long G4long
Definition: G4Types.hh:87
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetPolarization() const
void Initialize(const G4Track &) final
G4Step * UpdateStepForPostStep(G4Step *Step) final
G4double GetGlobalTime(G4double timeDelay=0.0) const
G4bool CheckIt(const G4Track &) final
G4Step * UpdateStepForAtRest(G4Step *Step) final
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
G4StepPoint * GetPostStepPoint() const
G4double GetGlobalTime() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4bool CheckIt(const G4Track &)
virtual void Initialize(const G4Track &)
static const G4int maxError
virtual void DumpInfo() const
const G4Track * theCurrentTrack
G4Step * UpdateStepInfo(G4Step *Step)
#define ns(x)
Definition: xmltok.c:1649