Geant4-11
G4AdjointPrimaryGenerator.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// G4AdjointPrimaryGenerator class implementation
27//
28// Author: L. Desorgher, SpaceIT GmbH - November 2009
29// Contract: ESA contract 21435/08/NL/AT
30// Customer: ESA/ESTEC
31// --------------------------------------------------------------------
32
35#include "G4Event.hh"
39#include "G4Navigator.hh"
41#include "G4VPhysicalVolume.hh"
42#include "G4Material.hh"
43#include "Randomize.hh"
44
45// --------------------------------------------------------------------
46//
48{
50 type_of_adjoint_source="Spherical";
52
57
59}
60
61// --------------------------------------------------------------------
62//
64{
66}
67
68// --------------------------------------------------------------------
69//
72 G4double E1, G4double E2)
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}
100
101// --------------------------------------------------------------------
102//
105 G4double E1, G4double E2)
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}
130
131// --------------------------------------------------------------------
132//
135{
137 center_spherical_source = center_pos;
138 type_of_adjoint_source = "Spherical";
146}
147
148// --------------------------------------------------------------------
149//
152{
154 type_of_adjoint_source ="ExternalSurfaceOfAVolume";
157}
158
159// --------------------------------------------------------------------
160//
163 G4ThreeVector direction,
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}
196
197// --------------------------------------------------------------------
198//
201{
202 G4double rand = G4UniformRand();
204 weight_corr=1.;
205 return distance;
206}
static const G4double pos
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double halfpi
Definition: G4SIunits.hh:57
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void DefinePhysicalVolume1(const G4String &aName)
static G4AdjointPosOnPhysVolGenerator * GetInstance()
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)
G4SingleParticleSource * theSingleParticleSource
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &v_name)
void GenerateAdjointPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
G4PhysicsFreeVector * theAccumulatedDepthVector
G4double SampleDistanceAlongBackRayAndComputeWeightCorrection(G4double &weight_corr)
void GenerateFwdPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
void ComputeAccumulatedDepthVectorAlongBackRay(G4ThreeVector glob_pos, G4ThreeVector direction, G4double ekin, G4ParticleDefinition *aPDef)
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)
G4double FindLinearEnergy(const G4double rand) const
void SetAngDistType(const G4String &)
void SetParticleMomentumDirection(const G4ParticleMomentum &aMomDirection)
void SetEnergyDisType(const G4String &)
void SetPosDisShape(const G4String &)
void SetCentreCoords(const G4ThreeVector &)
void SetPosDisType(const G4String &)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4SPSAngDistribution * GetAngDist() const
void GeneratePrimaryVertex(G4Event *evt)
G4SPSEneDistribution * GetEneDist() const
G4SPSPosDistribution * GetPosDist() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const