Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:22:04

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // This example is provided by the Geant4-DNA collaboration
0027 // Any report or published results obtained using the Geant4-DNA software
0028 // shall cite the following Geant4-DNA collaboration publication:
0029 // Med. Phys. 37 (2010) 4692-4708
0030 // and papers
0031 // M. Batmunkh et al. J Radiat Res Appl Sci 8 (2015) 498-507
0032 // O. Belov et al. Physica Medica 32 (2016) 1510-1520
0033 // The Geant4-DNA web site is available at http://geant4-dna.org
0034 //
0035 // -------------------------------------------------------------------
0036 // November 2016
0037 // -------------------------------------------------------------------
0038 //
0039 /// \file PrimaryGeneratorAction.cc
0040 /// \brief Implementation of the PrimaryGeneratorAction class
0041 
0042 #include "PrimaryGeneratorAction.hh"
0043 
0044 #include "G4SystemOfUnits.hh"
0045 //
0046 #include "CommandLineParser.hh"
0047 
0048 #include "G4Box.hh"
0049 #include "G4Event.hh"
0050 #include "G4LogicalVolume.hh"
0051 #include "G4LogicalVolumeStore.hh"
0052 #include "G4Orb.hh"
0053 #include "G4ParticleDefinition.hh"
0054 #include "G4ParticleGun.hh"
0055 #include "G4PhysicalConstants.hh"
0056 #include "G4PhysicalVolumeStore.hh"
0057 #include "G4Proton.hh"
0058 #include "G4VPhysicalVolume.hh"
0059 #include "Randomize.hh"
0060 
0061 using namespace G4DNAPARSER;
0062 
0063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0064 
0065 PrimaryGeneratorAction::PrimaryGeneratorAction() : G4VUserPrimaryGeneratorAction()
0066 {
0067   G4int n_particle = 1;
0068   fpParticleGun = new G4ParticleGun(n_particle);
0069   //
0070   G4ParticleDefinition* particle = G4Proton::Proton();
0071   fpParticleGun->SetParticleDefinition(particle);
0072   // default gun parameters
0073   fpParticleGun->SetParticleEnergy(10. * MeV);
0074   fpParticleGun->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.));
0075   fpParticleGun->SetParticlePosition(G4ThreeVector(0., 0., 0.));
0076 }
0077 
0078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0079 
0080 PrimaryGeneratorAction::~PrimaryGeneratorAction()
0081 {
0082   delete fpParticleGun;
0083 }
0084 
0085 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0086 
0087 void PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
0088 {
0089   // Initial kinetic energy of particles beam as Gaussion distribution!
0090 
0091   // G4double Ekin = 10.*MeV;
0092   // G4double deltaEkin = 9.0955e-5*MeV;
0093   // fpParticleGun->SetParticleEnergy(G4RandGauss::shoot(Ekin,deltaEkin));
0094 
0095   // In order to avoid dependence of PrimaryGeneratorAction
0096   // on DetectorConstruction class we get world volume
0097   // from G4LogicalVolumeStore
0098 
0099   // We included three options for particles direction:
0100 
0101   G4double mediumRadius = 0.;
0102   G4double boundingXHalfLength = 0.;
0103   G4double boundingYHalfLength = 0.;
0104   G4double boundingZHalfLength = 0.;
0105   G4LogicalVolume* mediumLV = G4LogicalVolumeStore::GetInstance()->GetVolume("Medium");
0106   G4LogicalVolume* boundingLV = G4LogicalVolumeStore::GetInstance()->GetVolume("BoundingSlice");
0107   G4Orb* mediumSphere = 0;
0108   G4Box* boundingSlice = 0;
0109   if (mediumLV && boundingLV) {
0110     mediumSphere = dynamic_cast<G4Orb*>(mediumLV->GetSolid());
0111     boundingSlice = dynamic_cast<G4Box*>(boundingLV->GetSolid());
0112   }
0113   if (mediumSphere && boundingSlice) {
0114     mediumRadius = mediumSphere->GetRadius();
0115     boundingXHalfLength = boundingSlice->GetXHalfLength();
0116     boundingYHalfLength = boundingSlice->GetYHalfLength();
0117     boundingZHalfLength = boundingSlice->GetZHalfLength();
0118 
0119     /// a) Partilces directed to "square" on the XY plane of bounding slice
0120     ///   (or YZ, XZ)
0121 
0122     if (CommandLineParser::GetParser()->GetCommandIfActive("-sXY")) {
0123       // INITIAL BEAM DIVERGENCE
0124       fpParticleGun->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.));
0125       // // INITIAL BEAM POSITION
0126       fpParticleGun->SetParticlePosition(G4ThreeVector(
0127         CLHEP::RandFlat::shoot(-boundingXHalfLength, boundingXHalfLength),
0128         CLHEP::RandFlat::shoot(-boundingYHalfLength, boundingYHalfLength), -mediumRadius));
0129     }
0130 
0131     /// b) Partilces directed to "disk" on the XY plane of
0132     ///    bounding slice (or YZ, XZ)
0133 
0134     else if (CommandLineParser::GetParser()->GetCommandIfActive("-dXY")) {
0135       fpParticleGun->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.));
0136       G4double x0, y0, z0;
0137       x0 = 100. * mm;
0138       y0 = 100. * mm;
0139       z0 = -mediumRadius;  // mediumRadius;
0140       while (!(std::sqrt(x0 * x0 + y0 * y0) <= mediumRadius)) {
0141         x0 = CLHEP::RandFlat::shoot(-mediumRadius, mediumRadius);
0142         y0 = CLHEP::RandFlat::shoot(-mediumRadius, mediumRadius);
0143       }
0144       fpParticleGun->SetParticlePosition(G4ThreeVector(x0, y0, z0));
0145     }
0146 
0147     /// c) Partilces directed towards the bounding slice (default option!)
0148     // Select a starting position on a sphere including the
0149     // target volume and neuron morphology
0150     else {
0151       G4double cosTheta = 2. * G4UniformRand() - 1;
0152       G4double sinTheta = std::sqrt(1. - cosTheta * cosTheta);
0153       G4double phi = twopi * G4UniformRand();
0154       G4ThreeVector positionStart(mediumRadius * sinTheta * std::cos(phi),
0155                                   mediumRadius * sinTheta * std::sin(phi), mediumRadius * cosTheta);
0156       fpParticleGun->SetParticlePosition(positionStart);
0157       // To compute the direction, select a point inside the target volume
0158       G4ThreeVector positionDir(boundingXHalfLength * (2. * G4UniformRand() - 1),
0159                                 boundingYHalfLength * (2. * G4UniformRand() - 1),
0160                                 boundingZHalfLength * (2. * G4UniformRand() - 1));
0161       fpParticleGun->SetParticleMomentumDirection((positionDir - positionStart).unit());
0162       // Surface area of sphere
0163       fGunArea = 4. * pi * mediumRadius * mediumRadius;
0164     }
0165   }
0166   else {
0167     G4cerr << "Bounding slice volume not found!" << G4endl;
0168     G4cerr << "Default particle kinematic used" << G4endl;
0169   }
0170 
0171   fpParticleGun->GeneratePrimaryVertex(anEvent);
0172 }