Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4RTPrimaryGeneratorAction Class Reference

#include <G4RTPrimaryGeneratorAction.hh>

Inheritance diagram for G4RTPrimaryGeneratorAction:
G4VUserPrimaryGeneratorAction

Public Member Functions

 G4RTPrimaryGeneratorAction ()
 
virtual ~G4RTPrimaryGeneratorAction ()
 
virtual void GeneratePrimaries (G4Event *anEvent)
 
void SetUp ()
 
- Public Member Functions inherited from G4VUserPrimaryGeneratorAction
 G4VUserPrimaryGeneratorAction ()
 
virtual ~G4VUserPrimaryGeneratorAction ()
 

Detailed Description

Definition at line 41 of file G4RTPrimaryGeneratorAction.hh.

Constructor & Destructor Documentation

G4RTPrimaryGeneratorAction::G4RTPrimaryGeneratorAction ( )

Definition at line 42 of file G4RTPrimaryGeneratorAction.cc.

References kInside.

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 }
CLHEP::Hep3Vector G4ThreeVector
G4RTPrimaryGeneratorAction::~G4RTPrimaryGeneratorAction ( )
virtual

Definition at line 68 of file G4RTPrimaryGeneratorAction.cc.

69 {;}

Member Function Documentation

void G4RTPrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Implements G4VUserPrimaryGeneratorAction.

Definition at line 71 of file G4RTPrimaryGeneratorAction.cc.

References G4Event::AddPrimaryVertex(), G4InuclParticleNames::gam, G4Event::GetEventID(), G4VPhysicalVolume::GetLogicalVolume(), G4ParticleDefinition::GetPDGMass(), G4LogicalVolume::GetSolid(), kInside, CLHEP::Hep3Vector::phi(), CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::rotateZ(), G4PrimaryParticle::SetKineticEnergy(), G4PrimaryParticle::SetMass(), G4PrimaryParticle::SetMomentumDirection(), G4PrimaryParticle::SetPolarization(), G4PrimaryVertex::SetPrimary(), G4InuclParticleNames::sp, CLHEP::Hep3Vector::theta(), CLHEP::Hep3Vector::unit(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

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 }
CLHEP::Hep3Vector G4ThreeVector
double x() const
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:143
int G4int
Definition: G4Types.hh:78
double z() const
void SetKineticEnergy(G4double eKin)
G4int GetEventID() const
Definition: G4Event.hh:140
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
Hep3Vector & rotateZ(double)
Definition: ThreeVector.cc:144
void SetMass(G4double mas)
double phi() const
double theta() const
G4LogicalVolume * GetLogicalVolume() const
G4double GetPDGMass() const
void SetMomentumDirection(const G4ThreeVector &p)
Hep3Vector unit() const
double y() const
void SetPrimary(G4PrimaryParticle *pp)
double G4double
Definition: G4Types.hh:76
void SetPolarization(const G4ThreeVector &pol)
G4VSolid * GetSolid() const
void G4RTPrimaryGeneratorAction::SetUp ( )

Definition at line 129 of file G4RTPrimaryGeneratorAction.cc.

References G4TheRayTracer::distortionOn, G4TheRayTracer::eyeDirection, G4TheRayTracer::eyePosition, FatalException, G4ParticleTable::FindParticle(), G4Exception(), G4VPhysicalVolume::GetLogicalVolume(), G4ParticleTable::GetParticleTable(), G4LogicalVolume::GetSolid(), G4TransportationManager::GetTransportationManager(), G4VSolid::Inside(), G4TheRayTracer::nColumn, G4TheRayTracer::nRow, and G4TheRayTracer::viewSpan.

Referenced by G4RTWorkerInitialization::WorkerRunStart().

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 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4ThreeVector eyeDirection
virtual EInside Inside(const G4ThreeVector &p) const =0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4LogicalVolume * GetLogicalVolume() const
static G4ParticleTable * GetParticleTable()
G4ThreeVector eyePosition
G4VSolid * GetSolid() const

The documentation for this class was generated from the following files: