Geant4_10
G4ParticleChangeForGamma.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: G4ParticleChangeForGamma.cc 68795 2013-04-05 13:24:46Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // ------------------------------------------------------------
34 // 15 April 2005 V.Ivanchenko for gamma EM processes
35 // --------------------------------------------------------------
36 //
37 // Modified::
38 // 28.08.06 V.Ivanchenko Add access to current track and polarizaion
39 //
40 // ------------------------------------------------------------
41 //
43 #include "G4SystemOfUnits.hh"
44 #include "G4Track.hh"
45 #include "G4Step.hh"
46 #include "G4DynamicParticle.hh"
47 #include "G4ExceptionSeverity.hh"
48 
50  : G4VParticleChange(), currentTrack(0), proposedKinEnergy(0.)
51 {
53  debugFlag = false;
54 #ifdef G4VERBOSE
55  if (verboseLevel>2) {
56  G4cout << "G4ParticleChangeForGamma::G4ParticleChangeForGamma() " << G4endl;
57  }
58 #endif
59 }
60 
62 {
63 #ifdef G4VERBOSE
64  if (verboseLevel>2) {
65  G4cout << "G4ParticleChangeForGamma::~G4ParticleChangeForGamma() " << G4endl;
66  }
67 #endif
68 }
69 
72 {
73  if (verboseLevel>1)
74  G4cout << "G4ParticleChangeForGamma:: copy constructor is called " << G4endl;
75 
76  currentTrack = right.currentTrack;
77  proposedKinEnergy = right.proposedKinEnergy;
78  proposedMomentumDirection = right.proposedMomentumDirection;
79  proposedPolarization = right.proposedPolarization;
80 }
81 
82 // assignment operator
85 {
86 #ifdef G4VERBOSE
87  if (verboseLevel>1)
88  G4cout << "G4ParticleChangeForGamma:: assignment operator is called " << G4endl;
89 #endif
90 
91  if (this != &right) {
92  if (theNumberOfSecondaries>0) {
93 #ifdef G4VERBOSE
94  if (verboseLevel>0) {
95  G4cout << "G4ParticleChangeForGamma: assignment operator Warning ";
96  G4cout << "theListOfSecondaries is not empty ";
97  }
98 #endif
100  if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
101  }
102  }
103  delete theListOfSecondaries;
107  G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
108  theListOfSecondaries->SetElement(index, newTrack); }
109 
114 
115  currentTrack = right.currentTrack;
116  proposedKinEnergy = right.proposedKinEnergy;
117  proposedMomentumDirection = right.proposedMomentumDirection;
118  proposedPolarization = right.proposedPolarization;
119  }
120  return *this;
121 }
122 
124 {
125  // create track
126  G4Track* aTrack = new G4Track(aParticle, currentTrack->GetGlobalTime(),
127  currentTrack->GetPosition());
128 
129  // Touchable handle is copied to keep the pointer
130  aTrack->SetTouchableHandle(currentTrack->GetTouchableHandle());
131 
132  // add a secondary
134 }
135 
136 //----------------------------------------------------------------
137 // method for updating G4Step
138 //
139 
141 {
143  pStep->SetStepLength( 0.0 );
144 
147  }
148 
149  return pStep;
150 }
151 
153 {
154  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
155  G4Track* pTrack = pStep->GetTrack();
156 
157  pPostStepPoint->SetKineticEnergy( proposedKinEnergy );
158  pPostStepPoint->SetMomentumDirection( proposedMomentumDirection );
159  pPostStepPoint->SetPolarization( proposedPolarization );
160 
161  // update velocity for scattering process and particles with mass
162  if(proposedKinEnergy > 0.0) {
163  if(pTrack->GetParticleDefinition()->GetPDGMass() > 0.0) {
164  pPostStepPoint->SetVelocity(pTrack->CalculateVelocity());
165  }
166  }
167 
169  pPostStepPoint->SetWeight( theParentWeight );
170  }
171 
174  return pStep;
175 }
176 
177 //----------------------------------------------------------------
178 // methods for printing messages
179 //
180 
182 {
183 // use base-class DumpInfo
185 
186  G4int oldprc = G4cout.precision(3);
187  G4cout << " Kinetic Energy (MeV): "
188  << std::setw(20) << proposedKinEnergy/MeV
189  << G4endl;
190  G4cout << " Momentum Direction: "
191  << std::setw(20) << proposedMomentumDirection
192  << G4endl;
193  G4cout << " Polarization: "
194  << std::setw(20) << proposedPolarization
195  << G4endl;
196  G4cout.precision(oldprc);
197 }
198 
200 {
201  G4bool itsOK = true;
202  G4bool exitWithError = false;
203 
204  G4double accuracy;
205 
206  // Energy should not be lager than initial value
207  accuracy = ( proposedKinEnergy - aTrack.GetKineticEnergy())/MeV;
208  if (accuracy > accuracyForWarning) {
209  itsOK = false;
210  exitWithError = (accuracy > accuracyForException);
211 #ifdef G4VERBOSE
212  G4cout << "G4ParticleChangeForGamma::CheckIt: ";
213  G4cout << "KinEnergy become larger than the initial value!"
214  << " Difference: " << accuracy << "[MeV] " <<G4endl;
215  G4cout << aTrack.GetDefinition()->GetParticleName()
216  << " E=" << aTrack.GetKineticEnergy()/MeV
217  << " pos=" << aTrack.GetPosition().x()/m
218  << ", " << aTrack.GetPosition().y()/m
219  << ", " << aTrack.GetPosition().z()/m
220  << G4endl;
221 #endif
222  }
223 
224  // dump out information of this particle change
225 #ifdef G4VERBOSE
226  if (!itsOK) DumpInfo();
227 #endif
228 
229  // Exit with error
230  if (exitWithError) {
231  G4Exception("G4ParticleChangeForGamma::CheckIt",
232  "TRACK004",EventMustBeAborted,
233  "energy was illegal");
234  }
235 
236  //correction
237  if (!itsOK) {
238  proposedKinEnergy = aTrack.GetKineticEnergy();
239  }
240 
241  itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
242  return itsOK;
243 }
G4ParticleDefinition * GetDefinition() const
void SetStepLength(G4double value)
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
void AddNonIonizingEnergyDeposit(G4double value)
Int_t index
Definition: macro.C:9
double x() const
G4ParticleChangeForGamma & operator=(const G4ParticleChangeForGamma &right)
const G4ThreeVector & GetPosition() const
G4TrackFastVector * theListOfSecondaries
void SetWeight(G4double aValue)
virtual G4bool CheckIt(const G4Track &)
void SetTouchableHandle(const G4TouchableHandle &apValue)
virtual void DumpInfo() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4Step * UpdateStepForPostStep(G4Step *Step)
double z() const
static const G4double accuracyForException
void AddSecondary(G4DynamicParticle *aParticle)
void SetMomentumDirection(const G4ThreeVector &aValue)
virtual G4bool CheckIt(const G4Track &)
G4double GetKineticEnergy() const
void SetPolarization(const G4ThreeVector &aValue)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4SteppingControl theSteppingControlFlag
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Step.hh:76
G4double GetGlobalTime() const
G4double CalculateVelocity() const
Definition: G4Track.cc:214
void AddSecondary(G4Track *aSecondary)
const G4TouchableHandle & GetTouchableHandle() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetPDGMass() const
void SetVelocity(G4double v)
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
G4StepPoint * GetPostStepPoint() const
void AddTotalEnergyDeposit(G4double value)
double y() const
#define G4endl
Definition: G4ios.hh:61
G4TrackStatus theStatusChange
double G4double
Definition: G4Types.hh:76
G4Step * UpdateStepForAtRest(G4Step *pStep)
static const G4double accuracyForWarning
void SetKineticEnergy(const G4double aValue)
G4Track * GetTrack() const
G4double theNonIonizingEnergyDeposit