Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
B5::DriftChamberHit Class Reference

Drift chamber hit. More...

#include <Doxymodules_basic.h>

Inheritance diagram for B5::DriftChamberHit:
G4VHit

Public Member Functions

 DriftChamberHit ()=default
 
 DriftChamberHit (G4int layerID)
 
 DriftChamberHit (const DriftChamberHit &right)=default
 
 ~DriftChamberHit () override=default
 
DriftChamberHitoperator= (const DriftChamberHit &right)=default
 
G4bool operator== (const DriftChamberHit &right) const
 
void * operator new (size_t)
 
void operator delete (void *aHit)
 
void Draw () override
 
const std::map< G4String, G4AttDef > * GetAttDefs () const override
 
std::vector< G4AttValue > * CreateAttValues () const override
 
void Print () override
 
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
 

Private Attributes

G4int fLayerID = -1
 
G4double fTime = 0.
 
G4ThreeVector fLocalPos
 
G4ThreeVector fWorldPos
 

Detailed Description

Drift chamber hit.

It records:

Definition at line 187 of file Doxymodules_basic.h.

Constructor & Destructor Documentation

◆ DriftChamberHit() [1/3]

B5::DriftChamberHit::DriftChamberHit ( )
default

◆ DriftChamberHit() [2/3]

B5::DriftChamberHit::DriftChamberHit ( G4int  layerID)

Definition at line 53 of file DriftChamberHit.cc.

54: fLayerID(layerID)
55{}

◆ DriftChamberHit() [3/3]

B5::DriftChamberHit::DriftChamberHit ( const DriftChamberHit right)
default

◆ ~DriftChamberHit()

B5::DriftChamberHit::~DriftChamberHit ( )
overridedefault

Member Function Documentation

◆ operator=()

DriftChamberHit & B5::DriftChamberHit::operator= ( const DriftChamberHit right)
default

◆ operator==()

G4bool B5::DriftChamberHit::operator== ( const DriftChamberHit right) const

Definition at line 59 of file DriftChamberHit.cc.

60{
61 return false;
62}

◆ operator new()

void * B5::DriftChamberHit::operator new ( size_t  )
inline

Definition at line 96 of file DriftChamberHit.hh.

97{
100 }
101 return (void*)DriftChamberHitAllocator->MallocSingle();
102}
G4ThreadLocal G4Allocator< DriftChamberHit > * DriftChamberHitAllocator

◆ operator delete()

void B5::DriftChamberHit::operator delete ( void *  aHit)
inline

Definition at line 104 of file DriftChamberHit.hh.

105{
106 DriftChamberHitAllocator->FreeSingle((DriftChamberHit*) aHit);
107}
DriftChamberHit()=default

◆ Draw()

void B5::DriftChamberHit::Draw ( )
override

Definition at line 66 of file DriftChamberHit.cc.

67{
68 auto visManager = G4VVisManager::GetConcreteInstance();
69 if (! visManager) return;
70
71 G4Circle circle(fWorldPos);
72 circle.SetScreenSize(2);
73 circle.SetFillStyle(G4Circle::filled);
74 G4VisAttributes attribs(G4Colour::Yellow());
75 circle.SetVisAttributes(attribs);
76 visManager->Draw(circle);
77}
G4ThreeVector fWorldPos

◆ GetAttDefs()

const std::map< G4String, G4AttDef > * B5::DriftChamberHit::GetAttDefs ( ) const
override

Definition at line 81 of file DriftChamberHit.cc.

82{
83 G4bool isNew;
84 auto store = G4AttDefStore::GetInstance("DriftChamberHit",isNew);
85
86 if (isNew) {
87 (*store)["HitType"]
88 = G4AttDef("HitType","Hit Type","Physics","","G4String");
89
90 (*store)["ID"]
91 = G4AttDef("ID","ID","Physics","","G4int");
92
93 (*store)["Time"]
94 = G4AttDef("Time","Time","Physics","G4BestUnit","G4double");
95
96 (*store)["Pos"]
97 = G4AttDef("Pos", "Position", "Physics","G4BestUnit","G4ThreeVector");
98 }
99
100 return store;
101}

◆ CreateAttValues()

std::vector< G4AttValue > * B5::DriftChamberHit::CreateAttValues ( ) const
override

Definition at line 105 of file DriftChamberHit.cc.

106{
107 auto values = new std::vector<G4AttValue>;
108
109 values
110 ->push_back(G4AttValue("HitType","DriftChamberHit",""));
111 values
112 ->push_back(G4AttValue("ID",G4UIcommand::ConvertToString(fLayerID),""));
113 values
114 ->push_back(G4AttValue("Time",G4BestUnit(fTime,"Time"),""));
115 values
116 ->push_back(G4AttValue("Pos",G4BestUnit(fWorldPos,"Length"),""));
117
118 return values;
119}

◆ Print()

void B5::DriftChamberHit::Print ( )
override

Definition at line 123 of file DriftChamberHit.cc.

124{
125 G4cout << " Layer[" << fLayerID << "] : time " << fTime/ns
126 << " (nsec) --- local (x,y) " << fLocalPos.x()
127 << ", " << fLocalPos.y() << G4endl;
128}
G4ThreeVector fLocalPos

◆ SetLayerID()

void B5::DriftChamberHit::SetLayerID ( G4int  z)
inline

Definition at line 73 of file DriftChamberHit.hh.

73{ fLayerID = z; }

◆ GetLayerID()

G4int B5::DriftChamberHit::GetLayerID ( ) const
inline

Definition at line 74 of file DriftChamberHit.hh.

74{ return fLayerID; }

◆ SetTime()

void B5::DriftChamberHit::SetTime ( G4double  t)
inline

Definition at line 76 of file DriftChamberHit.hh.

76{ fTime = t; }

◆ GetTime()

G4double B5::DriftChamberHit::GetTime ( ) const
inline

Definition at line 77 of file DriftChamberHit.hh.

77{ return fTime; }

◆ SetLocalPos()

void B5::DriftChamberHit::SetLocalPos ( G4ThreeVector  xyz)
inline

Definition at line 79 of file DriftChamberHit.hh.

79{ fLocalPos = xyz; }

◆ GetLocalPos()

G4ThreeVector B5::DriftChamberHit::GetLocalPos ( ) const
inline

Definition at line 80 of file DriftChamberHit.hh.

80{ return fLocalPos; }

◆ SetWorldPos()

void B5::DriftChamberHit::SetWorldPos ( G4ThreeVector  xyz)
inline

Definition at line 82 of file DriftChamberHit.hh.

82{ fWorldPos = xyz; }

◆ GetWorldPos()

G4ThreeVector B5::DriftChamberHit::GetWorldPos ( ) const
inline

Definition at line 83 of file DriftChamberHit.hh.

83{ return fWorldPos; }

Member Data Documentation

◆ fLayerID

G4int B5::DriftChamberHit::fLayerID = -1
private

Definition at line 86 of file DriftChamberHit.hh.

◆ fTime

G4double B5::DriftChamberHit::fTime = 0.
private

Definition at line 87 of file DriftChamberHit.hh.

◆ fLocalPos

G4ThreeVector B5::DriftChamberHit::fLocalPos
private

Definition at line 88 of file DriftChamberHit.hh.

◆ fWorldPos

G4ThreeVector B5::DriftChamberHit::fWorldPos
private

Definition at line 89 of file DriftChamberHit.hh.


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

Applications | User Support | Publications | Collaboration