Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/xray_TESdetector/src/XrayTESdetPrimaryGeneratorAction.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 /// \file XrayTESdetPrimaryGeneratorAction.cc
 28 /// \brief Implementation of the PrimaryGeneratorAction class
 29 //
 30 // Authors: P.Dondero (paolo.dondero@cern.ch), R.Stanzani (ronny.stanzani@cern.ch)
 31 //
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 34 
 35 
 36 #include "XrayTESdetPrimaryGeneratorAction.hh"
 37 #include "XrayTESdetDetectorConstruction.hh"
 38 
 39 #include "G4Event.hh"
 40 #include "G4ParticleGun.hh"
 41 #include "G4ThreeVector.hh"
 42 #include "G4ParticleTable.hh"
 43 #include "G4ParticleDefinition.hh"
 44 #include "Randomize.hh"
 45 #include "globals.hh"
 46 #include "G4NistManager.hh"
 47 #include "G4SystemOfUnits.hh"
 48 #include "G4PhysicalConstants.hh"
 49 
 50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 51 
 52 XrayTESdetPrimaryGeneratorAction::XrayTESdetPrimaryGeneratorAction()
 53 {
 54   G4int n_particle = 1;
 55   fParticleGun = new G4ParticleGun(n_particle);
 56   G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
 57   G4String particleName;
 58 
 59   // setting particle type
 60    fParticleGun->SetParticleDefinition(particleTable->FindParticle(particleName="proton"));
 61 }
 62 
 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 64 
 65 XrayTESdetPrimaryGeneratorAction::~XrayTESdetPrimaryGeneratorAction()
 66 {
 67   delete fParticleGun;
 68 }
 69 
 70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 71 
 72 void XrayTESdetPrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
 73 {
 74   G4double ax;
 75   G4double ay;
 76   G4double az;
 77 
 78   G4double bx;
 79   G4double by;
 80   G4double bz;
 81 
 82   G4double vx;
 83   G4double vy;
 84   G4double vz;
 85 
 86   G4double energy = 0;
 87   G4double u;
 88   G4double u1;
 89 
 90   G4double theta0;
 91   G4double phi0;
 92   G4double theta1;
 93   G4double phi1;
 94 
 95   G4double ene;
 96   G4double flux;
 97   G4double r;
 98 
 99   // Radius to shot particles
100   r = 27*cm;
101 
102   // Proton flux @L2 (Pamela2009), in p/cm2/s/sr [10 MeV - 10 GeV] = 0.407
103   G4double phi = 0.379;
104   G4double tr = 0.938;
105   G4double T;
106 
107   while (true)
108   {
109     u = G4UniformRand();
110     u1 = 0.268*G4UniformRand();      //max of the flux, in particles/cm2-s-GeV-sr
111     ene = 10.*MeV + u*100000.*MeV;
112     T = ene*0.001;
113     flux = 1.9*std::pow(((T + phi)*(T + phi + 2*tr)), -1.39)*T*(T + 2*tr)/((1 + 0.4866*std::pow(((T + phi)*(T + phi + 2*tr)),-1.255))*(T + phi)*(T + phi + 2*tr));
114     if (u1 <= flux)
115     {
116       energy = ene;
117       break;
118     }
119   }
120 
121   theta0 = std::acos(1. - 2.*G4UniformRand());
122   theta1 = std::acos(1. - 2.*G4UniformRand());
123   phi0 = twopi*G4UniformRand();
124   phi1 = twopi*G4UniformRand();
125   ax = r*std::sin(theta0)*std::cos(phi0);
126   ay = r*std::sin(theta0)*std::sin(phi0);
127   az = r*std::cos(theta0) + 80*mm;
128 
129   bx = r*std::sin(theta1)*std::cos(phi1);
130   by = r*std::sin(theta1)*std::sin(phi1);
131   bz = r*std::cos(theta1) + 80*mm;
132   vx = bx - ax;
133   vy = by - ay;
134   vz = bz - az;
135 
136   fParticleGun->SetParticleEnergy(energy);
137   fParticleGun->SetParticlePosition(G4ThreeVector(ax, ay, az));
138   fParticleGun->SetParticleMomentumDirection(G4ThreeVector(vx,vy,vz));
139   fParticleGun->GeneratePrimaryVertex(anEvent);
140 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143