LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
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 // ====================================================================
27 // MedicalBeam.cc
28 //
29 // 2005 Q
30 // ====================================================================
31 #include "MedicalBeam.hh"
32 #include "Randomize.hh"
33 #include "G4Event.hh"
34 #include "G4ParticleTable.hh"
35 #include "G4PrimaryVertex.hh"
36 
37 using namespace CLHEP;
38 
39 #include <cmath>
40 
41 // ====================================================================
42 //
43 // class description
44 //
45 // ====================================================================
46 
49  : particle(0),
50  kineticE(1.*MeV),
51  sourcePosition(G4ThreeVector()),
52  SSD(1.*m),
53  fieldShape(MedicalBeam::SQUARE),
54  fieldR(10.*cm)
56 {
57  fieldXY[0]= fieldXY[1]= 10.*cm;
58 }
59 
60 
64 {
65 }
66 
67 
71 {
72  // uniform distribution in a limitted solid angle
73  G4double dr;
75  dr= std::sqrt(sqr(fieldXY[0]/2.)+sqr(fieldXY[1]/2.));
76  } else {
77  dr= fieldR;
78  }
79 
80  G4double sin0= dr/SSD;
81  G4double cos0= std::sqrt(1.-sqr(sin0));
82 
83  G4double dcos, dsin, dphi, z;
84 
85  G4double x= DBL_MAX;
86  G4double y= DBL_MAX;
87 
88  G4double xmax, ymax;
90  xmax= fieldXY[0]/2./SSD;
91  ymax= fieldXY[1]/2./SSD;
92  } else {
93  xmax= ymax= DBL_MAX-1.;
94  }
95 
96  while(! (std::abs(x)< xmax && std::abs(y)< ymax) ) {
97  dcos= RandFlat::shoot(cos0, 1.);
98  dsin= std::sqrt(1.-sqr(dcos));
99  dphi= RandFlat::shoot(0., twopi);
100 
101  x= std::cos(dphi)*dsin*dcos;
102  y= std::sin(dphi)*dsin*dcos;
103  }
104  z= dcos;
105 
106  return G4ThreeVector(x,y,z);
107 }
108 
109 
111 void MedicalBeam::GeneratePrimaries(G4Event* anEvent)
113 {
114  if(particle==0) return;
115 
116  // create a new vertex
117  G4PrimaryVertex* vertex= new G4PrimaryVertex(sourcePosition, 0.*ns);
118 
119  // momentum
120  G4double mass= particle-> GetPDGMass();
121  G4double p= std::sqrt(sqr(mass+kineticE)-sqr(mass));
122  G4ThreeVector pmon= p*GenerateBeamDirection();
123  G4PrimaryParticle* primary= new G4PrimaryParticle(particle,
124  pmon.x(),
125  pmon.y(),
126  pmon.z());
127  // set primary to vertex
128  vertex-> SetPrimary(primary);
129 
130  // set vertex to event
131  anEvent-> AddPrimaryVertex(vertex);
132 }
133 
Float_t x
Definition: compare.C:6
G4ParticleDefinition * particle
Definition: MedicalBeam.hh:50
virtual void GeneratePrimaries(G4Event *anEvent)
Definition: MedicalBeam.cc:111
G4double fieldXY[2]
Definition: MedicalBeam.hh:56
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:276
constexpr auto abs(T v)
Returns the absolute value of the argument.
G4ThreeVector sourcePosition
Definition: MedicalBeam.hh:52
constexpr T sqr(T v)
G4double SSD
Definition: MedicalBeam.hh:54
FieldShape fieldShape
Definition: MedicalBeam.hh:55
G4double kineticE
Definition: MedicalBeam.hh:51
G4ThreeVector GenerateBeamDirection() const
Definition: MedicalBeam.cc:69
G4double fieldR
Definition: MedicalBeam.hh:57
vertex reconstruction