Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MSSteppingAction.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//
26// G4MSSteppingAction implementation
27//
28// Author: M.Asai, 5 May 2006
29// --------------------------------------------------------------------
30
31#include "G4MSSteppingAction.hh"
32
33#include "G4LogicalVolume.hh"
34#include "G4Material.hh"
35#include "G4Region.hh"
36#include "G4Step.hh"
37#include "G4VPhysicalVolume.hh"
38
39#include <map>
40
41// --------------------------------------------------------------------
43{
44 regionSensitive = rSens;
45 theRegion = reg;
46 length = 0.;
47 x0 = 0.;
48 lambda = 0.;
49 shape_mat_info_v.clear();
50}
51
52// --------------------------------------------------------------------
54{
55 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
56 G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
57
58 if (regionSensitive && (region != theRegion)) return;
59
60 G4double stlen = aStep->GetStepLength();
61 const G4Material* material = preStepPoint->GetMaterial();
62 length += stlen;
63 x0 += stlen / (material->GetRadlen());
64 lambda += stlen / (material->GetNuclearInterLength());
65
66 // store information per step (1 geantino step = 1 solid)
67 {
68 shape_mat_info_v.push_back({});
69
70 shape_mat_info_t & thisMaterialInfo = shape_mat_info_v.back();
71
72 // calculate average mass number and atomic number
73 {
74 const std::vector<const G4Element*> * ElementVector_ptr = material->GetElementVector();
75 for( auto & element : *ElementVector_ptr)
76 {
77 thisMaterialInfo.aveA += element->GetA();
78 thisMaterialInfo.aveZ += element->GetZ();
79 }
80 thisMaterialInfo.aveA /= ElementVector_ptr->size();
81 thisMaterialInfo.aveZ /= ElementVector_ptr->size();
82 }
83 thisMaterialInfo.density = material->GetDensity();
84 thisMaterialInfo.radiation_length = material->GetRadlen();
85 thisMaterialInfo.interaction_length = material->GetNuclearInterLength();
86 thisMaterialInfo.thickness = aStep->GetStepLength();
87 thisMaterialInfo.integrated_thickness = length;
88 thisMaterialInfo.lambda = stlen / (material->GetNuclearInterLength());
89 thisMaterialInfo.x0 = stlen / (material->GetRadlen());
90 thisMaterialInfo.integrated_lambda = lambda;
91 thisMaterialInfo.integrated_x0 = x0;
92 thisMaterialInfo.entry_point = preStepPoint->GetPosition();
93 thisMaterialInfo.exit_point = aStep->GetPostStepPoint()->GetPosition();
94 thisMaterialInfo.material_name = material->GetName();
95 }
96
97}
98
100{
101
102 G4int colwidth = 11;
103 G4int matname_colwidth = 15;
104
105 oss << G4endl<< G4endl;
106 oss << " Material Atomic properties (averaged) Radiation Interaction Integr. Lambda X0 Entry point Exit point\n name Mass Z density length length Thickness Thickness (cm) (cm)\n (g/mole) (g/cm3) (cm) (cm) (cm) (cm) \n";
107 oss << G4endl;
108
109
110 std::ios::fmtflags os_flags (oss.flags());
111 for( auto & matInfo : shape_mat_info_v)
112 {
113 oss << std::setw(matname_colwidth) << std::left << matInfo.GetName(matname_colwidth) << " ";
114
115 oss << std::setw(colwidth) << std::fixed << std::setprecision(2) << matInfo.aveA / (CLHEP::g/CLHEP::mole);
116 oss << std::setw(colwidth) << std::fixed << std::setprecision(2) << matInfo.aveZ;
117 oss << std::setw(colwidth) << std::scientific << std::setprecision(2) << matInfo.density / (CLHEP::g/CLHEP::cm3);
118 oss << std::setw(colwidth) << std::scientific << matInfo.radiation_length / CLHEP::cm;
119 oss << std::setw(colwidth) << std::scientific << matInfo.interaction_length / CLHEP::cm;
120 oss << std::setw(colwidth) << std::scientific << matInfo.thickness/CLHEP::cm;
121 oss << std::setw(colwidth) << std::scientific << matInfo.integrated_thickness/CLHEP::cm;
122 oss << std::setw(colwidth) << std::scientific << matInfo.lambda;
123 oss << std::setw(colwidth) << std::scientific << matInfo.x0;
124 oss << std::setw(colwidth) << " ";
125 oss << std::scientific << std::right << "(";
126 oss << std::scientific << std::left << matInfo.entry_point.x()/CLHEP::cm;
127 oss << ", " << matInfo.entry_point.y()/CLHEP::cm;
128 oss << ", " << matInfo.entry_point.z()/CLHEP::cm << ")";
129 oss << std::setw(colwidth) << " ";
130 oss << std::scientific << std::right << "(";
131 oss << std::scientific << std::left << matInfo.exit_point.x()/CLHEP::cm;
132 oss << ", " << matInfo.exit_point.y()/CLHEP::cm;
133 oss << ", " << matInfo.exit_point.z()/CLHEP::cm << ")";
134 oss << G4endl;
135 oss << G4endl;
136 }
137 oss.flags(os_flags); // Restore original stream format
138}
139
141{
142 // create database (db) of material name (key) and information
143 std::map<G4String, shape_mat_info_t> mat_db;
144 // accumulate information for each material name into mat_db
145 for(auto & matInfo : shape_mat_info_v)
146 {
147 G4String key = matInfo.material_name;
148 if( 0 == mat_db.count( key ) )
149 {
150 mat_db[key] = shape_mat_info_t{};
151 }
152
153 mat_db[key].x0 += matInfo.x0;
154 mat_db[key].thickness += matInfo.thickness;
155 mat_db[key].lambda += matInfo.lambda;
156 }
157
158 std::ios::fmtflags os_flags (oss.flags());
159 oss << std::scientific << std::setprecision(2) << '\t';
160 for(auto & [key,mat] : mat_db)
161 oss << '\t' << key
162 << '\t'<< mat.thickness/CLHEP::mm
163 << '\t'<< mat.x0
164 << '\t'<< mat.lambda;
165 oss.flags(os_flags); // Restore original stream format
166 return;
167}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4Region * GetRegion() const
void PrintEachMaterialVerbose(std::ostream &oss)
void UserSteppingAction(const G4Step *) override
void PrintIntegratedMaterialVerbose(std::ostream &oss)
void Initialize(G4bool rSens, G4Region *reg)
G4double GetDensity() const
const G4ElementVector * GetElementVector() const
G4double GetRadlen() const
const G4String & GetName() const
G4double GetNuclearInterLength() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4LogicalVolume * GetLogicalVolume() const