Geant4 11.2.2
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)
 
void SetProcess (G4VProcess *proc)
 

Detailed Description

Definition at line 55 of file GFlashHitMaker.hh.

Constructor & Destructor Documentation

◆ GFlashHitMaker()

GFlashHitMaker::GFlashHitMaker ( )

Definition at line 47 of file GFlashHitMaker.cc.

48{
49 fTouchableHandle = new G4TouchableHistory(); // talk to ?@@@
50 fpNavigator = new G4Navigator();
51 fNaviSetup = false;
52 fWorldWithSdName = "";
53 fpSpotS = new G4Step();
54 fpSpotP = new G4StepPoint();
55 // N.B. Pre and Post step points are common.
56 fpSpotS->SetPreStepPoint(fpSpotP);
57 fpSpotS->SetPostStepPoint(fpSpotP);
58}
void SetPostStepPoint(G4StepPoint *value)
void SetPreStepPoint(G4StepPoint *value)

◆ ~GFlashHitMaker()

GFlashHitMaker::~GFlashHitMaker ( )

Definition at line 60 of file GFlashHitMaker.cc.

61{
62 delete fpNavigator;
63 delete fpSpotP;
64 fpSpotS->ResetPreStepPoint();
65 fpSpotS->ResetPostStepPoint();
66 delete fpSpotS;
67}
G4StepPoint * ResetPreStepPoint(G4StepPoint *value=nullptr)
G4StepPoint * ResetPostStepPoint(G4StepPoint *value=nullptr)

Member Function Documentation

◆ make()

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

Definition at line 69 of file GFlashHitMaker.cc.

70{
71 // Locate the spot
72 if (!fNaviSetup)
73 {
74 // Choose the world volume that contains the sensitive detector based on its name (empty name for mass geometry)
75 G4VPhysicalVolume* worldWithSD = nullptr;
76 if(fWorldWithSdName.empty()) {
78 } else {
80 }
81 fpNavigator->SetWorldVolume(worldWithSD);
82 fpNavigator->
83 LocateGlobalPointAndUpdateTouchable(aSpot->GetPosition(),
84 fTouchableHandle(), false);
85 fNaviSetup = true;
86 }
87 else
88 {
89 fpNavigator->
90 LocateGlobalPointAndUpdateTouchable(aSpot->GetPosition(),
91 fTouchableHandle());
92 }
93
94 //--------------------------------------
95 // Produce Hits
96 // call sensitive part: taken/adapted from the stepping:
97 // Send G4Step information to Hit/Dig if the volume is sensitive
98 //--------------G4TouchableHistory----------------------------------------
99
100 G4VPhysicalVolume* pCurrentVolume = fTouchableHandle()->GetVolume();
101 G4VSensitiveDetector* pSensitive;
102 if( pCurrentVolume != 0 )
103 {
104 pSensitive = pCurrentVolume->GetLogicalVolume()->GetSensitiveDetector();
105 G4VGFlashSensitiveDetector * gflashSensitive =
106 dynamic_cast<G4VGFlashSensitiveDetector * > (pSensitive);
107 if( gflashSensitive )
108 {
109 // set spot information:
110 G4GFlashSpot theSpot(aSpot, aT, fTouchableHandle);
111 gflashSensitive->Hit(&theSpot);
112 }
113 else if( pSensitive )
114 {
115 fpSpotS->SetTotalEnergyDeposit(aSpot->GetEnergy());
116 fpSpotS->SetTrack(const_cast<G4Track*>(aT->GetPrimaryTrack()));
117 fpSpotP->SetWeight(aT->GetPrimaryTrack()->GetWeight());
118 fpSpotP->SetPosition(aSpot->GetPosition());
120 fpSpotP->SetLocalTime(aT->GetPrimaryTrack()->GetLocalTime());
122 fpSpotP->SetTouchableHandle(fTouchableHandle);
123 fpSpotP->SetProcessDefinedStep(fpProcess);
125 pSensitive->Hit(fpSpotS);
126 }
127 }
128 else
129 {
130 #ifdef GFLASH_DEBUG
131 G4cout << "GFlashHitMaker::Out of volume "<< G4endl;
132 #endif
133 }
134}
@ fUserDefinedLimit
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4Track * GetPrimaryTrack() const
G4VSensitiveDetector * GetSensitiveDetector() const
void SetWorldVolume(G4VPhysicalVolume *pWorld)
G4VPhysicalVolume * GetWorldVolume() const
void SetLocalTime(const G4double aValue)
void SetWeight(G4double aValue)
void SetStepStatus(const G4StepStatus aValue)
void SetProcessDefinedStep(const G4VProcess *aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetProperTime(const G4double aValue)
void SetGlobalTime(const G4double aValue)
void SetPosition(const G4ThreeVector &aValue)
void SetTotalEnergyDeposit(G4double value)
void SetTrack(G4Track *value)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4double GetWeight() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4double GetLocalTime() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4bool Hit(G4GFlashSpot *aSpot)
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)
G4ThreeVector GetPosition() const
G4double GetEnergy() const

◆ SetNameOfWorldWithSD()

void GFlashHitMaker::SetNameOfWorldWithSD ( const G4String & aName)
inline

Definition at line 63 of file GFlashHitMaker.hh.

63{fWorldWithSdName = aName;};

◆ SetProcess()

void GFlashHitMaker::SetProcess ( G4VProcess * proc)
inline

Definition at line 65 of file GFlashHitMaker.hh.

65{ fpProcess = proc; }

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