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

#include <G4PSTermination.hh>

+ Inheritance diagram for G4PSTermination:

Public Member Functions

 G4PSTermination (G4String name, G4int depth=0)
 
virtual ~G4PSTermination ()
 
void Weighted (G4bool flg=true)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
virtual void SetUnit (const G4String &unit)
 
- Public Member Functions inherited from G4VPrimitiveScorer
 G4VPrimitiveScorer (G4String name, G4int depth=0)
 
virtual ~G4VPrimitiveScorer ()
 
G4int GetCollectionID (G4int)
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
void SetUnit (const G4String &unit)
 
const G4StringGetUnit () const
 
G4double GetUnitValue () const
 
void SetMultiFunctionalDetector (G4MultiFunctionalDetector *d)
 
G4MultiFunctionalDetectorGetMultiFunctionalDetector () const
 
G4String GetName () const
 
void SetFilter (G4VSDFilter *f)
 
G4VSDFilterGetFilter () const
 
void SetVerboseLevel (G4int vl)
 
G4int GetVerboseLevel () const
 
void SetNijk (G4int i, G4int j, G4int k)
 

Protected Member Functions

virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
 
- Protected Member Functions inherited from G4VPrimitiveScorer
virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)=0
 
virtual G4int GetIndex (G4Step *)
 
void CheckAndSetUnit (const G4String &unit, const G4String &category)
 
G4VSolidComputeSolid (G4Step *aStep, G4int replicaIdx)
 
G4VSolidComputeCurrentSolid (G4Step *aStep)
 

Additional Inherited Members

- Protected Attributes inherited from G4VPrimitiveScorer
G4String primitiveName
 
G4MultiFunctionalDetectordetector
 
G4VSDFilterfilter
 
G4int verboseLevel
 
G4int indexDepth
 
G4String unitName
 
G4double unitValue
 
G4int fNi
 
G4int fNj
 
G4int fNk
 

Detailed Description

Definition at line 46 of file G4PSTermination.hh.

Constructor & Destructor Documentation

◆ G4PSTermination()

G4PSTermination::G4PSTermination ( G4String  name,
G4int  depth = 0 
)

Definition at line 42 of file G4PSTermination.cc.

43 :G4VPrimitiveScorer(name,depth),HCID(-1),EvtMap(0),weighted(false)
44{
45 SetUnit("");
46}
virtual void SetUnit(const G4String &unit)

◆ ~G4PSTermination()

G4PSTermination::~G4PSTermination ( )
virtual

Definition at line 48 of file G4PSTermination.cc.

49{;}

Member Function Documentation

◆ clear()

void G4PSTermination::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 72 of file G4PSTermination.cc.

72 {
73 EvtMap->clear();
74}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSTermination::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 76 of file G4PSTermination.cc.

77{;}

◆ EndOfEvent()

void G4PSTermination::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 69 of file G4PSTermination.cc.

70{;}

◆ Initialize()

void G4PSTermination::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 62 of file G4PSTermination.cc.

63{
65 if(HCID < 0) {HCID = GetCollectionID(0);}
66 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
67}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ PrintAll()

void G4PSTermination::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 79 of file G4PSTermination.cc.

80{
81 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
82 G4cout << " PrimitiveScorer " << GetName() << G4endl;
83 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
84 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
85 for(; itr != EvtMap->GetMap()->end(); itr++) {
86 G4cout << " copy no.: " << itr->first
87 << " terminations: " << *(itr->second)
88 << G4endl;
89 }
90}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Map_t * GetMap() const
Definition: G4THitsMap.hh:151
size_t entries() const
Definition: G4THitsMap.hh:155

◆ ProcessHits()

G4bool G4PSTermination::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 51 of file G4PSTermination.cc.

52{
53 if ( aStep->GetTrack()->GetTrackStatus() != fStopAndKill ) return FALSE;
54
55 G4int index = GetIndex(aStep);
56 G4double val = 1.0;
57 if(weighted) val *= aStep->GetPreStepPoint()->GetWeight();
58 EvtMap->add(index,val);
59 return TRUE;
60}
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4TrackStatus GetTrackStatus() const
virtual G4int GetIndex(G4Step *)
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

◆ SetUnit()

void G4PSTermination::SetUnit ( const G4String unit)
virtual

Definition at line 92 of file G4PSTermination.cc.

93{
94 if (unit == "" ){
95 unitName = unit;
96 unitValue = 1.0;
97 }else{
98 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
99 G4Exception("G4PSTermination::SetUnit","DetPS0017",JustWarning,msg);
100 }
101
102}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
const G4String & GetUnit() const

Referenced by G4PSTermination().

◆ Weighted()

void G4PSTermination::Weighted ( G4bool  flg = true)
inline

Definition at line 58 of file G4PSTermination.hh.

58{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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