Geant4  10.02.p02
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 // $Id: G4AdjointPrimaryGenerator.cc 81770 2014-06-05 08:31:31Z gcosmo $
27 //
29 // Class Name: G4AdjointCrossSurfChecker
30 // Author: L. Desorgher
31 // Organisation: SpaceIT GmbH
32 // Contract: ESA contract 21435/08/NL/AT
33 // Customer: ESA/ESTEC
35 
37 #include "G4PhysicalConstants.hh"
38 #include "G4Event.hh"
40 #include "G4ParticleDefinition.hh"
42 #include "G4Navigator.hh"
44 #include "G4VPhysicalVolume.hh"
45 #include "G4Material.hh"
46 #include "Randomize.hh"
47 /*
48 #include "G4AdjointCSManager.hh"
49 #include "G4MaterialCutsCouple.hh"
50 */
52 //
54  : radius_spherical_source(0.),fLinearNavigator(0),theAccumulatedDepthVector(0)
55 {
57  type_of_adjoint_source="Spherical";
59 
64 
66 
67 }
69 //
71 {
73 }
75 //
77 {
78  if (type_of_adjoint_source == "ExternalSurfaceOfAVolume") {
79 
80  //Generate position and direction relative to the external surface of sensitive volume
81  //-------------------------------------------------------------
82 
83  G4double costh_to_normal=1.;
84  G4ThreeVector pos =G4ThreeVector(0.,0.,0.);
85  G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
87  if (costh_to_normal <1.e-4) costh_to_normal =1.e-4;
88  //compute now the position along the ray backward direction
89 
92  }
93 
96 
99 
100 
101 
102 }
104 //
106 {
107  if (type_of_adjoint_source == "ExternalSurfaceOfAVolume") {
108 
109  //Generate position and direction relative to the external surface of sensitive volume
110  //-------------------------------------------------------------
111 
112  G4double costh_to_normal=1.;
113  G4ThreeVector pos =G4ThreeVector(0.,0.,0.);
114  G4ThreeVector direction = G4ThreeVector(0.,0.,1.);
116  if (costh_to_normal <1.e-4) costh_to_normal =1.e-4;
119  }
120 
123 
126 }
128 //
130 {
131  radius_spherical_source = radius;
132  center_spherical_source = center_pos;
133  type_of_adjoint_source ="Spherical";
141 }
143 //
145 {
147  type_of_adjoint_source ="ExternalSurfaceOfAVolume";
150 }
151 
153 //
155  G4ThreeVector glob_pos,
156  G4ThreeVector direction,
157  G4double,
162  G4ThreeVector position = glob_pos;
163  G4double safety=1.;
164  G4VPhysicalVolume* thePhysVolume =
166  G4double newStep =fLinearNavigator->ComputeStep(position,direction,1.e50,
167  safety);
170  //if (theAccumulatedCSDepthVector) delete theAccumulatedCSDepthVector;
171  //theAccumulatedCSDepthVector = new G4PhysicsOrderedFreeVector();
172 
173  G4double acc_depth=0.;
174  G4double acc_length=0.;
175  //G4double acc_cs_depth=0.;
176  //theAccumulatedCSDepthVector->InsertValues(acc_cs_depth, acc_length);
177  theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
178 
179  while (newStep > 0. && thePhysVolume) {
180  acc_length+=newStep;
181  /*
182  const G4MaterialCutsCouple* theMatCutsCouple=
183  thePhysVolume->GetLogicalVolume()->GetMaterialCutsCouple();
184 
185 
186  acc_cs_depth+=newStep*G4AdjointCSManager::GetAdjointCSManager()->GetTotalAdjointCS(aPartDef,
187  ekin,
188  theMatCutsCouple);
189  theAccumulatedCSDepthVector->InsertValues(acc_cs_depth, acc_length);*/
190 
191  acc_depth+=newStep*thePhysVolume->GetLogicalVolume()->GetMaterial()->GetDensity();
192  theAccumulatedDepthVector->InsertValues(acc_length,acc_depth);
193  position=position+newStep*direction;
194  thePhysVolume =
195  fLinearNavigator->LocateGlobalPointAndSetup(position,0,false);
196  newStep =fLinearNavigator->ComputeStep(position,direction,1.e50,
197  safety);
198  }
199 
200 
201 }
203 //
205 {G4double rand = G4UniformRand();
207  /*
208  G4double acc_cs_depth=theAccumulatedCSDepthVector->GetEnergy(distance);
209  weight_corr=std::exp(-acc_cs_depth);*/
210  weight_corr=1.;
211  return distance;
212 }
213 
214 
215 
216 
G4SPSAngDistribution * GetAngDist() const
G4SPSEneDistribution * GetEneDist() const
void SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector pos)
G4double FindLinearEnergy(G4double rand) const
static const double halfpi
Definition: G4SIunits.hh:76
G4SPSPosDistribution * GetPosDist() const
void SetPosDisType(G4String)
Important: This is a shared class between threads.
void ComputeAccumulatedDepthVectorAlongBackRay(G4ThreeVector glob_pos, G4ThreeVector direction, G4double ekin, G4ParticleDefinition *aPartDef)
CLHEP::Hep3Vector G4ThreeVector
void GenerateAdjointPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
G4Material * GetMaterial() const
void SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String &volume_name)
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:747
G4double GetDensity() const
Definition: G4Material.hh:180
G4Navigator * GetNavigatorForTracking() const
Andrea Dotti Feb 2015 Important: This is a shared class between threads.
G4double SampleDistanceAlongBackRayAndComputeWeightCorrection(G4double &weight_corr)
G4AdjointPosOnPhysVolGenerator * theG4AdjointPosOnPhysVolGenerator
void InsertValues(G4double energy, G4double value)
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
void SetCentreCoords(G4ThreeVector)
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4PhysicsOrderedFreeVector * theAccumulatedDepthVector
void GenerateFwdPrimaryVertex(G4Event *anEvt, G4ParticleDefinition *adj_part, G4double E1, G4double E2)
static G4TransportationManager * GetTransportationManager()
G4LogicalVolume * GetLogicalVolume() const
static const double pi
Definition: G4SIunits.hh:74
void GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(G4ThreeVector &p, G4ThreeVector &direction)
void DefinePhysicalVolume1(const G4String &aName)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:125
static G4AdjointPosOnPhysVolGenerator * GetInstance()
void GeneratePrimaryVertex(G4Event *evt)
double G4double
Definition: G4Types.hh:76
G4SingleParticleSource * theSingleParticleSource
static const G4double pos