ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimPlaGeantHit.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 04/2008
2 /******************************************************************
3  * Copyright (C) 2005-2016, Hector Alvarez-Pol *
4  * All rights reserved. *
5  * *
6  * License according to GNU LESSER GPL (see lgpl-3.0.txt). *
7  * For the list of contributors see CREDITS. *
8  ******************************************************************/
9 //////////////////////////////////////////////////////////////////
10 /// \class ActarSimPlaGeantHit
11 /// A Geant Hit in the plastic (scintillator) volume. It
12 /// represents the information of each step with energy deposited.
13 /////////////////////////////////////////////////////////////////
14 
15 #include "ActarSimPlaGeantHit.hh"
16 #include "G4UnitsTable.hh"
17 #include "G4VVisManager.hh"
18 #include "G4Circle.hh"
19 #include "G4Colour.hh"
20 #include "G4VisAttributes.hh"
21 
22 #include "G4PhysicalConstants.hh"
23 #include "G4SystemOfUnits.hh"
24 
25 G4Allocator<ActarSimPlaGeantHit> ActarSimPlaGeantHitAllocator;
26 
27 //////////////////////////////////////////////////////////////////
28 /// Constructor
30 }
31 
32 //////////////////////////////////////////////////////////////////
33 /// Destructor
35 }
36 
37 //////////////////////////////////////////////////////////////////
38 /// Copy constructor
40  edep = right.edep;
41  eBeforePla = right.eBeforePla;
42  eAfterPla = right.eAfterPla;
43  pos = right.pos;
44  prePos = right.prePos;
45  localPos = right.localPos;
46  localPrePos = right.localPrePos;
47  detName = right.detName;
48  postDetName = right.postDetName;
49  preDetName = right.preDetName;
50  detID = right.detID;
51  toF = right.toF;
52  trackID = right.trackID;
53  parentID = right.parentID;
54  particleID = right.particleID;
56  particleMass = right.particleMass;
57 }
58 
59 //////////////////////////////////////////////////////////////////
60 /// Operator =
62  edep = right.edep;
63  eBeforePla = right.eBeforePla;
64  eAfterPla = right.eAfterPla;
65  pos = right.pos;
66  prePos = right.prePos;
67  localPos = right.localPos;
68  localPrePos = right.localPrePos;
69  detName = right.detName;
70  postDetName = right.postDetName;
71  preDetName = right.preDetName;
72  detID = right.detID;
73  toF = right.toF;
74  trackID = right.trackID;
75  parentID = right.parentID;
76  particleID = right.particleID;
78  particleMass = right.particleMass;
79  return *this;
80 }
81 
82 //////////////////////////////////////////////////////////////////
83 /// Operator ==
85  return (this==&right) ? 1 : 0;
86 }
87 
88 //////////////////////////////////////////////////////////////////
89 /// Draws the Hit. A clear red point on the Hit position
91  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
92  if(pVVisManager) {
93  G4Circle circle(pos);
94  circle.SetScreenSize(4);
95  circle.SetFillStyle(G4Circle::filled);
96  G4Colour colour(1.,0.,0.);
97  G4VisAttributes attribs(colour);
98  circle.SetVisAttributes(attribs);
99  pVVisManager->Draw(circle);
100  }
101 }
102 
103 //////////////////////////////////////////////////////////////////
104 /// Prints full information about the calGeantHit
106  G4cout << "##################################################################"
107  << G4endl
108  << "############### ActarSimSciGeantHit::Print() ###################" << G4endl
109  << "trackID: " << trackID
110  << ", parentID: " << parentID
111  << ", particleID: " << particleID
112  << ", particleCharge: " << particleCharge << G4endl;
113  G4cout << "detID: " << detID
114  << ", detName: " << detName
115  << ", postDetName: " << postDetName
116  << ", preDetName: " << preDetName
117  << G4endl;
118  G4cout << "edep: " << edep / MeV << " MeV"
119  << ", pos: " << pos << " mm" << G4endl
120  << ", prePos: " << prePos << " mm" << G4endl;
121  G4cout << "toF: " << toF / ns << " ns" << ", localPos: " << localPos << " mm"
122  << ", localPrePos: " << localPrePos << " mm" << G4endl;
123  G4cout << "##################################################################"
124  << G4endl;
125 }
void Print()
Prints full information about the calGeantHit.
G4ThreeVector localPrePos
Local (for the given detName and detID) coordinates of interaction (postStep)
G4double eBeforePla
Energy before entering plastic.
G4ThreeVector pos
PostStep position of the step.
G4String detName
Name of the volume where the interaction takes place.
G4ThreeVector localPos
Local (for the given detName and detID) coordinates of interaction (postStep)
G4ThreeVector prePos
PreStep position of the step.
G4int detID
ID (copy) of the detector where the interaction takes place.
G4double particleCharge
Particle charge.
const ActarSimPlaGeantHit & operator=(const ActarSimPlaGeantHit &)
Operator =.
G4Allocator< ActarSimPlaGeantHit > ActarSimPlaGeantHitAllocator
ActarSimPlaGeantHit()
Constructor.
void Draw()
Draws the Hit. A clear red point on the Hit position.
G4double toF
ToF of the interaction (postStep)
G4int parentID
Parent ID.
~ActarSimPlaGeantHit()
Destructor.
G4double particleMass
Particle mass.
G4int operator==(const ActarSimPlaGeantHit &) const
Operator ==.
G4int particleID
Particle ID.
G4String preDetName
Name of the volume at the following step.
G4String postDetName
Name of the volume at the previous step.
G4double edep
Energy deposited in the step.
G4double eAfterPla
Energy when exiting plastic.