Geant4 10.7.0
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"
32#include "G4GeometryManager.hh"
34#include "G4VPhysicalVolume.hh"
35#include "G4Event.hh"
36#include "G4PrimaryVertex.hh"
37#include "G4PrimaryParticle.hh"
38
39#include "G4TheMTRayTracer.hh"
40
42{
43 G4ThreeVector zero;
44 particle_definition = 0;
45 particle_energy = 1.0*CLHEP::GeV;
46 particle_time = 0.0;
47 particle_polarization = zero;
48
49 pWorld = 0;
50 whereisit = kInside;
51
52 nRow = 0;
53 nColumn = 0;
54
55 eyePosition = zero;
56 eyeDirection = zero;
57 up = G4ThreeVector(0,1,0);
58 headAngle = 0.0;
59 viewSpan = 0.0;
60 stepAngle = 0.0;
61 viewSpanX = 0.0;
62 viewSpanY = 0.0;
63
64 distortionOn = false;
65}
66
68{;}
69
71{
72 // Note: We don't use G4ParticleGun here, as instantiating a G4ParticleGun
73 // object causes creation of UI commands and corresponding UI messenger
74 // that interfare with normal G4ParticleGun UI commands.
75
76 // evId = iRow * nColumn + iColumn
77 G4int evId = anEvent->GetEventID();
78 G4int iRow = evId / nColumn;
79 G4int iColumn = evId % nColumn;
80 G4double angleX = -(viewSpanX/2. - G4double(iColumn)*stepAngle);
81 G4double angleY = viewSpanY/2. - G4double(iRow)*stepAngle;
82 G4ThreeVector rayDirection;
83 if(distortionOn)
84 { rayDirection = G4ThreeVector(-std::tan(angleX)/std::cos(angleY),std::tan(angleY)/std::cos(angleX),1.0); }
85 else
86 { rayDirection = G4ThreeVector(-std::tan(angleX),std::tan(angleY),1.0); }
87 G4double cp = std::cos(eyeDirection.phi());
88 G4double sp = std::sqrt(1.-cp*cp);
89 G4double ct = std::cos(eyeDirection.theta());
90 G4double st = std::sqrt(1.-ct*ct);
91 G4double gam = std::atan2(ct*cp*up.x()+ct*sp*up.y()-st*up.z(), -sp*up.x()+cp*up.y());
92 rayDirection.rotateZ(-gam);
93 rayDirection.rotateZ(headAngle);
94 rayDirection.rotateUz(eyeDirection);
95
96 G4ThreeVector rayPosition(eyePosition);
97 if (whereisit != kInside) {
98 // Eye position is outside the world, so move it inside.
99 G4double outsideDistance = pWorld->GetLogicalVolume()->GetSolid()->
100 DistanceToIn(rayPosition,rayDirection);
101 if(outsideDistance != kInfinity)
102 { rayPosition = rayPosition + (outsideDistance+0.001)*rayDirection; }
103 else
104 {
105 // Ray does not intercept world at all.
106 // Return without primary particle.
107 return;
108 }
109 }
110
111 // create a new vertex
112 G4PrimaryVertex* vertex = new G4PrimaryVertex(rayPosition,particle_time);
113
114 // create new primaries and set them to the vertex
115 G4double mass = particle_definition->GetPDGMass();
116 G4PrimaryParticle* particle = new G4PrimaryParticle(particle_definition);
117 particle->SetKineticEnergy( particle_energy );
118 particle->SetMass( mass );
119 particle->SetMomentumDirection( rayDirection.unit() );
120 particle->SetPolarization(particle_polarization.x(),
121 particle_polarization.y(),
122 particle_polarization.z());
123 vertex->SetPrimary( particle );
124
125 anEvent->AddPrimaryVertex( vertex );
126}
127
129{
131 particle_definition = particleTable->FindParticle("geantino");
132 if(!particle_definition)
133 {
134 G4String msg;
135 msg = " G4RayTracer uses geantino to trace the ray, but your physics list does not\n";
136 msg += "define G4Geantino. Please add G4Geantino in your physics list.";
137 G4Exception("G4RTPrimaryGeneratorAction::SetUp","VisRayTracer00101",FatalException,msg);
138 }
139
140 G4TheMTRayTracer* rt = G4TheMTRayTracer::theInstance;
141 nRow = rt->nRow;
142 nColumn = rt->nColumn;
143 eyePosition = rt->eyePosition;
144 eyeDirection = rt->eyeDirection;
145 viewSpan = rt->viewSpan;
146 stepAngle = viewSpan/100.;
147 viewSpanX = stepAngle*nColumn;
148 viewSpanY = stepAngle*nRow;
149 distortionOn = rt->distortionOn;
150
152 GetNavigatorForTracking()->GetWorldVolume();
153 whereisit = pWorld->GetLogicalVolume()->GetSolid()->Inside(eyePosition);
154}
155
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
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)
Definition: ThreeVector.cc:107
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
G4int GetEventID() const
Definition: G4Event.hh:118
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:121
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