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

#include <GFlashHitMaker.hh>

Public Member Functions

 GFlashHitMaker ()
 
 ~GFlashHitMaker ()
 
void make (GFlashEnergySpot *aSpot, const G4FastTrack *aT)
 
void SetNameOfWorldWithSD (const G4String &aName)
 

Detailed Description

Definition at line 51 of file GFlashHitMaker.hh.

Constructor & Destructor Documentation

◆ GFlashHitMaker()

GFlashHitMaker::GFlashHitMaker ( )

Definition at line 45 of file GFlashHitMaker.cc.

46{
47 fTouchableHandle = new G4TouchableHistory(); // talk to ?@@@
48 fpNavigator = new G4Navigator();
49 fNaviSetup = false;
50 fWorldWithSdName = "";
51}

◆ ~GFlashHitMaker()

GFlashHitMaker::~GFlashHitMaker ( )

Definition at line 53 of file GFlashHitMaker.cc.

54{
55 delete fpNavigator;
56}

Member Function Documentation

◆ make()

void GFlashHitMaker::make ( GFlashEnergySpot aSpot,
const G4FastTrack aT 
)

Navigator

Definition at line 58 of file GFlashHitMaker.cc.

59{
60 // Locate the spot
61 if (!fNaviSetup)
62 {
63 // Choose the world volume that contains the sensitive detector based on its name (empty name for mass geometry)
64 G4VPhysicalVolume* worldWithSD = nullptr;
65 if(fWorldWithSdName.empty()) {
67 } else {
69 }
70 fpNavigator->SetWorldVolume(worldWithSD);
71 fpNavigator->
72 LocateGlobalPointAndUpdateTouchable(aSpot->GetPosition(),
73 fTouchableHandle(), false);
74 fNaviSetup = true;
75 }
76 else
77 {
78 fpNavigator->
79 LocateGlobalPointAndUpdateTouchable(aSpot->GetPosition(),
80 fTouchableHandle());
81 }
82
83 //--------------------------------------
84 // Fills attribute of the G4Step needed
85 // by our sensitive detector:
86 //-------------------------------------
87 // set spot information:
88 G4GFlashSpot theSpot(aSpot, aT, fTouchableHandle);
89 ///Navigator
90 //--------------------------------------
91 // Produce Hits
92 // call sensitive part: taken/adapted from the stepping:
93 // Send G4Step information to Hit/Dig if the volume is sensitive
94 //--------------G4TouchableHistory----------------------------------------
95
96 G4VPhysicalVolume* pCurrentVolume = fTouchableHandle()->GetVolume();
97 G4VSensitiveDetector* pSensitive;
98 if( pCurrentVolume != 0 )
99 {
100 pSensitive = pCurrentVolume->GetLogicalVolume()->GetSensitiveDetector();
101 G4VGFlashSensitiveDetector * gflashSensitive =
102 dynamic_cast<G4VGFlashSensitiveDetector * > (pSensitive);
103 if( gflashSensitive )
104 {
105 gflashSensitive->Hit(&theSpot);
106 }
107 else if (( pSensitive ) &&
108 ( pCurrentVolume->GetLogicalVolume()->GetFastSimulationManager() )
109 ) // Using gflash without implementing the
110 // gflashSensitive detector interface -> not allowed!
111
112 {
113 G4cerr << "ERROR - GFlashHitMaker::make()" << G4endl
114 << " It is required to implement the "<< G4endl
115 << " G4VGFlashSensitiveDetector interface in "<< G4endl
116 << " addition to the usual SensitiveDetector class."
117 << G4endl;
118 G4Exception("GFlashHitMaker::make()", "InvalidSetup", FatalException,
119 "G4VGFlashSensitiveDetector interface not implemented.");
120 }
121 }
122 else
123 {
124 #ifdef GFLASH_DEBUG
125 G4cout << "GFlashHitMaker::Out of volume "<< G4endl;
126 #endif
127 }
128}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VSensitiveDetector * GetSensitiveDetector() const
G4FastSimulationManager * GetFastSimulationManager() const
void SetWorldVolume(G4VPhysicalVolume *pWorld)
G4VPhysicalVolume * GetWorldVolume() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4bool Hit(G4GFlashSpot *aSpot)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41
G4ThreeVector GetPosition() const

◆ SetNameOfWorldWithSD()

void GFlashHitMaker::SetNameOfWorldWithSD ( const G4String aName)
inline

Definition at line 59 of file GFlashHitMaker.hh.

59{fWorldWithSdName = aName;};

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