source: trunk/environments/g4py/site-modules/primaries/MedicalBeam/MedicalBeam.cc

Last change on this file was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 4.3 KB
Line 
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// $Id: MedicalBeam.cc,v 1.6 2006/06/29 15:30:10 gunter Exp $
27// $Name: geant4-09-04-beta-01 $
28// ====================================================================
29//   MedicalBeam.cc
30//
31//                                         2005 Q
32// ====================================================================
33#include "MedicalBeam.hh"
34#include "Randomize.hh"
35#include "G4Event.hh"
36#include "G4ParticleTable.hh"
37#include "G4PrimaryVertex.hh"
38
39using namespace CLHEP;
40
41#include <cmath>
42
43// ====================================================================
44//
45// class description
46//
47// ====================================================================
48
49///////////////////////////////////
50MedicalBeam::MedicalBeam()
51  : particle(0), 
52    kineticE(1.*MeV),
53    sourcePosition(G4ThreeVector()),
54    SSD(1.*m), 
55    fieldShape(MedicalBeam::SQUARE),
56    fieldR(10.*cm)
57///////////////////////////////////
58{
59  fieldXY[0]= fieldXY[1]= 10.*cm;
60}
61
62
63///////////////////////////
64MedicalBeam::~MedicalBeam()
65///////////////////////////
66{
67}
68
69
70////////////////////////////////////////////////////////
71G4ThreeVector MedicalBeam::GenerateBeamDirection() const
72////////////////////////////////////////////////////////
73{
74  // uniform distribution in a limitted solid angle
75  G4double dr;
76  if(fieldShape == MedicalBeam::SQUARE) {
77    dr= std::sqrt(sqr(fieldXY[0]/2.)+sqr(fieldXY[1]/2.));
78  } else {
79    dr= fieldR;
80  }
81
82  G4double sin0= dr/SSD;
83  G4double cos0= std::sqrt(1.-sqr(sin0));
84
85  G4double dcos, dsin, dphi, z;
86
87  G4double x= DBL_MAX;
88  G4double y= DBL_MAX;
89
90  G4double xmax, ymax;
91  if(fieldShape == MedicalBeam::SQUARE) {
92    xmax= fieldXY[0]/2./SSD;
93    ymax= fieldXY[1]/2./SSD;
94  } else {
95    xmax= ymax= DBL_MAX-1.;
96  }
97
98  while(! (std::abs(x)< xmax && std::abs(y)< ymax) ) {
99    dcos= RandFlat::shoot(cos0, 1.);
100    dsin= std::sqrt(1.-sqr(dcos));
101    dphi= RandFlat::shoot(0., twopi);
102
103    x= std::cos(dphi)*dsin*dcos;
104    y= std::sin(dphi)*dsin*dcos;
105  }
106  z= dcos;
107
108  return G4ThreeVector(x,y,z);
109}
110
111
112/////////////////////////////////////////////////////
113void MedicalBeam::GeneratePrimaries(G4Event* anEvent)
114/////////////////////////////////////////////////////
115{
116  if(particle==0) return;
117
118  // create a new vertex
119  G4PrimaryVertex* vertex= new G4PrimaryVertex(sourcePosition, 0.*ns);
120
121  // momentum
122  G4double mass= particle-> GetPDGMass();
123  G4double p= std::sqrt(sqr(mass+kineticE)-sqr(mass));
124  G4ThreeVector pmon= p*GenerateBeamDirection();
125  G4PrimaryParticle* primary= new G4PrimaryParticle(particle, 
126                                                    pmon.x(), 
127                                                    pmon.y(), 
128                                                    pmon.z());
129  // set primary to vertex
130  vertex-> SetPrimary(primary);
131
132  // set vertex to event
133  anEvent-> AddPrimaryVertex(vertex);
134}
135
Note: See TracBrowser for help on using the repository browser.