Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-10-14 08:09:02

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 //
0027 /// @file MedicalBeam.cc
0028 /// @brief Define beam profile as primary generator
0029 
0030 #include "MedicalBeam.hh"
0031 
0032 #include "G4Event.hh"
0033 #include "G4ParticleTable.hh"
0034 #include "G4PrimaryVertex.hh"
0035 #include "Randomize.hh"
0036 
0037 #include <cmath>
0038 
0039 using namespace CLHEP;
0040 
0041 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0042 MedicalBeam::MedicalBeam()
0043   : fparticle(0),
0044     fkineticE(1. * MeV),
0045     fsourcePosition(G4ThreeVector()),
0046     fSSD(1. * m),
0047     ffieldShape(MedicalBeam::kSQUARE),
0048     ffieldR(10. * cm)
0049 {
0050   G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
0051   fparticle = particleTable->FindParticle("proton");
0052 
0053   fkineticE = 200. * MeV;
0054   fsourcePosition = G4ThreeVector(0., 0., -125. * cm);
0055   fSSD = 100. * cm;
0056   ffieldXY[0] = ffieldXY[1] = 5. * cm;
0057 }
0058 
0059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0060 MedicalBeam::~MedicalBeam() {}
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0063 G4ThreeVector MedicalBeam::GenerateBeamDirection() const
0064 {
0065   // uniform distribution in a limitted solid angle
0066   G4double dr;
0067   if (ffieldShape == MedicalBeam::kSQUARE) {
0068     dr = std::sqrt(sqr(ffieldXY[0] / 2.) + sqr(ffieldXY[1] / 2.));
0069   }
0070   else {
0071     dr = ffieldR;
0072   }
0073 
0074   G4double sin0 = dr / fSSD;
0075   G4double cos0 = std::sqrt(1. - sqr(sin0));
0076 
0077   G4double dcos = 0.;
0078   G4double dsin, dphi, z;
0079 
0080   G4double x = DBL_MAX;
0081   G4double y = DBL_MAX;
0082 
0083   G4double xmax, ymax;
0084   if (ffieldShape == MedicalBeam::kSQUARE) {
0085     xmax = ffieldXY[0] / 2. / fSSD;
0086     ymax = ffieldXY[1] / 2. / fSSD;
0087   }
0088   else {
0089     xmax = ymax = DBL_MAX - 1.;
0090   }
0091 
0092   while (!(std::abs(x) < xmax && std::abs(y) < ymax)) {
0093     dcos = G4RandFlat::shoot(cos0, 1.);
0094     dsin = std::sqrt(1. - sqr(dcos));
0095     dphi = G4RandFlat::shoot(0., twopi);
0096 
0097     x = std::cos(dphi) * dsin * dcos;
0098     y = std::sin(dphi) * dsin * dcos;
0099   }
0100   z = dcos;
0101 
0102   return G4ThreeVector(x, y, z);
0103 }
0104 
0105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0106 void MedicalBeam::GeneratePrimaries(G4Event* anEvent)
0107 {
0108   if (fparticle == NULL) return;
0109 
0110   // create a new vertex
0111   G4PrimaryVertex* vertex = new G4PrimaryVertex(fsourcePosition, 0. * ns);
0112 
0113   // momentum
0114   G4double mass = fparticle->GetPDGMass();
0115   G4double p = std::sqrt(sqr(mass + fkineticE) - sqr(mass));
0116   G4ThreeVector pmon = p * GenerateBeamDirection();
0117   G4PrimaryParticle* primary = new G4PrimaryParticle(fparticle, pmon.x(), pmon.y(), pmon.z());
0118 
0119   // set primary to vertex
0120   vertex->SetPrimary(primary);
0121 
0122   // set vertex to event
0123   anEvent->AddPrimaryVertex(vertex);
0124 }