Geant4  10.02.p03
G4PSTermination Class Reference

#include <G4PSTermination.hh>

Inheritance diagram for G4PSTermination:
Collaboration 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)
 
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 G4int GetIndex (G4Step *)
 
void CheckAndSetUnit (const G4String &unit, const G4String &category)
 

Private Attributes

G4int HCID
 
G4THitsMap< G4double > * EvtMap
 
G4bool weighted
 

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 47 of file G4PSTermination.hh.

Constructor & Destructor Documentation

◆ G4PSTermination()

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

Definition at line 43 of file G4PSTermination.cc.

44  :G4VPrimitiveScorer(name,depth),HCID(-1),EvtMap(0),weighted(false)
45 {
46  SetUnit("");
47 }
G4THitsMap< G4double > * EvtMap
virtual void SetUnit(const G4String &unit)
G4VPrimitiveScorer(G4String name, G4int depth=0)
Here is the call graph for this function:

◆ ~G4PSTermination()

G4PSTermination::~G4PSTermination ( )
virtual

Definition at line 49 of file G4PSTermination.cc.

50 {;}

Member Function Documentation

◆ clear()

void G4PSTermination::clear ( void  )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 73 of file G4PSTermination.cc.

73  {
74  EvtMap->clear();
75 }
void clear()
Definition: G4THitsMap.hh:209
G4THitsMap< G4double > * EvtMap
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DrawAll()

void G4PSTermination::DrawAll ( void  )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 77 of file G4PSTermination.cc.

78 {;}
Here is the caller graph for this function:

◆ EndOfEvent()

void G4PSTermination::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 70 of file G4PSTermination.cc.

71 {;}
Here is the caller graph for this function:

◆ Initialize()

void G4PSTermination::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 63 of file G4PSTermination.cc.

64 {
66  if(HCID < 0) {HCID = GetCollectionID(0);}
68 }
G4String GetName() const
G4int GetCollectionID(G4int)
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4THitsMap< G4double > * EvtMap
G4MultiFunctionalDetector * detector
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PrintAll()

void G4PSTermination::PrintAll ( void  )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 80 of file G4PSTermination.cc.

81 {
82  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
83  G4cout << " PrimitiveScorer " << GetName() << G4endl;
84  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
85  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
86  for(; itr != EvtMap->GetMap()->end(); itr++) {
87  G4cout << " copy no.: " << itr->first
88  << " terminations: " << *(itr->second)
89  << G4endl;
90  }
91 }
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
G4String GetName() const
G4GLOB_DLL std::ostream G4cout
G4int entries() const
Definition: G4THitsMap.hh:79
G4THitsMap< G4double > * EvtMap
#define G4endl
Definition: G4ios.hh:61
G4MultiFunctionalDetector * detector
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 52 of file G4PSTermination.cc.

53 {
54  if ( aStep->GetTrack()->GetTrackStatus() != fStopAndKill ) return FALSE;
55 
56  G4int index = GetIndex(aStep);
57  G4double val = 1.0;
58  if(weighted) val *= aStep->GetPreStepPoint()->GetWeight();
59  EvtMap->add(index,val);
60  return TRUE;
61 }
Int_t index
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:138
int G4int
Definition: G4Types.hh:78
#define FALSE
Definition: globals.hh:52
virtual G4int GetIndex(G4Step *)
#define TRUE
Definition: globals.hh:55
G4THitsMap< G4double > * EvtMap
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SetUnit()

void G4PSTermination::SetUnit ( const G4String unit)
virtual

Definition at line 93 of file G4PSTermination.cc.

94 {
95  if (unit == "" ){
96  unitName = unit;
97  unitValue = 1.0;
98  }else{
99  G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
100  G4Exception("G4PSTermination::SetUnit","DetPS0017",JustWarning,msg);
101  }
102 
103 }
G4String GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetUnit() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Weighted()

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

Definition at line 59 of file G4PSTermination.hh.

59 { weighted = flg; }
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ EvtMap

G4THitsMap<G4double>* G4PSTermination::EvtMap
private

Definition at line 75 of file G4PSTermination.hh.

◆ HCID

G4int G4PSTermination::HCID
private

Definition at line 74 of file G4PSTermination.hh.

◆ weighted

G4bool G4PSTermination::weighted
private

Definition at line 76 of file G4PSTermination.hh.


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