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 eventgenerator/particleGun/src/TrackingAction.cc 27 /// \brief Implementation of the TrackingAction class 28 // 29 // 30 // 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 33 34 #include "TrackingAction.hh" 35 36 #include "HistoManager.hh" 37 #include "PrimaryGeneratorAction.hh" 38 #include "PrimaryGeneratorAction2.hh" 39 #include "PrimaryGeneratorAction3.hh" 40 #include "PrimaryGeneratorAction4.hh" 41 #include "RunAction.hh" 42 43 #include "G4PhysicalConstants.hh" 44 #include "G4Track.hh" 45 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 47 48 TrackingAction::TrackingAction(PrimaryGeneratorAction* prim) : fPrimary(prim) 49 { 50 // parameters for generator action #3 51 fNewUz = fPrimary->GetAction3()->GetNewUz(); 52 } 53 54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 55 56 void TrackingAction::PreUserTrackingAction(const G4Track* track) 57 { 58 G4int selectedGeneratorAction = fPrimary->GetSelectedAction(); 59 G4AnalysisManager* analysis = G4AnalysisManager::Instance(); 60 G4int id = 0; 61 62 if (selectedGeneratorAction == 2) { 63 // energy spectrum 64 // 65 id = 1; 66 analysis->FillH1(id, track->GetKineticEnergy()); 67 } 68 69 else if (selectedGeneratorAction == 0) { 70 // particle direction: cos(alpha) 71 // 72 id = 5; 73 G4ThreeVector um = track->GetMomentumDirection(); 74 G4double cosalpha = um.z(); 75 analysis->FillH1(id, cosalpha); 76 77 // particle direction: psi 78 // 79 id = 6; 80 G4double psi = std::atan2(um.y(), um.x()); 81 if (psi < 0.) psi += twopi; 82 analysis->FillH1(id, psi); 83 } 84 else if (selectedGeneratorAction == 3) { 85 // particle direction in local frame: cos(alpha) 86 // 87 id = 5; 88 G4ThreeVector um = track->GetMomentumDirection(); 89 G4double cosalpha = fNewUz * um; 90 analysis->FillH1(id, cosalpha); 91 92 // particle direction in local frame: psi 93 // 94 id = 6; 95 // complete local frame 96 G4ThreeVector u1(1., 0., 0.); 97 u1.rotateUz(fNewUz); 98 G4ThreeVector u2(0., 1., 0.); 99 u2.rotateUz(fNewUz); 100 // 101 G4double psi = std::atan2(u2 * um, u1 * um); 102 if (psi < 0.) psi += twopi; 103 analysis->FillH1(id, psi); 104 } 105 106 else if (selectedGeneratorAction == 4) { 107 G4ThreeVector vertex = track->GetVertexPosition(); 108 G4double r = vertex.mag(); 109 if (r <= 0.0) return; 110 // local frame : new uz = ur 111 G4ThreeVector ur = vertex / r; 112 113 // vertex position. radial distribution dN/dv = f(r) 114 // 115 id = 2; 116 G4double dr = analysis->GetH1Width(id); 117 G4double dv = 2 * twopi * r * r * dr; 118 if (dv > 0.) analysis->FillH1(id, r, 1. / dv); 119 120 // vertex position: cos(theta) 121 // 122 id = 3; 123 G4double costheta = ur.z(); 124 analysis->FillH1(id, costheta); 125 126 // vertex position: phi 127 // 128 id = 4; 129 G4double phi = std::atan2(ur.y(), ur.x()); 130 if (phi < 0.) phi += twopi; 131 analysis->FillH1(id, phi); 132 133 // particle direction in local frame: cos(alpha) 134 // 135 id = 5; 136 G4ThreeVector um = track->GetMomentumDirection(); 137 G4double cosalpha = ur * um; 138 analysis->FillH1(id, cosalpha); 139 140 // particle direction in local frame: psi 141 // 142 id = 6; 143 // complete local frame 144 G4ThreeVector u1(1., 0., 0.); 145 u1.rotateUz(ur); 146 G4ThreeVector u2(0., 1., 0.); 147 u2.rotateUz(ur); 148 // 149 G4double psi = std::atan2(u2 * um, u1 * um); 150 if (psi < 0.) psi += twopi; 151 analysis->FillH1(id, psi); 152 } 153 } 154 155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 156 157 void TrackingAction::PostUserTrackingAction(const G4Track*) {} 158 159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 160