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 // G4VRestDiscreteProcess class implementation 27 // 28 // Authors: 29 // - 2 December 1995, G.Cosmo - First implementation, based on object model 30 // - 8 January 1997, H.Kurashige - New Physics scheme 31 // -------------------------------------------------------------------- 32 33 #include "G4VRestDiscreteProcess.hh" 34 #include "G4SystemOfUnits.hh" 35 36 #include "Randomize.hh" 37 #include "G4Step.hh" 38 #include "G4Track.hh" 39 #include "G4MaterialTable.hh" 40 #include "G4VParticleChange.hh" 41 42 // -------------------------------------------------------------------- 43 G4VRestDiscreteProcess::G4VRestDiscreteProcess() 44 : G4VProcess("No Name Discrete Process") 45 { 46 G4Exception("G4VRestDiscreteProcess::G4VRestDiscreteProcess", "ProcMan102", 47 JustWarning, "Default constructor is called"); 48 } 49 50 // -------------------------------------------------------------------- 51 G4VRestDiscreteProcess:: 52 G4VRestDiscreteProcess(const G4String& aName, G4ProcessType aType) 53 : G4VProcess(aName, aType) 54 { 55 enableAlongStepDoIt = false; 56 } 57 58 // -------------------------------------------------------------------- 59 G4VRestDiscreteProcess::~G4VRestDiscreteProcess() 60 { 61 } 62 63 // -------------------------------------------------------------------- 64 G4VRestDiscreteProcess:: 65 G4VRestDiscreteProcess(G4VRestDiscreteProcess& right) 66 : G4VProcess(right) 67 { 68 } 69 70 // -------------------------------------------------------------------- 71 G4double G4VRestDiscreteProcess:: 72 PostStepGetPhysicalInteractionLength( const G4Track& track, 73 G4double previousStepSize, 74 G4ForceCondition* condition ) 75 { 76 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) 77 { 78 // beginning of tracking (or just after DoIt() of this process) 79 ResetNumberOfInteractionLengthLeft(); 80 } 81 else if ( previousStepSize > 0.0) 82 { 83 // subtract NumberOfInteractionLengthLeft 84 SubtractNumberOfInteractionLengthLeft(previousStepSize); 85 } 86 else 87 { 88 // zero step 89 // DO NOTHING 90 } 91 92 // condition is set to "Not Forced" 93 *condition = NotForced; 94 95 // get mean free path 96 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition); 97 98 G4double value; 99 if (currentInteractionLength < DBL_MAX) 100 { 101 value = theNumberOfInteractionLengthLeft * currentInteractionLength; 102 } 103 else 104 { 105 value = DBL_MAX; 106 } 107 #ifdef G4VERBOSE 108 if (verboseLevel>1) 109 { 110 G4double length = (value < DBL_MAX) ? value/cm : value; 111 G4cout << "G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength() - "; 112 G4cout << "[ " << GetProcessName() << "]" << G4endl; 113 track.GetDynamicParticle()->DumpInfo(); 114 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl; 115 G4cout << "InteractionLength= " << length <<"[cm] " << G4endl; 116 } 117 #endif 118 return value; 119 } 120 121 // -------------------------------------------------------------------- 122 G4VParticleChange* 123 G4VRestDiscreteProcess::PostStepDoIt( const G4Track& , 124 const G4Step& ) 125 { 126 ClearNumberOfInteractionLengthLeft(); 127 128 return pParticleChange; 129 } 130 131 // -------------------------------------------------------------------- 132 G4double G4VRestDiscreteProcess:: 133 AtRestGetPhysicalInteractionLength( const G4Track& track, 134 G4ForceCondition* condition ) 135 { 136 // beginning of tracking 137 ResetNumberOfInteractionLengthLeft(); 138 139 // condition is set to "Not Forced" 140 *condition = NotForced; 141 142 // get mean life time 143 currentInteractionLength = GetMeanLifeTime(track, condition); 144 145 G4double time = (currentInteractionLength < DBL_MAX) ? 146 theNumberOfInteractionLengthLeft * currentInteractionLength : DBL_MAX; 147 148 #ifdef G4VERBOSE 149 if ((currentInteractionLength <0.0) || (verboseLevel>2)) 150 { 151 G4double t = (currentInteractionLength < DBL_MAX) ? 152 currentInteractionLength/ns : DBL_MAX; 153 G4cout << "G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength() - "; 154 G4cout << "[ " << GetProcessName() << "]" << G4endl; 155 track.GetDynamicParticle()->DumpInfo(); 156 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl; 157 G4cout << "MeanLifeTime = " << t << " [ns]" 158 << G4endl; 159 } 160 #endif 161 162 return time; 163 } 164 165 // -------------------------------------------------------------------- 166 G4VParticleChange* 167 G4VRestDiscreteProcess::AtRestDoIt( const G4Track&, 168 const G4Step& ) 169 { 170 ClearNumberOfInteractionLengthLeft(); 171 172 return pParticleChange; 173 } 174