Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RTPrimaryGeneratorAction.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28
31#include "G4ParticleTable.hh"
33#include "G4Event.hh"
34#include "G4PrimaryVertex.hh"
35#include "G4PrimaryParticle.hh"
36
37#include "G4TheMTRayTracer.hh"
38
40{
41 G4ThreeVector zero;
42 particle_definition = 0;
43 particle_energy = 1.0*CLHEP::GeV;
44 particle_time = 0.0;
45 particle_polarization = zero;
46
47 pWorld = 0;
48 whereisit = kInside;
49
50 nRow = 0;
51 nColumn = 0;
52
53 eyePosition = zero;
54 eyeDirection = zero;
55 up = G4ThreeVector(0,1,0);
56 headAngle = 0.0;
57 viewSpan = 0.0;
58 stepAngle = 0.0;
59 viewSpanX = 0.0;
60 viewSpanY = 0.0;
61
62 distortionOn = false;
63}
64
67
69{
70 // Note: We don't use G4ParticleGun here, as instantiating a G4ParticleGun
71 // object causes creation of UI commands and corresponding UI messenger
72 // that interfare with normal G4ParticleGun UI commands.
73
74 // evId = iRow * nColumn + iColumn
75 G4int evId = anEvent->GetEventID();
76 G4int iRow = evId / nColumn;
77 G4int iColumn = evId % nColumn;
78 G4double angleX = -(viewSpanX/2. - G4double(iColumn)*stepAngle);
79 G4double angleY = viewSpanY/2. - G4double(iRow)*stepAngle;
80 G4ThreeVector rayDirection;
81 if(distortionOn)
82 { rayDirection = G4ThreeVector(-std::tan(angleX)/std::cos(angleY),std::tan(angleY)/std::cos(angleX),1.0); }
83 else
84 { rayDirection = G4ThreeVector(-std::tan(angleX),std::tan(angleY),1.0); }
85 G4double cp = std::cos(eyeDirection.phi());
86 G4double sp = std::sqrt(1.-cp*cp);
87 G4double ct = std::cos(eyeDirection.theta());
88 G4double st = std::sqrt(1.-ct*ct);
89 G4double gam = std::atan2(ct*cp*up.x()+ct*sp*up.y()-st*up.z(), -sp*up.x()+cp*up.y());
90 rayDirection.rotateZ(-gam);
91 rayDirection.rotateZ(headAngle);
92 rayDirection.rotateUz(eyeDirection);
93
94 G4ThreeVector rayPosition(eyePosition);
95 if (whereisit != kInside) {
96 // Eye position is outside the world, so move it inside.
97 G4double outsideDistance = pWorld->GetLogicalVolume()->GetSolid()->
98 DistanceToIn(rayPosition,rayDirection);
99 if(outsideDistance != kInfinity)
100 { rayPosition = rayPosition + (outsideDistance+0.001)*rayDirection; }
101 else
102 {
103 // Ray does not intercept world at all.
104 // Return without primary particle.
105 return;
106 }
107 }
108
109 // create a new vertex
110 G4PrimaryVertex* vertex = new G4PrimaryVertex(rayPosition,particle_time);
111
112 // create new primaries and set them to the vertex
113 G4double mass = particle_definition->GetPDGMass();
114 G4PrimaryParticle* particle = new G4PrimaryParticle(particle_definition);
115 particle->SetKineticEnergy( particle_energy );
116 particle->SetMass( mass );
117 particle->SetMomentumDirection( rayDirection.unit() );
118 particle->SetPolarization(particle_polarization.x(),
119 particle_polarization.y(),
120 particle_polarization.z());
121 vertex->SetPrimary( particle );
122
123 anEvent->AddPrimaryVertex( vertex );
124}
125
127{
129 particle_definition = particleTable->FindParticle("geantino");
130 if(!particle_definition)
131 {
132 G4String msg;
133 msg = " G4RayTracer uses geantino to trace the ray, but your physics list does not\n";
134 msg += "define G4Geantino. Please add G4Geantino in your physics list.";
135 G4Exception("G4RTPrimaryGeneratorAction::SetUp","VisRayTracer00101",FatalException,msg);
136 }
137
138 G4TheMTRayTracer* rt = G4TheMTRayTracer::theInstance;
139 nRow = rt->nRow;
140 nColumn = rt->nColumn;
141 eyePosition = rt->eyePosition;
142 eyeDirection = rt->eyeDirection;
143 viewSpan = rt->viewSpan;
144 stepAngle = viewSpan/100.;
145 viewSpanX = stepAngle*nColumn;
146 viewSpanY = stepAngle*nRow;
147 distortionOn = rt->distortionOn;
148
150 GetNavigatorForTracking()->GetWorldVolume();
151 whereisit = pWorld->GetLogicalVolume()->GetSolid()->Inside(eyePosition);
152}
153
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
double z() const
Hep3Vector unit() const
double phi() const
double theta() const
double x() const
Hep3Vector & rotateZ(double)
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
G4int GetEventID() const
Definition G4Event.hh:123
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition G4Event.hh:126
G4VSolid * GetSolid() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
void SetPolarization(const G4ThreeVector &pol)
void SetKineticEnergy(G4double eKin)
void SetMomentumDirection(const G4ThreeVector &p)
void SetMass(G4double mas)
void SetPrimary(G4PrimaryParticle *pp)
virtual void GeneratePrimaries(G4Event *anEvent)
G4ThreeVector eyeDirection
G4ThreeVector eyePosition
static G4TransportationManager * GetTransportationManager()
G4LogicalVolume * GetLogicalVolume() const
virtual EInside Inside(const G4ThreeVector &p) const =0
@ kInside
Definition geomdefs.hh:70