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 "EventAction.hh" 28 29 #include "SensitiveDetectorHit.hh" 30 31 #include "G4AnalysisManager.hh" 32 #include "G4Event.hh" 33 #include "G4EventManager.hh" 34 #include "G4HCofThisEvent.hh" 35 #include "G4RunManager.hh" 36 #include "G4SDManager.hh" 37 #include "G4SystemOfUnits.hh" 38 #include "G4Trajectory.hh" 39 #include "G4TrajectoryContainer.hh" 40 #include "G4UImanager.hh" 41 #include "G4VHitsCollection.hh" 42 #include "G4VVisManager.hh" 43 #include "G4ios.hh" 44 45 EventAction::EventAction() : fSDHT_ID(-1) {} 46 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 48 49 EventAction::~EventAction() 50 { 51 ; 52 } 53 54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 55 56 void EventAction::BeginOfEventAction(const G4Event*) 57 { 58 ; 59 } 60 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 62 63 void EventAction::EndOfEventAction(const G4Event* evt) 64 { 65 G4SDManager* SDman = G4SDManager::GetSDMpointer(); 66 67 G4ThreeVector ssd[3]; 68 ssd[0] = G4ThreeVector(0., 0., 0.); 69 ssd[1] = G4ThreeVector(0., 0., 0.); 70 ssd[2] = G4ThreeVector(0., 0., 0.); 71 72 if (fSDHT_ID == -1) { 73 G4String sdName; 74 if (SDman->FindSensitiveDetector(sdName = "telescope", 0)) { 75 fSDHT_ID = SDman->GetCollectionID(sdName = "telescope/collection"); 76 } 77 } 78 79 SensitiveDetectorHitsCollection* sdht = 0; 80 G4HCofThisEvent* hce = evt->GetHCofThisEvent(); 81 82 if (hce) { 83 if (fSDHT_ID != -1) { 84 G4VHitsCollection* aHCSD = hce->GetHC(fSDHT_ID); 85 sdht = (SensitiveDetectorHitsCollection*)(aHCSD); 86 } 87 } 88 89 int bTotalHits = 0; 90 if (sdht) { 91 int n_hit_sd = sdht->entries(); 92 for (int i2 = 0; i2 < 3; i2++) { 93 for (int i1 = 0; i1 < n_hit_sd; i1++) { 94 SensitiveDetectorHit* aHit = (*sdht)[i1]; 95 if (aHit->GetLayerID() == i2) { 96 ssd[i2] = aHit->GetWorldPos(); 97 bTotalHits++; 98 } 99 } 100 } 101 } 102 103 if (bTotalHits > 2) { 104 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 105 G4double angXin = (ssd[1].x() - ssd[0].x()) / (ssd[1].z() - ssd[0].z()); 106 G4double angYin = (ssd[1].y() - ssd[0].y()) / (ssd[1].z() - ssd[0].z()); 107 108 analysisManager->FillNtupleDColumn(0, angXin * 1.E6 * CLHEP::rad); 109 analysisManager->FillNtupleDColumn(1, angYin * 1.E6 * CLHEP::rad); 110 111 double posXin = ssd[1].x() - angXin * ssd[1].z(); 112 double posYin = ssd[1].y() - angYin * ssd[1].z(); 113 114 analysisManager->FillNtupleDColumn(2, posXin / CLHEP::mm); 115 analysisManager->FillNtupleDColumn(3, posYin / CLHEP::mm); 116 117 G4double angXout = (ssd[2].x() - posXin) / (ssd[2].z()); 118 G4double angYout = (ssd[2].y() - posYin) / (ssd[2].z()); 119 analysisManager->FillNtupleDColumn(4, angXout * 1.E6 * CLHEP::rad); 120 analysisManager->FillNtupleDColumn(5, angYout * 1.E6 * CLHEP::rad); 121 122 analysisManager->AddNtupleRow(); 123 } 124 } 125 126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 127