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 electromagnetic/TestEm5/src/StackingAction.cc 27 /// \brief Implementation of the StackingAction class 28 // 29 // 30 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 33 #include "StackingAction.hh" 34 35 #include "EventAction.hh" 36 #include "HistoManager.hh" 37 #include "Run.hh" 38 #include "StackingMessenger.hh" 39 40 #include "G4EmSecondaryParticleType.hh" 41 #include "G4RunManager.hh" 42 #include "G4Track.hh" 43 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 45 46 StackingAction::StackingAction(EventAction* event) : fEventAction(event) 47 { 48 fStackMessenger = new StackingMessenger(this); 49 } 50 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 52 53 StackingAction::~StackingAction() 54 { 55 delete fStackMessenger; 56 } 57 58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 59 60 G4ClassificationOfNewTrack StackingAction::ClassifyNewTrack(const G4Track* aTrack) 61 { 62 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 63 64 // keep primary particle 65 if (aTrack->GetParentID() == 0) { 66 return fUrgent; 67 } 68 69 G4int procID = aTrack->GetCreatorProcess()->GetProcessSubType(); 70 G4int modelID = aTrack->GetCreatorModelID(); 71 72 // count secondary particles 73 Run* run = static_cast<Run*>(G4RunManager::GetRunManager()->GetNonConstCurrentRun()); 74 run->CountParticles(aTrack->GetDefinition()); 75 /* 76 G4cout << "###StackingAction: new " 77 << aTrack->GetDefinition()->GetParticleName() 78 << " E(MeV)= " << aTrack->GetKineticEnergy() 79 << " " << aTrack->GetMomentumDirection() << G4endl; 80 */ 81 // 82 // energy spectrum of secondaries 83 // 84 G4double energy = aTrack->GetKineticEnergy(); 85 G4double charge = aTrack->GetDefinition()->GetPDGCharge(); 86 87 if (charge != 0.) { 88 analysisManager->FillH1(2, energy); 89 analysisManager->FillH1(4, energy); 90 if (procID >= 51 && procID <= 65) { 91 analysisManager->FillH1(58, energy); 92 analysisManager->FillH1(60, energy); 93 } 94 else if (_AugerElectron == modelID) { 95 analysisManager->FillH1(50, energy); 96 analysisManager->FillH1(52, energy); 97 } 98 else if (_ePIXE == modelID) { 99 analysisManager->FillH1(54, energy); 100 analysisManager->FillH1(56, energy); 101 } 102 } 103 104 if (aTrack->GetDefinition() == G4Gamma::Gamma()) { 105 analysisManager->FillH1(3, energy); 106 analysisManager->FillH1(5, energy); 107 if (procID >= 51 && procID <= 65) { 108 analysisManager->FillH1(59, energy); 109 analysisManager->FillH1(61, energy); 110 } 111 else if (_Fluorescence == modelID) { 112 analysisManager->FillH1(51, energy); 113 analysisManager->FillH1(53, energy); 114 } 115 else if (_GammaPIXE == modelID) { 116 analysisManager->FillH1(55, energy); 117 analysisManager->FillH1(57, energy); 118 } 119 } 120 121 // stack or delete secondaries 122 G4ClassificationOfNewTrack status = fUrgent; 123 if (0 < fKillSecondary) { 124 if (fKillSecondary == 1) { 125 fEventAction->AddEnergy(energy); 126 } 127 status = fKill; 128 } 129 130 return status; 131 } 132 133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 134