Geant4-11
Public Member Functions | Private Attributes
G4AdjointPrimaryGenerator Class Reference

#include <G4AdjointPrimaryGenerator.hh>

Public Member Functions

void ComputeAccumulatedDepthVectorAlongBackRay (G4ThreeVector glob_pos, G4ThreeVector direction, G4double ekin, G4ParticleDefinition *aPDef)
 
 G4AdjointPrimaryGenerator ()
 
 G4AdjointPrimaryGenerator (const G4AdjointPrimaryGenerator &)=delete
 
void GenerateAdjointPrimaryVertex (G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
 
void GenerateFwdPrimaryVertex (G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
 
G4AdjointPrimaryGeneratoroperator= (const G4AdjointPrimaryGenerator &)=delete
 
G4double SampleDistanceAlongBackRayAndComputeWeightCorrection (G4double &weight_corr)
 
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume (const G4String &v_name)
 
void SetSphericalAdjointPrimarySource (G4double radius, G4ThreeVector pos)
 
 ~G4AdjointPrimaryGenerator ()
 

Private Attributes

G4ThreeVector center_spherical_source
 
G4NavigatorfLinearNavigator = nullptr
 
G4double radius_spherical_source = 0.0
 
G4PhysicsFreeVectortheAccumulatedDepthVector = nullptr
 
G4AdjointPosOnPhysVolGeneratortheG4AdjointPosOnPhysVolGenerator = nullptr
 
G4SingleParticleSourcetheSingleParticleSource = nullptr
 
G4String type_of_adjoint_source
 

Detailed Description

Definition at line 58 of file G4AdjointPrimaryGenerator.hh.

Constructor & Destructor Documentation

◆ G4AdjointPrimaryGenerator() [1/2]

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( )

Definition at line 47 of file G4AdjointPrimaryGenerator.cc.

48{
50 type_of_adjoint_source="Spherical";
52
57
59}
CLHEP::Hep3Vector G4ThreeVector
static G4AdjointPosOnPhysVolGenerator * GetInstance()
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
G4SingleParticleSource * theSingleParticleSource
void SetAngDistType(const G4String &)
void SetEnergyDisType(const G4String &)
void SetPosDisType(const G4String &)
G4SPSAngDistribution * GetAngDist() const
G4SPSEneDistribution * GetEneDist() const
G4SPSPosDistribution * GetPosDist() const

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

◆ ~G4AdjointPrimaryGenerator()

G4AdjointPrimaryGenerator::~G4AdjointPrimaryGenerator ( )

Definition at line 63 of file G4AdjointPrimaryGenerator.cc.

64{
66}

References theSingleParticleSource.

◆ G4AdjointPrimaryGenerator() [2/2]

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( const G4AdjointPrimaryGenerator )
delete

Member Function Documentation

◆ ComputeAccumulatedDepthVectorAlongBackRay()

void G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay ( G4ThreeVector  glob_pos,
G4ThreeVector  direction,
G4double  ekin,
G4ParticleDefinition aPDef 
)

Definition at line 161 of file G4AdjointPrimaryGenerator.cc.

165{
166 if (fLinearNavigator == nullptr)
167 {
170 }
171 G4ThreeVector position = glob_pos;
172 G4double safety=1.;
173 G4VPhysicalVolume* thePhysVolume =
175 G4double newStep = fLinearNavigator->ComputeStep(position,direction,1.e50,
176 safety);
179
180 G4double acc_depth=0.;
181 G4double acc_length=0.;
182 theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
183
184 while (newStep > 0. && thePhysVolume != nullptr)
185 {
186 acc_length+=newStep;
187 acc_depth+=newStep*thePhysVolume->GetLogicalVolume()
188 ->GetMaterial()->GetDensity();
189 theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
190 position=position+newStep*direction;
191 thePhysVolume = fLinearNavigator
192 ->LocateGlobalPointAndSetup(position,nullptr,false);
193 newStep = fLinearNavigator->ComputeStep(position,direction,1.e50,safety);
194 }
195}
double G4double
Definition: G4Types.hh:83
G4PhysicsFreeVector * theAccumulatedDepthVector
G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:176
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:770
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:132
void InsertValues(const G4double energy, const G4double value)
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const

References G4Navigator::ComputeStep(), fLinearNavigator, G4Material::GetDensity(), G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetMaterial(), G4TransportationManager::GetNavigatorForTracking(), G4TransportationManager::GetTransportationManager(), G4PhysicsFreeVector::InsertValues(), G4Navigator::LocateGlobalPointAndSetup(), and theAccumulatedDepthVector.

◆ GenerateAdjointPrimaryVertex()

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

Definition at line 70 of file G4AdjointPrimaryGenerator.cc.

73{
74 if (type_of_adjoint_source == "ExternalSurfaceOfAVolume")
75 {
76 // Generate position and direction relative to the external surface
77 // of sensitive volume
78
79 G4double costh_to_normal=1.;
81 G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
84 costh_to_normal);
85 if (costh_to_normal <1.e-4) { costh_to_normal = 1.e-4; }
86
87 // compute now the position along the ray backward direction
88 //
90 ->SetParticleMomentumDirection(-direction);
92 }
93
96
99}
static const G4double pos
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void SetParticleMomentumDirection(const G4ParticleMomentum &aMomDirection)
void SetCentreCoords(const G4ThreeVector &)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
void GeneratePrimaryVertex(G4Event *evt)

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

◆ GenerateFwdPrimaryVertex()

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

Definition at line 103 of file G4AdjointPrimaryGenerator.cc.

106{
107 if (type_of_adjoint_source == "ExternalSurfaceOfAVolume")
108 {
109 // Generate position and direction relative to the external surface
110 // of sensitive volume
111
112 G4double costh_to_normal=1.;
114 G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
117 costh_to_normal);
118 if (costh_to_normal <1.e-4) { costh_to_normal =1.e-4; }
120 ->SetParticleMomentumDirection(direction);
122 }
123
126
129}

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

