9#include <arich/simulation/SensitiveDetector.h>
10#include <arich/dataobjects/ARICHSimHit.h>
13#include <G4OpBoundaryProcess.hh>
14#include <G4ProcessManager.hh>
17#include <simulation/kernel/UserInfo.h>
21#include <framework/datastore/StoreArray.h>
22#include <framework/datastore/RelationArray.h>
23#include <Math/Vector2D.h>
36 Simulation::SensitiveDetectorBase(
"ARICH",
Const::ARICH)
43 hits.registerInDataStore();
44 particles.registerRelationTo(hits);
52 G4Track& track = *aStep->GetTrack();
53 if (track.GetDefinition()->GetParticleName() !=
"opticalphoton")
return false;
56 const G4double globalTime = track.GetGlobalTime();
58 const G4StepPoint& postStep = *aStep->GetPostStepPoint();
59 const G4ThreeVector& postworldPosition = postStep.GetPosition();
61 const G4StepPoint& preStep = *aStep->GetPreStepPoint();
62 const G4ThreeVector& preworldPosition = preStep.GetPosition();
66 if (postworldPosition.z() - preworldPosition.z() > 0) dir = +1;
69 if (dir == 1 && postStep.GetStepStatus() != fGeomBoundary) {
return false;}
70 if (dir == -1 && preStep.GetStepStatus() != fGeomBoundary) {
return false;}
72 if ((track.GetNextVolume())->GetName() ==
"ARICH.HAPDWall") {
73 track.SetTrackStatus(fStopAndKill);
78 G4OpBoundaryProcessStatus theStatus = Undefined;
80 G4ProcessManager* OpManager = G4OpticalPhoton::OpticalPhoton()->GetProcessManager();
83 G4int MAXofPostStepLoops =
84 OpManager->GetPostStepProcessVector()->entries();
85 G4ProcessVector* fPostStepDoItVector =
86 OpManager->GetPostStepProcessVector(typeDoIt);
87 for (G4int i = 0; i < MAXofPostStepLoops; i++) {
88 G4VProcess* fCurrentProcess = (*fPostStepDoItVector)[i];
89 G4OpBoundaryProcess* opProcess =
dynamic_cast<G4OpBoundaryProcess*
>(fCurrentProcess);
90 if (opProcess) { theStatus = opProcess->GetStatus();
break;}
95 if (theStatus == 3 && dir < 0)
return 0;
96 if (theStatus == 4 && dir < 0) {
if (gRandom->Uniform() < 0.5) track.SetTrackStatus(fStopAndKill);
return 0; }
102 if (info) applyQE = info->getStatus() < 2;
106 double energy = track.GetKineticEnergy() / CLHEP::eV;
109 double fraction = 1.;
110 if (info) fraction = info->getFraction();
112 if (theStatus == 3) qeffi *=
m_simPar->getQEScaling();
113 if (gRandom->Uniform() * fraction > qeffi) {
115 if (theStatus == 3 && gRandom->Uniform() <
m_simPar->getWindowAbsorbtion()) track.SetTrackStatus(fStopAndKill);
121 const G4int moduleID = dir > 0 ? postStep.GetTouchableHandle()->GetReplicaNumber(1) :
122 preStep.GetTouchableHandle()->GetReplicaNumber(2);
124 if (moduleID <= 0)
return false;
128 G4ThreeVector localPosition = dir > 0 ? postStep.GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(
130 preStep.GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(preworldPosition);
133 const G4double energy = track.GetKineticEnergy() / CLHEP::eV;
139 ROOT::Math::XYVector locpos(localPosition.x() / CLHEP::cm, localPosition.y() / CLHEP::cm);
149 track.SetTrackStatus(fStopAndKill);
Class ARICHSimHit - Geant4 simulated hit for ARICH.
This class provides a set of constants for the framework.
Low-level class to create/modify relations between StoreArrays.
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
UserInfo class which is used to attach additional information to Geant4 particles and tracks.
Accessor to arrays stored in the data store.
T * appendNew()
Construct a new T object at the end of the array.
DBObjPtr< ARICHSimulationPar > m_simPar
simulation parameters from the DB
bool step(G4Step *aStep, G4TouchableHistory *) override
Process each step and calculate variables defined in PXDSimHit.
SensitiveDetector()
Constructor.
Abstract base class for different kinds of events.