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

#include <B5DriftChamberHit.hh>

Inheritance diagram for B5DriftChamberHit:
G4VHit

Public Member Functions

 B5DriftChamberHit ()
 
 B5DriftChamberHit (G4int z)
 
 B5DriftChamberHit (const B5DriftChamberHit &right)
 
virtual ~B5DriftChamberHit ()
 
const B5DriftChamberHitoperator= (const B5DriftChamberHit &right)
 
int operator== (const B5DriftChamberHit &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 SetLayerID (G4int z)
 
G4int GetLayerID () const
 
void SetTime (G4double t)
 
G4double GetTime () 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

Drift chamber hit

It records:

Definition at line 52 of file B5DriftChamberHit.hh.

Constructor & Destructor Documentation

B5DriftChamberHit::B5DriftChamberHit ( )

Definition at line 51 of file B5DriftChamberHit.cc.

52 : G4VHit(), fLayerID(-1), fTime(0.), fLocalPos(0), fWorldPos(0)
53 {}
G4VHit()
Definition: G4VHit.cc:34
B5DriftChamberHit::B5DriftChamberHit ( G4int  z)

Definition at line 57 of file B5DriftChamberHit.cc.

58 : G4VHit(), fLayerID(z), fTime(0.), fLocalPos(0), fWorldPos(0)
59 {}
G4VHit()
Definition: G4VHit.cc:34
G4double z
Definition: TRTMaterials.hh:39
B5DriftChamberHit::B5DriftChamberHit ( const B5DriftChamberHit right)

Definition at line 68 of file B5DriftChamberHit.cc.

69 : G4VHit() {
70  fLayerID = right.fLayerID;
71  fWorldPos = right.fWorldPos;
72  fLocalPos = right.fLocalPos;
73  fTime = right.fTime;
74 }
G4VHit()
Definition: G4VHit.cc:34
B5DriftChamberHit::~B5DriftChamberHit ( )
virtual

Definition at line 63 of file B5DriftChamberHit.cc.

64 {}

Member Function Documentation

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

Reimplemented from G4VHit.

Definition at line 137 of file B5DriftChamberHit.cc.

References G4UIcommand::ConvertToString(), and G4BestUnit.

138 {
139  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
140 
141  values
142  ->push_back(G4AttValue("HitType","DriftChamberHit",""));
143  values
144  ->push_back(G4AttValue("ID",G4UIcommand::ConvertToString(fLayerID),""));
145  values
146  ->push_back(G4AttValue("Time",G4BestUnit(fTime,"Time"),""));
147  values
148  ->push_back(G4AttValue("Pos",G4BestUnit(fWorldPos,"Length"),""));
149 
150  return values;
151 }
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void B5DriftChamberHit::Draw ( )
virtual

Reimplemented from G4VHit.

Definition at line 96 of file B5DriftChamberHit.cc.

References G4VVisManager::Draw(), G4VMarker::filled, G4VVisManager::GetConcreteInstance(), G4VMarker::SetFillStyle(), G4VMarker::SetScreenSize(), and G4Visible::SetVisAttributes().

97 {
99  if(pVVisManager)
100  {
101  G4Circle circle(fWorldPos);
102  circle.SetScreenSize(2);
103  circle.SetFillStyle(G4Circle::filled);
104  G4Colour colour(1.,1.,0.);
105  G4VisAttributes attribs(colour);
106  circle.SetVisAttributes(attribs);
107  pVVisManager->Draw(circle);
108  }
109 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
static G4VVisManager * GetConcreteInstance()
const std::map< G4String, G4AttDef > * B5DriftChamberHit::GetAttDefs ( ) const
virtual

Reimplemented from G4VHit.

Definition at line 113 of file B5DriftChamberHit.cc.

References G4AttDefStore::GetInstance().

114 {
115  G4bool isNew;
116  std::map<G4String,G4AttDef>* store
117  = G4AttDefStore::GetInstance("B5DriftChamberHit",isNew);
118 
119  if (isNew) {
120  (*store)["HitType"]
121  = G4AttDef("HitType","Hit Type","Physics","","G4String");
122 
123  (*store)["ID"]
124  = G4AttDef("ID","ID","Physics","","G4int");
125 
126  (*store)["Time"]
127  = G4AttDef("Time","Time","Physics","G4BestUnit","G4double");
128 
129  (*store)["Pos"]
130  = G4AttDef("Pos", "Position", "Physics","G4BestUnit","G4ThreeVector");
131  }
132  return store;
133 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
G4int B5DriftChamberHit::GetLayerID ( ) const
inline

Definition at line 72 of file B5DriftChamberHit.hh.

Referenced by B5EventAction::EndOfEventAction().

72 { return fLayerID; }
G4ThreeVector B5DriftChamberHit::GetLocalPos ( ) const
inline

Definition at line 78 of file B5DriftChamberHit.hh.

Referenced by B5EventAction::EndOfEventAction().

78 { return fLocalPos; }
G4double B5DriftChamberHit::GetTime ( ) const
inline

Definition at line 75 of file B5DriftChamberHit.hh.

75 { return fTime; }
G4ThreeVector B5DriftChamberHit::GetWorldPos ( ) const
inline

Definition at line 81 of file B5DriftChamberHit.hh.

81 { return fWorldPos; }
void B5DriftChamberHit::operator delete ( void aHit)
inline

Definition at line 101 of file B5DriftChamberHit.hh.

References B5DriftChamberHitAllocator.

102 {
103  B5DriftChamberHitAllocator->FreeSingle((B5DriftChamberHit*) aHit);
104 }
G4ThreadLocal G4Allocator< B5DriftChamberHit > * B5DriftChamberHitAllocator
void * B5DriftChamberHit::operator new ( size_t  )
inline

Definition at line 94 of file B5DriftChamberHit.hh.

References B5DriftChamberHitAllocator.

95 {
98  return (void*)B5DriftChamberHitAllocator->MallocSingle();
99 }
G4ThreadLocal G4Allocator< B5DriftChamberHit > * B5DriftChamberHitAllocator
const B5DriftChamberHit & B5DriftChamberHit::operator= ( const B5DriftChamberHit right)

Definition at line 78 of file B5DriftChamberHit.cc.

79 {
80  fLayerID = right.fLayerID;
81  fWorldPos = right.fWorldPos;
82  fLocalPos = right.fLocalPos;
83  fTime = right.fTime;
84  return *this;
85 }
int B5DriftChamberHit::operator== ( const B5DriftChamberHit right) const

Definition at line 89 of file B5DriftChamberHit.cc.

90 {
91  return 0;
92 }
void B5DriftChamberHit::Print ( void  )
virtual

Reimplemented from G4VHit.

Definition at line 155 of file B5DriftChamberHit.cc.

References G4cout, G4endl, ns, CLHEP::Hep3Vector::x(), and CLHEP::Hep3Vector::y().

Referenced by B5EventAction::EndOfEventAction().

156 {
157  G4cout << " Layer[" << fLayerID << "] : time " << fTime/ns
158  << " (nsec) --- local (x,y) " << fLocalPos.x()
159  << ", " << fLocalPos.y() << G4endl;
160 }
double x() const
G4GLOB_DLL std::ostream G4cout
double y() const
#define G4endl
Definition: G4ios.hh:61
#define ns
Definition: xmlparse.cc:597
void B5DriftChamberHit::SetLayerID ( G4int  z)
inline

Definition at line 71 of file B5DriftChamberHit.hh.

References z.

71 { fLayerID = z; }
G4double z
Definition: TRTMaterials.hh:39
void B5DriftChamberHit::SetLocalPos ( G4ThreeVector  xyz)
inline

Definition at line 77 of file B5DriftChamberHit.hh.

Referenced by B5DriftChamberSD::ProcessHits().

77 { fLocalPos = xyz; }
void B5DriftChamberHit::SetTime ( G4double  t)
inline

Definition at line 74 of file B5DriftChamberHit.hh.

Referenced by B5DriftChamberSD::ProcessHits().

74 { fTime = t; }
void B5DriftChamberHit::SetWorldPos ( G4ThreeVector  xyz)
inline

Definition at line 80 of file B5DriftChamberHit.hh.

Referenced by B5DriftChamberSD::ProcessHits().

80 { fWorldPos = xyz; }

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