Referenced by G4AdjointPrimaryGeneratorAction::GeneratePrimaries().

◆ operator=()

G4AdjointPrimaryGenerator & G4AdjointPrimaryGenerator::operator= ( const G4AdjointPrimaryGenerator )
delete

◆ SampleDistanceAlongBackRayAndComputeWeightCorrection()

G4double G4AdjointPrimaryGenerator::SampleDistanceAlongBackRayAndComputeWeightCorrection ( G4double weight_corr)

Definition at line 199 of file G4AdjointPrimaryGenerator.cc.

201{
202 G4double rand = G4UniformRand();
204 weight_corr=1.;
205 return distance;
206}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double FindLinearEnergy(const G4double rand) const

References G4PhysicsVector::FindLinearEnergy(), G4UniformRand, and theAccumulatedDepthVector.

◆ SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume()

void G4AdjointPrimaryGenerator::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume ( const G4String v_name)

◆ SetSphericalAdjointPrimarySource()

void G4AdjointPrimaryGenerator::SetSphericalAdjointPrimarySource ( G4double  radius,
G4ThreeVector  pos 
)

Definition at line 133 of file G4AdjointPrimaryGenerator.cc.

References center_spherical_source, G4SingleParticleSource::GetAngDist(), G4SingleParticleSource::GetPosDist(), halfpi, pi, radius_spherical_source, G4SPSAngDistribution::SetAngDistType(), G4SPSPosDistribution::SetCentreCoords(), G4SPSAngDistribution::SetMaxTheta(), G4SPSAngDistribution::SetMinTheta(), G4SPSPosDistribution::SetPosDisShape(), G4SPSPosDistribution::SetPosDisType(), G4SPSPosDistribution::SetRadius(), theSingleParticleSource, and type_of_adjoint_source.

Referenced by G4AdjointPrimaryGeneratorAction::SetSphericalAdjointPrimarySource().

Field Documentation

◆ center_spherical_source

G4ThreeVector G4AdjointPrimaryGenerator::center_spherical_source
private

◆ fLinearNavigator

G4Navigator* G4AdjointPrimaryGenerator::fLinearNavigator = nullptr
private

◆ radius_spherical_source

G4double G4AdjointPrimaryGenerator::radius_spherical_source = 0.0
private

Definition at line 98 of file G4AdjointPrimaryGenerator.hh.

Referenced by SetSphericalAdjointPrimarySource().

◆ theAccumulatedDepthVector

G4PhysicsFreeVector* G4AdjointPrimaryGenerator::theAccumulatedDepthVector = nullptr
private

◆ theG4AdjointPosOnPhysVolGenerator

G4AdjointPosOnPhysVolGenerator* G4AdjointPrimaryGenerator::theG4AdjointPosOnPhysVolGenerator = nullptr
private

◆ theSingleParticleSource

G4SingleParticleSource* G4AdjointPrimaryGenerator::theSingleParticleSource = nullptr
private

◆ type_of_adjoint_source

G4String G4AdjointPrimaryGenerator::type_of_adjoint_source
private

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