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

#include <G4PSNofCollision.hh>

+ Inheritance diagram for G4PSNofCollision:

Public Member Functions

 G4PSNofCollision (G4String name, G4int depth=0)
 
virtual ~G4PSNofCollision ()
 
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 G4PSNofCollision.hh.

Constructor & Destructor Documentation

◆ G4PSNofCollision()

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

Definition at line 41 of file G4PSNofCollision.cc.

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

◆ ~G4PSNofCollision()

G4PSNofCollision::~G4PSNofCollision ( )
virtual

Definition at line 47 of file G4PSNofCollision.cc.

48{;}

Member Function Documentation

◆ clear()

void G4PSNofCollision::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 71 of file G4PSNofCollision.cc.

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

◆ DrawAll()

void G4PSNofCollision::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 75 of file G4PSNofCollision.cc.

76{;}

◆ EndOfEvent()

void G4PSNofCollision::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 68 of file G4PSNofCollision.cc.

69{;}

◆ Initialize()

void G4PSNofCollision::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 61 of file G4PSNofCollision.cc.

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

◆ PrintAll()

void G4PSNofCollision::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 78 of file G4PSNofCollision.cc.

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

◆ ProcessHits()

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

Implements G4VPrimitiveScorer.

Definition at line 50 of file G4PSNofCollision.cc.

51{
52 if ( aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ) return TRUE;
53
54 G4int index = GetIndex(aStep);
55 G4double val = 1.0;
56 if(weighted) val *= aStep->GetPreStepPoint()->GetWeight();
57 EvtMap->add(index,val);
58 return TRUE;
59}
@ fGeomBoundary
Definition: G4StepStatus.hh:43
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define TRUE
Definition: Globals.hh:27
G4StepStatus GetStepStatus() const
G4double GetWeight() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
virtual G4int GetIndex(G4Step *)
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

◆ SetUnit()

void G4PSNofCollision::SetUnit ( const G4String unit)
virtual

Definition at line 93 of file G4PSNofCollision.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("G4PSNofCollision::SetUnit","DetPS0009",JustWarning,msg);
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 G4PSNofCollision().

◆ Weighted()

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

Definition at line 58 of file G4PSNofCollision.hh.

58{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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