CGEM BOSS 6.6.5.g
BESIII Offline Software System
Loading...
Searching...
No Matches
BesMdcHit.cc
Go to the documentation of this file.
1//---------------------------------------------------------------------------//
2// BOOST --- BESIII Object_Oriented Simulation Tool //
3//---------------------------------------------------------------------------//
4//Description: Track type hit class for BESIII MDC
5//Author: Yuan Ye([email protected])
6//Created: 4 Dec, 2003
7//Modified:
8//Comment:
9//---------------------------------------------------------------------------//
10
11#include "BesMdcHit.hh"
12#include "G4UnitsTable.hh"
13#include "G4VVisManager.hh"
14#include "G4Circle.hh"
15#include "G4Colour.hh"
16#include "G4VisAttributes.hh"
17
18G4Allocator<BesMdcHit> BesMdcHitAllocator;
19
20//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
21
23
24//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
25
27
28//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
29
31 : G4VHit()
32{
33 trackID = right.trackID;
34 layerNo = right.layerNo;
35 cellNo = right.cellNo;
36 edep = right.edep;
37 pos = right.pos;
38 driftD = right.driftD;
39 driftT = right.driftT;
40 globalT = right.globalT;
41 theta = right.theta;
42 enterAngle= right.enterAngle;
43 posFlag = right.posFlag;
44 m_p3 = right.m_p3;
45 m_flightLength = right.m_flightLength;
46 m_pdg_code = right.m_pdg_code;
47 m_isSecondary = right.m_isSecondary;
48 m_creatorProcess = right.m_creatorProcess;
49 m_digiId = right.m_digiId;
50}
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
55{
56 trackID = right.trackID;
57 layerNo = right.layerNo;
58 cellNo = right.cellNo;
59 edep = right.edep;
60 pos = right.pos;
61 driftD = right.driftD;
62 driftT = right.driftT;
63 globalT = right.globalT;
64 theta = right.theta;
65 enterAngle= right.enterAngle;
66 posFlag = right.posFlag;
67 m_p3 = right.m_p3;
68 m_flightLength = right.m_flightLength;
69 m_pdg_code = right.m_pdg_code;
70 m_isSecondary = right.m_isSecondary;
71 m_creatorProcess = right.m_creatorProcess;
72 m_digiId = right.m_digiId;
73 return *this;
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77
78int BesMdcHit::operator==(const BesMdcHit& right) const
79{
80 return (this==&right) ? 1 : 0;
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86{
87 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
88 if(pVVisManager)
89 {
90 G4Circle circle(pos);
91 circle.SetScreenSize(2.);
92 circle.SetFillStyle(G4Circle::filled);
93 G4Colour colour(1.,0.,0.);
94 G4VisAttributes attribs(colour);
95 circle.SetVisAttributes(attribs);
96 pVVisManager->Draw(circle);
97 }
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
103{
104 G4cout << " trackID: " << trackID
105 << " layerNo: " << layerNo
106 << " cellNo: " << cellNo
107 << " energy deposit: " << G4BestUnit(edep,"Energy")
108 << " position: " << G4BestUnit(pos,"Length")
109 << " driftD: " << G4BestUnit(driftD,"Length")
110 << " driftT: " << G4BestUnit(driftT,"Time")
111 << " globalT: " << G4BestUnit(globalT,"Time")
112 << " theta: " << G4BestUnit(theta,"Angle")
113 << " enterAngle: " << G4BestUnit(enterAngle,"Angle")
114 << " posFlag: " << posFlag
115 <<G4endl;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119
G4Allocator< BesMdcHit > BesMdcHitAllocator
Definition: BesMdcHit.cc:18
void Print()
Definition: BesMdcHit.cc:102
void Draw()
Definition: BesMdcHit.cc:85
~BesMdcHit()
Definition: BesMdcHit.cc:26
int operator==(const BesMdcHit &) const
Definition: BesMdcHit.cc:78
const BesMdcHit & operator=(const BesMdcHit &)
Definition: BesMdcHit.cc:54