Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B5HadCalorimeterHit Class Reference

#include <B5HadCalorimeterHit.hh>

Inheritance diagram for B5HadCalorimeterHit:
Collaboration diagram for B5HadCalorimeterHit:

Public Member Functions

 B5HadCalorimeterHit ()
 
 B5HadCalorimeterHit (G4int iCol, G4int iRow)
 
 B5HadCalorimeterHit (const B5HadCalorimeterHit &right)
 
virtual ~B5HadCalorimeterHit ()
 
const B5HadCalorimeterHitoperator= (const B5HadCalorimeterHit &right)
 
int operator== (const B5HadCalorimeterHit &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 SetColumnID (G4int z)
 
G4int GetColumnID () const
 
void SetRowID (G4int z)
 
G4int GetRowID () const
 
void SetEdep (G4double de)
 
void AddEdep (G4double de)
 
G4double GetEdep () const
 
void SetPos (G4ThreeVector xyz)
 
G4ThreeVector GetPos () const
 
void SetRot (G4RotationMatrix rmat)
 
G4RotationMatrix GetRot () const
 
- Public Member Functions inherited from G4VHit
 G4VHit ()
 
virtual ~G4VHit ()
 
G4int operator== (const G4VHit &right) const
 

Detailed Description

Hadron Calorimeter hit

It records:

  • the cell column ID and row ID
  • the energy deposit
  • the cell position and rotation

Definition at line 52 of file B5HadCalorimeterHit.hh.

Constructor & Destructor Documentation

B5HadCalorimeterHit::B5HadCalorimeterHit ( )

Definition at line 52 of file B5HadCalorimeterHit.cc.

53 : G4VHit(),
54  fColumnID(-1), fRowID(-1), fEdep(0.), fPos(0)
55 {}
G4VHit()
Definition: G4VHit.cc:34
B5HadCalorimeterHit::B5HadCalorimeterHit ( G4int  iCol,
G4int  iRow 
)

Definition at line 59 of file B5HadCalorimeterHit.cc.

60 : G4VHit(),
61  fColumnID(columnID), fRowID(rowID), fEdep(0.), fPos(0)
62 {}
G4VHit()
Definition: G4VHit.cc:34
B5HadCalorimeterHit::B5HadCalorimeterHit ( const B5HadCalorimeterHit right)

Definition at line 71 of file B5HadCalorimeterHit.cc.

72 : G4VHit(),
73  fColumnID(right.fColumnID),
74  fRowID(right.fRowID),
75  fEdep(right.fEdep),
76  fPos(right.fPos),
77  fRot(right.fRot)
78 {}
G4VHit()
Definition: G4VHit.cc:34
B5HadCalorimeterHit::~B5HadCalorimeterHit ( )
virtual

Definition at line 66 of file B5HadCalorimeterHit.cc.

67 {}

Member Function Documentation

void B5HadCalorimeterHit::AddEdep ( G4double  de)
inline

Definition at line 78 of file B5HadCalorimeterHit.hh.

78 { fEdep += de; }
std::vector< G4AttValue > * B5HadCalorimeterHit::CreateAttValues ( ) const
virtual

Reimplemented from G4VHit.

Definition at line 147 of file B5HadCalorimeterHit.cc.

148 {
149  auto values = new std::vector<G4AttValue>;
150 
151  values
152  ->push_back(G4AttValue("HitType","HadCalorimeterHit",""));
153  values
154  ->push_back(G4AttValue("Column",G4UIcommand::ConvertToString(fColumnID),
155  ""));
156  values
157  ->push_back(G4AttValue("Row",G4UIcommand::ConvertToString(fRowID),""));
158  values
159  ->push_back(G4AttValue("Energy",G4BestUnit(fEdep,"Energy"),""));
160  values
161  ->push_back(G4AttValue("Pos",G4BestUnit(fPos,"Length"),""));
162 
163  return values;
164 }
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:372
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1

Here is the call graph for this function:

void B5HadCalorimeterHit::Draw ( )
virtual

Reimplemented from G4VHit.

Definition at line 102 of file B5HadCalorimeterHit.cc.

103 {
104  auto visManager = G4VVisManager::GetConcreteInstance();
105  if (! visManager || (fEdep==0.)) return;
106 
107  // Draw a calorimeter cell with depth propotional to the energy deposition
108  G4Transform3D trans(fRot.inverse(),fPos);
109  G4VisAttributes attribs;
110  G4Colour colour(1.,0.,0.);
111  attribs.SetColour(colour);
112  attribs.SetForceSolid(true);
113  G4Box box("dummy",15.*cm,15.*cm,1.*m*fEdep/(0.1*GeV));
114  visManager->Draw(box,attribs,trans);
115 }
static G4VVisManager * GetConcreteInstance()
Definition: G4Box.hh:64
HepRotation inverse() const
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double GeV
Definition: G4SIunits.hh:217

Here is the call graph for this function:

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

Reimplemented from G4VHit.

Definition at line 119 of file B5HadCalorimeterHit.cc.

120 {
121  G4bool isNew;
122  auto store = G4AttDefStore::GetInstance("B5HadCalorimeterHit",isNew);
123 
124  if (isNew) {
125  (*store)["HitType"]
126  = G4AttDef("HitType","Hit Type","Physics","","G4String");
127 
128  (*store)["Column"]
129  = G4AttDef("Column","Column ID","Physics","","G4int");
130 
131  (*store)["Row"]
132  = G4AttDef("Row","Row ID","Physics","","G4int");
133 
134  (*store)["Energy"]
135  = G4AttDef("Energy","Energy Deposited","Physics","G4BestUnit",
136  "G4double");
137 
138  (*store)["Pos"]
139  = G4AttDef("Pos", "Position", "Physics","G4BestUnit",
140  "G4ThreeVector");
141  }
142  return store;
143 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Here is the call graph for this function:

G4int B5HadCalorimeterHit::GetColumnID ( ) const
inline

Definition at line 72 of file B5HadCalorimeterHit.hh.

72 { return fColumnID; }
G4double B5HadCalorimeterHit::GetEdep ( ) const
inline

Definition at line 79 of file B5HadCalorimeterHit.hh.

79 { return fEdep; }
G4ThreeVector B5HadCalorimeterHit::GetPos ( ) const
inline

Definition at line 82 of file B5HadCalorimeterHit.hh.

82 { return fPos; }
G4RotationMatrix B5HadCalorimeterHit::GetRot ( ) const
inline

Definition at line 85 of file B5HadCalorimeterHit.hh.

85 { return fRot; }
G4int B5HadCalorimeterHit::GetRowID ( ) const
inline

Definition at line 75 of file B5HadCalorimeterHit.hh.

75 { return fRowID; }
void B5HadCalorimeterHit::operator delete ( void aHit)
inline

Definition at line 107 of file B5HadCalorimeterHit.hh.

108 {
110 }
G4ThreadLocal G4Allocator< B5HadCalorimeterHit > * B5HadCalorimeterHitAllocator
void * B5HadCalorimeterHit::operator new ( size_t  )
inline

Definition at line 99 of file B5HadCalorimeterHit.hh.

100 {
103  }
104  return (void*)B5HadCalorimeterHitAllocator->MallocSingle();
105 }
G4ThreadLocal G4Allocator< B5HadCalorimeterHit > * B5HadCalorimeterHitAllocator
const B5HadCalorimeterHit & B5HadCalorimeterHit::operator= ( const B5HadCalorimeterHit right)

Definition at line 82 of file B5HadCalorimeterHit.cc.

84 {
85  fColumnID = right.fColumnID;
86  fRowID = right.fRowID;
87  fEdep = right.fEdep;
88  fPos = right.fPos;
89  fRot = right.fRot;
90  return *this;
91 }
int B5HadCalorimeterHit::operator== ( const B5HadCalorimeterHit right) const

Definition at line 95 of file B5HadCalorimeterHit.cc.

96 {
97  return (fColumnID==right.fColumnID&&fRowID==right.fRowID);
98 }
void B5HadCalorimeterHit::Print ( void  )
virtual

Reimplemented from G4VHit.

Definition at line 168 of file B5HadCalorimeterHit.cc.

169 {
170  G4cout << " Cell[" << fRowID << ", " << fColumnID << "] "
171  << fEdep/MeV << " (MeV) " << fPos << G4endl;
172 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void B5HadCalorimeterHit::SetColumnID ( G4int  z)
inline

Definition at line 71 of file B5HadCalorimeterHit.hh.

71 { fColumnID = z; }
tuple z
Definition: test.py:28
void B5HadCalorimeterHit::SetEdep ( G4double  de)
inline

Definition at line 77 of file B5HadCalorimeterHit.hh.

77 { fEdep = de; }
void B5HadCalorimeterHit::SetPos ( G4ThreeVector  xyz)
inline

Definition at line 81 of file B5HadCalorimeterHit.hh.

81 { fPos = xyz; }
void B5HadCalorimeterHit::SetRot ( G4RotationMatrix  rmat)
inline

Definition at line 84 of file B5HadCalorimeterHit.hh.

84 { fRot = rmat; }
void B5HadCalorimeterHit::SetRowID ( G4int  z)
inline

Definition at line 74 of file B5HadCalorimeterHit.hh.

74 { fRowID = z; }
tuple z
Definition: test.py:28

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