Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VGFlashSensitiveDetector Class Referenceabstract

#include <G4VGFlashSensitiveDetector.hh>

Inheritance diagram for G4VGFlashSensitiveDetector:

Public Member Functions

 G4VGFlashSensitiveDetector ()
 
 G4VGFlashSensitiveDetector (const G4VGFlashSensitiveDetector &)
 
virtual ~G4VGFlashSensitiveDetector ()
 
G4int operator== (const G4VGFlashSensitiveDetector &right) const
 
G4int operator!= (const G4VGFlashSensitiveDetector &right) const
 
G4bool Hit (G4GFlashSpot *aSpot)
 

Protected Member Functions

virtual G4bool ProcessHits (G4GFlashSpot *aSpot, G4TouchableHistory *ROhist)=0
 

Detailed Description

Definition at line 53 of file G4VGFlashSensitiveDetector.hh.

Constructor & Destructor Documentation

G4VGFlashSensitiveDetector::G4VGFlashSensitiveDetector ( )
inline

Definition at line 58 of file G4VGFlashSensitiveDetector.hh.

58 {}
G4VGFlashSensitiveDetector::G4VGFlashSensitiveDetector ( const G4VGFlashSensitiveDetector )
inline

Definition at line 59 of file G4VGFlashSensitiveDetector.hh.

59 {}
virtual G4VGFlashSensitiveDetector::~G4VGFlashSensitiveDetector ( )
inlinevirtual

Definition at line 67 of file G4VGFlashSensitiveDetector.hh.

67 {}

Member Function Documentation

G4bool G4VGFlashSensitiveDetector::Hit ( G4GFlashSpot aSpot)
inline

Definition at line 76 of file G4VGFlashSensitiveDetector.hh.

77  {
78  // This is the public method invoked by GFlashHitMaker for generating
79  // hits. The actual user's implementation for generating hits must be
80  // implemented in GenerateHits() virtual protected method.
81 
82  G4bool result = true;
84  = dynamic_cast<G4VSensitiveDetector *>(this);
85  if(!This)
86  {
87  G4Exception("G4VGFlashSensitiveDetector::Hit()",
88  "InvalidSetup", FatalException,
89  "Needs also to inherit from G4VSensitiveDetector!");
90  return false;
91  }
92  if(This->isActive())
93  {
94  G4VReadOutGeometry * ROgeometry = 0;
95  G4TouchableHistory* ROhis = 0;
96 
97  if(This) ROgeometry = This->GetROgeometry();
98  if(ROgeometry)
99  {
100  // fake pre-step point for touchable from read-out geometry.
101  G4Step fakeStep;
102  G4StepPoint * tmpPoint = fakeStep.GetPreStepPoint();
103  tmpPoint->SetTouchableHandle(aSpot->GetTouchableHandle());
104  tmpPoint->SetPosition(aSpot->GetPosition());
105  tmpPoint->SetMomentumDirection(aSpot->GetOriginatorTrack()
107  result = ROgeometry->CheckROVolume(&fakeStep, ROhis);
108  }
109  if(result) result = ProcessHits(aSpot, ROhis);
110  }
111  else
112  {
113  result = false;
114  }
115  return result;
116  }
G4double G4ParticleHPJENDLHEData::G4double result
G4TouchableHandle GetTouchableHandle() const
Definition: G4GFlashSpot.hh:62
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
void SetPosition(const G4ThreeVector &aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
G4StepPoint * GetPreStepPoint() const
const G4FastTrack * GetOriginatorTrack() const
Definition: G4GFlashSpot.hh:60
virtual G4bool ProcessHits(G4GFlashSpot *aSpot, G4TouchableHistory *ROhist)=0
bool G4bool
Definition: G4Types.hh:79
G4VReadOutGeometry * GetROgeometry() const
Definition: G4Step.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool CheckROVolume(G4Step *, G4TouchableHistory *&)
const G4ThreeVector & GetMomentumDirection() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4ThreeVector GetPosition() const
Definition: G4GFlashSpot.hh:64

Here is the call graph for this function:

Here is the caller graph for this function:

G4int G4VGFlashSensitiveDetector::operator!= ( const G4VGFlashSensitiveDetector right) const
inline

Definition at line 71 of file G4VGFlashSensitiveDetector.hh.

72  {return this != &right;}
G4int G4VGFlashSensitiveDetector::operator== ( const G4VGFlashSensitiveDetector right) const
inline

Definition at line 69 of file G4VGFlashSensitiveDetector.hh.

70  {return this == &right;}
virtual G4bool G4VGFlashSensitiveDetector::ProcessHits ( G4GFlashSpot aSpot,
G4TouchableHistory ROhist 
)
protectedpure virtual

Implemented in ExGflashSensitiveDetector.

Here is the caller graph for this function:


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