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 // ------------------------------------------------------------ 28 // GEANT 4 class implementation file 29 // ------------------------------------------------------------ 30 // 31 32 #include "G4ErrorTrackLengthTarget.hh" 33 34 #include "G4ParticleTable.hh" 35 #include "G4ParticleDefinition.hh" 36 #include "G4VProcess.hh" 37 #include "G4ProcessVector.hh" 38 #include "G4ProcessManager.hh" 39 40 #ifdef G4VERBOSE 41 # include "G4ErrorPropagatorData.hh" //for verbosity checking 42 #endif 43 44 //---------------------------------------------------------------------------- 45 G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget(const G4double maxTrkLength) 46 : G4VDiscreteProcess("G4ErrorTrackLengthTarget") 47 , theMaximumTrackLength(maxTrkLength) 48 { 49 theType = G4ErrorTarget_TrkL; 50 51 G4ParticleTable::G4PTblDicIterator* theParticleIterator = 52 G4ParticleTable::GetParticleTable()->GetIterator(); 53 54 // loop over all particles in G4ParticleTable 55 56 theParticleIterator->reset(); 57 while((*theParticleIterator)()) // Loop checking, 06.08.2015, G.Cosmo 58 { 59 G4ParticleDefinition* particle = theParticleIterator->value(); 60 G4ProcessManager* pmanager = particle->GetProcessManager(); 61 if(!particle->IsShortLived()) 62 { 63 // Add transportation process for all particles other than "shortlived" 64 if(pmanager == 0) 65 { 66 // Error !! no process manager 67 G4String particleName = particle->GetParticleName(); 68 G4Exception("G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget", 69 "No process manager", RunMustBeAborted, particleName); 70 } 71 else 72 { 73 G4ProcessVector* procvec = pmanager->GetProcessList(); 74 std::size_t isiz = (G4int)procvec->size(); 75 76 for(G4int ii = 0; ii < (G4int)isiz; ++ii) 77 { 78 if(((*procvec)[ii])->GetProcessName() == "G4ErrorTrackLengthTarget") 79 { 80 pmanager->RemoveProcess((*procvec)[ii]); 81 } 82 } 83 pmanager->AddDiscreteProcess(this, 4); 84 isiz = procvec->size(); 85 } 86 } 87 else 88 { 89 // shortlived particle case 90 } 91 } 92 } 93 94 //----------------------------------------------------------------------- 95 G4double G4ErrorTrackLengthTarget::PostStepGetPhysicalInteractionLength( 96 const G4Track& track, G4double, G4ForceCondition* condition) 97 { 98 *condition = NotForced; 99 return GetMeanFreePath(track, 0., condition); 100 } 101 102 //----------------------------------------------------------------------- 103 G4double G4ErrorTrackLengthTarget::GetMeanFreePath(const class G4Track& track, 104 G4double, 105 enum G4ForceCondition*) 106 { 107 #ifdef G4VERBOSE 108 if(G4ErrorPropagatorData::verbose() >= 3) 109 { 110 G4cout << " G4ErrorTrackLengthTarget::GetMeanFreePath " 111 << theMaximumTrackLength - track.GetTrackLength() << G4endl; 112 } 113 #endif 114 115 return theMaximumTrackLength - track.GetTrackLength(); 116 } 117 118 G4VParticleChange* G4ErrorTrackLengthTarget::PostStepDoIt(const G4Track& aTrack, 119 const G4Step&) 120 { 121 theParticleChange.Initialize(aTrack); 122 return &theParticleChange; 123 } 124 125 //----------------------------------------------------------------------- 126 void G4ErrorTrackLengthTarget::Dump(const G4String& msg) const 127 { 128 G4cout << msg << "G4ErrorTrackLengthTarget: max track length = " 129 << theMaximumTrackLength << G4endl; 130 } 131