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

#include <G4AdjointPrimaryGenerator.hh>

Public Member Functions

 G4AdjointPrimaryGenerator ()
 
 ~G4AdjointPrimaryGenerator ()
 
void GenerateAdjointPrimaryVertex (G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
 
void GenerateFwdPrimaryVertex (G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
 
void SetSphericalAdjointPrimarySource (G4double radius, G4ThreeVector pos)
 
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume (const G4String &volume_name)
 

Detailed Description

Definition at line 63 of file G4AdjointPrimaryGenerator.hh.

Constructor & Destructor Documentation

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( )

Definition at line 45 of file G4AdjointPrimaryGenerator.cc.

References G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetEneDist(), G4AdjointPosOnPhysVolGenerator::GetInstance(), G4SingleParticleSource::GetPosDist(), G4SPSEneDistribution::SetAlpha(), G4SPSAngDistribution::SetAngDistType(), G4SPSEneDistribution::SetEnergyDisType(), and G4SPSPosDistribution::SetPosDisType().

46  : radius_spherical_source(0.)
47 {
48  center_spherical_source = G4ThreeVector(0.,0.,0.);
49  type_of_adjoint_source="Spherical";
50  theSingleParticleSource = new G4SingleParticleSource();
51 
52  theSingleParticleSource->GetEneDist()->SetEnergyDisType("Pow");
53  theSingleParticleSource->GetEneDist()->SetAlpha(-1.);
54  theSingleParticleSource->GetPosDist()->SetPosDisType("Point");
55  theSingleParticleSource->GetAngDist()->SetAngDistType("planar");
56 
57  theG4AdjointPosOnPhysVolGenerator = G4AdjointPosOnPhysVolGenerator::GetInstance();
58 }
G4SPSEneDistribution * GetEneDist()
CLHEP::Hep3Vector G4ThreeVector
G4SPSPosDistribution * GetPosDist()
void SetEnergyDisType(G4String)
G4SPSAngDistribution * GetAngDist()
static G4AdjointPosOnPhysVolGenerator * GetInstance()
G4AdjointPrimaryGenerator::~G4AdjointPrimaryGenerator ( )

Definition at line 61 of file G4AdjointPrimaryGenerator.cc.

62 {
63  delete theSingleParticleSource;
64 }

Member Function Documentation

void G4AdjointPrimaryGenerator::GenerateAdjointPrimaryVertex ( G4Event anEvt,
G4ParticleDefinition adj_part,
G4double  E1,
G4double  E2 
)

Definition at line 67 of file G4AdjointPrimaryGenerator.cc.

References G4AdjointPosOnPhysVolGenerator::GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(), G4SingleParticleSource::GeneratePrimaryVertex(), G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetEneDist(), G4SingleParticleSource::GetPosDist(), G4SPSPosDistribution::SetCentreCoords(), G4SPSEneDistribution::SetEmax(), G4SPSEneDistribution::SetEmin(), G4SingleParticleSource::SetParticleDefinition(), and G4SPSAngDistribution::SetParticleMomentumDirection().

68 {
69  if (type_of_adjoint_source == "ExternalSurfaceOfAVolume") {
70 
71  //Generate position and direction relative to the external surface of sensitive volume
72  //-------------------------------------------------------------
73 
74  G4double costh_to_normal=1.;
75  G4ThreeVector pos =G4ThreeVector(0.,0.,0.);
76  G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
77  theG4AdjointPosOnPhysVolGenerator->GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(pos, direction,costh_to_normal);
78  if (costh_to_normal <1.e-4) costh_to_normal =1.e-4;
79  theSingleParticleSource->GetAngDist()->SetParticleMomentumDirection(-direction);
80  theSingleParticleSource->GetPosDist()->SetCentreCoords(pos);
81  }
82 
83  theSingleParticleSource->GetEneDist()->SetEmin(E1);
84  theSingleParticleSource->GetEneDist()->SetEmax(E2);
85 
86  theSingleParticleSource->SetParticleDefinition(adj_part);
87  theSingleParticleSource->GeneratePrimaryVertex(anEvent);
88 }
G4SPSEneDistribution * GetEneDist()
CLHEP::Hep3Vector G4ThreeVector
G4SPSPosDistribution * GetPosDist()
void SetCentreCoords(G4ThreeVector)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4SPSAngDistribution * GetAngDist()
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
void GeneratePrimaryVertex(G4Event *evt)
double G4double
Definition: G4Types.hh:76
void G4AdjointPrimaryGenerator::GenerateFwdPrimaryVertex ( G4Event anEvt,
G4ParticleDefinition adj_part,
G4double  E1,
G4double  E2 
)

Definition at line 91 of file G4AdjointPrimaryGenerator.cc.

References G4AdjointPosOnPhysVolGenerator::GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(), G4SingleParticleSource::GeneratePrimaryVertex(), G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetEneDist(), G4SingleParticleSource::GetPosDist(), G4SPSPosDistribution::SetCentreCoords(), G4SPSEneDistribution::SetEmax(), G4SPSEneDistribution::SetEmin(), G4SingleParticleSource::SetParticleDefinition(), and G4SPSAngDistribution::SetParticleMomentumDirection().

Referenced by G4AdjointPrimaryGeneratorAction::GeneratePrimaries().

92 {
93  if (type_of_adjoint_source == "ExternalSurfaceOfAVolume") {
94 
95  //Generate position and direction relative to the external surface of sensitive volume
96  //-------------------------------------------------------------
97 
98  G4double costh_to_normal=1.;
99  G4ThreeVector pos =G4ThreeVector(0.,0.,0.);
100  G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
101  theG4AdjointPosOnPhysVolGenerator->GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(pos, direction,costh_to_normal);
102  if (costh_to_normal <1.e-4) costh_to_normal =1.e-4;
103  theSingleParticleSource->GetAngDist()->SetParticleMomentumDirection(direction);
104  theSingleParticleSource->GetPosDist()->SetCentreCoords(pos);
105  }
106 
107  theSingleParticleSource->GetEneDist()->SetEmin(E1);
108  theSingleParticleSource->GetEneDist()->SetEmax(E2);
109 
110  theSingleParticleSource->SetParticleDefinition(fwd_part);
111  theSingleParticleSource->GeneratePrimaryVertex(anEvent);
112 }
G4SPSEneDistribution * GetEneDist()
CLHEP::Hep3Vector G4ThreeVector
G4SPSPosDistribution * GetPosDist()
void SetCentreCoords(G4ThreeVector)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4SPSAngDistribution * GetAngDist()
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
void GeneratePrimaryVertex(G4Event *evt)
double G4double
Definition: G4Types.hh:76
void G4AdjointPrimaryGenerator::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume ( const G4String volume_name)

Definition at line 130 of file G4AdjointPrimaryGenerator.cc.

References G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume1(), G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetPosDist(), G4SPSAngDistribution::SetAngDistType(), and G4SPSPosDistribution::SetPosDisType().

Referenced by G4AdjointPrimaryGeneratorAction::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume().

131 {
132  theG4AdjointPosOnPhysVolGenerator->DefinePhysicalVolume1(volume_name);
133  type_of_adjoint_source ="ExternalSurfaceOfAVolume";
134  theSingleParticleSource->GetPosDist()->SetPosDisType("Point");
135  theSingleParticleSource->GetAngDist()->SetAngDistType("planar");
136 }
G4SPSPosDistribution * GetPosDist()
G4SPSAngDistribution * GetAngDist()
void DefinePhysicalVolume1(const G4String &aName)
void G4AdjointPrimaryGenerator::SetSphericalAdjointPrimarySource ( G4double  radius,
G4ThreeVector  pos 
)

Definition at line 115 of file G4AdjointPrimaryGenerator.cc.

References G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetPosDist(), python.hepunit::halfpi, python.hepunit::pi, G4SPSAngDistribution::SetAngDistType(), G4SPSPosDistribution::SetCentreCoords(), G4SPSAngDistribution::SetMaxTheta(), G4SPSAngDistribution::SetMinTheta(), G4SPSPosDistribution::SetPosDisShape(), G4SPSPosDistribution::SetPosDisType(), and G4SPSPosDistribution::SetRadius().

Referenced by G4AdjointPrimaryGeneratorAction::SetSphericalAdjointPrimarySource().

116 {
117  radius_spherical_source = radius;
118  center_spherical_source = center_pos;
119  type_of_adjoint_source ="Spherical";
120  theSingleParticleSource->GetPosDist()->SetPosDisType("Surface");
121  theSingleParticleSource->GetPosDist()->SetPosDisShape("Sphere");
122  theSingleParticleSource->GetPosDist()->SetCentreCoords(center_pos);
123  theSingleParticleSource->GetPosDist()->SetRadius(radius);
124  theSingleParticleSource->GetAngDist()->SetAngDistType("cos");
125  theSingleParticleSource->GetAngDist()->SetMaxTheta(pi);
126  theSingleParticleSource->GetAngDist()->SetMinTheta(halfpi);
127 }
G4SPSPosDistribution * GetPosDist()
void SetCentreCoords(G4ThreeVector)
G4SPSAngDistribution * GetAngDist()

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