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 27 #include "G4ErrorEnergyLoss.hh" 28 #include "G4EnergyLossForExtrapolator.hh" 29 #include "G4ErrorPropagatorData.hh" 30 31 //------------------------------------------------------------------- 32 G4ErrorEnergyLoss::G4ErrorEnergyLoss(const G4String& processName, 33 G4ProcessType type) 34 : G4VContinuousProcess(processName, type) 35 { 36 if (verboseLevel>2) { 37 G4cout << GetProcessName() << " is created " << G4endl; 38 } 39 40 theELossForExtrapolator = new G4EnergyLossForExtrapolator; 41 theStepLimit = 1.*CLHEP::mm; 42 } 43 44 //------------------------------------------------------------------- 45 void G4ErrorEnergyLoss::InstantiateEforExtrapolator() 46 {} 47 48 //------------------------------------------------------------------- 49 G4ErrorEnergyLoss::~G4ErrorEnergyLoss() 50 { 51 delete theELossForExtrapolator; 52 } 53 54 //------------------------------------------------------------------- 55 56 G4bool G4ErrorEnergyLoss::IsApplicable(const G4ParticleDefinition& aParticleType) 57 { 58 return (aParticleType.GetPDGCharge() != 0); 59 } 60 61 //------------------------------------------------------------------- 62 G4VParticleChange* 63 G4ErrorEnergyLoss::AlongStepDoIt(const G4Track& aTrack, const G4Step& aStep) 64 { 65 aParticleChange.Initialize(aTrack); 66 67 G4ErrorPropagatorData* g4edata = G4ErrorPropagatorData::GetErrorPropagatorData(); 68 69 G4double kinEnergyStart = aTrack.GetKineticEnergy(); 70 G4double step_length = aStep.GetStepLength(); 71 72 const G4Material* aMaterial = aTrack.GetMaterial(); 73 const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition(); 74 G4double kinEnergyEnd = kinEnergyStart; 75 76 // backward - energy increased 77 if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) { 78 kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart, 79 step_length, 80 aMaterial, 81 aParticleDef ); 82 G4double kinEnergyHalfStep = (kinEnergyStart + kinEnergyEnd) * 0.5; 83 84 #ifdef G4VERBOSE 85 if(G4ErrorPropagatorData::verbose() >= 3 ) 86 G4cout << " G4ErrorEnergyLoss FWD end " << kinEnergyEnd 87 << " halfstep " << kinEnergyHalfStep << G4endl; 88 #endif 89 90 //--- rescale to energy lost at 1/2 step 91 kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyHalfStep, 92 step_length, 93 aMaterial, 94 aParticleDef ); 95 kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd ); 96 97 // forward - energy decreased 98 } else { 99 kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyStart, 100 step_length, 101 aMaterial, 102 aParticleDef ); 103 G4double kinEnergyHalfStep = (kinEnergyStart + kinEnergyEnd) * 0.5; 104 #ifdef G4VERBOSE 105 if(G4ErrorPropagatorData::verbose() >= 3 ) 106 G4cout << " G4ErrorEnergyLoss BCKD end " << kinEnergyEnd 107 << " halfstep " << kinEnergyHalfStep << G4endl; 108 #endif 109 110 //--- rescale to energy lost at 1/2 step 111 kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyHalfStep, 112 step_length, 113 aMaterial, 114 aParticleDef ); 115 kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd ); 116 } 117 118 G4double edepo = kinEnergyEnd - kinEnergyStart; 119 120 #ifdef G4VERBOSE 121 if( G4ErrorPropagatorData::verbose() >= 2 ) 122 G4cout << "AlongStepDoIt Estart= " << kinEnergyStart << " Eend " << kinEnergyEnd 123 << " Ediff " << kinEnergyStart-kinEnergyEnd << " step= " << step_length 124 << " mate= " << aMaterial->GetName() 125 << " particle= " << aParticleDef->GetParticleName() << G4endl; 126 #endif 127 128 aParticleChange.ClearDebugFlag(); 129 aParticleChange.ProposeLocalEnergyDeposit( edepo ); 130 aParticleChange.SetNumberOfSecondaries(0); 131 132 aParticleChange.ProposeEnergy( kinEnergyEnd ); 133 134 return &aParticleChange; 135 } 136 137 138 //------------------------------------------------------------------- 139 G4double G4ErrorEnergyLoss::GetContinuousStepLimit(const G4Track& aTrack, 140 G4double, G4double, G4double& ) 141 { 142 G4double ekin = aTrack.GetKineticEnergy(); 143 const G4Material* mat = aTrack.GetMaterial(); 144 const G4ParticleDefinition* part = 145 aTrack.GetDynamicParticle()->GetDefinition(); 146 G4double range = theELossForExtrapolator->ComputeRange(ekin, part, mat); 147 G4double delta = std::max(range*theFractionLimit, theStepLimit); 148 #ifdef G4VERBOSE 149 if(G4ErrorPropagatorData::verbose() >= 2 ) { 150 G4cout << " G4ErrorEnergyLoss: limiting Step " << delta 151 << " energy(GeV) " << ekin / CLHEP::GeV 152 << " for " << part->GetParticleName() << G4endl; 153 } 154 #endif 155 return delta; 156 } 157