Geant4  10.02.p03
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
 

Private Attributes

B5HadCalorimeterHitsCollectionfHitsCollection
 
G4int fHCID
 

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::B5HadCalorimeterSD ( G4String  name)

Definition at line 43 of file B5HadCalorimeterSD.cc.

45 {
46  collectionName.insert("HadCalorimeterColl");
47 }
B5HadCalorimeterHitsCollection * fHitsCollection
G4VSensitiveDetector(G4String name)
G4CollectionNameVector collectionName
Here is the call graph for this function:

◆ ~B5HadCalorimeterSD()

B5HadCalorimeterSD::~B5HadCalorimeterSD ( )
virtual

Definition at line 51 of file B5HadCalorimeterSD.cc.

52 {}

Member Function Documentation

◆ Initialize()

void B5HadCalorimeterSD::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 56 of file B5HadCalorimeterSD.cc.

57 {
60  if (fHCID<0)
62  hce->AddHitsCollection(fHCID,fHitsCollection);
63 
64  // fill calorimeter hits with zero energy deposition
65  for (G4int iColumn=0;iColumn<10;iColumn++)
66  for (G4int iRow=0;iRow<2;iRow++)
67  {
69  fHitsCollection->insert(hit);
70  }
71 }
B5HadCalorimeterHitsCollection * fHitsCollection
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
int G4int
Definition: G4Types.hh:78
G4THitsCollection< B5HadCalorimeterHit > B5HadCalorimeterHitsCollection
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
G4CollectionNameVector collectionName
Here is the call graph for this function:

◆ ProcessHits()

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

Implements G4VSensitiveDetector.

Definition at line 75 of file B5HadCalorimeterSD.cc.

76 {
77  G4double edep = step->GetTotalEnergyDeposit();
78  if (edep==0.) return true;
79 
80  G4TouchableHistory* touchable
81  = (G4TouchableHistory*)(step->GetPreStepPoint()->GetTouchable());
82 
83  G4int rowNo = touchable->GetCopyNumber(2);
84  G4int columnNo = touchable->GetCopyNumber(3);
85 
86  G4int hitID = 2*columnNo+rowNo;
87  B5HadCalorimeterHit* hit = (*fHitsCollection)[hitID];
88 
89  // check if it is first touch
90  if (hit->GetColumnID()<0)
91  {
92  hit->SetColumnID(columnNo);
93  hit->SetRowID(rowNo);
94  G4int depth = touchable->GetHistory()->GetDepth();
95  G4AffineTransform transform
96  = touchable->GetHistory()->GetTransform(depth-2);
97  transform.Invert();
98  hit->SetRot(transform.NetRotation());
99  hit->SetPos(transform.NetTranslation());
100  }
101  // add energy deposition
102  hit->AddEdep(edep);
103 
104  return true;
105 }
const G4AffineTransform & GetTransform(G4int n) const
const G4NavigationHistory * GetHistory() const
int G4int
Definition: G4Types.hh:78
Double_t edep
void AddEdep(G4double de)
G4AffineTransform & Invert()
G4RotationMatrix NetRotation() const
G4ThreeVector NetTranslation() const
void SetPos(G4ThreeVector xyz)
void SetRot(G4RotationMatrix rmat)
G4int GetCopyNumber(G4int depth=0) const
double G4double
Definition: G4Types.hh:76
G4int GetDepth() const
Here is the call graph for this function:

Member Data Documentation

◆ fHCID

G4int B5HadCalorimeterSD::fHCID
private

Definition at line 55 of file B5HadCalorimeterSD.hh.

◆ fHitsCollection

B5HadCalorimeterHitsCollection* B5HadCalorimeterSD::fHitsCollection
private

Definition at line 54 of file B5HadCalorimeterSD.hh.


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