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

#include <G4PSNofStep.hh>

+ Inheritance diagram for G4PSNofStep:

Public Member Functions

 G4PSNofStep (G4String name, G4int depth=0)
 
virtual ~G4PSNofStep ()
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
virtual void SetUnit (const G4String &unit)
 
void SetBoundaryFlag (G4bool flg=true)
 
- Public Member Functions inherited from G4VPrimitivePlotter
 G4VPrimitivePlotter (G4String name, G4int depth=0)
 
virtual ~G4VPrimitivePlotter ()
 
void Plot (G4int copyNo, G4int histID)
 
G4int GetNumberOfHist () const
 
- 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 G4VPrimitivePlotter
std::map< G4int, G4inthitIDMap
 
- 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 G4PSNofStep.hh.

Constructor & Destructor Documentation

◆ G4PSNofStep()

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

Definition at line 43 of file G4PSNofStep.cc.

44 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0),boundaryFlag(false)
45{
46 SetUnit("");
47}
virtual void SetUnit(const G4String &unit)
Definition: G4PSNofStep.cc:109

◆ ~G4PSNofStep()

G4PSNofStep::~G4PSNofStep ( )
virtual

Definition at line 49 of file G4PSNofStep.cc.

50{;}

Member Function Documentation

◆ clear()

void G4PSNofStep::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 88 of file G4PSNofStep.cc.

88 {
89 EvtMap->clear();
90}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSNofStep::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 92 of file G4PSNofStep.cc.

93{;}

◆ EndOfEvent()

void G4PSNofStep::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 85 of file G4PSNofStep.cc.

86{;}

◆ Initialize()

void G4PSNofStep::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 78 of file G4PSNofStep.cc.

79{
81 if(HCID < 0) {HCID = GetCollectionID(0);}
82 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
83}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ PrintAll()

void G4PSNofStep::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 95 of file G4PSNofStep.cc.

96{
97 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
98 G4cout << " PrimitiveScorer " << GetName() << G4endl;
99 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
100 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
101 for(; itr != EvtMap->GetMap()->end(); itr++) {
102 G4cout << " copy no.: " << itr->first
103 << " num of step: " << *(itr->second)
104 << " [steps] "
105 << G4endl;
106 }
107}
#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 G4PSNofStep::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 52 of file G4PSNofStep.cc.

53{
54 if ( boundaryFlag ) {
55 if ( aStep->GetStepLength() == 0. ) return FALSE;
56 }
57 G4int index = GetIndex(aStep);
58 G4double val = 1.0;
59 EvtMap->add(index,val);
60
61 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
62 {
63 auto filler = G4VScoreHistFiller::Instance();
64 if(!filler)
65 {
66 G4Exception("G4PSNofStep::ProcessHits","SCORER0123",JustWarning,
67 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
68 }
69 else
70 {
71 filler->FillH1(hitIDMap[index],aStep->GetStepLength(),val);
72 }
73 }
74
75 return TRUE;
76}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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 GetStepLength() const
std::map< G4int, G4int > hitIDMap
virtual G4int GetIndex(G4Step *)
static G4VScoreHistFiller * Instance()
size_t add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:177

◆ SetBoundaryFlag()

void G4PSNofStep::SetBoundaryFlag ( G4bool  flg = true)
inline

Definition at line 71 of file G4PSNofStep.hh.

72 {boundaryFlag = flg;}

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetUnit()

void G4PSNofStep::SetUnit ( const G4String unit)
virtual

Definition at line 109 of file G4PSNofStep.cc.

110{
111 if (unit == "" ){
112 unitName = unit;
113 unitValue = 1.0;
114 }else{
115 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
116 G4Exception("G4PSNofStep::SetUnit","DetPS0011",JustWarning,msg);
117 }
118
119}
const G4String & GetUnit() const

Referenced by G4PSNofStep().


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