Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 &v_name)
void ComputeAccumulatedDepthVectorAlongBackRay (G4ThreeVector glob_pos, G4ThreeVector direction, G4double ekin, G4ParticleDefinition *aPDef)
G4double SampleDistanceAlongBackRayAndComputeWeightCorrection (G4double &weight_corr)

Static Public Member Functions

static G4AdjointPrimaryGeneratorGetInstance ()

Detailed Description

Definition at line 58 of file G4AdjointPrimaryGenerator.hh.

Constructor & Destructor Documentation

◆ G4AdjointPrimaryGenerator()

G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator ( )

Definition at line 60 of file G4AdjointPrimaryGenerator.cc.

61{
62 center_spherical_source = G4ThreeVector(0.,0.,0.);
63 type_of_adjoint_source="Spherical";
64 theSingleParticleSource = new G4SingleParticleSource();
65 theSingleParticleSource->GetEneDist()->SetEnergyDisType("Pow");
66 theSingleParticleSource->GetEneDist()->SetAlpha(-1.);
67 theSingleParticleSource->GetPosDist()->SetPosDisType("Point");
68 theSingleParticleSource->GetAngDist()->SetAngDistType("planar");
69
70 theG4AdjointPosOnPhysVolGenerator = G4AdjointPosOnPhysVolGenerator::GetInstance();
71}
CLHEP::Hep3Vector G4ThreeVector
static G4AdjointPosOnPhysVolGenerator * GetInstance()

Referenced by GetInstance().

◆ ~G4AdjointPrimaryGenerator()

G4AdjointPrimaryGenerator::~G4AdjointPrimaryGenerator ( )

Definition at line 75 of file G4AdjointPrimaryGenerator.cc.

76{
77 delete theSingleParticleSource;
78}

Member Function Documentation

◆ ComputeAccumulatedDepthVectorAlongBackRay()

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

Definition at line 173 of file G4AdjointPrimaryGenerator.cc.

177{
178 if (fLinearNavigator == nullptr)
179 {
182 }
183 G4ThreeVector position = glob_pos;
184 G4double safety=1.;
185 G4VPhysicalVolume* thePhysVolume =
186 fLinearNavigator->LocateGlobalPointAndSetup(position);
187 G4double newStep = fLinearNavigator->ComputeStep(position,direction,1.e50,
188 safety);
189 delete theAccumulatedDepthVector;
190 theAccumulatedDepthVector = new G4PhysicsFreeVector();
191
192 G4double acc_depth=0.;
193 G4double acc_length=0.;
194 theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
195
196 while (newStep > 0. && thePhysVolume != nullptr)
197 {
198 acc_length+=newStep;
199 acc_depth+=newStep*thePhysVolume->GetLogicalVolume()
200 ->GetMaterial()->GetDensity();
201 theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
202 position=position+newStep*direction;
203 thePhysVolume = fLinearNavigator
204 ->LocateGlobalPointAndSetup(position,nullptr,false);
205 newStep = fLinearNavigator->ComputeStep(position,direction,1.e50,safety);
206 }
207}
double G4double
Definition G4Types.hh:83
G4Material * GetMaterial() const
G4double GetDensity() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const

◆ GenerateAdjointPrimaryVertex()

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

Definition at line 82 of file G4AdjointPrimaryGenerator.cc.

85{
86 if (type_of_adjoint_source == "ExternalSurfaceOfAVolume")
87 {
88 // Generate position and direction relative to the external surface
89 // of sensitive volume
90
91 G4double costh_to_normal=1.;
93 G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
94 theG4AdjointPosOnPhysVolGenerator
95 ->GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(pos, direction,
96 costh_to_normal);
97 if (costh_to_normal <1.e-4) { costh_to_normal = 1.e-4; }
98
99 // compute now the position along the ray backward direction
100 //
101 theSingleParticleSource->GetAngDist()
102 ->SetParticleMomentumDirection(-direction);
103 theSingleParticleSource->GetPosDist()->SetCentreCoords(pos);
104 }
105
106 theSingleParticleSource->GetEneDist()->SetEmin(E1);
107 theSingleParticleSource->GetEneDist()->SetEmax(E2);
108
109 theSingleParticleSource->SetParticleDefinition(adj_part);
110 theSingleParticleSource->GeneratePrimaryVertex(anEvent);
111}

◆ GenerateFwdPrimaryVertex()

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

Definition at line 115 of file G4AdjointPrimaryGenerator.cc.

118{
119 if (type_of_adjoint_source == "ExternalSurfaceOfAVolume")
120 {
121 // Generate position and direction relative to the external surface
122 // of sensitive volume
123
124 G4double costh_to_normal=1.;
126 G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
127 theG4AdjointPosOnPhysVolGenerator
128 ->GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(pos, direction,
129 costh_to_normal);
130 if (costh_to_normal <1.e-4) { costh_to_normal =1.e-4; }
131 theSingleParticleSource->GetAngDist()
132 ->SetParticleMomentumDirection(direction);
133 theSingleParticleSource->GetPosDist()->SetCentreCoords(pos);
134 }
135
136 theSingleParticleSource->GetEneDist()->SetEmin(E1);
137 theSingleParticleSource->GetEneDist()->SetEmax(E2);
138
139 theSingleParticleSource->SetParticleDefinition(fwd_part);
140 theSingleParticleSource->GeneratePrimaryVertex(anEvent);
141}

Referenced by G4AdjointPrimaryGeneratorAction::GeneratePrimaries().

◆ GetInstance()

◆ SampleDistanceAlongBackRayAndComputeWeightCorrection()

G4double G4AdjointPrimaryGenerator::SampleDistanceAlongBackRayAndComputeWeightCorrection ( G4double & weight_corr)

Definition at line 211 of file G4AdjointPrimaryGenerator.cc.

213{
214 G4double rand = G4UniformRand();
215 G4double distance = theAccumulatedDepthVector->FindLinearEnergy(rand);
216 weight_corr=1.;
217 return distance;
218}
#define G4UniformRand()
Definition Randomize.hh:52

◆ SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume()

void G4AdjointPrimaryGenerator::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume ( const G4String & v_name)

Definition at line 162 of file G4AdjointPrimaryGenerator.cc.

164{
165 theG4AdjointPosOnPhysVolGenerator->DefinePhysicalVolume1(volume_name);
166 type_of_adjoint_source ="ExternalSurfaceOfAVolume";
167 theSingleParticleSource->GetPosDist()->SetPosDisType("Point");
168 theSingleParticleSource->GetAngDist()->SetAngDistType("planar");
169}

Referenced by G4AdjointPrimaryGeneratorAction::SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume().

◆ SetSphericalAdjointPrimarySource()

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

Definition at line 145 of file G4AdjointPrimaryGenerator.cc.

147{
148 radius_spherical_source = radius;
149 center_spherical_source = center_pos;
150 type_of_adjoint_source = "Spherical";
151 theSingleParticleSource->GetPosDist()->SetPosDisType("Surface");
152 theSingleParticleSource->GetPosDist()->SetPosDisShape("Sphere");
153 theSingleParticleSource->GetPosDist()->SetCentreCoords(center_pos);
154 theSingleParticleSource->GetPosDist()->SetRadius(radius);
155 theSingleParticleSource->GetAngDist()->SetAngDistType("cos");
156 theSingleParticleSource->GetAngDist()->SetMaxTheta(pi);
157 theSingleParticleSource->GetAngDist()->SetMinTheta(halfpi);
158}

Referenced by G4AdjointPrimaryGeneratorAction::SetSphericalAdjointPrimarySource().


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