Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:21:09

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 /// \file eventgenerator/particleGun/src/PrimaryGeneratorAction2.cc
0027 /// \brief Implementation of the PrimaryGeneratorAction2 class
0028 //
0029 //
0030 //
0031 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0032 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0033 
0034 #include "PrimaryGeneratorAction2.hh"
0035 
0036 #include "PrimaryGeneratorAction.hh"
0037 
0038 #include "G4Event.hh"
0039 #include "G4ParticleDefinition.hh"
0040 #include "G4ParticleGun.hh"
0041 #include "G4ParticleTable.hh"
0042 #include "G4PhysicalConstants.hh"
0043 #include "G4SystemOfUnits.hh"
0044 #include "Randomize.hh"
0045 
0046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0047 
0048 PrimaryGeneratorAction2::PrimaryGeneratorAction2(G4ParticleGun* gun) : fParticleGun(gun)
0049 {
0050   // energy distribution
0051   //
0052   InitFunction();
0053 }
0054 
0055 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0056 
0057 void PrimaryGeneratorAction2::GeneratePrimaries(G4Event* anEvent)
0058 {
0059   // uniform solid angle
0060   G4double cosAlpha = 1. - 2 * G4UniformRand();  // cosAlpha uniform in [-1,+1]
0061   G4double sinAlpha = std::sqrt(1. - cosAlpha * cosAlpha);
0062   G4double psi = twopi * G4UniformRand();  // psi uniform in [0, 2*pi]
0063   G4ThreeVector dir(sinAlpha * std::cos(psi), sinAlpha * std::sin(psi), cosAlpha);
0064 
0065   fParticleGun->SetParticleMomentumDirection(dir);
0066 
0067   // set energy from a tabulated distribution
0068   //
0069   // G4double energy = RejectAccept();
0070   G4double energy = InverseCumul();
0071   fParticleGun->SetParticleEnergy(energy);
0072 
0073   // create vertex
0074   //
0075   fParticleGun->GeneratePrimaryVertex(anEvent);
0076 }
0077 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0078 
0079 void PrimaryGeneratorAction2::InitFunction()
0080 {
0081   // tabulated function
0082   // Y is assumed positive, linear per segment, continuous
0083   //
0084   fNPoints = 16;
0085   const G4double xx[] = {37 * keV,  39 * keV,  45 * keV,  51 * keV, 57 * keV, 69 * keV,
0086                          71 * keV,  75 * keV,  83 * keV,  91 * keV, 97 * keV, 107 * keV,
0087                          125 * keV, 145 * keV, 159 * keV, 160 * keV};
0088 
0089   const G4double yy[] = {0.000,  0.077,  0.380,  2.044,  5.535, 15.077, 12.443, 14.766,
0090                          17.644, 18.518, 17.772, 14.776, 8.372, 3.217,  0.194,  0.000};
0091 
0092   // copy arrays in std::vector and compute fMax
0093   //
0094   fX.resize(fNPoints);
0095   fY.resize(fNPoints);
0096   fYmax = 0.;
0097   for (G4int j = 0; j < fNPoints; j++) {
0098     fX[j] = xx[j];
0099     fY[j] = yy[j];
0100     if (fYmax < fY[j]) fYmax = fY[j];
0101   };
0102 
0103   // compute slopes
0104   //
0105   fSlp.resize(fNPoints);
0106   for (G4int j = 0; j < fNPoints - 1; j++) {
0107     fSlp[j] = (fY[j + 1] - fY[j]) / (fX[j + 1] - fX[j]);
0108   };
0109 
0110   // compute cumulative function
0111   //
0112   fYC.resize(fNPoints);
0113   fYC[0] = 0.;
0114   for (G4int j = 1; j < fNPoints; j++) {
0115     fYC[j] = fYC[j - 1] + 0.5 * (fY[j] + fY[j - 1]) * (fX[j] - fX[j - 1]);
0116   };
0117 }
0118 
0119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0120 
0121 G4double PrimaryGeneratorAction2::RejectAccept()
0122 {
0123   // tabulated function
0124   // Y is assumed positive, linear per segment, continuous
0125   // (see Particle Data Group: pdg.lbl.gov --> Monte Carlo techniques)
0126   //
0127   G4double Xrndm = 0., Yrndm = 0., Yinter = -1.;
0128 
0129   while (Yrndm > Yinter) {
0130     // choose a point randomly
0131     Xrndm = fX[0] + G4UniformRand() * (fX[fNPoints - 1] - fX[0]);
0132     Yrndm = G4UniformRand() * fYmax;
0133     // find bin
0134     G4int j = fNPoints - 2;
0135     while ((fX[j] > Xrndm) && (j > 0))
0136       j--;
0137     // compute Y(x_rndm) by linear interpolation
0138     Yinter = fY[j] + fSlp[j] * (Xrndm - fX[j]);
0139   };
0140   return Xrndm;
0141 }
0142 
0143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0144 
0145 G4double PrimaryGeneratorAction2::InverseCumul()
0146 {
0147   // tabulated function
0148   // Y is assumed positive, linear per segment, continuous
0149   // --> cumulative function is second order polynomial
0150   // (see Particle Data Group: pdg.lbl.gov --> Monte Carlo techniques)
0151 
0152   // choose y randomly
0153   G4double Yrndm = G4UniformRand() * fYC[fNPoints - 1];
0154   // find bin
0155   G4int j = fNPoints - 2;
0156   while ((fYC[j] > Yrndm) && (j > 0))
0157     j--;
0158   // y_rndm --> x_rndm :  fYC(x) is second order polynomial
0159   G4double Xrndm = fX[j];
0160   G4double a = fSlp[j];
0161   if (a != 0.) {
0162     G4double b = fY[j] / a, c = 2 * (Yrndm - fYC[j]) / a;
0163     G4double delta = b * b + c;
0164     G4int sign = 1;
0165     if (a < 0.) sign = -1;
0166     Xrndm += sign * std::sqrt(delta) - b;
0167   }
0168   else if (fY[j] > 0.) {
0169     Xrndm += (Yrndm - fYC[j]) / fY[j];
0170   };
0171   return Xrndm;
0172 }
0173 
0174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......