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

Helper class for hit creation. More...

#include <G4FastSimHitMaker.hh>

Public Member Functions

 G4FastSimHitMaker ()
 
 ~G4FastSimHitMaker ()
 
void make (const G4FastHit &aHit, const G4FastTrack &aTrack)
 
void SetNameOfWorldWithSD (const G4String &aName)
 

Detailed Description

Helper class for hit creation.

Helper class that can be employed in the fast simulation models. It allows to deposit energy at given position (G4FastHit), provided it is located within the sensitive detector that derives from G4VFastSimSensitiveDetector base class. An extended example extended/parameterisations/Par03 demonstrates how to use G4FastSimHitMaker to create multiple deposits from the fast simulation model.

Definition at line 47 of file G4FastSimHitMaker.hh.

Constructor & Destructor Documentation

◆ G4FastSimHitMaker()

G4FastSimHitMaker::G4FastSimHitMaker ( )

Definition at line 35 of file G4FastSimHitMaker.cc.

36{
37 fTouchableHandle = new G4TouchableHistory();
38 fpNavigator = new G4Navigator();
39 fNaviSetup = false;
40 fWorldWithSdName = "";
41}

◆ ~G4FastSimHitMaker()

G4FastSimHitMaker::~G4FastSimHitMaker ( )

Definition at line 45 of file G4FastSimHitMaker.cc.

45{ delete fpNavigator; }

Member Function Documentation

◆ make()

void G4FastSimHitMaker::make ( const G4FastHit aHit,
const G4FastTrack aTrack 
)

Deposit energy at given position.

Parameters
[in]aHitCreated hit (energy and position)
[in]aTrackFast track with access to particle's track and properties in envelope's local coordinates

Definition at line 49 of file G4FastSimHitMaker.cc.

50{
51 // do not make empty deposit
52 if(aHit.GetEnergy() <= 0)
53 return;
54 // Locate the spot
55 if(!fNaviSetup)
56 {
57 // Choose the world volume that contains the sensitive detector based on its
58 // name (empty name for mass geometry)
59 G4VPhysicalVolume* worldWithSD = nullptr;
60 if(fWorldWithSdName.empty())
61 {
65 }
66 else
67 {
68 worldWithSD =
70 fWorldWithSdName);
71 }
72 fpNavigator->SetWorldVolume(worldWithSD);
73 // use track global position
75 aTrack.GetPrimaryTrack()->GetPosition(), fTouchableHandle(), false);
76 fNaviSetup = true;
77 }
78 else
79 {
80 // for further deposits use hit (local) position and local->global
81 // transformation
84 aHit.GetPosition()),
85 fTouchableHandle());
86 }
87 G4VPhysicalVolume* currentVolume = fTouchableHandle()->GetVolume();
88
89 G4VSensitiveDetector* sensitive;
90 if(currentVolume != 0)
91 {
92 sensitive = currentVolume->GetLogicalVolume()->GetSensitiveDetector();
93 G4VFastSimSensitiveDetector* fastSimSensitive =
94 dynamic_cast<G4VFastSimSensitiveDetector*>(sensitive);
95 if(fastSimSensitive)
96 {
97 fastSimSensitive->Hit(&aHit, &aTrack, &fTouchableHandle);
98 }
99 else if(sensitive &&
100 currentVolume->GetLogicalVolume()->GetFastSimulationManager())
101 {
102 G4cerr << "ERROR - G4FastSimHitMaker::make()" << G4endl
103 << " It is required to derive from the " << G4endl
104 << " G4VFastSimSensitiveDetector in " << G4endl
105 << " addition to the usual G4VSensitiveDetector class."
106 << G4endl;
107 G4Exception("G4FastSimHitMaker::make()", "InvalidSetup", FatalException,
108 "G4VFastSimSensitiveDetector interface not implemented.");
109 }
110 }
111}
@ 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
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector GetPosition() const
Get position.
Definition: G4FastHit.hh:65
G4double GetEnergy() const
Get energy.
Definition: G4FastHit.hh:58
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
const G4AffineTransform * GetInverseAffineTransformation() const
Definition: G4FastTrack.hh:236
G4VSensitiveDetector * GetSensitiveDetector() const
G4FastSimulationManager * GetFastSimulationManager() const
void LocateGlobalPointAndUpdateTouchable(const G4ThreeVector &position, const G4ThreeVector &direction, G4VTouchable *touchableToUpdate, const G4bool RelativeSearch=true)
void SetWorldVolume(G4VPhysicalVolume *pWorld)
G4VPhysicalVolume * GetWorldVolume() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
Base class for the sensitive detector used within the fast simulation.
G4bool Hit(const G4FastHit *aHit, const G4FastTrack *aTrack, G4TouchableHandle *aTouchable)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41

◆ SetNameOfWorldWithSD()

void G4FastSimHitMaker::SetNameOfWorldWithSD ( const G4String aName)
inline

If sensitive detector class is in the parallel world, it must be specified, otherwise no sensitive detector will be found (mass geometry will be checked).

Parameters
[in]aNameName of the parallel world

Definition at line 62 of file G4FastSimHitMaker.hh.

63 {
64 fWorldWithSdName = aName;
65 };

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