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

Sensitive detector. More...

#include <Doxymodules_parameterisations.h>

Inheritance diagram for Par04SensitiveDetector:
G4VSensitiveDetector G4VFastSimSensitiveDetector

Public Member Functions

 Par04SensitiveDetector (G4String aName)
 
 Par04SensitiveDetector (G4String aName, G4ThreeVector aNbOfCells, G4ThreeVector aNSizeOfCells)
 
virtual ~Par04SensitiveDetector ()
 
virtual void Initialize (G4HCofThisEvent *HCE) final
 Create hit collection.
 
virtual G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *aROhist) final
 Process energy deposit from the full simulation.
 
virtual G4bool ProcessHits (const G4FastHit *aHit, const G4FastTrack *aTrack, G4TouchableHistory *aROhist) final
 Process energy deposit from the fast simulation.
 
Par04HitRetrieveAndSetupHit (G4ThreeVector aPosition)
 Process energy deposit - common part for full and fast simulation It is invoked from ProcessHits() methods, and sets basic hit properties (position, etc.), common for hit from fast and full simulation.
 
virtual void EndOfEvent (G4HCofThisEvent *aHC) final
 Rewrite hits map to a vector.
 
- Public Member Functions inherited from G4VFastSimSensitiveDetector
G4bool Hit (const G4FastHit *aHit, const G4FastTrack *aTrack, G4TouchableHandle *aTouchable)
 

Private Attributes

Par04HitsCollectionfHitsCollection = nullptr
 Collection of hits.
 
std::unordered_map< G4int, std::unique_ptr< Par04Hit > > fHitsMap
 Map of hits to be used in runtime.
 
G4int fHitCollectionID = -1
 ID of collection of hits.
 
G4ThreeVector fMeshNbOfCells = { 10, 10, 10 }
 Number of mesh readout cells in cylindrical coordinates.
 
G4ThreeVector fMeshSizeOfCells = { 1 * m, 2 * CLHEP::pi / 10., 1 * m }
 Size of mesh readout cells in cylindrical coordinates.
 
G4ThreeVector fEntrancePosition = { -1, -1, -1 }
 Retrieved once per event: position of entering particle.
 
G4ThreeVector fEntranceDirection = { -1, -1, -1 }
 Retrieved once per event: direction of entering particle.
 

Detailed Description

Sensitive detector.

Describes how to store the energy deposited within the detector. It derives from two classes: G4VSensitiveDetector and G4VFastSimSensitiveDetector. Addition of G4VFastSimSensitiveDetector is necessary in order to handle the energy deposits from the fast simulation.

Two ProcessHits() methods are introduced to handle energy deposited from full (detailed) simulation, and from fast simulation. The common part is handled by RetrieveAdnSetupHit() method.

Definition at line 92 of file Doxymodules_parameterisations.h.

Constructor & Destructor Documentation

◆ Par04SensitiveDetector() [1/2]

Par04SensitiveDetector::Par04SensitiveDetector ( G4String  aName)

Definition at line 52 of file Par04SensitiveDetector.cc.

54{
55 collectionName.insert("hits");
56}

◆ Par04SensitiveDetector() [2/2]

Par04SensitiveDetector::Par04SensitiveDetector ( G4String  aName,
G4ThreeVector  aNbOfCells,
G4ThreeVector  aNSizeOfCells 
)

Definition at line 59 of file Par04SensitiveDetector.cc.

62 , fMeshNbOfCells(aNb)
63 , fMeshSizeOfCells(aSize)
64{
65 collectionName.insert("hits");
66}
G4ThreeVector fMeshSizeOfCells
Size of mesh readout cells in cylindrical coordinates.
G4ThreeVector fMeshNbOfCells
Number of mesh readout cells in cylindrical coordinates.

◆ ~Par04SensitiveDetector()

Par04SensitiveDetector::~Par04SensitiveDetector ( )
virtualdefault

Member Function Documentation

◆ Initialize()

void Par04SensitiveDetector::Initialize ( G4HCofThisEvent HCE)
finalvirtual

Create hit collection.

Definition at line 74 of file Par04SensitiveDetector.cc.

