Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CellScorer Class Reference

#include <G4CellScorer.hh>

+ Inheritance diagram for G4CellScorer:

Public Member Functions

 G4CellScorer ()
 
virtual ~G4CellScorer ()
 
virtual void ScoreAnExitingStep (const G4Step &aStep, const G4GeometryCell &gCell)
 
virtual void ScoreAnEnteringStep (const G4Step &aStep, const G4GeometryCell &gCell)
 
virtual void ScoreAnInVolumeStep (const G4Step &aStep, const G4GeometryCell &gCell)
 
const G4CellScoreComposerGetCellScoreComposer () const
 
const G4CellScoreValuesGetCellScoreValues () const
 
- Public Member Functions inherited from G4VCellScorer
 G4VCellScorer ()
 
virtual ~G4VCellScorer ()
 
virtual void ScoreAnExitingStep (const G4Step &aStep, const G4GeometryCell &gCell)=0
 
virtual void ScoreAnEnteringStep (const G4Step &aStep, const G4GeometryCell &gCell)=0
 
virtual void ScoreAnInVolumeStep (const G4Step &aStep, const G4GeometryCell &gCell)=0
 

Detailed Description

Definition at line 51 of file G4CellScorer.hh.

Constructor & Destructor Documentation

◆ G4CellScorer()

G4CellScorer::G4CellScorer ( )

Definition at line 42 of file G4CellScorer.cc.

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}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ ~G4CellScorer()

G4CellScorer::~G4CellScorer ( )
virtual

Definition at line 53 of file G4CellScorer.cc.

54{}

Member Function Documentation

◆ GetCellScoreComposer()

const G4CellScoreComposer & G4CellScorer::GetCellScoreComposer ( ) const

Definition at line 96 of file G4CellScorer.cc.

96 {
97 return fCellScoreComposer;
98}

◆ GetCellScoreValues()

const G4CellScoreValues & G4CellScorer::GetCellScoreValues ( ) const

Definition at line 100 of file G4CellScorer.cc.

100 {
101 return fCellScoreComposer.GetStandardCellScoreValues();
102}
const G4CellScoreValues & GetStandardCellScoreValues() const

◆ ScoreAnEnteringStep()

void G4CellScorer::ScoreAnEnteringStep ( const G4Step aStep,
const G4GeometryCell gCell 
)
virtual

Implements G4VCellScorer.

Definition at line 62 of file G4CellScorer.cc.

63 {
64 // population counting
65 ScorePopulation(post_gCell, aStep);
66 // entering tracks
67 fCellScoreComposer.TrackEnters();
68}

◆ ScoreAnExitingStep()

void G4CellScorer::ScoreAnExitingStep ( const G4Step aStep,
const G4GeometryCell gCell 
)
virtual

Implements G4VCellScorer.

Definition at line 56 of file G4CellScorer.cc.

57 {
58 fCellScoreComposer.EstimatorCalculation(aStep);
59 ScorePopulation(pre_gCell, aStep);
60}
void EstimatorCalculation(const G4Step &step)

◆ ScoreAnInVolumeStep()

void G4CellScorer::ScoreAnInVolumeStep ( const G4Step aStep,
const G4GeometryCell gCell 
)
virtual

Implements G4VCellScorer.

Definition at line 70 of file G4CellScorer.cc.

71 {
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}
@ fGeomBoundary
Definition: G4StepStatus.hh:54
void SetCollisionWeight(G4double weight)
G4StepStatus GetStepStatus() const
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
G4double GetWeight() const

The documentation for this class was generated from the following files: