Geant4
10.01
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
MedicalBeam.cc
Go to the documentation of this file.
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 78126 2013-12-03 17:43:56Z gcosmo $
27
//
30
31
#include <cmath>
32
#include "
G4Event.hh
"
33
#include "
G4ParticleTable.hh
"
34
#include "
G4PrimaryVertex.hh
"
35
#include "
Randomize.hh
"
36
#include "MedicalBeam.hh"
37
38
using namespace
CLHEP
;
39
40
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
MedicalBeam::MedicalBeam
()
42
: fparticle(0), fkineticE(1.*
MeV
), fsourcePosition(
G4ThreeVector
()),
43
fSSD(1.*
m
), ffieldShape(
MedicalBeam
::kSQUARE), ffieldR(10.*
cm
)
44
{
45
G4ParticleTable
* particleTable =
G4ParticleTable::GetParticleTable
();
46
fparticle
= particleTable-> FindParticle(
"proton"
);
47
48
fkineticE
= 200.*
MeV
;
49
fsourcePosition
=
G4ThreeVector
(0.,0.,-125.*
cm
);
50
fSSD
= 100.*
cm
;
51
ffieldXY
[0] =
ffieldXY
[1] = 5.*
cm
;
52
}
53
54
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55
MedicalBeam::~MedicalBeam
()
56
{
57
}
58
59
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
G4ThreeVector
MedicalBeam::GenerateBeamDirection
()
const
61
{
62
// uniform distribution in a limitted solid angle
63
G4double
dr;
64
if
(
ffieldShape
==
MedicalBeam::kSQUARE
) {
65
dr = std::sqrt(
sqr
(
ffieldXY
[0]/2.) +
sqr
(
ffieldXY
[1]/2.));
66
}
else
{
67
dr =
ffieldR
;
68
}
69
70
G4double
sin0 = dr/
fSSD
;
71
G4double
cos0 = std::sqrt(1.-
sqr
(sin0));
72
73
G4double
dcos = 0.;
74
G4double
dsin, dphi,
z
;
75
76
G4double
x =
DBL_MAX
;
77
G4double
y =
DBL_MAX
;
78
79
G4double
xmax, ymax;
80
if
(
ffieldShape
==
MedicalBeam::kSQUARE
) {
81
xmax =
ffieldXY
[0]/2./
fSSD
;
82
ymax =
ffieldXY
[1]/2./
fSSD
;
83
}
else
{
84
xmax = ymax =
DBL_MAX
-1.;
85
}
86
87
while
(! (std::abs(x) < xmax && std::abs(y) < ymax) ) {
88
dcos =
RandFlat::shoot
(cos0, 1.);
89
dsin = std::sqrt(1.-
sqr
(dcos));
90
dphi =
RandFlat::shoot
(0., twopi);
91
92
x = std::cos(dphi)*dsin*dcos;
93
y = std::sin(dphi)*dsin*dcos;
94
}
95
z = dcos;
96
97
return
G4ThreeVector
(x,y,z);
98
}
99
100
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
void
MedicalBeam::GeneratePrimaries
(
G4Event
* anEvent)
102
{
103
if
(
fparticle
== NULL )
return
;
104
105
// create a new vertex
106
G4PrimaryVertex
* vertex =
new
G4PrimaryVertex
(
fsourcePosition
, 0.*
ns
);
107
108
// momentum
109
G4double
mass =
fparticle
-> GetPDGMass();
110
G4double
p = std::sqrt(
sqr
(mass+
fkineticE
)-
sqr
(mass));
111
G4ThreeVector
pmon = p*
GenerateBeamDirection
();
112
G4PrimaryParticle
* primary =
113
new
G4PrimaryParticle
(
fparticle
, pmon.x(), pmon.y(), pmon.z());
114
115
// set primary to vertex
116
vertex-> SetPrimary(primary);
117
118
// set vertex to event
119
anEvent-> AddPrimaryVertex(vertex);
120
}
cm
static const double cm
Definition:
G4SIunits.hh:106
G4INCL::DeJongSpin::shoot
ThreeVector shoot(const G4int Ap, const G4int Af)
Definition:
G4INCLDeJongSpin.cc:68
MeV
static const double MeV
Definition:
G4SIunits.hh:193
G4ThreeVector
CLHEP::Hep3Vector G4ThreeVector
Definition:
G4ThreeVector.hh:42
MedicalBeam::ffieldXY
G4double ffieldXY[2]
Definition:
MedicalBeam.hh:83
MedicalBeam::kSQUARE
Definition:
MedicalBeam.hh:42
z
G4double z
Definition:
TRTMaterials.hh:39
MedicalBeam::GenerateBeamDirection
G4ThreeVector GenerateBeamDirection() const
Definition:
MedicalBeam.cc:70
G4PrimaryVertex
Definition:
G4PrimaryVertex.hh:50
G4ParticleTable.hh
MedicalBeam::MedicalBeam
MedicalBeam()
Definition:
MedicalBeam.cc:49
G4ParticleTable
Definition:
G4ParticleTable.hh:65
G4PrimaryVertex.hh
Randomize.hh
MedicalBeam::ffieldR
G4double ffieldR
Definition:
MedicalBeam.hh:84
MedicalBeam::fSSD
G4double fSSD
Definition:
MedicalBeam.hh:81
MedicalBeam::~MedicalBeam
~MedicalBeam()
Definition:
MedicalBeam.cc:63
MedicalBeam::fkineticE
G4double fkineticE
Definition:
MedicalBeam.hh:78
G4ParticleTable::GetParticleTable
static G4ParticleTable * GetParticleTable()
Definition:
G4ParticleTable.cc:96
MedicalBeam
Definition:
MedicalBeam.hh:46
CLHEP
Definition:
AxisAngle.cc:16
MedicalBeam::fsourcePosition
G4ThreeVector fsourcePosition
Definition:
MedicalBeam.hh:79
G4PrimaryParticle
Definition:
G4PrimaryParticle.hh:68
m
static const double m
Definition:
G4SIunits.hh:110
sqr
T sqr(const T &x)
Definition:
templates.hh:145
G4double
double G4double
Definition:
G4Types.hh:76
G4Event.hh
MedicalBeam::GeneratePrimaries
virtual void GeneratePrimaries(G4Event *anEvent)
Definition:
MedicalBeam.cc:112
DBL_MAX
#define DBL_MAX
Definition:
templates.hh:83
MedicalBeam::ffieldShape
FieldShape ffieldShape
Definition:
MedicalBeam.hh:82
ns
#define ns
Definition:
xmlparse.cc:597
G4Event
Definition:
G4Event.hh:52
MedicalBeam::fparticle
G4ParticleDefinition * fparticle
Definition:
MedicalBeam.hh:77
source
examples
extended
parallel
MPI
examples
exMPI02
src
MedicalBeam.cc
Generated on Fri Dec 12 2014 22:47:01 for Geant4 by
1.8.8