Geant4 Cross Reference |
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 33 #include "G4AdjointPrimaryGenerator.hh" 34 #include "G4PhysicalConstants.hh" 35 #include "G4Event.hh" 36 #include "G4SingleParticleSource.hh" 37 #include "G4ParticleDefinition.hh" 38 #include "G4AdjointPosOnPhysVolGenerator.hh" 39 #include "G4Navigator.hh" 40 #include "G4TransportationManager.hh" 41 #include "G4VPhysicalVolume.hh" 42 #include "G4Material.hh" 43 #include "Randomize.hh" 44 45 // -------------------------------------------------------------------- 46 // 47 G4AdjointPrimaryGenerator::G4AdjointPrimaryGenerator() 48 { 49 center_spherical_source = G4ThreeVector(0.,0.,0.); 50 type_of_adjoint_source="Spherical"; 51 theSingleParticleSource = new G4SingleParticleSource(); 52 53 theSingleParticleSource->GetEneDist()->SetEnergyDisType("Pow"); 54 theSingleParticleSource->GetEneDist()->SetAlpha(-1.); 55 theSingleParticleSource->GetPosDist()->SetPosDisType("Point"); 56 theSingleParticleSource->GetAngDist()->SetAngDistType("planar"); 57 58 theG4AdjointPosOnPhysVolGenerator = G4AdjointPosOnPhysVolGenerator::GetInstance(); 59 } 60 61 // -------------------------------------------------------------------- 62 // 63 G4AdjointPrimaryGenerator::~G4AdjointPrimaryGenerator() 64 { 65 delete theSingleParticleSource; 66 } 67 68 // -------------------------------------------------------------------- 69 // 70 void G4AdjointPrimaryGenerator:: 71 GenerateAdjointPrimaryVertex(G4Event* anEvent, G4ParticleDefinition* adj_part, 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.; 80 G4ThreeVector pos =G4ThreeVector(0.,0.,0.); 81 G4ThreeVector direction = G4ThreeVector(0.,0.,1.); 82 theG4AdjointPosOnPhysVolGenerator 83 ->GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(pos, direction, 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 // 89 theSingleParticleSource->GetAngDist() 90 ->SetParticleMomentumDirection(-direction); 91 theSingleParticleSource->GetPosDist()->SetCentreCoords(pos); 92 } 93 94 theSingleParticleSource->GetEneDist()->SetEmin(E1); 95 theSingleParticleSource->GetEneDist()->SetEmax(E2); 96 97 theSingleParticleSource->SetParticleDefinition(adj_part); 98 theSingleParticleSource->GeneratePrimaryVertex(anEvent); 99 } 100 101 // -------------------------------------------------------------------- 102 // 103 void G4AdjointPrimaryGenerator:: 104 GenerateFwdPrimaryVertex(G4Event* anEvent,G4ParticleDefinition* fwd_part, 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.; 113 G4ThreeVector pos =G4ThreeVector(0.,0.,0.); 114 G4ThreeVector direction = G4ThreeVector(0.,0.,1.); 115 theG4AdjointPosOnPhysVolGenerator 116 ->GenerateAPositionOnTheExtSurfaceOfThePhysicalVolume(pos, direction, 117 costh_to_normal); 118 if (costh_to_normal <1.e-4) { costh_to_normal =1.e-4; } 119 theSingleParticleSource->GetAngDist() 120 ->SetParticleMomentumDirection(direction); 121 theSingleParticleSource->GetPosDist()->SetCentreCoords(pos); 122 } 123 124 theSingleParticleSource->GetEneDist()->SetEmin(E1); 125 theSingleParticleSource->GetEneDist()->SetEmax(E2); 126 127 theSingleParticleSource->SetParticleDefinition(fwd_part); 128 theSingleParticleSource->GeneratePrimaryVertex(anEvent); 129 } 130 131 // -------------------------------------------------------------------- 132 // 133 void G4AdjointPrimaryGenerator:: 134 SetSphericalAdjointPrimarySource(G4double radius, G4ThreeVector center_pos) 135 { 136 radius_spherical_source = radius; 137 center_spherical_source = center_pos; 138 type_of_adjoint_source = "Spherical"; 139 theSingleParticleSource->GetPosDist()->SetPosDisType("Surface"); 140 theSingleParticleSource->GetPosDist()->SetPosDisShape("Sphere"); 141 theSingleParticleSource->GetPosDist()->SetCentreCoords(center_pos); 142 theSingleParticleSource->GetPosDist()->SetRadius(radius); 143 theSingleParticleSource->GetAngDist()->SetAngDistType("cos"); 144 theSingleParticleSource->GetAngDist()->SetMaxTheta(pi); 145 theSingleParticleSource->GetAngDist()->SetMinTheta(halfpi); 146 } 147 148 // -------------------------------------------------------------------- 149 // 150 void G4AdjointPrimaryGenerator:: 151 SetAdjointPrimarySourceOnAnExtSurfaceOfAVolume(const G4String& volume_name) 152 { 153 theG4AdjointPosOnPhysVolGenerator->DefinePhysicalVolume1(volume_name); 154 type_of_adjoint_source ="ExternalSurfaceOfAVolume"; 155 theSingleParticleSource->GetPosDist()->SetPosDisType("Point"); 156 theSingleParticleSource->GetAngDist()->SetAngDistType("planar"); 157 } 158 159 // -------------------------------------------------------------------- 160 // 161 void G4AdjointPrimaryGenerator:: 162 ComputeAccumulatedDepthVectorAlongBackRay(G4ThreeVector glob_pos, 163 G4ThreeVector direction, 164 G4double, G4ParticleDefinition*) 165 { 166 if (fLinearNavigator == nullptr) 167 { 168 fLinearNavigator = G4TransportationManager::GetTransportationManager() 169 ->GetNavigatorForTracking(); 170 } 171 G4ThreeVector position = glob_pos; 172 G4double safety=1.; 173 G4VPhysicalVolume* thePhysVolume = 174 fLinearNavigator->LocateGlobalPointAndSetup(position); 175 G4double newStep = fLinearNavigator->ComputeStep(position,direction,1.e50, 176 safety); 177 delete theAccumulatedDepthVector; 178 theAccumulatedDepthVector = new G4PhysicsFreeVector(); 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 // 199 G4double G4AdjointPrimaryGenerator:: 200 SampleDistanceAlongBackRayAndComputeWeightCorrection(G4double& weight_corr) 201 { 202 G4double rand = G4UniformRand(); 203 G4double distance = theAccumulatedDepthVector->FindLinearEnergy(rand); 204 weight_corr=1.; 205 return distance; 206 } 207