Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
XAluminumElectrodeHit Class Reference

#include <XAluminumElectrodeHit.hh>

Inheritance diagram for XAluminumElectrodeHit:
G4VHit

Public Member Functions

 XAluminumElectrodeHit ()
 
virtual ~XAluminumElectrodeHit ()
 
 XAluminumElectrodeHit (const XAluminumElectrodeHit &right)
 
const XAluminumElectrodeHitoperator= (const XAluminumElectrodeHit &right)
 
int operator== (const XAluminumElectrodeHit &right) const
 
voidoperator new (size_t)
 
void operator delete (void *aHit)
 
virtual void Draw ()
 
virtual const std::map
< G4String, G4AttDef > * 
GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
virtual void Print ()
 
void SetTime (G4double t)
 
G4double GetTime () const
 
void SetEDep (G4double e)
 
G4double GetEDep () const
 
void SetLocalPos (G4ThreeVector xyz)
 
G4ThreeVector GetLocalPos () const
 
void SetWorldPos (G4ThreeVector xyz)
 
G4ThreeVector GetWorldPos () const
 
- Public Member Functions inherited from G4VHit
 G4VHit ()
 
virtual ~G4VHit ()
 
G4int operator== (const G4VHit &right) const
 

Detailed Description

Definition at line 42 of file XAluminumElectrodeHit.hh.

Constructor & Destructor Documentation

XAluminumElectrodeHit::XAluminumElectrodeHit ( )

Definition at line 50 of file XAluminumElectrodeHit.cc.

51 {
52  fTime = 0.;
53  fEdep = 0.;
54 }
XAluminumElectrodeHit::~XAluminumElectrodeHit ( )
virtual

Definition at line 58 of file XAluminumElectrodeHit.cc.

59 {;}
XAluminumElectrodeHit::XAluminumElectrodeHit ( const XAluminumElectrodeHit right)

Definition at line 63 of file XAluminumElectrodeHit.cc.

64 : G4VHit() {
65  fTime = right.fTime;
66  fEdep = right.fEdep;
67  fWorldPos = right.fWorldPos;
68  fLocalPos = right.fLocalPos;
69 }
G4VHit()
Definition: G4VHit.cc:34

Member Function Documentation

std::vector< G4AttValue > * XAluminumElectrodeHit::CreateAttValues ( ) const
virtual

Reimplemented from G4VHit.

Definition at line 134 of file XAluminumElectrodeHit.cc.

References G4BestUnit.

