FDsimG4OpticalHit.cc
Go to the documentation of this file.
1 #include "FDsimG4OpticalHit.hh"
2 #include "FDsimG4Colours.hh"
3 #include "G4VVisManager.hh"
4 #include "G4Colour.hh"
5 #include "G4Circle.hh"
6 #include "G4VisAttributes.hh"
7 
8 using namespace TelescopeSimulatorLX;
9 
10 G4Allocator<FDsimG4OpticalHit> FDsimG4OpticalHitAllocator;
11 
12 
14  fPhotEne(0.0),
15  fPhotWvl(0.0),
16  fPhotTime(0.0),
17  fPhotPos(0.0),
18  fPhotDir(0.0),
19  fPixelPos(0.0),
20  fHitID(0),
21  fWeight(0)
22 {
23 }
24 
25 
27 {
28 }
29 
30 
32 {
33  fPhotEne = right.fPhotEne;
34  fPhotWvl = right.fPhotWvl;
35  fPhotTime = right.fPhotTime;
36  fPhotPos = right.fPhotPos;
37  fPhotDir = right.fPhotDir;
38  fPixelPos = right.fPixelPos;
39  fHitID = right.fHitID;
40  fWeight = right.fWeight;
41 }
42 
43 
45 {
46  fPhotEne = right.fPhotEne;
47  fPhotWvl = right.fPhotWvl;
48  fPhotTime = right.fPhotTime;
49  fPhotPos = right.fPhotPos;
50  fPhotDir = right.fPhotDir;
51  fPixelPos = right.fPixelPos;
52  fHitID = right.fHitID;
53  fWeight = right.fWeight;
54 
55  return *this;
56 }
57 
58 
60 {
61  return (this==&right) ? 1 : 0;
62 }
63 
64 
66 {
67  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
68  if(pVVisManager)
69  {
70  // define a circle in a 3D space
71  G4double tMin = 0*ns;
72  G4double tMax = 20000*ns;
73  G4double iRed = 0;
74  G4double iBlue = ((tMax-fPhotTime)/(tMax-tMin));
75  G4double iGreen = 0;
76 
77  G4Circle circle(fPhotWorldPos);
78  // circle.SetScreenSize(1.0);
79 
80  circle.SetWorldSize(10*mm);
81  // circle.SetScreenDiameter(2.0);
82  circle.SetFillStyle(G4Circle::filled);
83  // make the circle red
84  G4Colour colour(iRed,iGreen,iBlue);
85  // G4Colour colour(0,200,0);
86  // G4VisAttributes attribs(FDsimG4Colours::yellow);
87 
88  G4VisAttributes attribs(colour);
89  circle.SetVisAttributes(attribs);
90  pVVisManager->Draw(circle);
91  }
92 }
93 
94 
96 {
97  G4cerr << "*****************************************" << G4endl;
98  G4cerr << " Optical HIT " << G4endl;
99  G4cerr << " Energy (eV) : " << fPhotEne/eV << G4endl;
100  G4cerr << " Wavelength (nm) : " << fPhotWvl/nanometer << G4endl;
101  G4cerr << " Position (mm) : " << fPhotWorldPos/mm << G4endl;
102  G4cerr << " Direction(mm) : " << fPhotDir << G4endl;
103  G4cerr << " Time (ns) : " << fPhotTime/ns << G4endl;
104  G4cerr << " Pixel Position (mm) : " << fPixelPos/mm << G4endl;
105  G4cerr << " Weight : " << fWeight << G4endl;
106  G4cerr << "*****************************************" << G4endl;
107 }
const double eV
Definition: GalacticUnits.h:35
constexpr double mm
Definition: AugerUnits.h:113
const FDsimG4OpticalHit & operator=(const FDsimG4OpticalHit &)
constexpr double nanometer
Definition: AugerUnits.h:102
G4Allocator< FDsimG4OpticalHit > FDsimG4OpticalHitAllocator
int operator==(const FDsimG4OpticalHit &) const
const double ns

, generated on Tue Sep 26 2023.