Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastSimHitMaker.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26
27#include "G4FastSimHitMaker.hh"
28
29#include "G4Step.hh"
30#include "G4StepPoint.hh"
31#include "G4TouchableHandle.hh"
35
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}
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
52{
53 delete fpNavigator;
54 delete fpSpotP;
55 fpSpotS->ResetPreStepPoint();
56 fpSpotS->ResetPostStepPoint();
57 delete fpSpotS;
58}
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62void G4FastSimHitMaker::make(const G4FastHit& aHit, const G4FastTrack& aTrack)
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
83 fTouchableHandle(), false);
84 fNaviSetup = true;
85 }
86 else {
87 // for further deposits use hit (local) position and local->global
88 // transformation
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);
114 sensitive->Hit(fpSpotS);
115 }
116 }
117}
@ fUserDefinedLimit
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
Minimal hit created in the fast simulation.
Definition G4FastHit.hh:48
G4ThreeVector GetPosition() const
Get position.
Definition G4FastHit.hh:63
G4double GetEnergy() const
Get energy.
Definition G4FastHit.hh:59
void make(const G4FastHit &aHit, const G4FastTrack &aTrack)
const G4Track * GetPrimaryTrack() const
const G4AffineTransform * GetInverseAffineTransformation() const
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
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)
G4StepPoint * ResetPreStepPoint(G4StepPoint *value=nullptr)
void SetPostStepPoint(G4StepPoint *value)
void SetPreStepPoint(G4StepPoint *value)
void SetTotalEnergyDeposit(G4double value)
G4StepPoint * ResetPostStepPoint(G4StepPoint *value=nullptr)
void SetTrack(G4Track *value)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) 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
Base class for the sensitive detector used within the fast simulation.
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)