ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimGasGeantHit.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 04/2006
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 ActarSimGasGeantHit
11 /// A Geant Hit in the calorimeter volume. It represents the
12 /// information of each step with energy deposited in the
13 /// calorimeter volume.
14 /////////////////////////////////////////////////////////////////
15 
16 #include "ActarSimGasGeantHit.hh"
17 #include "G4UnitsTable.hh"
18 #include "G4VVisManager.hh"
19 #include "G4Circle.hh"
20 #include "G4Colour.hh"
21 #include "G4VisAttributes.hh"
22 
23 #include "G4PhysicalConstants.hh"
24 #include "G4SystemOfUnits.hh"
25 
26 G4Allocator<ActarSimGasGeantHit> ActarSimGasGeantHitAllocator;
27 
28 //////////////////////////////////////////////////////////////////
29 /// Constructor
31 }
32 
33 //////////////////////////////////////////////////////////////////
34 /// Destructor
36 }
37 
38 //////////////////////////////////////////////////////////////////
39 /// Copy constructor
41  trackID = right.trackID;
42  parentID = right.parentID;
43  edep = right.edep;
45  particleMass = right.particleMass;
46  particleID = right.particleID;
47  prePos = right.prePos;
48  postPos = right.postPos;
49  detName = right.detName;
50  detID = right.detID;
51  preToF = right.preToF;
52  postToF = right.postToF;
53  stepLength = right.stepLength;
54  stepEnergy = right.stepEnergy;
55 }
56 
57 //////////////////////////////////////////////////////////////////
58 /// Operator =
60  trackID = right.trackID;
61  parentID = right.parentID;
62  edep = right.edep;
64  particleMass = right.particleMass;
65  particleID = right.particleID;
66  prePos = right.prePos;
67  postPos = right.postPos;
68  detName = right.detName;
69  detID = right.detID;
70  preToF = right.preToF;
71  postToF = right.postToF;
72  stepLength = right.stepLength;
73  stepEnergy = right.stepEnergy;
74 
75  return *this;
76 }
77 
78 //////////////////////////////////////////////////////////////////
79 /// Operator ==
81  return (this==&right) ? 1 : 0;
82 }
83 
84 //////////////////////////////////////////////////////////////////
85 /// Draws the Hit. A clear red point on the Hit position
87  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
88  if(pVVisManager) {
89  G4Circle circle(prePos);
90  circle.SetScreenSize(4);
91  circle.SetFillStyle(G4Circle::filled);
92  G4Colour colour(1.,0.,0.);
93  G4VisAttributes attribs(colour);
94  circle.SetVisAttributes(attribs);
95  pVVisManager->Draw(circle);
96  }
97 }
98 
99 //////////////////////////////////////////////////////////////////
100 /// Prints full information about the calGeantHit
102  G4cout << "##################################################################"
103  << G4endl
104  << "############ ActarSimGasGeantHit::Print() ################" << G4endl
105  << "trackID: " << trackID
106  << "parentID: " << parentID
107  << ", detID: " << detID
108  << ", detName: " << detName << G4endl;
109  G4cout << "edep: " << edep / MeV << " MeV"
110  << ", stepEnergy: " << stepEnergy
111  << ", prePos: " << prePos
112  << ", postPos: " << postPos
113  << ", stepLength: " << stepLength / mm << " mm"
114  << ", preToF: " << preToF / ns << " ns"
115  << ", posToF: " << postToF / ns << " ns" << G4endl;
116  G4cout << "##################################################################"
117  << G4endl;
118 }
const ActarSimGasGeantHit & operator=(const ActarSimGasGeantHit &)
Operator =.
G4ThreeVector postPos
Position after step.
G4double particleMass
Mass of the particle.
void Print()
Prints full information about the calGeantHit.
void Draw()
Draws the Hit. A clear red point on the Hit position.
G4Allocator< ActarSimGasGeantHit > ActarSimGasGeantHitAllocator
G4int operator==(const ActarSimGasGeantHit &) const
Operator ==.
G4double particleCharge
Charge of the particle.
G4int detID
Detector ID.
G4double postToF
Time after step.
G4int parentID
ID for the parent track.
~ActarSimGasGeantHit()
Destructor.
ActarSimGasGeantHit()
Constructor.
G4int trackID
ID for each track.
G4ThreeVector prePos
Position before step.
G4double stepLength
Length of the step.
G4int particleID
Particle ID according to the GDP-coding.
G4double stepEnergy
Particle energy before step.
G4String detName
Detector name where energy is deposited.
G4double edep
Energy deposited.
G4double preToF
Time before step.