Geant4
v4-10.4-release
메인 페이지
관련된 페이지
모듈
네임스페이스
클래스
파일들
파일 목록
파일 멤버
모두
클래스
네임스페이스들
파일들
함수
변수
타입정의
열거형 타입
열거형 멤버
Friends
매크로
그룹들
페이지들
examples
extended
parallel
MPI
examples
exMPI03
src
examples/extended/parallel/MPI/examples/exMPI03/src/MedicalBeam.cc
이 파일의 문서화 페이지로 가기
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 76431 2013-11-10 20:28:49Z kmura $
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 =
G4RandFlat::shoot
(cos0, 1.);
89
dsin = std::sqrt(1.-
sqr
(dcos));
90
dphi =
G4RandFlat::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
}
x
Float_t x
Definition:
compare.C:6
MedicalBeam::GeneratePrimaries
virtual void GeneratePrimaries(G4Event *anEvent)
Definition:
environments/g4py/site-modules/primaries/MedicalBeam/MedicalBeam.cc:112
G4ThreeVector
CLHEP::Hep3Vector G4ThreeVector
Definition:
G4ThreeVector.hh:42
G4ParticleTable::GetParticleTable
static G4ParticleTable * GetParticleTable()
Definition:
G4ParticleTable.cc:105
MeV
static constexpr double MeV
Definition:
G4SIunits.hh:214
G4PrimaryVertex.hh
Randomize.hh
MedicalBeam::fSSD
G4double fSSD
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:81
y
Float_t y
Definition:
compare.C:6
p
const char * p
Definition:
xmltok.h:285
z
Double_t z
Definition:
advanced/microbeam/plot.C:277
MedicalBeam::fkineticE
G4double fkineticE
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:78
CLHEP::Hep3Vector::z
double z() const
G4PrimaryParticle
Definition:
G4PrimaryParticle.hh:68
ymax
Double_t ymax
Definition:
extended/medical/dna/microyz/plot.C:25
G4Event.hh
CLHEP::cm
static constexpr double cm
Definition:
SystemOfUnits.h:99
G4double
double G4double
Definition:
G4Types.hh:76
MedicalBeam::kSQUARE
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:42
G4Event
Definition:
G4Event.hh:52
CLHEP::MeV
static constexpr double MeV
Definition:
SystemOfUnits.h:174
twopi
static constexpr double twopi
Definition:
G4SIunits.hh:76
MedicalBeam::MedicalBeam
MedicalBeam()
Definition:
environments/g4py/site-modules/primaries/MedicalBeam/MedicalBeam.cc:49
MedicalBeam::ffieldXY
G4double ffieldXY[2]
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:83
MedicalBeam::fparticle
G4ParticleDefinition * fparticle
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:77
G4ParticleTable.hh
G4INCL::DeJongSpin::shoot
ThreeVector shoot(const G4int Ap, const G4int Af)
Definition:
G4INCLDeJongSpin.cc:68
CLHEP::m
static constexpr double m
Definition:
SystemOfUnits.h:109
CLHEP::Hep3Vector
Definition:
ThreeVector.h:41
MedicalBeam::~MedicalBeam
~MedicalBeam()
Definition:
environments/g4py/site-modules/primaries/MedicalBeam/MedicalBeam.cc:63
cm
static constexpr double cm
Definition:
G4SIunits.hh:119
CLHEP::Hep3Vector::x
double x() const
MedicalBeam::ffieldR
G4double ffieldR
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:84
sqr
T sqr(const T &x)
Definition:
templates.hh:145
CLHEP::Hep3Vector::y
double y() const
G4ParticleTable
Definition:
G4ParticleTable.hh:66
MedicalBeam::GenerateBeamDirection
G4ThreeVector GenerateBeamDirection() const
Definition:
environments/g4py/site-modules/primaries/MedicalBeam/MedicalBeam.cc:70
DBL_MAX
#define DBL_MAX
Definition:
templates.hh:83
ns
#define ns
Definition:
xmlparse.cc:614
MedicalBeam::ffieldShape
FieldShape ffieldShape
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:82
MedicalBeam::fsourcePosition
G4ThreeVector fsourcePosition
Definition:
examples/extended/parallel/MPI/examples/exMPI02/include/MedicalBeam.hh:79
MedicalBeam
Definition:
environments/g4py/site-modules/primaries/MedicalBeam/MedicalBeam.hh:46
G4PrimaryVertex
Definition:
G4PrimaryVertex.hh:50
다음에 의해 생성됨 :
1.8.5