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

#include <G4PSNofSecondary.hh>

+ Inheritance diagram for G4PSNofSecondary:

Public Member Functions

 G4PSNofSecondary (G4String name, G4int depth=0)
 
void SetParticle (const G4String &particleName)
 
void Weighted (G4bool flg=true)
 
virtual ~G4PSNofSecondary ()
 
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 *)
 
- 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 54 of file G4PSNofSecondary.hh.

Constructor & Destructor Documentation

◆ G4PSNofSecondary()

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

Definition at line 43 of file G4PSNofSecondary.cc.

44 :G4VPrimitivePlotter(name,depth),HCID(-1),EvtMap(0),particleDef(0),
45 weighted(true)
46{;}

◆ ~G4PSNofSecondary()

G4PSNofSecondary::~G4PSNofSecondary ( )
virtual

Definition at line 48 of file G4PSNofSecondary.cc.

49{;}

Member Function Documentation

◆ clear()

void G4PSNofSecondary::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 107 of file G4PSNofSecondary.cc.

107 {
108 EvtMap->clear();
109}
void clear()
Definition: G4THitsMap.hh:537

◆ DrawAll()

void G4PSNofSecondary::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 111 of file G4PSNofSecondary.cc.

112{;}

◆ EndOfEvent()

void G4PSNofSecondary::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 104 of file G4PSNofSecondary.cc.

105{;}

◆ Initialize()

void G4PSNofSecondary::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 97 of file G4PSNofSecondary.cc.

98{
100 if(HCID < 0) {HCID = GetCollectionID(0);}
101 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
102}
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4String GetName() const
G4MultiFunctionalDetector * detector
G4int GetCollectionID(G4int)

◆ PrintAll()

void G4PSNofSecondary::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 114 of file G4PSNofSecondary.cc.

115{
116 G4cout << " PrimitiveScorer " << GetName() << G4endl;
117 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
118 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
119 for(; itr != EvtMap->GetMap()->end(); itr++) {
120 G4cout << " copy no.: " << itr->first
121 << " num of secondaries: " << *(itr->second)/GetUnitValue()
122 << G4endl;
123 }
124}
#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 G4PSNofSecondary::ProcessHits ( G4Step aStep,
G4TouchableHistory  
)
protectedvirtual

Implements G4VPrimitiveScorer.

Definition at line 51 of file G4PSNofSecondary.cc.

52{
53 //- check for newly produced particle. e.g. first step.
54 if ( aStep->GetTrack()->GetCurrentStepNumber() != 1) return FALSE;
55 //- check for this is not a primary particle. e.g. ParentID > 0 .
56 if ( aStep->GetTrack()->GetParentID() == 0 ) return FALSE;
57 //- check the particle if the partifle definition is given.
58 if ( particleDef && particleDef != aStep->GetTrack()->GetDefinition() )
59 return FALSE;
60 //
61 //- This is a newly produced secondary particle.
62 G4int index = GetIndex(aStep);
63 G4double weight = 1.0;
64 if ( weighted ) weight *= aStep->GetPreStepPoint()->GetWeight();
65 EvtMap->add(index,weight);
66
67 if(hitIDMap.size()>0 && hitIDMap.find(index)!=hitIDMap.end())
68 {
69 auto filler = G4VScoreHistFiller::Instance();
70 if(!filler)
71 {
72 G4Exception("G4PSVolumeFlux::ProcessHits","SCORER0123",JustWarning,
73 "G4TScoreHistFiller is not instantiated!! Histogram is not filled.");
74 }
75 else
76 {
77 filler->FillH1(hitIDMap[index],aStep->GetPreStepPoint()->GetKineticEnergy(),weight);
78 }
79 }
80
81 return TRUE;
82}
@ 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 GetKineticEnergy() const
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4int GetCurrentStepNumber() const
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() 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

◆ SetParticle()

void G4PSNofSecondary::SetParticle ( const G4String particleName)

Definition at line 84 of file G4PSNofSecondary.cc.

84 {
87 if(!pd) {
88 G4String msg = "Particle <";
89 msg += particleName;
90 msg += "> not found.";
91 G4Exception("G4PSNofSecondary::SetParticle",
92 "DetPS0101",FatalException,msg);
93 }
94 particleDef = pd;
95}
@ FatalException
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

◆ SetUnit()

void G4PSNofSecondary::SetUnit ( const G4String unit)
virtual

Definition at line 126 of file G4PSNofSecondary.cc.

127{
128 if (unit == "" ){
129 unitName = unit;
130 unitValue = 1.0;
131 }else{
132 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
133 G4Exception("G4PSNofSecondary::SetUnit","DetPS0010",JustWarning,msg);
134 }
135}
const G4String & GetUnit() const

◆ Weighted()

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

Definition at line 63 of file G4PSNofSecondary.hh.

63{ weighted = flg; }

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