source: trunk/examples/extended/parallel/MPI/exMPI02/src/MedicalBeam.cc@ 1036

Last change on this file since 1036 was 807, checked in by garnier, 17 years ago

update

File size: 4.5 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.1 2007/11/16 14:29:33 kmura Exp $
27// $Name: $
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
60 G4ParticleTable* particleTable= G4ParticleTable::GetParticleTable();
61 particle= particleTable-> FindParticle("proton");
62
63 kineticE= 200.*MeV;
64 sourcePosition= G4ThreeVector(0.,0.,-125.*cm);
65 SSD= 100.*cm;
66 fieldXY[0]= fieldXY[1]= 5.*cm;
67}
68
69
70///////////////////////////
71MedicalBeam::~MedicalBeam()
72///////////////////////////
73{
74}
75
76
77////////////////////////////////////////////////////////
78G4ThreeVector MedicalBeam::GenerateBeamDirection() const
79////////////////////////////////////////////////////////
80{
81 // uniform distribution in a limitted solid angle
82 G4double dr;
83 if(fieldShape == MedicalBeam::SQUARE) {
84 dr= std::sqrt(sqr(fieldXY[0]/2.)+sqr(fieldXY[1]/2.));
85 } else {
86 dr= fieldR;
87 }
88
89 G4double sin0= dr/SSD;
90 G4double cos0= std::sqrt(1.-sqr(sin0));
91
92 G4double dcos=0.;
93 G4double dsin, dphi, z;
94
95 G4double x= DBL_MAX;
96 G4double y= DBL_MAX;
97
98 G4double xmax, ymax;
99 if(fieldShape == MedicalBeam::SQUARE) {
100 xmax= fieldXY[0]/2./SSD;
101 ymax= fieldXY[1]/2./SSD;
102 } else {
103 xmax= ymax= DBL_MAX-1.;
104 }
105
106 while(! (std::abs(x)< xmax && std::abs(y)< ymax) ) {
107 dcos= RandFlat::shoot(cos0, 1.);
108 dsin= std::sqrt(1.-sqr(dcos));
109 dphi= RandFlat::shoot(0., twopi);
110
111 x= std::cos(dphi)*dsin*dcos;
112 y= std::sin(dphi)*dsin*dcos;
113 }
114 z= dcos;
115
116 return G4ThreeVector(x,y,z);
117}
118
119
120/////////////////////////////////////////////////////
121void MedicalBeam::GeneratePrimaries(G4Event* anEvent)
122/////////////////////////////////////////////////////
123{
124 if(particle==0) return;
125
126 // create a new vertex
127 G4PrimaryVertex* vertex= new G4PrimaryVertex(sourcePosition, 0.*ns);
128
129 // momentum
130 G4double mass= particle-> GetPDGMass();
131 G4double p= std::sqrt(sqr(mass+kineticE)-sqr(mass));
132 G4ThreeVector pmon= p*GenerateBeamDirection();
133 G4PrimaryParticle* primary= new G4PrimaryParticle(particle,
134 pmon.x(),
135 pmon.y(),
136 pmon.z());
137 // set primary to vertex
138 vertex-> SetPrimary(primary);
139
140 // set vertex to event
141 anEvent-> AddPrimaryVertex(vertex);
142}
143
Note: See TracBrowser for help on using the repository browser.