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

#include <RE01CalorimeterSD.hh>

Inheritance diagram for RE01CalorimeterSD:
Collaboration diagram for RE01CalorimeterSD:

Public Member Functions

 RE01CalorimeterSD (G4String name)
 
virtual ~RE01CalorimeterSD ()
 
virtual void Initialize (G4HCofThisEvent *HCE)
 
- 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
 

Protected Member Functions

virtual G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist)
 
- Protected Member Functions inherited from G4VSensitiveDetector
virtual G4int GetCollectionID (G4int i)
 

Additional Inherited Members

- Protected Attributes inherited from G4VSensitiveDetector
G4CollectionNameVector collectionName
 
G4String SensitiveDetectorName
 
G4String thePathName
 
G4String fullPathName
 
G4int verboseLevel
 
G4bool active
 
G4VReadOutGeometryROgeometry
 
G4VSDFilterfilter
 

Detailed Description

Definition at line 42 of file RE01CalorimeterSD.hh.

Constructor & Destructor Documentation

RE01CalorimeterSD::RE01CalorimeterSD ( G4String  name)

Definition at line 44 of file RE01CalorimeterSD.cc.

45  :G4VSensitiveDetector(name), fCalCollection(0),
46  fNumberOfCellsInZ(20),fNumberOfCellsInPhi(48)
47 {
48  // Initialize data member.
49  for(G4int j=0;j<fNumberOfCellsInZ;j++)
50  for(G4int k=0;k<fNumberOfCellsInPhi;k++)
51  {
52  fCellID[j][k] = -1;
53  }
54  //
55  G4String HCname;
56  collectionName.insert(HCname="calCollection");
57 }
int G4int
Definition: G4Types.hh:78
G4VSensitiveDetector(G4String name)
G4CollectionNameVector collectionName

Here is the call graph for this function:

RE01CalorimeterSD::~RE01CalorimeterSD ( )
virtual

Definition at line 60 of file RE01CalorimeterSD.cc.

61 {;}

Member Function Documentation

void RE01CalorimeterSD::Initialize ( G4HCofThisEvent HCE)
virtual

Reimplemented from G4VSensitiveDetector.

Definition at line 64 of file RE01CalorimeterSD.cc.

65 {
66  fCalCollection = new RE01CalorimeterHitsCollection
68  for(G4int j=0;j<fNumberOfCellsInZ;j++)
69  for(G4int k=0;k<fNumberOfCellsInPhi;k++)
70  {
71  fCellID[j][k] = -1;
72  }
73 
74  static G4int HCID = -1;
75  if(HCID<0)
76  { HCID = GetCollectionID(0); }
77  HCE->AddHitsCollection( HCID, fCalCollection );
78 }
G4THitsCollection< RE01CalorimeterHit > RE01CalorimeterHitsCollection
int G4int
Definition: G4Types.hh:78
virtual G4int GetCollectionID(G4int i)
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4CollectionNameVector collectionName

Here is the call graph for this function:

G4bool RE01CalorimeterSD::ProcessHits ( G4Step aStep,
G4TouchableHistory ROhist 
)
protectedvirtual

Implements G4VSensitiveDetector.

Definition at line 81 of file RE01CalorimeterSD.cc.

82 {
83 //***** RE05CalorimeterSD has been migrated to Geant4 version 10 that does not
84 //***** support Readout Geometry in multi-threaded mode. Now RE05CalorimeterSD
85 //***** is assigned to a dedicaed parallel world. The pointer "aStep" points to
86 //***** a G4Step object for the parallel world.
87 
88  G4double edep = aStep->GetTotalEnergyDeposit();
89  if(edep==0.) return false;
90 
91  const G4VTouchable* ROhist = aStep->GetPreStepPoint()->GetTouchable();
92  G4int copyIDinZ = ROhist->GetReplicaNumber();
93  G4int copyIDinPhi = ROhist->GetReplicaNumber(1);
94 
95  if(fCellID[copyIDinZ][copyIDinPhi]==-1)
96  {
98  (ROhist->GetVolume()->GetLogicalVolume(),copyIDinZ,copyIDinPhi);
99  calHit->SetEdep( edep );
100  G4AffineTransform aTrans = ROhist->GetHistory()->GetTopTransform();
101  aTrans.Invert();
102  calHit->SetPos(aTrans.NetTranslation());
103  calHit->SetRot(aTrans.NetRotation());
104  calHit->SetTrackInformation(aStep->GetTrack());
105  G4int icell = fCalCollection->insert( calHit );
106  fCellID[copyIDinZ][copyIDinPhi] = icell - 1;
107  if(verboseLevel>0)
108  { G4cout << " New Calorimeter Hit on CellID "
109  << copyIDinZ << " " << copyIDinPhi << G4endl; }
110  }
111  else
112  {
113  (*fCalCollection)[fCellID[copyIDinZ][copyIDinPhi]]->AddEdep(edep);
114  (*fCalCollection)[fCellID[copyIDinZ][copyIDinPhi]]
115  ->SetTrackInformation(aStep->GetTrack());
116  if(verboseLevel>0)
117  { G4cout << " Energy added to CellID "
118  << copyIDinZ << " " << copyIDinPhi << G4endl; }
119  }
120 
121  return true;
122 }
G4ThreeVector NetTranslation() const
virtual const G4NavigationHistory * GetHistory() const
Definition: G4VTouchable.cc:86
G4int insert(T *aHit)
const G4VTouchable * GetTouchable() const
void SetTrackInformation(const G4Track *aTrack)
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
G4AffineTransform & Invert()
G4GLOB_DLL std::ostream G4cout
void SetPos(G4ThreeVector xyz)
void SetRot(G4RotationMatrix rmat)
void SetEdep(G4double de)
G4double GetTotalEnergyDeposit() const
G4RotationMatrix NetRotation() const
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
#define G4endl
Definition: G4ios.hh:61
const G4AffineTransform & GetTopTransform() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const

Here is the call graph for this function:


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