135 {
136  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
137 
138  values->push_back(G4AttValue("HitType","XAluminumElectrodeHit",""));
139 
140  values->push_back
141  (G4AttValue("Time",G4BestUnit(fTime,"Time"),""));
142 
143  values->push_back
144  (G4AttValue("EDep",G4BestUnit(fEdep,"Energy"),""));
145 
146  values->push_back
147  (G4AttValue("Pos",G4BestUnit(fWorldPos,"Length"),""));
148 
149  return values;
150 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void XAluminumElectrodeHit::Draw ( )
virtual

Reimplemented from G4VHit.

Definition at line 91 of file XAluminumElectrodeHit.cc.

References G4VVisManager::Draw(), G4VMarker::filled, G4VVisManager::GetConcreteInstance(), python.hepunit::ms, G4VisAttributes::SetEndTime(), G4VMarker::SetFillStyle(), G4VMarker::SetScreenSize(), G4VisAttributes::SetStartTime(), and G4Visible::SetVisAttributes().

92 {
94  if(pVVisManager)
95  {
96  G4Circle circle(fWorldPos);
97  circle.SetScreenSize(15);
98  circle.SetFillStyle(G4Circle::filled);
99  G4Colour colour(0.65,0.65,0.);
100  G4VisAttributes attribs(colour);
101  attribs.SetStartTime(fTime);
102  attribs.SetEndTime(fTime+1*ms);
103  circle.SetVisAttributes(attribs);
104  pVVisManager->Draw(circle);
105  }
106 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
static G4VVisManager * GetConcreteInstance()
const std::map< G4String, G4AttDef > * XAluminumElectrodeHit::GetAttDefs ( ) const
virtual

Reimplemented from G4VHit.

Definition at line 110 of file XAluminumElectrodeHit.cc.

References G4AttDefStore::GetInstance().

111 {
112  G4bool isNew;
113  std::map<G4String,G4AttDef>* store
114  = G4AttDefStore::GetInstance("XAluminumElectrodeHit",isNew);
115  if (isNew) {
116  G4String HitType("HitType");
117  (*store)[HitType] = G4AttDef(HitType,"Hit Type","Physics","","G4String");
118 
119  G4String Time("Time");
120  (*store)[Time] = G4AttDef(Time,"Time","Physics","G4BestUnit","G4double");
121 
122  G4String EDep("EDep");
123  (*store)[EDep] = G4AttDef(Time,"EDep","Physics","G4BestUnit","G4double");
124 
125  G4String Pos("Pos");
126  (*store)[Pos] = G4AttDef(Pos, "Position",
127  "Physics","G4BestUnit","G4ThreeVector");
128  }
129  return store;
130 }
ush Pos
Definition: deflate.h:89
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
G4double XAluminumElectrodeHit::GetEDep ( ) const
inline

Definition at line 70 of file XAluminumElectrodeHit.hh.

Referenced by XAluminumElectrodeSensitivity::ProcessHits().

70 { return fEdep; }
G4ThreeVector XAluminumElectrodeHit::GetLocalPos ( ) const
inline

Definition at line 72 of file XAluminumElectrodeHit.hh.

72 { return fLocalPos; }
G4double XAluminumElectrodeHit::GetTime ( ) const
inline

Definition at line 68 of file XAluminumElectrodeHit.hh.

68 { return fTime; }
G4ThreeVector XAluminumElectrodeHit::GetWorldPos ( ) const
inline

Definition at line 74 of file XAluminumElectrodeHit.hh.

74 { return fWorldPos; }
void XAluminumElectrodeHit::operator delete ( void aHit)
inline

Definition at line 88 of file XAluminumElectrodeHit.hh.

References G4Allocator< Type >::FreeSingle().

89 {
91 }
G4Allocator< XAluminumElectrodeHit > XAluminumElectrodeHitAllocator
void * XAluminumElectrodeHit::operator new ( size_t  )
inline

Definition at line 81 of file XAluminumElectrodeHit.hh.

References G4Allocator< Type >::MallocSingle().

82 {
83  void* aHit;
84  aHit = (void*)XAluminumElectrodeHitAllocator.MallocSingle();
85  return aHit;
86 }
G4Allocator< XAluminumElectrodeHit > XAluminumElectrodeHitAllocator
const XAluminumElectrodeHit & XAluminumElectrodeHit::operator= ( const XAluminumElectrodeHit right)

Definition at line 73 of file XAluminumElectrodeHit.cc.

74 {
75  fTime = right.fTime;
76  fEdep = right.fEdep;
77  fWorldPos = right.fWorldPos;
78  fLocalPos = right.fLocalPos;
79  return *this;
80 }
int XAluminumElectrodeHit::operator== ( const XAluminumElectrodeHit right) const

Definition at line 84 of file XAluminumElectrodeHit.cc.

85 {
86  return 0;
87 }
void XAluminumElectrodeHit::Print ( void  )
virtual

Reimplemented from G4VHit.

Definition at line 154 of file XAluminumElectrodeHit.cc.

References python.hepunit::eV, G4cout, G4endl, and ns.

155 {
156  G4cout << " time " << fTime/ns << " (nsec) : at " << fLocalPos
157  << " -- fEdep = " << fEdep/eV << " [eV]" << G4endl;
158 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
#define ns
Definition: xmlparse.cc:597
void XAluminumElectrodeHit::SetEDep ( G4double  e)
inline

Definition at line 69 of file XAluminumElectrodeHit.hh.

Referenced by XAluminumElectrodeSensitivity::ProcessHits().

69 { fEdep = e; }
void XAluminumElectrodeHit::SetLocalPos ( G4ThreeVector  xyz)
inline

Definition at line 71 of file XAluminumElectrodeHit.hh.

Referenced by XAluminumElectrodeSensitivity::ProcessHits().

71 { fLocalPos = xyz; }
void XAluminumElectrodeHit::SetTime ( G4double  t)
inline

Definition at line 67 of file XAluminumElectrodeHit.hh.

Referenced by XAluminumElectrodeSensitivity::ProcessHits().

67 { fTime = t; }
void XAluminumElectrodeHit::SetWorldPos ( G4ThreeVector  xyz)
inline

Definition at line 73 of file XAluminumElectrodeHit.hh.

Referenced by XAluminumElectrodeSensitivity::ProcessHits().

73 { fWorldPos = xyz; }

The documentation for this class was generated from the following files: