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 /// \file SteppingAction.cc 27 /// \brief Implementation of the SteppingAction class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 33 #include "SteppingAction.hh" 34 35 #include "DetectorConstruction.hh" 36 #include "EventAction.hh" 37 #include "HistoManager.hh" 38 #include "Run.hh" 39 40 #include "G4RunManager.hh" 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 44 SteppingAction::SteppingAction(DetectorConstruction* det, EventAction* evt) 45 : fDetector(det), fEventAct(evt) 46 {} 47 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 49 50 void SteppingAction::UserSteppingAction(const G4Step* aStep) 51 { 52 // count processes 53 // 54 const G4StepPoint* prePoint = aStep->GetPreStepPoint(); 55 const G4StepPoint* endPoint = aStep->GetPostStepPoint(); 56 const G4VProcess* process = endPoint->GetProcessDefinedStep(); 57 Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun()); 58 run->CountProcesses(process); 59 60 // if World, return 61 // 62 G4VPhysicalVolume* volume = prePoint->GetTouchableHandle()->GetVolume(); 63 // if sum of absorbers do not fill exactly a layer: check material, not volume. 64 const G4Material* mat = volume->GetLogicalVolume()->GetMaterial(); 65 if (mat == fDetector->GetWorldMaterial()) return; 66 67 const G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition(); 68 69 // here we are in an absorber. Locate it 70 // 71 G4int absorNum = prePoint->GetTouchableHandle()->GetCopyNumber(0); 72 G4int layerNum = prePoint->GetTouchableHandle()->GetCopyNumber(1); 73 74 // collect energy deposit (taking into account track weight) 75 G4double edep = aStep->GetTotalEnergyDeposit(); 76 ////edep *= (aStep->GetTrack()->GetWeight()); 77 78 // collect step length of charged particles 79 G4double stepl = 0.; 80 if (particle->GetPDGCharge() != 0.) stepl = aStep->GetStepLength(); 81 82 // sum up per event 83 fEventAct->SumEnergy(absorNum, edep, stepl); 84 85 // longitudinal profile of edep per absorber 86 if (edep > 0.) { 87 G4AnalysisManager::Instance()->FillH1(kMaxAbsor + absorNum, G4double(layerNum + 1), edep); 88 } 89 90 // energy flow 91 // 92 // unique identificator of layer+absorber 93 G4int Idnow = (fDetector->GetNbOfAbsor()) * layerNum + absorNum; 94 G4int plane; 95 // 96 // leaving the absorber ? 97 if (endPoint->GetStepStatus() == fGeomBoundary) { 98 G4ThreeVector position = endPoint->GetPosition(); 99 G4ThreeVector direction = endPoint->GetMomentumDirection(); 100 G4double Eflow = endPoint->GetKineticEnergy(); 101 if (direction.x() >= 0.) 102 run->SumEnergyFlow(plane = Idnow + 1, Eflow); 103 else 104 run->SumEnergyFlow(plane = Idnow, -Eflow); 105 } 106 107 //// example of Birk attenuation 108 /// G4double destep = aStep->GetTotalEnergyDeposit(); 109 /// G4double response = BirksAttenuation(aStep); 110 /// G4cout << " Destep: " << destep/keV << " keV" 111 /// << " response after Birks: " << response/keV << " keV" << G4endl; 112 } 113 114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 115 116 G4double SteppingAction::BirksAttenuation(const G4Step* aStep) 117 { 118 // Example of Birk attenuation law in organic scintillators. 119 // adapted from Geant3 PHYS337. See MIN 80 (1970) 239-244 120 // 121 const G4Material* material = aStep->GetTrack()->GetMaterial(); 122 G4double birk1 = material->GetIonisation()->GetBirksConstant(); 123 G4double destep = aStep->GetTotalEnergyDeposit(); 124 G4double stepl = aStep->GetStepLength(); 125 G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge(); 126 // 127 G4double response = destep; 128 if (birk1 * destep * stepl * charge != 0.) { 129 response = destep / (1. + birk1 * destep / stepl); 130 } 131 return response; 132 } 133 134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 135