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