Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CellScorer.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$
28//
29// ----------------------------------------------------------------------
30// GEANT 4 class source file
31//
32// G4CellScorer.cc
33//
34// ----------------------------------------------------------------------
35
36#include "G4CellScorer.hh"
37#include "G4GeometryCell.hh"
38#include "G4Track.hh"
39#include "G4EventManager.hh"
40#include "G4Event.hh"
41
43{
44 G4cout << "--------------------------------------------------------" << G4endl
45 << "WARNING: Class <G4CellScorer> is now obsolete |" << G4endl
46 << " and will be removed starting from next Geant4 |" << G4endl
47 << " major release. Please, consider switching to |" << G4endl
48 << " general purpose scoring functionality. |" << G4endl
49 << "--------------------------------------------------------"
50 << G4endl;
51}
52
54{}
55
57 const G4GeometryCell &pre_gCell){
58 fCellScoreComposer.EstimatorCalculation(aStep);
59 ScorePopulation(pre_gCell, aStep);
60}
61
63 const G4GeometryCell &post_gCell){
64 // population counting
65 ScorePopulation(post_gCell, aStep);
66 // entering tracks
67 fCellScoreComposer.TrackEnters();
68}
69
71 const G4GeometryCell &post_gCell){
72 // population counting
73 ScorePopulation(post_gCell, aStep);
74 // collisions counting
76 fCellScoreComposer.SetCollisionWeight(aStep.GetTrack()->GetWeight());
77 }
78
79 // counting for step length estimators
80 fCellScoreComposer.EstimatorCalculation(aStep);
81
82}
83
84void G4CellScorer::ScorePopulation(const G4GeometryCell &,
85 const G4Step &aStep) {
86 // check for new event
88 GetConstCurrentEvent()->
89 GetEventID());
90 // increase population
91 if (fTrackLogger.FirstEnterance(aStep.GetTrack()->GetTrackID())) {
92 fCellScoreComposer.NewTrackPopedUp();
93 }
94}
95
97 return fCellScoreComposer;
98}
99
101 return fCellScoreComposer.GetStandardCellScoreValues();
102}
103
@ fGeomBoundary
Definition: G4StepStatus.hh:54
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetCollisionWeight(G4double weight)
const G4CellScoreValues & GetStandardCellScoreValues() const
void EstimatorCalculation(const G4Step &step)
virtual void ScoreAnEnteringStep(const G4Step &aStep, const G4GeometryCell &gCell)
Definition: G4CellScorer.cc:62
virtual void ScoreAnInVolumeStep(const G4Step &aStep, const G4GeometryCell &gCell)
Definition: G4CellScorer.cc:70
virtual void ScoreAnExitingStep(const G4Step &aStep, const G4GeometryCell &gCell)
Definition: G4CellScorer.cc:56
const G4CellScoreValues & GetCellScoreValues() const
virtual ~G4CellScorer()
Definition: G4CellScorer.cc:53
const G4CellScoreComposer & GetCellScoreComposer() const
Definition: G4CellScorer.cc:96
static G4EventManager * GetEventManager()
G4StepStatus GetStepStatus() const
Definition: G4Step.hh:78
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
void SetEventID(G4int id)
G4bool FirstEnterance(G4int trid)
G4int GetTrackID() const
G4double GetWeight() const