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

Hadron calorimeter sensitive detector. More...

#include <B5HadCalorimeterSD.hh>

Inheritance diagram for B5HadCalorimeterSD:
Collaboration diagram for B5HadCalorimeterSD:

Public Member Functions

 B5HadCalorimeterSD (G4String name)
 
virtual ~B5HadCalorimeterSD ()
 
virtual void Initialize (G4HCofThisEvent *HCE)
 
virtual G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist)
 
- Public Member Functions inherited from G4VSensitiveDetector
 G4VSensitiveDetector (G4String name)
 
 G4VSensitiveDetector (const G4VSensitiveDetector &right)
 
virtual ~G4VSensitiveDetector ()
 
G4VSensitiveDetectoroperator= (const G4VSensitiveDetector &right)
 
G4int operator== (const G4VSensitiveDetector &right) const
 
G4int operator!= (const G4VSensitiveDetector &right) const
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 
G4bool Hit (G4Step *aStep)
 
void SetROgeometry (G4VReadOutGeometry *value)
 
void SetFilter (G4VSDFilter *value)
 
G4int GetNumberOfCollections () const
 
G4String GetCollectionName (G4int id) const
 
void SetVerboseLevel (G4int vl)
 
void Activate (G4bool activeFlag)
 
G4bool isActive () const
 
G4String GetName () const
 
G4String GetPathName () const
 
G4String GetFullPathName () const
 
G4VReadOutGeometryGetROgeometry () const
 
G4VSDFilterGetFilter () const
 
virtual G4VSensitiveDetectorClone () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VSensitiveDetector
virtual G4int GetCollectionID (G4int i)
 
- Protected Attributes inherited from G4VSensitiveDetector
G4CollectionNameVector collectionName
 
G4String SensitiveDetectorName
 
G4String thePathName
 
G4String fullPathName
 
G4int verboseLevel
 
G4bool active
 
G4VReadOutGeometryROgeometry
 
G4VSDFilterfilter
 

Detailed Description

Hadron calorimeter sensitive detector.

Definition at line 44 of file B5HadCalorimeterSD.hh.

Constructor & Destructor Documentation

B5HadCalorimeterSD::B5HadCalorimeterSD ( G4String  name)

Definition at line 44 of file B5HadCalorimeterSD.cc.

45 : G4VSensitiveDetector(name),
46  fHitsCollection(nullptr), fHCID(-1)
47 {
48  collectionName.insert("HadCalorimeterColl");
49 }
G4VSensitiveDetector(G4String name)
G4CollectionNameVector collectionName

Here is the call graph for this function:

B5HadCalorimeterSD::~B5HadCalorimeterSD ( )
virtual

Definition at line 53 of file B5HadCalorimeterSD.cc.

54 {}

Member Function Documentation

void B5HadCalorimeterSD::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 58 of file B5HadCalorimeterSD.cc.

59 {
60  fHitsCollection
62  if (fHCID<0) {
63  fHCID = G4SDManager::GetSDMpointer()->GetCollectionID(fHitsCollection);
64  }
65  hce->AddHitsCollection(fHCID,fHitsCollection);
66 
67  // fill calorimeter hits with zero energy deposition
68  for (auto column=0;column<kNofHadColumns;column++) {
69  for (auto row=0;row<kNofHadRows;row++) {
70  fHitsCollection->insert(new B5HadCalorimeterHit());
71  }
72  }
73 }
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
G4int insert(T *aHit)
constexpr G4int kNofHadRows
Definition: B5Constants.hh:43
constexpr G4int kNofHadColumns
Definition: B5Constants.hh:42
G4THitsCollection< B5HadCalorimeterHit > B5HadCalorimeterHitsCollection
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
G4CollectionNameVector collectionName

Here is the call graph for this function:

G4bool B5HadCalorimeterSD::ProcessHits ( G4Step aStep,
G4TouchableHistory ROhist 
)
virtual

Implements G4VSensitiveDetector.

Definition at line 77 of file B5HadCalorimeterSD.cc.

78 {
79  auto edep = step->GetTotalEnergyDeposit();
80  if (edep==0.) return true;
81 
82  auto touchable = step->GetPreStepPoint()->GetTouchable();
83  auto rowNo = touchable->GetCopyNumber(2);
84  auto columnNo = touchable->GetCopyNumber(3);
85  auto hitID = kNofHadRows*columnNo+rowNo;
86  auto hit = (*fHitsCollection)[hitID];
87 
88  // check if it is first touch
89  if (hit->GetColumnID()<0) {
90  hit->SetColumnID(columnNo);
91  hit->SetRowID(rowNo);
92  auto depth = touchable->GetHistory()->GetDepth();
93  auto transform = touchable->GetHistory()->GetTransform(depth-2);
94  transform.Invert();
95  hit->SetRot(transform.NetRotation());
96  hit->SetPos(transform.NetTranslation());
97  }
98  // add energy deposition
99  hit->AddEdep(edep);
100 
101  return true;
102 }
constexpr G4int kNofHadRows
Definition: B5Constants.hh:43

Here is the call graph for this function:


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