Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SteppingAction.cc
Go to the documentation of this file.
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 //
28 //
29 // $Id$
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33 
34 #include "SteppingAction.hh"
35 #include "PrimaryGeneratorAction.hh"
36 #include "RunAction.hh"
37 #include "HistoManager.hh"
38 
39 #include "G4ParticleTypes.hh"
40 #include "G4RunManager.hh"
41 #include "G4HadronicProcess.hh"
42 
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44 
46 :fPrimary(prim),fRunAction(RuAct)
47 { }
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
52 { }
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
57 {
58  // count processes
59  //
60  const G4StepPoint* endPoint = aStep->GetPostStepPoint();
61  const G4VProcess* process = endPoint->GetProcessDefinedStep();
62  fRunAction->CountProcesses(process);
63 
64  // check that an real interaction occured (eg. not a transportation)
65  G4StepStatus stepStatus = endPoint->GetStepStatus();
66  G4bool transmit = (stepStatus==fGeomBoundary || stepStatus==fWorldBoundary);
67  if (transmit) return;
68 
69  //real processes : sum track length
70  //
71  G4double stepLength = aStep->GetStepLength();
72  fRunAction->SumTrack(stepLength);
73 
74  //energy-momentum balance initialisation
75  //
76  const G4StepPoint* prePoint = aStep->GetPreStepPoint();
77  G4double Q = - prePoint->GetKineticEnergy();
78  G4ThreeVector Pbalance = - prePoint->GetMomentum();
79 
80  //initialisation of the nuclear channel identification
81  //
82  G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition();
83  G4String partName = particle->GetParticleName();
84  G4String nuclearChannel = partName;
85  G4HadronicProcess* hproc = (G4HadronicProcess*) process;
86  const G4Isotope* target = hproc->GetTargetIsotope();
87  G4String targetName = "XXXX";
88  if (target) targetName = target->GetName();
89  nuclearChannel += " + " + targetName + " --> ";
90 
91  //scattered primary particle (if any)
92  //
93  G4AnalysisManager* analysis = G4AnalysisManager::Instance();
94  G4int ih = 1;
95  if (aStep->GetTrack()->GetTrackStatus() == fAlive) {
96  G4double energy = endPoint->GetKineticEnergy();
97  analysis->FillH1(ih,energy);
98  //
99  G4ThreeVector momentum = endPoint->GetMomentum();
100  Q += energy;
101  Pbalance += momentum;
102  //
103  nuclearChannel += partName + " + ";
104  }
105 
106  //secondaries
107  //
108  const G4TrackVector* secondary = fpSteppingManager->GetSecondary();
109  for (size_t lp=0; lp<(*secondary).size(); lp++) {
110  particle = (*secondary)[lp]->GetDefinition();
111  G4String name = particle->GetParticleName();
112  G4String type = particle->GetParticleType();
113  G4double charge = particle->GetPDGCharge();
114  G4double energy = (*secondary)[lp]->GetKineticEnergy();
115  fRunAction->ParticleCount(name,energy);
116  //energy spectrum
117  if (charge > 3.) ih = 2;
118  else if (particle == G4Gamma::Gamma()) ih = 3;
119  else if (particle == G4Neutron::Neutron()) ih = 4;
120  else if (particle == G4Proton::Proton()) ih = 5;
121  else if (particle == G4Deuteron::Deuteron()) ih = 6;
122  else if (particle == G4Alpha::Alpha()) ih = 7;
123  else if (type == "nucleus") ih = 8;
124  else if (type == "meson") ih = 9;
125  else if (type == "baryon") ih = 10;
126  analysis->FillH1(ih,energy);
127  //energy-momentum balance
128  G4ThreeVector momentum = (*secondary)[lp]->GetMomentum();
129  Q += energy;
130  Pbalance += momentum;
131  //particle flag
132  fParticleFlag[particle]++;
133  }
134 
135  //energy-momentum balance
136  G4double Pbal = Pbalance.mag();
137  fRunAction->Balance(Pbal);
138  ih = 11;
139  analysis->FillH1(ih,Q);
140  ih = 12;
141  analysis->FillH1(ih,Pbal);
142 
143  // nuclear channel
144  const G4int kMax = 16;
145  const G4String conver[] = {"0","","2 ","3 ","4 ","5 ","6 ","7 ","8 ","9 ",
146  "10 ","11 ","12 ","13 ","14 ","15 ","16 "};
147  std::map<G4ParticleDefinition*,G4int>::iterator ip;
148  for (ip = fParticleFlag.begin(); ip != fParticleFlag.end(); ip++) {
149  particle = ip->first;
150  G4String name = particle->GetParticleName();
151  G4int nb = ip->second;
152  if (nb > kMax) nb = kMax;
153  G4String Nb = conver[nb];
154  if (particle == G4Gamma::Gamma()) {
155  fRunAction->CountGamma(nb);
156  Nb = "N ";
157  }
158  if (ip != fParticleFlag.begin()) nuclearChannel += " + ";
159  nuclearChannel += Nb + name;
160  }
161 
163  fRunAction->CountNuclearChannel(nuclearChannel, Q);
164 
165  fParticleFlag.clear();
166 
167  // kill event after first interaction
168  //
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173 
174