Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
UltraPMTSD.cc
이 파일의 문서화 페이지로 가기
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 // --------------------------------------------------------------
28 // GEANT 4 - ULTRA experiment example
29 // --------------------------------------------------------------
30 //
31 // Code developed by:
32 // B. Tome, M.C. Espirito-Santo, A. Trindade, P. Rodrigues
33 //
34 // ****************************************************
35 // * UltraPMTSD.cc
36 // ****************************************************
37 //
38 // Class used to define the Ultra photomultiplier as a sensitive detector.
39 // Hits in this sensitive detector are defined in the UltraOpticalHit class
40 //
41 #include "UltraPMTSD.hh"
42 
43 #include "G4Material.hh"
44 #include "G4Step.hh"
45 #include "G4VTouchable.hh"
46 #include "G4TouchableHistory.hh"
47 #include "G4SDManager.hh"
48 #include "G4OpticalPhoton.hh"
49 #include "G4ParticleDefinition.hh"
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 
54 {
55 
56  collectionName.insert("OpticalHitsCollection");
57 
58 }
59 
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 
64 {;}
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67 
69 {
70 
71 
72  static int HCID1 = -1;
73 
74 
75  // SensitiveDetectorName and collectionName are data members of G4VSensitiveDetector
76 
79 
80  if(HCID1<0)
81  { HCID1 = GetCollectionID(0); }
83 
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87 
89 {
90 
91  // Get Material
92 
93  G4String thisVolume = aStep->GetTrack()->GetVolume()->GetName() ;
94  const G4ParticleDefinition* particle = aStep->GetTrack()->GetDefinition();
95 
96 
97  if (thisVolume != "PMT1" && thisVolume != "PMT2")
98  return false;
99  if (particle != G4OpticalPhoton::Definition() )
100  return false;
101 
102  if(particle == G4OpticalPhoton::Definition())
104 
106  G4ThreeVector HitPosition = aStep->GetPreStepPoint()->GetPosition() ;
107 
108  UltraOpticalHit* OpticalHit = new UltraOpticalHit ;
109  OpticalHit->SetEnergy(kineticEnergy);
110  OpticalHit->SetPosition(HitPosition);
111 
112 
113  OpticalHitsCollection->insert(OpticalHit);
114 
115 
116 #ifdef ULTRA_VERBOSE
117  G4cout << "*******************************" << G4endl;
118  G4cout << " PMT HIT " << G4endl;
119  G4cout << " Volume: " << thisVolume << G4endl;
120  G4cout << " Photon energy (eV) : " << kineticEnergy/eV << G4endl;
121  G4cout << " POSITION (mm) : "
122  << HitPosition.x()/mm << " " << HitPosition.y()/mm << " " << HitPosition.z()/mm << G4endl;
123  G4cout << "*******************************" << G4endl;
124 #endif
125 
126 
127  return true;
128 }
129 
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132 
134 {
135  static G4int HCID = -1;
136  if(HCID<0)
137  {
139  }
141 
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145 
147 {;}
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150 
152 {;}
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
157 {;}
158 
159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
virtual G4int GetCollectionID(G4int i)
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
G4double GetKineticEnergy() const
const XML_Char * name
Definition: expat.h:151
void EndOfEvent(G4HCofThisEvent *)
Definition: UltraPMTSD.cc:133
G4StepPoint * GetPreStepPoint() const
static constexpr double mm
Definition: G4SIunits.hh:115
void PrintAll()
Definition: UltraPMTSD.cc:156
#define G4endl
Definition: G4ios.hh:61
void SetEnergy(G4double fEn)
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
double z() const
G4THitsCollection< UltraOpticalHit > UltraOpticalHitsCollection
UltraOpticalHitsCollection * OpticalHitsCollection
Definition: UltraPMTSD.hh:72
G4ParticleDefinition * GetDefinition() const
void clear()
Definition: UltraPMTSD.cc:146
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetPosition() const
static G4OpticalPhoton * Definition()
G4Track * GetTrack() const
static constexpr double eV
Definition: G4SIunits.hh:215
void SetPosition(G4ThreeVector xyz)
UltraPMTSD(G4String)
Definition: UltraPMTSD.cc:53
Definition: G4Step.hh:76
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
int G4int
Definition: G4Types.hh:78
G4CollectionNameVector collectionName
G4GLOB_DLL std::ostream G4cout
double x() const
void Initialize(G4HCofThisEvent *)
Definition: UltraPMTSD.cc:68
void DrawAll()
Definition: UltraPMTSD.cc:151
void SetTrackStatus(const G4TrackStatus aTrackStatus)
double y() const
G4bool ProcessHits(G4Step *astep, G4TouchableHistory *ROHist)
Definition: UltraPMTSD.cc:88
const G4String & GetName() const
G4VPhysicalVolume * GetVolume() const