Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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{}
44
46
47 G4StepPoint *p = aStep.GetPreStepPoint();
48 if (!p) {
49 G4Exception("G4CellScoreComposer::EstimatorCalculation","Det0191",FatalException," no pointer to pre PreStepPoint!");
50 }
51 G4double sl = aStep.GetStepLength();
52 G4double slw = sl * p->GetWeight();
53 G4double slwe = slw * p->GetKineticEnergy();
54
55 G4double v = p->GetVelocity();
56 if (!(v>0.)) {
57 v = 10e-9;
58 }
59
60 fSCScoreValues.fSumSL += sl;
61 fSCScoreValues.fSumSLW += slw;
62 fSCScoreValues.fSumSLW_v += slw / v;
63 fSCScoreValues.fSumSLWE += slwe;
64 fSCScoreValues.fSumSLWE_v += slwe / v;
65
66}
68 fSCScoreValues.fSumTracksEntering++;
69}
71 fSCScoreValues.fSumPopulation++;
72}
73
75 fSCScoreValues.fSumCollisions++;
76 fSCScoreValues.fSumCollisionsWeight+=weight;
77}
78
79
82 if (fSCScoreValues.fSumSLW > 0.) {
83 //divide by SumSLW or SumSLW_v ?
84 fSCScoreValues.fNumberWeightedEnergy =
85 fSCScoreValues.fSumSLWE_v / fSCScoreValues.fSumSLW_v;
86
87 fSCScoreValues.fFluxWeightedEnergy =
88 fSCScoreValues.fSumSLWE / fSCScoreValues.fSumSLW;
89
90 fSCScoreValues.fAverageTrackWeight =
91 fSCScoreValues.fSumSLW / fSCScoreValues.fSumSL;
92 }
93 return fSCScoreValues;
94}
95
97 fSCScoreValues.fImportance = importance;
98}
99
100std::ostream& operator<<(std::ostream &out,
101 const G4CellScoreComposer &ps) {
103 out << "Tracks entering: " << scores.fSumTracksEntering << G4endl;
104 out << "Population: " << scores.fSumPopulation << G4endl;
105 out << "Collisions: " << scores.fSumCollisions << G4endl;
106 out << "Collisions*Wgt: " << scores.fSumCollisionsWeight << G4endl;
107 out << "NumWGTedEnergy: " << scores.fNumberWeightedEnergy << G4endl;
108 out << "FluxWGTedEnergy: " << scores.fFluxWeightedEnergy << G4endl;
109 out << "Aver.TrackWGT*I: " << scores.fAverageTrackWeight*
110 scores.fImportance << G4endl;
111 return out;
112}
113
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:35
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