G4PSPassageTrackLength Class Reference

#include <G4PSPassageTrackLength.hh>

Inheritance diagram for G4PSPassageTrackLength:

G4VPrimitiveScorer G4PSPassageTrackLength3D

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)

Protected Member Functions

virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
G4bool IsPassed (G4Step *)

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.

References SetUnit().

00047   :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
00048    weighted(false)
00049 {
00050     SetUnit("mm");
00051 }

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

Definition at line 53 of file G4PSPassageTrackLength.cc.

References SetUnit().

00056   :G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fTrackLength(0.),
00057    weighted(false)
00058 {
00059     SetUnit(unit);
00060 }

G4PSPassageTrackLength::~G4PSPassageTrackLength (  )  [virtual]

Definition at line 63 of file G4PSPassageTrackLength.cc.

00064 {;}


Member Function Documentation

void G4PSPassageTrackLength::clear (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 119 of file G4PSPassageTrackLength.cc.

00119                                   {
00120   EvtMap->clear();
00121 }

void G4PSPassageTrackLength::DrawAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 123 of file G4PSPassageTrackLength.cc.

00124 {;}

void G4PSPassageTrackLength::EndOfEvent ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 116 of file G4PSPassageTrackLength.cc.

00117 {;}

void G4PSPassageTrackLength::Initialize ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 107 of file G4PSPassageTrackLength.cc.

References G4HCofThisEvent::AddHitsCollection(), G4VPrimitiveScorer::detector, G4VPrimitiveScorer::GetCollectionID(), G4VPrimitiveScorer::GetName(), and G4VSensitiveDetector::GetName().

00108 {
00109   fCurrentTrkID = -1;
00110 
00111   EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
00112   if ( HCID < 0 ) HCID = GetCollectionID(0);
00113   HCE->AddHitsCollection(HCID,EvtMap);
00114 }

G4bool G4PSPassageTrackLength::IsPassed ( G4Step  )  [protected]

Definition at line 77 of file G4PSPassageTrackLength.cc.

References FALSE, fGeomBoundary, G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4Step::GetStepLength(), G4StepPoint::GetStepStatus(), G4Step::GetTrack(), G4Track::GetTrackID(), G4StepPoint::GetWeight(), and TRUE.

Referenced by ProcessHits().

00077                                                     {
00078   G4bool Passed = FALSE;
00079 
00080   G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
00081   G4bool IsExit  = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
00082 
00083   G4int  trkid  = aStep->GetTrack()->GetTrackID();
00084   G4double trklength  = aStep->GetStepLength();
00085   if(weighted) trklength *= aStep->GetPreStepPoint()->GetWeight();
00086 
00087   if ( IsEnter &&IsExit ){         // Passed at one step
00088     fTrackLength = trklength;      // Track length is absolutely given.
00089     Passed = TRUE;                 
00090   }else if ( IsEnter ){            // Enter a new geometry
00091     fCurrentTrkID = trkid;         // Resetting the current track.
00092     fTrackLength  = trklength;     
00093   }else if ( IsExit ){             // Exit a current geometry
00094     if ( fCurrentTrkID == trkid ) {
00095       fTrackLength  += trklength;  // Adding the track length to current one,
00096       Passed = TRUE;               // if the track is same as entered.
00097     }
00098   }else{                           // Inside geometry
00099     if ( fCurrentTrkID == trkid ){ // Adding the track length to current one ,
00100       fTrackLength  += trklength;  // if the track is same as entered.
00101     }
00102   }
00103 
00104   return Passed;
00105 }

void G4PSPassageTrackLength::PrintAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 126 of file G4PSPassageTrackLength.cc.

References G4VPrimitiveScorer::detector, G4cout, G4endl, G4VPrimitiveScorer::GetName(), G4VSensitiveDetector::GetName(), G4VPrimitiveScorer::GetUnit(), and G4VPrimitiveScorer::GetUnitValue().

00127 {
00128   G4cout << " MultiFunctionalDet  " << detector->GetName() << G4endl;
00129   G4cout << " PrimitiveSenstivity " << GetName() << G4endl;
00130   G4cout << " Number of entries " << EvtMap->entries() << G4endl;
00131   std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
00132   for(; itr != EvtMap->GetMap()->end(); itr++) {
00133     G4cout << "  copy no.: " << itr->first
00134            << "  track length : " 
00135            << *(itr->second)/GetUnitValue()
00136            << " [" << GetUnit()<< "]"
00137            << G4endl;
00138   }
00139 }

G4bool G4PSPassageTrackLength::ProcessHits ( G4Step ,
G4TouchableHistory  
) [protected, virtual]

Implements G4VPrimitiveScorer.

Definition at line 66 of file G4PSPassageTrackLength.cc.

References G4VPrimitiveScorer::GetIndex(), IsPassed(), and TRUE.

00067 {
00068 
00069   if ( IsPassed(aStep) ) {
00070     G4int index = GetIndex(aStep);
00071     EvtMap->add(index,fTrackLength);
00072   }
00073 
00074   return TRUE;
00075 }

void G4PSPassageTrackLength::SetUnit ( const G4String unit  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 141 of file G4PSPassageTrackLength.cc.

References G4VPrimitiveScorer::CheckAndSetUnit().

Referenced by G4PSPassageTrackLength(), and G4PSPassageTrackLength3D::G4PSPassageTrackLength3D().

00142 {
00143     CheckAndSetUnit(unit,"Length");
00144 }

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

Definition at line 56 of file G4PSPassageTrackLength.hh.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

00056 { weighted = flg; }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:03 2013 for Geant4 by  doxygen 1.4.7