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

#include <G4PSPassageCellCurrent.hh>

+ Inheritance diagram for G4PSPassageCellCurrent:

Public Member Functions

 G4PSPassageCellCurrent (G4String name, G4int depth=0)
 
virtual ~G4PSPassageCellCurrent ()
 
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 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 *)
 
virtual G4bool IsPassed (G4Step *)
 
- 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 50 of file G4PSPassageCellCurrent.hh.

Constructor & Destructor Documentation

◆ G4PSPassageCellCurrent()

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

Definition at line 50 of file G4PSPassageCellCurrent.cc.

51 :G4VPrimitivePlotter(name,depth),HCID(-1),fCurrentTrkID(-1),fCurrent(0),
52 EvtMap(0),weighted(true)
53{
54 SetUnit("");
55}
virtual void SetUnit(const G4String &unit)

◆ ~G4PSPassageCellCurrent()

G4PSPassageCellCurrent::~G4PSPassageCellCurrent ( )
virtual

Definition at line 57 of file G4PSPassageCellCurrent.cc.

58{;}

Member Function Documentation

◆ clear()

void G4PSPassageCellCurrent::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 126 of file G4PSPassageCellCurrent.cc.

126 {
127 EvtMap->clear();
128}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSPassageCellCurrent::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 130 of file G4PSPassageCellCurrent.cc.

131{;}

◆ EndOfEvent()

void G4PSPassageCellCurrent::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 122 of file G4PSPassageCellCurrent.cc.

123{
124}

◆ Initialize()

void G4PSPassageCellCurrent::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 111 of file G4PSPassageCellCurrent.cc.

112{
113 fCurrentTrkID = -1;
114
115 EvtMap = new G4THitsMap<G4double>(detector->GetName(),
116 GetName());
117 if ( HCID < 0 ) HCID = GetCollectionID(0);
118 HCE->AddHitsCollection(HCID,EvtMap);
119
120}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ IsPassed()

G4bool G4PSPassageCellCurrent::IsPassed ( G4Step aStep)
protectedvirtual

Definition at line 88 of file G4PSPassageCellCurrent.cc.

88 {
89 G4bool Passed = FALSE;
90
91 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
92 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
93
94 G4int trkid = aStep->GetTrack()->GetTrackID();
95
96 if ( IsEnter &&IsExit ){ // Passed at one step
97 Passed = TRUE;
98 }else if ( IsEnter ){ // Enter a new geometry
99 fCurrentTrkID = trkid; // Resetting the current track.
100 }else if ( IsExit ){ // Exit a current geometry
101 if ( fCurrentTrkID == trkid ) {
102 Passed = TRUE; // if the track is same as entered.
103 }
104 }else{ // Inside geometry
105 if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
106 }
107 }
108 return Passed;
109}
@ fGeomBoundary
Definition: G4StepStatus.hh:43
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
G4StepStatus GetStepStatus() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const

Referenced by ProcessHits().

◆ PrintAll()

void G4PSPassageCellCurrent::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 133 of file G4PSPassageCellCurrent.cc.

134{
135 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
136 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
137 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
138 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
139 for(; itr != EvtMap->GetMap()->end(); itr++) {
140 G4cout << " copy no.: " << itr->first
141 << " cell current : " << *(itr->second)
142 << " [tracks] "
143 << G4endl;
144 }
145}
#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 G4PSPassageCellCurrent::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 60 of file G4PSPassageCellCurrent.cc.

61{
62
63 if ( IsPassed(aStep) ) {
64 fCurrent = 1.;
65 if(weighted) fCurrent = aStep->GetPreStepPoint()->GetWeight();
66 G4int index = GetIndex(aStep);
67 EvtMap->add(index,fCurrent);
68
69 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
70 {
71 auto filler = G4VScoreHistFiller::Instance();
72 if(!filler)
73 {
74 G4Exception("G4PSVolumeFlux::ProcessHits","SCORER0123",JustWarning,
75 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
76 }
77 else
78 {
79 filler->FillH1(hitIDMap[index],aStep->GetPreStepPoint()->GetKineticEnergy(),fCurrent);
80 }
81 }
82
83 }
84
85 return TRUE;
86}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
virtual G4bool IsPassed(G4Step *)
G4double GetKineticEnergy() const
G4double GetWeight() 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

◆ SetUnit()

void G4PSPassageCellCurrent::SetUnit ( const G4String unit)
virtual

Definition at line 147 of file G4PSPassageCellCurrent.cc.

148{
149 if (unit == "" ){
150 unitName = unit;
151 unitValue = 1.0;
152 }else{
153 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
154 G4Exception("G4PSPassageCellCurrent::SetUnit","DetPS0012",JustWarning,msg);
155 }
156}
const G4String & GetUnit() const

Referenced by G4PSPassageCellCurrent().

◆ Weighted()

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

Definition at line 57 of file G4PSPassageCellCurrent.hh.

57{ weighted = flg; }

Referenced by G4ScoreQuantityMessenger::SetNewValue().


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