Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PSPassageTrackLength Class Reference

#include <G4PSPassageTrackLength.hh>

Inheritance diagram for G4PSPassageTrackLength:
Collaboration diagram for G4PSPassageTrackLength:

Public Member Functions

 G4PSPassageTrackLength (G4String name, G4int depth=0)
 
 G4PSPassageTrackLength (G4String name, const G4String &unit, G4int depth=0)
 
virtual ~G4PSPassageTrackLength ()
 
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 *)
 
G4bool IsPassed (G4Step *)
 
- Protected Member Functions inherited from G4VPrimitiveScorer
virtual G4int GetIndex (G4Step *)
 
void CheckAndSetUnit (const G4String &unit, const G4String &category)
 

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 G4PSPassageTrackLength.hh.

Constructor & Destructor Documentation

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

Definition at line 46 of file G4PSPassageTrackLength.cc.

47  :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
48  EvtMap(0),weighted(false)
49 {
50  SetUnit("mm");
51 }
virtual void SetUnit(const G4String &unit)
G4VPrimitiveScorer(G4String name, G4int depth=0)

Here is the call graph for this function:

G4PSPassageTrackLength::G4PSPassageTrackLength ( G4String  name,
const G4String unit,
G4int  depth = 0 
)

Definition at line 53 of file G4PSPassageTrackLength.cc.

56  :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
57  EvtMap(0),weighted(false)
58 {
59  SetUnit(unit);
60 }
virtual void SetUnit(const G4String &unit)
G4VPrimitiveScorer(G4String name, G4int depth=0)

Here is the call graph for this function:

G4PSPassageTrackLength::~G4PSPassageTrackLength ( )
virtual

Definition at line 63 of file G4PSPassageTrackLength.cc.

64 {;}

Member Function Documentation

void G4PSPassageTrackLength::clear ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 119 of file G4PSPassageTrackLength.cc.

119  {
120  EvtMap->clear();
121 }
void clear()
Definition: G4THitsMap.hh:267

Here is the call graph for this function:

void G4PSPassageTrackLength::DrawAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 123 of file G4PSPassageTrackLength.cc.

124 {;}
void G4PSPassageTrackLength::EndOfEvent ( G4HCofThisEvent )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 116 of file G4PSPassageTrackLength.cc.

117 {;}
void G4PSPassageTrackLength::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 107 of file G4PSPassageTrackLength.cc.

108 {
109  fCurrentTrkID = -1;
110 
111  EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
112  if ( HCID < 0 ) HCID = GetCollectionID(0);
113  HCE->AddHitsCollection(HCID,EvtMap);
114 }
G4String GetName() const
G4int GetCollectionID(G4int)
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4MultiFunctionalDetector * detector

Here is the call graph for this function:

G4bool G4PSPassageTrackLength::IsPassed ( G4Step aStep)
protected

Definition at line 77 of file G4PSPassageTrackLength.cc.

77  {
78  G4bool Passed = FALSE;
79 
80  G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
81  G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
82 
83  G4int trkid = aStep->GetTrack()->GetTrackID();
84  G4double trklength = aStep->GetStepLength();
85  if(weighted) trklength *= aStep->GetPreStepPoint()->GetWeight();
86 
87  if ( IsEnter &&IsExit ){ // Passed at one step
88  fTrackLength = trklength; // Track length is absolutely given.
89  Passed = TRUE;
90  }else if ( IsEnter ){ // Enter a new geometry
91  fCurrentTrkID = trkid; // Resetting the current track.
92  fTrackLength = trklength;
93  }else if ( IsExit ){ // Exit a current geometry
94  if ( fCurrentTrkID == trkid ) {
95  fTrackLength += trklength; // Adding the track length to current one,
96  Passed = TRUE; // if the track is same as entered.
97  }
98  }else{ // Inside geometry
99  if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
100  fTrackLength += trklength; // if the track is same as entered.
101  }
102  }
103 
104  return Passed;
105 }
G4double GetWeight() const
G4double GetStepLength() const
G4StepStatus GetStepStatus() const
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
bool G4bool
Definition: G4Types.hh:79
#define FALSE
Definition: globals.hh:52
#define TRUE
Definition: globals.hh:55
G4int GetTrackID() const
G4StepPoint * GetPostStepPoint() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PSPassageTrackLength::PrintAll ( )
virtual

Reimplemented from G4VPrimitiveScorer.

Definition at line 126 of file G4PSPassageTrackLength.cc.

127 {
128  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
129  G4cout << " PrimitiveSenstivity " << GetName() << G4endl;
130  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
131  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
132  for(; itr != EvtMap->GetMap()->end(); itr++) {
133  G4cout << " copy no.: " << itr->first
134  << " track length : "
135  << *(itr->second)/GetUnitValue()
136  << " [" << GetUnit()<< "]"
137  << G4endl;
138  }
139 }
G4String GetName() const
G4GLOB_DLL std::ostream G4cout
G4double GetUnitValue() const
G4int entries() const
Definition: G4THitsMap.hh:200
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:99
#define G4endl
Definition: G4ios.hh:61
G4MultiFunctionalDetector * detector
const G4String & GetUnit() const

Here is the call graph for this function:

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

Implements G4VPrimitiveScorer.

Definition at line 66 of file G4PSPassageTrackLength.cc.

67 {
68 
69  if ( IsPassed(aStep) ) {
70  G4int index = GetIndex(aStep);
71  EvtMap->add(index,fTrackLength);
72  }
73 
74  return TRUE;
75 }
int G4int
Definition: G4Types.hh:78
virtual G4int GetIndex(G4Step *)
#define TRUE
Definition: globals.hh:55
G4int add(const G4int &key, T *&aHit) const
Definition: G4THitsMap.hh:106

Here is the call graph for this function:

void G4PSPassageTrackLength::SetUnit ( const G4String unit)
virtual

Definition at line 141 of file G4PSPassageTrackLength.cc.

142 {
143  CheckAndSetUnit(unit,"Length");
144 }
void CheckAndSetUnit(const G4String &unit, const G4String &category)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 56 of file G4PSPassageTrackLength.hh.

56 { weighted = flg; }

Here is the caller graph for this function:


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