/
G4EicDircOpBoundaryProcess.cc
97 lines (81 loc) · 3.95 KB
/
G4EicDircOpBoundaryProcess.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
#include "G4EicDircOpBoundaryProcess.h"
#include <Geant4/G4Step.hh>
#include <Geant4/G4TouchableHistory.hh>
#include <Geant4/G4Track.hh>
#include <set>
G4VParticleChange* G4EicDircOpBoundaryProcess::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
{
G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
G4VParticleChange* pParticleChange = G4OpBoundaryProcess::PostStepDoIt(aTrack, aStep);
// static std::set<std::string> prevol;
// static std::set<std::string> postvol;
// std::string prevolnam = pPreStepPoint->GetPhysicalVolume()->GetName();
// if (prevol.find(prevolnam) == prevol.end())
// {
// std::cout << "PreStep in " <<prevolnam << std::endl;
// prevol.insert(prevolnam);
// }
// std::string postvolnam = pPostStepPoint->GetPhysicalVolume()->GetName();
// if (postvol.find(postvolnam) == postvol.end())
// {
// std::cout << "PostStep in " <<postvolnam << std::endl;
// postvol.insert(postvolnam);
// }
// int parentId = aTrack.GetParentID();
// std::cout<<"parentId "<<parentId <<std::endl;
// if(parentId==1) pParticleChange->ProposeTrackStatus(fStopAndKill);
/*double endofbar = 0.5*(4200+4*0.05); //1250/2.;
// ideal focusing
// if(PrtManager::Instance()->GetLens() == 10)
{
G4ThreeVector theGlobalPoint1 = pPostStepPoint->GetPosition();
G4TouchableHistory* touchable = (G4TouchableHistory*)(pPostStepPoint->GetTouchable());
G4ThreeVector lpoint = touchable->GetHistory()->GetTransform( 1 ).TransformPoint(theGlobalPoint1);
if(lpoint.getZ() < endofbar+0.0001 && lpoint.getZ() > endofbar-0.0001){
G4ThreeVector ww = pPreStepPoint->GetTouchableHandle()->GetHistory()->
GetTopTransform().Inverse().TransformPoint(G4ThreeVector(0,0,endofbar));
if(aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName()!="wGlue")
pParticleChange->ProposeTrackStatus(fStopAndKill);
else
aParticleChange.ProposePosition(ww.getX(), ww.getY(),lpoint.getZ()-0.0005);
return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
}
}
{
pParticleChange->ProposeTrackStatus(fStopAndKill);
if(pPreStepPoint->GetPosition().z() < endofbar) pParticleChange->ProposeTrackStatus(fStopAndKill);
}
if(aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName()=="wExpVol" && pPostStepPoint->GetPosition().z()<pPreStepPoint->GetPosition().z()){
pParticleChange->ProposeTrackStatus(fStopAndKill);
}*/
if (aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName() == "wLens3" && pPostStepPoint->GetPosition().z() > pPreStepPoint->GetPosition().z())
{
pParticleChange->ProposeTrackStatus(fStopAndKill);
}
// kill photons outside bar and prizm
if (GetStatus() == FresnelRefraction && aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName() == "wDirc")
{
pParticleChange->ProposeTrackStatus(fStopAndKill);
}
if ((aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName() == "wLens1" || aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName() == "wLens2") && aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName() == "wDirc")
{
pParticleChange->ProposeTrackStatus(fStopAndKill);
}
// // black edge of the lens3
// if((aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName()=="wLens3"
// && aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName()=="wDirc")
// || (aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName()=="wLens3"
// && aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName()=="wLens3")){
// pParticleChange->ProposeTrackStatus(fStopAndKill);
// }
if (aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName() == "wLens1" && aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName() == "wLens1")
{
pParticleChange->ProposeTrackStatus(fStopAndKill);
}
if (aStep.GetPreStepPoint()->GetPhysicalVolume()->GetName() == "wLens2" && aStep.GetPostStepPoint()->GetPhysicalVolume()->GetName() == "wLens2")
{
pParticleChange->ProposeTrackStatus(fStopAndKill);
}
return pParticleChange;
}