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

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 50 of file G4FastSimHitMaker.hh.

Constructor & Destructor Documentation

◆ G4FastSimHitMaker()

G4FastSimHitMaker::G4FastSimHitMaker ( )

Definition at line 36 of file G4FastSimHitMaker.cc.

37{
38 fTouchableHandle = new G4TouchableHistory();
39 fpNavigator = new G4Navigator();
40 fNaviSetup = false;
41 fWorldWithSdName = "";
42 fpSpotS = new G4Step();
43 fpSpotP = new G4StepPoint();
44 // N.B. Pre and Post step points are common.
45 fpSpotS->SetPreStepPoint(fpSpotP);
46 fpSpotS->SetPostStepPoint(fpSpotP);
47}

◆ ~G4FastSimHitMaker()

G4FastSimHitMaker::~G4FastSimHitMaker ( )

Definition at line 51 of file G4FastSimHitMaker.cc.

52{
53 delete fpNavigator;
54 delete fpSpotP;
55 fpSpotS->ResetPreStepPoint();
56 fpSpotS->ResetPostStepPoint();
57 delete fpSpotS;
58}

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 62 of file G4FastSimHitMaker.cc.

63{
64 // do not make empty deposit
65 if (aHit.GetEnergy() <= 0) return;
66 // Locate the spot
67 if (!fNaviSetup) {
68 // Choose the world volume that contains the sensitive detector based on its
69 // name (empty name for mass geometry)
70 G4VPhysicalVolume* worldWithSD = nullptr;
71 if (fWorldWithSdName.empty()) {
75 }
76 else {
77 worldWithSD =
79 }
80 fpNavigator->SetWorldVolume(worldWithSD);
81 // use track global position
82 fpNavigator->LocateGlobalPointAndUpdateTouchable(aTrack.GetPrimaryTrack()->GetPosition(),
83 fTouchableHandle(), false);
84 fNaviSetup = true;
85 }
86 else {
87 // for further deposits use hit (local) position and local->global
88 // transformation
89 fpNavigator->LocateGlobalPointAndUpdateTouchable(
91 fTouchableHandle());
92 }
93 G4VPhysicalVolume* currentVolume = fTouchableHandle()->GetVolume();
94
95 if (currentVolume != nullptr) {
96 G4VSensitiveDetector* sensitive = currentVolume->GetLogicalVolume()->GetSensitiveDetector();
97 auto fastSimSensitive = dynamic_cast<G4VFastSimSensitiveDetector*>(sensitive);
98 if (fastSimSensitive != nullptr) {
99 fastSimSensitive->Hit(&aHit, &aTrack, &fTouchableHandle);
100 }
101 else if ((sensitive != nullptr)
102 && (currentVolume->GetLogicalVolume()->GetFastSimulationManager() != nullptr))
103 {
104 fpSpotS->SetTotalEnergyDeposit(aHit.GetEnergy());
105 fpSpotS->SetTrack(const_cast<G4Track*>(aTrack.GetPrimaryTrack()));
106 fpSpotP->SetWeight(aTrack.GetPrimaryTrack()->GetWeight());
107 fpSpotP->SetPosition(aHit.GetPosition());
108 fpSpotP->SetGlobalTime(aTrack.GetPrimaryTrack()->GetGlobalTime());
109 fpSpotP->SetLocalTime(aTrack.GetPrimaryTrack()->GetLocalTime());
110 fpSpotP->SetProperTime(aTrack.GetPrimaryTrack()->GetProperTime());
111 fpSpotP->SetTouchableHandle(fTouchableHandle);
112 fpSpotP->SetProcessDefinedStep(fpProcess);
113 fpSpotP->SetStepStatus(fUserDefinedLimit);
114 sensitive->Hit(fpSpotS);
115 }
116 }
117}
@ fUserDefinedLimit
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector GetPosition() const
Get position.
Definition G4FastHit.hh:63
G4double GetEnergy() const
Get energy.
Definition G4FastHit.hh:59
const G4Track * GetPrimaryTrack() const
const G4AffineTransform * GetInverseAffineTransformation() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4FastSimulationManager * GetFastSimulationManager() const
G4VPhysicalVolume * GetWorldVolume() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4double GetProperTime() const
G4double GetLocalTime() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)

◆ 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 65 of file G4FastSimHitMaker.hh.

65{ fWorldWithSdName = aName; };

◆ SetProcess()

void G4FastSimHitMaker::SetProcess ( G4VProcess * proc)
inline

Definition at line 66 of file G4FastSimHitMaker.hh.

66{ fpProcess = proc; }

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