Geant4-11
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
34
36{
39 fNaviSetup = false;
41}
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
49void G4FastSimHitMaker::make(const G4FastHit& aHit, const G4FastTrack& aTrack)
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 =
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()),
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
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
G4TouchableHandle fTouchableHandle
Touchable.
G4Navigator * fpNavigator
Navigator.
G4bool fNaviSetup
Flag specifying if navigator has been already set up.
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
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:34