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 XraySPOSteppingAction.cc 27 /// \brief Implementation of the SteppingAction class 28 // 29 // Authors: P.Dondero (paolo.dondero@cern.ch), R.Stanzani (ronny.stanzani@cern.ch) 30 // 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 33 34 35 #include "XraySPOSteppingAction.hh" 36 #include "XraySPODetectorConstruction.hh" 37 #include "XraySPOHistoManager.hh" 38 #include "G4SteppingManager.hh" 39 #include "G4RunManager.hh" 40 #include "Randomize.hh" 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 44 void XraySPOSteppingAction::UserSteppingAction(const G4Step* step) 45 { 46 G4int parentID = step->GetTrack()->GetParentID(); 47 48 // Only primaries 49 if (parentID == 0) 50 { 51 G4RunManager *rm = G4RunManager::GetRunManager(); 52 G4int eventID = rm->GetCurrentEvent()->GetEventID(); 53 54 G4double len_fraction = 0.1; // approx 1/20 of the lenght of the pore. Particles scatters an avg of 5 times, so by selecting all events with a distance >10mm we have all the scatters inside. 55 G4bool proceed = false; 56 57 if (fPrevEventID != eventID) 58 { 59 fNumReflections = -2; 60 } 61 62 G4ThreeVector direction = step->GetTrack()->GetMomentumDirection(); 63 auto theta = (G4double)direction.theta(); 64 auto phi = (G4double)direction.phi(); 65 66 G4double x = step->GetPreStepPoint()->GetPosition().x(); 67 G4double y = step->GetPreStepPoint()->GetPosition().y(); 68 G4double z = step->GetPreStepPoint()->GetPosition().z(); 69 70 const G4String& vol_name = step->GetTrack()->GetVolume()->GetName(); 71 G4int trackid = step->GetTrack()->GetTrackID(); 72 73 G4StepPoint * postStep = step->GetPostStepPoint(); 74 G4StepPoint* post_point = step->GetPostStepPoint(); 75 76 const G4String& proc_name = postStep->GetProcessDefinedStep()->GetProcessName(); 77 78 // Count the number of reflections 79 G4String s = ""; 80 G4String& post_phys_name = s; 81 if (fPrevEventID == eventID && step->GetStepLength() >= len_fraction) 82 { 83 if (post_point != nullptr) 84 { 85 G4VPhysicalVolume* post_phys = post_point->GetPhysicalVolume(); 86 if (post_phys != nullptr) 87 { 88 post_phys_name = post_phys->GetName(); 89 } 90 91 if (post_phys_name == "cone_1" || post_phys_name == "cone_2" || post_phys_name == "cone_3" || post_phys_name == "cone_4" || post_phys_name == "cone_5" || post_phys_name == "cone_6" || post_phys_name == "cone_7" || post_phys_name == "cone_8") 92 { 93 if (fNumReflections < 0) 94 { 95 fNumReflections = 1; 96 } 97 else 98 { 99 fNumReflections += 1; 100 } 101 proceed = true; 102 } 103 } 104 } 105 106 if (vol_name == "pDummyEntrance" || vol_name == "pDummyExit" || vol_name == "pDummySphere") 107 { 108 proceed = true; 109 } 110 111 if (proceed) 112 { 113 // G4cout << "Save with " << fNumReflections << " reflections" << G4endl; 114 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 115 analysisManager->FillNtupleIColumn(0, eventID); 116 analysisManager->FillNtupleSColumn(1, vol_name); 117 analysisManager->FillNtupleIColumn(2, trackid); 118 analysisManager->FillNtupleDColumn(3, x); 119 analysisManager->FillNtupleDColumn(4, y); 120 analysisManager->FillNtupleDColumn(5, z); 121 analysisManager->FillNtupleDColumn(6, theta); 122 analysisManager->FillNtupleDColumn(7, phi); 123 analysisManager->FillNtupleSColumn(8, proc_name); 124 analysisManager->FillNtupleIColumn(9, parentID); 125 analysisManager->FillNtupleIColumn(10, fNumReflections); 126 analysisManager->AddNtupleRow(); 127 } 128 fPrevx = x; 129 fPrevy = y; 130 fPrevz = z; 131 fPrevEventID = eventID; 132 } 133 } 134 135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 136