G4PSMinKinEAtGeneration Class Reference

#include <G4PSMinKinEAtGeneration.hh>

Inheritance diagram for G4PSMinKinEAtGeneration:

G4VPrimitiveScorer G4PSMinKinEAtGeneration3D

Public Member Functions

 G4PSMinKinEAtGeneration (G4String name, G4int depth=0)
 G4PSMinKinEAtGeneration (G4String name, const G4String &unit, G4int depth=0)
virtual ~G4PSMinKinEAtGeneration ()
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 *)

Detailed Description

Definition at line 48 of file G4PSMinKinEAtGeneration.hh.


Constructor & Destructor Documentation

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

Definition at line 42 of file G4PSMinKinEAtGeneration.cc.

References SetUnit().

00043   :G4VPrimitiveScorer(name,depth),HCID(-1)
00044 {
00045     SetUnit("MeV");
00046 }

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

Definition at line 48 of file G4PSMinKinEAtGeneration.cc.

References SetUnit().

00051   :G4VPrimitiveScorer(name,depth),HCID(-1)
00052 {
00053     SetUnit(unit);
00054 }

G4PSMinKinEAtGeneration::~G4PSMinKinEAtGeneration (  )  [virtual]

Definition at line 56 of file G4PSMinKinEAtGeneration.cc.

00057 {;}


Member Function Documentation

void G4PSMinKinEAtGeneration::clear (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 109 of file G4PSMinKinEAtGeneration.cc.

00109                                    {
00110   EvtMap->clear();
00111 }

void G4PSMinKinEAtGeneration::DrawAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 113 of file G4PSMinKinEAtGeneration.cc.

00114 {;}

void G4PSMinKinEAtGeneration::EndOfEvent ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 106 of file G4PSMinKinEAtGeneration.cc.

00107 {;}

void G4PSMinKinEAtGeneration::Initialize ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 99 of file G4PSMinKinEAtGeneration.cc.

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

00100 {
00101   EvtMap = new G4THitsMap<G4double>(detector->GetName(),GetName());
00102   if(HCID < 0) {HCID = GetCollectionID(0);}
00103   HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
00104 }

void G4PSMinKinEAtGeneration::PrintAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 116 of file G4PSMinKinEAtGeneration.cc.

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

00117 {
00118   G4cout << " PrimitiveScorer " << GetName() << G4endl;
00119   G4cout << " Number of entries " << EvtMap->entries() << G4endl;
00120   std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
00121   for(; itr != EvtMap->GetMap()->end(); itr++) {
00122     G4cout << "  copy no.: " << itr->first
00123            << "  energy: " << *(itr->second)/GetUnitValue()
00124            << " ["<<GetUnit()<<"]"
00125            << G4endl;
00126   }
00127 }

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

Implements G4VPrimitiveScorer.

Definition at line 59 of file G4PSMinKinEAtGeneration.cc.

References FALSE, G4Track::GetCurrentStepNumber(), G4VPrimitiveScorer::GetIndex(), G4StepPoint::GetKineticEnergy(), G4Track::GetParentID(), G4Step::GetPreStepPoint(), G4Step::GetTrack(), and TRUE.

00060 {
00061   // ===============
00062   //  First Step, 
00063   //  Confirm this is a newly produced secondary.
00064   //  
00065   //- check for newly produced particle. e.g. Step number is 1.
00066   if ( aStep->GetTrack()->GetCurrentStepNumber() != 1) return FALSE;
00067   //- check for this is not a primary particle. e.g. ParentID != 0 .
00068   if ( aStep->GetTrack()->GetParentID() == 0 ) return FALSE;
00069 
00070   //===============================================
00071   //- This is a newly produced secondary particle.
00072   //===============================================
00073 
00074   // ===============
00075   //  Second Step,
00076   //  Confirm this track has lower energy than previous one.
00077   //
00078 
00079   // -Kinetic energy of this particle at the starting point.
00080   G4double kinetic = aStep->GetPreStepPoint()->GetKineticEnergy();
00081 
00082   // -Stored value in the current HitsMap.
00083   G4int  index = GetIndex(aStep);
00084   G4double*   mapValue=  ((*EvtMap)[index]);
00085 
00086   // -
00087   // If mapValue exits (e.g not NULL ), compare it with 
00088   // current track's kinetic energy.
00089   if ( mapValue && ( kinetic > *mapValue ) ) return FALSE;
00090 
00091   // -
00092   // Current Track is a newly produced secondary and has lower
00093   // kinetic energy than previous one in this geometry.
00094   //
00095   EvtMap->set(index, kinetic);
00096   return TRUE;
00097 }

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

Reimplemented from G4VPrimitiveScorer.

Definition at line 129 of file G4PSMinKinEAtGeneration.cc.

References G4VPrimitiveScorer::CheckAndSetUnit().

Referenced by G4PSMinKinEAtGeneration(), G4PSMinKinEAtGeneration3D::G4PSMinKinEAtGeneration3D(), and G4ScoreQuantityMessenger::SetNewValue().

00130 {
00131     CheckAndSetUnit(unit,"Energy");
00132 }


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