Geant4 11.1.1
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
31#include "G4TouchableHandle.hh"
32#include "G4Step.hh"
33#include "G4StepPoint.hh"
34
36
38{
39 fTouchableHandle = new G4TouchableHistory();
40 fpNavigator = new G4Navigator();
41 fNaviSetup = false;
42 fWorldWithSdName = "";
43 fpSpotS = new G4Step();
44 fpSpotP = new G4StepPoint();
45 // N.B. Pre and Post step points are common.
46 fpSpotS->SetPreStepPoint(fpSpotP);
47 fpSpotS->SetPostStepPoint(fpSpotP);
48}
49
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51
53{
54 delete fpNavigator;
55 delete fpSpotP;
56 fpSpotS->ResetPreStepPoint();
57 fpSpotS->ResetPostStepPoint();
58 delete fpSpotS;
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
63void G4FastSimHitMaker::make(const G4FastHit& aHit, const G4FastTrack& aTrack)
64{
65 // do not make empty deposit
66 if(aHit.GetEnergy() <= 0)
67 return;
68 // Locate the spot
69 if(!fNaviSetup)
70 {
71 // Choose the world volume that contains the sensitive detector based on its
72 // name (empty name for mass geometry)
73 G4VPhysicalVolume* worldWithSD = nullptr;
74 if(fWorldWithSdName.empty())
75 {
79 }
80 else
81 {
82 worldWithSD =
84 fWorldWithSdName);
85 }
86 fpNavigator->SetWorldVolume(worldWithSD);
87 // use track global position
89 aTrack.GetPrimaryTrack()->GetPosition(), fTouchableHandle(), false);
90 fNaviSetup = true;
91 }
92 else
93 {
94 // for further deposits use hit (local) position and local->global
95 // transformation
98 aHit.GetPosition()),
99 fTouchableHandle());
100 }
101 G4VPhysicalVolume* currentVolume = fTouchableHandle()->GetVolume();
102
103 if(currentVolume != 0)
104 {
105 G4VSensitiveDetector* sensitive = currentVolume->GetLogicalVolume()->GetSensitiveDetector();
106 G4VFastSimSensitiveDetector* fastSimSensitive =
107 dynamic_cast<G4VFastSimSensitiveDetector*>(sensitive);
108 if(fastSimSensitive)
109 {
110 fastSimSensitive->Hit(&aHit, &aTrack, &fTouchableHandle);
111 }
112 else if(sensitive &&
113 currentVolume->GetLogicalVolume()->GetFastSimulationManager())
114 {
115 fpSpotS->SetTotalEnergyDeposit(aHit.GetEnergy());
116 fpSpotS->SetTrack(const_cast<G4Track*>(aTrack.GetPrimaryTrack()));
117 fpSpotP->SetWeight(aTrack.GetPrimaryTrack()->GetWeight());
118 fpSpotP->SetPosition(aHit.GetPosition());
119 fpSpotP->SetGlobalTime(aTrack.GetPrimaryTrack()->GetGlobalTime());
120 fpSpotP->SetLocalTime(aTrack.GetPrimaryTrack()->GetLocalTime());
121 fpSpotP->SetProperTime(aTrack.GetPrimaryTrack()->GetProperTime());
122 fpSpotP->SetTouchableHandle(fTouchableHandle);
123 fpSpotP->SetProcessDefinedStep(fpProcess);
125 sensitive->Hit(fpSpotS);
126 }
127 }
128}
@ fUserDefinedLimit
Definition: G4StepStatus.hh:51
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:65
G4double GetEnergy() const
Get energy.
Definition: G4FastHit.hh:58
void make(const G4FastHit &aHit, const G4FastTrack &aTrack)
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
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)
Definition: G4Step.hh:62
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)
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.
G4bool Hit(const G4FastHit *aHit, const G4FastTrack *aTrack, G4TouchableHandle *aTouchable)
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34