Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
XAluminumElectrodeHit Class Reference

#include <Doxymodules_exoticphysics.h>

Inheritance diagram for XAluminumElectrodeHit:
G4VHit

Public Member Functions

 XAluminumElectrodeHit ()
 
virtual ~XAluminumElectrodeHit ()
 
 XAluminumElectrodeHit (const XAluminumElectrodeHit &right)
 
const XAluminumElectrodeHitoperator= (const XAluminumElectrodeHit &right)
 
G4bool operator== (const XAluminumElectrodeHit &right) const
 
void * operator 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 ()
 

Public Attributes

G4double fTime
 
G4double fEdep
 
G4ThreeVector fLocalPos
 
G4ThreeVector fWorldPos
 

Detailed Description

Definition at line 46 of file Doxymodules_exoticphysics.h.

Constructor & Destructor Documentation

◆ XAluminumElectrodeHit() [1/2]

XAluminumElectrodeHit::XAluminumElectrodeHit ( )

Definition at line 50 of file XAluminumElectrodeHit.cc.

51{
52 fTime = 0.;
53 fEdep = 0.;
54}

◆ ~XAluminumElectrodeHit()

XAluminumElectrodeHit::~XAluminumElectrodeHit ( )
virtual

Definition at line 58 of file XAluminumElectrodeHit.cc.

59{;}

◆ XAluminumElectrodeHit() [2/2]

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}

Member Function Documentation

◆ operator=()

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}

◆ operator==()

G4bool XAluminumElectrodeHit::operator== ( const XAluminumElectrodeHit right) const

Definition at line 84 of file XAluminumElectrodeHit.cc.

85{
86 return false;
87}

◆ operator new()

void * XAluminumElectrodeHit::operator new ( size_t  )
inline

Definition at line 71 of file XAluminumElectrodeHit.hh.

72{
73 if (!XAluminumElectrodeHitAllocator) // Singleton
75
76 return (void*)XAluminumElectrodeHitAllocator->MallocSingle();
77}
G4ThreadLocal G4Allocator< XAluminumElectrodeHit > * XAluminumElectrodeHitAllocator

◆ operator delete()

void XAluminumElectrodeHit::operator delete ( void *  aHit)
inline

◆ Draw()

void XAluminumElectrodeHit::Draw ( )
virtual

Definition at line 91 of file XAluminumElectrodeHit.cc.

92{
93 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
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}

◆ GetAttDefs()

const std::map< G4String, G4AttDef > * XAluminumElectrodeHit::GetAttDefs ( ) const
virtual

Definition at line 110 of file XAluminumElectrodeHit.cc.

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}

◆ CreateAttValues()

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

Definition at line 134 of file XAluminumElectrodeHit.cc.

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}

◆ Print()

void XAluminumElectrodeHit::Print ( )
virtual

Definition at line 154 of file XAluminumElectrodeHit.cc.

155{
156 G4cout << " time " << fTime/ns << " (nsec) : at " << fLocalPos
157 << " -- fEdep = " << fEdep/eV << " [eV]" << G4endl;
158}

Member Data Documentation

◆ fTime

G4double XAluminumElectrodeHit::fTime

Definition at line 46 of file XAluminumElectrodeHit.hh.

◆ fEdep

G4double XAluminumElectrodeHit::fEdep

Definition at line 47 of file XAluminumElectrodeHit.hh.

◆ fLocalPos

G4ThreeVector XAluminumElectrodeHit::fLocalPos

Definition at line 48 of file XAluminumElectrodeHit.hh.

◆ fWorldPos

G4ThreeVector XAluminumElectrodeHit::fWorldPos

Definition at line 49 of file XAluminumElectrodeHit.hh.


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

Applications | User Support | Publications | Collaboration