75{
76 fHitsCollection = new Par04HitsCollection(SensitiveDetectorName, collectionName[0]);
77 if(fHitCollectionID < 0)
78 {
79 fHitCollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(fHitsCollection);
80 }
81 aHCE->AddHitsCollection(fHitCollectionID, fHitsCollection);
82
83 // reset entrance position
84 fEntrancePosition.set(-1, -1, -1);
85 fEntranceDirection.set(-1, -1, -1);
86}
G4THitsCollection< Par04Hit > Par04HitsCollection
Definition Par04Hit.hh:145
G4ThreeVector fEntrancePosition
Retrieved once per event: position of entering particle.
G4int fHitCollectionID
ID of collection of hits.
G4ThreeVector fEntranceDirection
Retrieved once per event: direction of entering particle.
Par04HitsCollection * fHitsCollection
Collection of hits.

◆ ProcessHits() [1/2]

G4bool Par04SensitiveDetector::ProcessHits ( G4Step aStep,
G4TouchableHistory aROhist 
)
finalvirtual

Process energy deposit from the full simulation.

Definition at line 90 of file Par04SensitiveDetector.cc.

91{
92 G4double edep = aStep->GetTotalEnergyDeposit();
93 if(edep == 0.)
94 return true;
95
96 auto hit = RetrieveAndSetupHit(aStep->GetPostStepPoint()->GetPosition());
97 if(hit == nullptr)
98 return true;
99
100 // Add energy deposit from G4Step
101 hit->AddEdep(edep);
102 // Increment the counter
103 hit->AddNdep(1);
104
105 // Fill time information from G4Step
106 // If it's already filled, choose hit with earliest global time
107 if(hit->GetTime() == -1 || hit->GetTime() > aStep->GetTrack()->GetGlobalTime())
108 hit->SetTime(aStep->GetTrack()->GetGlobalTime());
109
110 // Set hit type to full simulation (only if hit is not already marked as fast
111 // sim)
112 if(hit->GetType() != 1)
113 hit->SetType(0);
114
115 return true;
116}
Par04Hit * RetrieveAndSetupHit(G4ThreeVector aPosition)
Process energy deposit - common part for full and fast simulation It is invoked from ProcessHits() me...

◆ ProcessHits() [2/2]

G4bool Par04SensitiveDetector::ProcessHits ( const G4FastHit aHit,
const G4FastTrack aTrack,
G4TouchableHistory aROhist 
)
finalvirtual

Process energy deposit from the fast simulation.

Implements G4VFastSimSensitiveDetector.

Definition at line 120 of file Par04SensitiveDetector.cc.

122{
123 G4double edep = aHit->GetEnergy();
124 if(edep == 0.)
125 return true;
126
127 auto hit = RetrieveAndSetupHit(aHit->GetPosition());
128 if(hit == nullptr)
129 return true;
130
131 // Add energy deposit from G4FastHit
132 hit->AddEdep(edep);
133 // Increment the counter
134 hit->AddNdep(1);
135
136 // Fill time information from G4FastTrack
137 // If it's already filled, choose hit with earliest global time
138 if(hit->GetTime() == -1 || hit->GetTime() > aTrack->GetPrimaryTrack()->GetGlobalTime())
139 {
140 hit->SetTime(aTrack->GetPrimaryTrack()->GetGlobalTime());
141 }
142
143 // Set hit type to fast simulation (even if hit was already marked as full
144 // sim, overwrite it)
145 hit->SetType(1);
146
147 return true;
148}
G4double GetEnergy() const
G4ThreeVector GetPosition() const

◆ RetrieveAndSetupHit()

Par04Hit * Par04SensitiveDetector::RetrieveAndSetupHit ( G4ThreeVector  aPosition)

Process energy deposit - common part for full and fast simulation It is invoked from ProcessHits() methods, and sets basic hit properties (position, etc.), common for hit from fast and full simulation.

Definition at line 152 of file Par04SensitiveDetector.cc.

