Geant4_10
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 //
27 // $Id: G4ParticleChangeForDecay.cc 68795 2013-04-05 13:24:46Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 //
34 //
35 // ------------------------------------------------------------
36 // Implemented for the new scheme 23 Mar. 1998 H.Kurahige
37 // Remove modification of energy/momentum 20 Jul, 1998 H.Kurashige
38 // --------------------------------------------------------------
39 
41 #include "G4SystemOfUnits.hh"
42 #include "G4Track.hh"
43 #include "G4Step.hh"
44 #include "G4TrackFastVector.hh"
45 #include "G4DynamicParticle.hh"
46 #include "G4ExceptionSeverity.hh"
47 
49  : G4VParticleChange(),
50  theGlobalTime0(0.),theLocalTime0(0.),theTimeChange(0.)
51 {
52 #ifdef G4VERBOSE
53  if (verboseLevel>2) {
54  G4cout << "G4ParticleChangeForDecay::G4ParticleChangeForDecay() " << G4endl;
55  }
56 #endif
57 }
58 
60 {
61 #ifdef G4VERBOSE
62  if (verboseLevel>2) {
63  G4cout << "G4ParticleChangeForDecay::~G4ParticleChangeForDecay() " << G4endl;
64  }
65 #endif
66 }
67 
69 {
74 }
75 
76 
78 {
79  if (this != &right){
80  if (theNumberOfSecondaries>0) {
81 #ifdef G4VERBOSE
82  if (verboseLevel>0) {
83  G4cout << "G4ParticleChangeForDecay: assignment operator Warning ";
84  G4cout << "theListOfSecondaries is not empty ";
85  }
86 #endif
88  if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
89  }
90  }
91  delete theListOfSecondaries;
95  G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
97 
102 
107  }
108  return *this;
109 }
110 
112 {
113  return ((G4VParticleChange *)this == (G4VParticleChange *) &right);
114 }
115 
117 {
118  return ((G4VParticleChange *)this != (G4VParticleChange *) &right);
119 }
120 
121 //----------------------------------------------------------------
122 // methods for Initialization
123 //
125 {
126  // use base class's method at first
128 
129  const G4DynamicParticle* pParticle = track.GetDynamicParticle();
130 
131  // set TimeChange equal to local time of the parent track
132  theTimeChange = track.GetLocalTime();
133 
134  // set initial Local/Global time of the parent track
135  theLocalTime0 = track.GetLocalTime();
136  theGlobalTime0 = track.GetGlobalTime();
137 
138  // set the Polarization equal to those of the parent track
139  thePolarizationChange = pParticle->GetPolarization();
140 }
141 
142 //----------------------------------------------------------------
143 // methods for updating G4Step
144 //
145 
147 {
150  }
151 
152  // Update the G4Step specific attributes
153  return UpdateStepInfo(pStep);
154 }
155 
156 
158 {
159  // A physics process always calculates the final state of the particle
160 
161  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
162 
163  // update polarization
164  pPostStepPoint->SetPolarization( thePolarizationChange );
165 
166  // update time
167  pPostStepPoint->SetGlobalTime( GetGlobalTime() );
168  pPostStepPoint->SetLocalTime( theTimeChange );
169  pPostStepPoint->AddProperTime (theTimeChange-theLocalTime0);
170 
171 #ifdef G4VERBOSE
172  G4Track* aTrack = pStep->GetTrack();
173  if (debugFlag) CheckIt(*aTrack);
174 #endif
175 
176  if (isParentWeightProposed )pPostStepPoint->SetWeight( theParentWeight );
177 
178  // Update the G4Step specific attributes
179  return UpdateStepInfo(pStep);
180 }
181 
183 {
184 // Show header
186 
187  G4int oldprc = G4cout.precision(3);
188  G4cout << " Time (ns) : "
189  << std::setw(20) << theTimeChange/ns << G4endl;
190  G4cout.precision(oldprc);
191 }
192 
194 {
195  G4bool exitWithError = false;
196 
197  G4double accuracy;
198 
199  // local time should not go back
200  G4bool itsOK =true;
201  accuracy = -1.0*(theTimeChange - theLocalTime0)/ns;
202  if (accuracy > accuracyForWarning) {
203  itsOK = false;
204  exitWithError = (accuracy > accuracyForException);
205 #ifdef G4VERBOSE
206  G4cout << " G4ParticleChangeForDecay::CheckIt : ";
207  G4cout << "the local time goes back !!"
208  << " Difference: " << accuracy << "[ns] " <<G4endl;
209  G4cout << aTrack.GetDefinition()->GetParticleName()
210  << " E=" << aTrack.GetKineticEnergy()/MeV
211  << " pos=" << aTrack.GetPosition().x()/m
212  << ", " << aTrack.GetPosition().y()/m
213  << ", " << aTrack.GetPosition().z()/m
214  <<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  G4Exception("G4ParticleChangeForDecay::CheckIt",
226  "TRACK005",EventMustBeAborted,
227  "time was illegal");
228  }
229 
230  // correction
231  if (!itsOK) {
232  theTimeChange = aTrack.GetLocalTime();
233  }
234 
235  itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
236  return itsOK;
237 }
238 
239 
240 
241 
242 
G4ParticleDefinition * GetDefinition() const
virtual void Initialize(const G4Track &)
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4bool operator==(const G4ParticleChangeForDecay &right) const
G4double GetLocalTime() const
Int_t index
Definition: macro.C:9
double x() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetGlobalTime(G4double timeDelay=0.0) const
virtual G4bool CheckIt(const G4Track &)
const G4ThreeVector & GetPosition() const
G4TrackFastVector * theListOfSecondaries
void SetWeight(G4double aValue)
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
virtual void DumpInfo() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
double z() const
static const G4double accuracyForException
void SetGlobalTime(const G4double aValue)
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
void SetLocalTime(const G4double aValue)
virtual G4bool CheckIt(const G4Track &)
virtual void Initialize(const G4Track &)
G4double GetKineticEnergy() const
void SetPolarization(const G4ThreeVector &aValue)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4SteppingControl theSteppingControlFlag
Definition: G4Step.hh:76
G4double GetGlobalTime() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4Step * UpdateStepInfo(G4Step *Step)
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
G4StepPoint * GetPostStepPoint() const
G4ParticleChangeForDecay & operator=(const G4ParticleChangeForDecay &right)
double y() const
const G4ThreeVector & GetPolarization() const
G4bool operator!=(const G4ParticleChangeForDecay &right) const
#define G4endl
Definition: G4ios.hh:61
G4TrackStatus theStatusChange
double G4double
Definition: G4Types.hh:76
static const G4double accuracyForWarning
G4Track * GetTrack() const
#define ns
Definition: xmlparse.cc:597
void AddProperTime(const G4double aValue)