Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id: G4CellScoreComposer.cc 67992 2013-03-13 10:59:57Z gcosmo $
28 //
29 // ----------------------------------------------------------------------
30 // GEANT 4 class source file
31 //
32 // G4CellSCoreComposer.cc
33 //
34 // ----------------------------------------------------------------------
35 
36 #include "G4CellScoreComposer.hh"
37 #include "G4Step.hh"
38 
40  fSCScoreValues()
41 {}
42 
44 {}
45 
47 
48  G4StepPoint *p = 0;
49  p = aStep.GetPreStepPoint();
50  if (!p) {
51  G4Exception("G4CellScoreComposer::EstimatorCalculation","Det0191",FatalException," no pointer to pre PreStepPoint!");
52  }
53  G4double sl = aStep.GetStepLength();
54  G4double slw = sl * p->GetWeight();
55  G4double slwe = slw * p->GetKineticEnergy();
56 
57  G4double v = p->GetVelocity();
58  if (!(v>0.)) {
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 
68 }
70  fSCScoreValues.fSumTracksEntering++;
71 }
73  fSCScoreValues.fSumPopulation++;
74 }
75 
77  fSCScoreValues.fSumCollisions++;
78  fSCScoreValues.fSumCollisionsWeight+=weight;
79 }
80 
81 
84  if (fSCScoreValues.fSumSLW > 0.) {
85  //divide by SumSLW or SumSLW_v ?
86  fSCScoreValues.fNumberWeightedEnergy =
87  fSCScoreValues.fSumSLWE_v / fSCScoreValues.fSumSLW_v;
88 
89  fSCScoreValues.fFluxWeightedEnergy =
90  fSCScoreValues.fSumSLWE / fSCScoreValues.fSumSLW;
91 
92  fSCScoreValues.fAverageTrackWeight =
93  fSCScoreValues.fSumSLW / fSCScoreValues.fSumSL;
94  }
95  return fSCScoreValues;
96 }
97 
99  fSCScoreValues.fImportance = importance;
100 }
101 
102 std::ostream& operator<<(std::ostream &out,
103  const G4CellScoreComposer &ps) {
105  out << "Tracks entering: " << scores.fSumTracksEntering << G4endl;
106  out << "Population: " << scores.fSumPopulation << G4endl;
107  out << "Collisions: " << scores.fSumCollisions << G4endl;
108  out << "Collisions*Wgt: " << scores.fSumCollisionsWeight << G4endl;
109  out << "NumWGTedEnergy: " << scores.fNumberWeightedEnergy << G4endl;
110  out << "FluxWGTedEnergy: " << scores.fFluxWeightedEnergy << G4endl;
111  out << "Aver.TrackWGT*I: " << scores.fAverageTrackWeight*
112  scores.fImportance << G4endl;
113  return out;
114 }
115 
G4double GetWeight() const
G4double GetStepLength() const
const char * p
Definition: xmltok.h:285
G4double fNumberWeightedEnergy
G4double GetVelocity() const
G4StepPoint * GetPreStepPoint() const
void EstimatorCalculation(const G4Step &step)
Definition: G4Step.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
#define G4endl
Definition: G4ios.hh:61
void SetImportnace(G4double importance)
void SetCollisionWeight(G4double weight)
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
const G4CellScoreValues & GetStandardCellScoreValues() const