153{
154 if(fEntrancePosition.x() == -1)
155 {
156 auto info = dynamic_cast<Par04EventInformation*>(
157 G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetUserInformation());
158 if(info == nullptr)
159 return nullptr;
161 fEntranceDirection = info->GetDirection();
162 }
163
164 auto delta = aGlobalPosition - fEntrancePosition;
165
166 // Calculate rotation matrix along the particle momentum direction
167 // It will rotate the shower axes to match the incoming particle direction
168 G4RotationMatrix rotMatrix = G4RotationMatrix();
169 double particleTheta = fEntranceDirection.theta();
170 double particlePhi = fEntranceDirection.phi();
171 rotMatrix.rotateZ(-particlePhi);
172 rotMatrix.rotateY(-particleTheta);
173 G4RotationMatrix rotMatrixInv = CLHEP::inverseOf(rotMatrix);
174
175 delta = rotMatrix * delta;
176
177 G4int rhoNo = std::floor(delta.perp() / fMeshSizeOfCells.x());
178 G4int phiNo = std::floor((CLHEP::pi + delta.phi()) / fMeshSizeOfCells.y());
179 G4int zNo = std::floor(delta.z() / fMeshSizeOfCells.z());
180
181 std::size_t hitID =
182 fMeshNbOfCells.x() * fMeshNbOfCells.z() * phiNo + fMeshNbOfCells.z() * rhoNo + zNo;
183
184 if(zNo >= fMeshNbOfCells.z() || rhoNo >= fMeshNbOfCells.x() || zNo < 0)
185 {
186 return nullptr;
187 }
188
189 auto hit = fHitsMap[hitID].get();
190 if (hit==nullptr) {
191 fHitsMap[hitID] = std::unique_ptr<Par04Hit>(new Par04Hit());
192 hit = fHitsMap[hitID].get();
193 hit->SetPhiId(phiNo);
194 hit->SetRhoId(rhoNo);
195 hit->SetZid(zNo);
196 hit->SetRot(rotMatrixInv);
197 hit->SetPos(fEntrancePosition +
198 rotMatrixInv * G4ThreeVector(0, 0, (zNo + 0.5) * fMeshSizeOfCells.z()));
199
200 }
201 return hit;
202}
G4ThreeVector GetPosition() const
Get particle position.
Hit class to store energy deposited in the sensitive detector.
std::unordered_map< G4int, std::unique_ptr< Par04Hit > > fHitsMap
Map of hits to be used in runtime.

◆ EndOfEvent()

void Par04SensitiveDetector::EndOfEvent ( G4HCofThisEvent aHC)
finalvirtual

Rewrite hits map to a vector.

Definition at line 206 of file Par04SensitiveDetector.cc.

207{
208 for(const auto& hits: fHitsMap){
209 fHitsCollection->insert(new Par04Hit(*hits.second.get()));
210 }
211 fHitsMap.clear();
212}

Member Data Documentation

◆ fHitsCollection

Par04HitsCollection* Par04SensitiveDetector::fHitsCollection = nullptr
private

Collection of hits.

Definition at line 82 of file Par04SensitiveDetector.hh.

◆ fHitsMap

std::unordered_map<G4int, std::unique_ptr<Par04Hit> > Par04SensitiveDetector::fHitsMap
private

Map of hits to be used in runtime.

Definition at line 84 of file Par04SensitiveDetector.hh.

◆ fHitCollectionID

G4int Par04SensitiveDetector::fHitCollectionID = -1
private

ID of collection of hits.

Definition at line 86 of file Par04SensitiveDetector.hh.

◆ fMeshNbOfCells

G4ThreeVector Par04SensitiveDetector::fMeshNbOfCells = { 10, 10, 10 }
private

Number of mesh readout cells in cylindrical coordinates.

Definition at line 88 of file Par04SensitiveDetector.hh.

88{ 10, 10, 10 };

◆ fMeshSizeOfCells

G4ThreeVector Par04SensitiveDetector::fMeshSizeOfCells = { 1 * m, 2 * CLHEP::pi / 10., 1 * m }
private

Size of mesh readout cells in cylindrical coordinates.

Definition at line 90 of file Par04SensitiveDetector.hh.

90{ 1 * m, 2 * CLHEP::pi / 10., 1 * m };

◆ fEntrancePosition

G4ThreeVector Par04SensitiveDetector::fEntrancePosition = { -1, -1, -1 }
private

Retrieved once per event: position of entering particle.

Definition at line 92 of file Par04SensitiveDetector.hh.

92{ -1, -1, -1 };

◆ fEntranceDirection

G4ThreeVector Par04SensitiveDetector::fEntranceDirection = { -1, -1, -1 }
private

Retrieved once per event: direction of entering particle.

Definition at line 94 of file Par04SensitiveDetector.hh.

94{ -1, -1, -1 };

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

Applications | User Support | Publications | Collaboration