Geant4 Cross Reference |
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 // G4ParticleChangeForGamma class implementation 27 // 28 // Author: Hisaya Kurashige, 23 March 1998 29 // Revision: Vladimir Ivantchenko, 15 April 2005 30 // -------------------------------------------------------------------- 31 32 #include "G4ParticleChangeForGamma.hh" 33 #include "G4SystemOfUnits.hh" 34 #include "G4Track.hh" 35 #include "G4Step.hh" 36 #include "G4DynamicParticle.hh" 37 #include "G4ExceptionSeverity.hh" 38 39 // -------------------------------------------------------------------- 40 G4ParticleChangeForGamma::G4ParticleChangeForGamma() 41 { 42 // Disable flag to avoid check of each secondary at each step 43 debugFlag = false; 44 SetNumberOfSecondaries(2); 45 } 46 47 // -------------------------------------------------------------------- 48 void G4ParticleChangeForGamma::AddSecondary(G4DynamicParticle* aParticle) 49 { 50 // create track 51 G4Track* aTrack = new G4Track(aParticle, theCurrentTrack->GetGlobalTime(), 52 theCurrentTrack->GetPosition()); 53 54 // touchable handle is copied to keep the pointer 55 aTrack->SetTouchableHandle(theCurrentTrack->GetTouchableHandle()); 56 57 // add a secondary 58 G4VParticleChange::AddSecondary(aTrack); 59 } 60 61 // -------------------------------------------------------------------- 62 G4Step* G4ParticleChangeForGamma::UpdateStepForAtRest(G4Step* pStep) 63 { 64 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit); 65 pStep->SetStepLength(0.0); 66 67 if(isParentWeightProposed) 68 { 69 pStep->GetPostStepPoint()->SetWeight(theParentWeight); 70 } 71 72 return pStep; 73 } 74 75 // -------------------------------------------------------------------- 76 G4Step* G4ParticleChangeForGamma::UpdateStepForPostStep(G4Step* pStep) 77 { 78 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint(); 79 80 pPostStepPoint->SetMomentumDirection(proposedMomentumDirection); 81 pPostStepPoint->SetPolarization(proposedPolarization); 82 83 // update velocity for scattering process and particles with mass 84 if(proposedKinEnergy > 0.0) 85 { 86 pPostStepPoint->SetKineticEnergy(proposedKinEnergy); 87 G4double mass = theCurrentTrack->GetDefinition()->GetPDGMass(); 88 G4double v = CLHEP::c_light; 89 if(mass > 0.0) { 90 v *= std::sqrt(proposedKinEnergy*(proposedKinEnergy + 2*mass))/ 91 (proposedKinEnergy + mass); 92 } 93 pPostStepPoint->SetVelocity(v); 94 } 95 else 96 { 97 pPostStepPoint->SetKineticEnergy(0.0); 98 pPostStepPoint->SetVelocity(0.0); 99 } 100 101 if(isParentWeightProposed) 102 { 103 pPostStepPoint->SetWeight(theParentWeight); 104 } 105 106 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit); 107 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit); 108 return pStep; 109 } 110 111 // -------------------------------------------------------------------- 112 void G4ParticleChangeForGamma::DumpInfo() const 113 { 114 // use base-class DumpInfo 115 G4VParticleChange::DumpInfo(); 116 117 G4long oldprc = G4cout.precision(8); 118 G4cout << " -----------------------------------------------" << G4endl; 119 G4cout << " G4ParticleChangeForGamma proposes: " << G4endl; 120 G4cout << " Kinetic Energy (MeV): " << std::setw(20) 121 << proposedKinEnergy / MeV << G4endl; 122 G4cout << " Momentum Direction: " << std::setw(20) 123 << proposedMomentumDirection << G4endl; 124 G4cout << " Polarization: " << std::setw(20) << proposedPolarization 125 << G4endl; 126 G4cout.precision(oldprc); 127 } 128 129