OMSim
Geant4 for IceCube optical module studies
Loading...
Searching...
No Matches
OMSimSensitiveDetector.hh
1#pragma once
2
3#include "G4VSensitiveDetector.hh"
4#include "G4ThreeVector.hh"
5#include "OMSimPMTResponse.hh"
6#include "OMSimOpBoundaryProcess.hh"
7#include <vector>
8
9class G4Step;
10class G4TouchableHistory;
11
12
18enum class DetectorType {
19 PMT,
20 VolumePhotonDetector,
21 BoundaryPhotonDetector,
22 PerfectPMT
23};
24
29struct PhotonInfo {
30 G4int eventID;
31 G4double globalTime;
32 G4double localTime;
33 G4double trackLength;
34 G4double kineticEnergy;
35 G4double wavelength;
36 G4ThreeVector globalPosition;
37 G4ThreeVector localPosition;
38 G4ThreeVector momentumDirection;
39 G4ThreeVector deltaPosition;
40 G4int pmtNumber;
41 G4int detectorID;
43};
44
45
46
51class OMSimSensitiveDetector : public G4VSensitiveDetector
52{
53public:
54 OMSimSensitiveDetector(G4String pName, DetectorType pDetectorType);
56
57 G4bool ProcessHits(G4Step *pStep, G4TouchableHistory *pTouchableHistory) override;
58 void setPMTResponse(OMSimPMTResponse *pResponse);
59
60private:
61 bool m_QEcut;
62 OMSimPMTResponse *m_PMTResponse;
63 DetectorType m_detectorType;
65
66 G4bool checkVolumeAbsorption(G4Step *pStep);
67 G4bool checkBoundaryAbsorption(G4Step *pStep);
68 PhotonInfo getPhotonInfo(G4Step *pStep);
69 G4bool handlePMT(G4Step *pStep, G4TouchableHistory *pTouchableHistory);
70 G4bool handleGeneralPhotonDetector(G4Step *pStep, G4TouchableHistory *pTouchableHistory);
71 bool isPhotonDetected(double p_efficiency);
72 void storePhotonHit(PhotonInfo &pInfo);
74 void killParticle(G4Track *pTrack);
75};
Simulation of PMT response.
Definition OMSimOpBoundaryProcess.hh:159
Class to simulate PMT response.
Definition OMSimPMTResponse.hh:20
Represents a sensitive detector.
Definition OMSimSensitiveDetector.hh:52
G4bool handlePMT(G4Step *pStep, G4TouchableHistory *pTouchableHistory)
Handles hits for PMT detectors.
Definition OMSimSensitiveDetector.cc:207
G4bool ProcessHits(G4Step *pStep, G4TouchableHistory *pTouchableHistory) override
Processes hits for optical photons in the detector.
Definition OMSimSensitiveDetector.cc:90
void storePhotonHit(PhotonInfo &pInfo)
Stores photon hit information into the HitManager.
Definition OMSimSensitiveDetector.cc:252
void setPMTResponse(OMSimPMTResponse *pResponse)
Sets the PMT response model.
Definition OMSimSensitiveDetector.cc:74
PhotonInfo getPhotonInfo(G4Step *pStep)
Retrieves photon information from a given step.
Definition OMSimSensitiveDetector.cc:130
G4bool handleGeneralPhotonDetector(G4Step *pStep, G4TouchableHistory *pTouchableHistory)
Handles hits for general photon detectors.
Definition OMSimSensitiveDetector.cc:239
~OMSimSensitiveDetector()
Destructor for OMSimSensitiveDetector.
Definition OMSimSensitiveDetector.cc:40
static thread_local G4OpBoundaryProcess * m_boundaryProcess
Definition OMSimSensitiveDetector.hh:64
bool isPhotonDetected(double p_efficiency)
Monte carlo if the photon was detected based on the detection probability.
Definition OMSimSensitiveDetector.cc:196
void fetchBoundaryProcess()
Fetches the boundary process for detecting boundary absorptions.
Definition OMSimSensitiveDetector.cc:52
G4bool checkVolumeAbsorption(G4Step *pStep)
Checks if the photon was absorbed in the volume.
Definition OMSimSensitiveDetector.cc:162
void killParticle(G4Track *pTrack)
Stop the particle from propagating further. Necessary for 100% efficient detectors.
Definition OMSimSensitiveDetector.cc:273
G4bool checkBoundaryAbsorption(G4Step *pStep)
Checks if the photon was detected at a boundary.
Definition OMSimSensitiveDetector.cc:172
Represents the output pulse for a detected photon.
Definition OMSimPMTResponse.hh:29
Contains information about a detected photon which will be appended in HitManager.
Definition OMSimSensitiveDetector.hh:29
G4double globalTime
Global time of the photon hit.
Definition OMSimSensitiveDetector.hh:31
OMSimPMTResponse::PMTPulse PMTResponse
PMT response to the photon hit.
Definition OMSimSensitiveDetector.hh:42
G4ThreeVector deltaPosition
Change in position of the photon hit.
Definition OMSimSensitiveDetector.hh:39
G4double kineticEnergy
Kinetic energy of the photon.
Definition OMSimSensitiveDetector.hh:34
G4double wavelength
Wavelength of the photon.
Definition OMSimSensitiveDetector.hh:35
G4int eventID
Event ID of the photon hit.
Definition OMSimSensitiveDetector.hh:30
G4ThreeVector globalPosition
Global position of the photon hit.
Definition OMSimSensitiveDetector.hh:36
G4ThreeVector momentumDirection
Direction of the photon's momentum.
Definition OMSimSensitiveDetector.hh:38
G4int pmtNumber
PMT number associated with the photon hit.
Definition OMSimSensitiveDetector.hh:40
G4double trackLength
Length of the photon's track.
Definition OMSimSensitiveDetector.hh:33
G4double localTime
Local time of the photon hit.
Definition OMSimSensitiveDetector.hh:32
G4int detectorID
ID of the detector registering the photon.
Definition OMSimSensitiveDetector.hh:41
G4ThreeVector localPosition
Local position of the photon hit.
Definition OMSimSensitiveDetector.hh:37