Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CellScoreComposer.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//
27//
28// ----------------------------------------------------------------------
29// GEANT 4 class source file
30//
31// G4CellSCoreComposer.cc
32//
33// ----------------------------------------------------------------------
34
36#include "G4Step.hh"
37
39 : fSCScoreValues()
40{}
41
43
45{
46 G4StepPoint* p = aStep.GetPreStepPoint();
47 if(!p)
48 {
49 G4Exception("G4CellScoreComposer::EstimatorCalculation", "Det0191",
50 FatalException, " no pointer to pre PreStepPoint!");
51 }
52 G4double sl = aStep.GetStepLength();
53 G4double slw = sl * p->GetWeight();
54 G4double slwe = slw * p->GetKineticEnergy();
55
56 G4double v = p->GetVelocity();
57 if(!(v > 0.))
58 {
59 v = 10e-9;
60 }
61
62 fSCScoreValues.fSumSL += sl;
63 fSCScoreValues.fSumSLW += slw;
64 fSCScoreValues.fSumSLW_v += slw / v;
65 fSCScoreValues.fSumSLWE += slwe;
66 fSCScoreValues.fSumSLWE_v += slwe / v;
67}
70
72{
73 fSCScoreValues.fSumCollisions++;
74 fSCScoreValues.fSumCollisionsWeight += weight;
75}
76
78{
79 if(fSCScoreValues.fSumSLW > 0.)
80 {
81 // divide by SumSLW or SumSLW_v ?
82 fSCScoreValues.fNumberWeightedEnergy =
83 fSCScoreValues.fSumSLWE_v / fSCScoreValues.fSumSLW_v;
84
85 fSCScoreValues.fFluxWeightedEnergy =
86 fSCScoreValues.fSumSLWE / fSCScoreValues.fSumSLW;
87
88 fSCScoreValues.fAverageTrackWeight =
89 fSCScoreValues.fSumSLW / fSCScoreValues.fSumSL;
90 }
91 return fSCScoreValues;
92}
93
95{
96 fSCScoreValues.fImportance = importance;
97}
98
99std::ostream& operator<<(std::ostream& out, const G4CellScoreComposer& ps)
100{
102 out << "Tracks entering: " << scores.fSumTracksEntering << G4endl;
103 out << "Population: " << scores.fSumPopulation << G4endl;
104 out << "Collisions: " << scores.fSumCollisions << G4endl;
105 out << "Collisions*Wgt: " << scores.fSumCollisionsWeight << G4endl;
106 out << "NumWGTedEnergy: " << scores.fNumberWeightedEnergy << G4endl;
107 out << "FluxWGTedEnergy: " << scores.fFluxWeightedEnergy << G4endl;
108 out << "Aver.TrackWGT*I: " << scores.fAverageTrackWeight * scores.fImportance
109 << G4endl;
110 return out;
111}
std::ostream & operator<<(std::ostream &out, const G4CellScoreComposer &ps)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
void SetImportnace(G4double importance)
void SetCollisionWeight(G4double weight)
const G4CellScoreValues & GetStandardCellScoreValues() const
void EstimatorCalculation(const G4Step &step)
G4double fNumberWeightedEnergy
G4double GetVelocity() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const