Fork me on GitHub

source: git/external/Hector/H_Beam.h@ 5b822e5

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5b822e5 was 5b822e5, checked in by pavel <pavel@…>, 11 years ago

add Hector module

  • Property mode set to 100644
File size: 7.3 KB
Line 
1#ifndef _H_Beam_
2#define _H_Beam_
3
4 /* * * * * * * * * * * * * * * * * * * * * * * * * * * *
5 * *
6* --<--<-- A fast simulator --<--<-- *
7* / --<--<-- of particle --<--<-- *
8* ----HECTOR----< *
9* \ -->-->-- transport through -->-->-- *
10* -->-->-- generic beamlines -->-->-- *
11* *
12* JINST 2:P09005 (2007) *
13* X Rouby, J de Favereau, K Piotrzkowski (CP3) *
14* http://www.fynu.ucl.ac.be/hector.html *
15* *
16* Center for Cosmology, Particle Physics and Phenomenology *
17* Universite catholique de Louvain *
18* Louvain-la-Neuve, Belgium *
19 * *
20 * * * * * * * * * * * * * * * * * * * * * * * * * * * */
21
22/// \file H_Beam.h
23/// \brief Describes a set a particles as a beam
24///
25
26// c++ #includes
27//#include <iostream>
28#include <vector>
29#include <sstream>
30#include <cmath>
31
32// ROOT #includes
33#include "TH2F.h"
34#include "TGraphErrors.h"
35#include "TMultiGraph.h"
36#include "TMath.h"
37
38// local #includes
39#include "H_BeamParticle.h"
40#include "H_AbstractBeamLine.h"
41#include "H_Parameters.h"
42#include "H_OpticalElement.h"
43using namespace std;
44
45/// \brief Beam made from a STL vector of H_BeamParticle.
46///
47/// Example of usage :
48/// <PRE>
49/// H_Beam mybeam;
50/// mybeam.createBeamParticles(100);
51/// H_AbstractBeamLine * beamline = new H_AbstractBeamLine(200);
52/// ... (<-- fills the beam)
53/// mybeam.computePath(beamline)
54/// </PRE>
55class H_Beam {
56
57 public:
58 /// Constructors, destructor and operator
59 //@{
60 H_Beam();
61 H_Beam(const H_Beam&);
62 H_Beam& operator=(const H_Beam&);
63 ~H_Beam();
64 //@}
65 /// Fills the beam with particles in given position/angle/energy intervals (flat distribution)
66 void particleGun(const unsigned int Number_of_particles, const float E_min=BE, const float E_max=BE, const float fs_min=0, const float fs_max=0, const float fx_min=0, const float fx_max=0, const float fy_min=0, const float fy_max=0, const float tx_min=-TMath::Pi()/2., const float tx_max=TMath::Pi()/2., const float ty_min=-TMath::Pi()/2., const float ty_max=TMath::Pi()/2., const float p_mass=MP, const double p_charge=QP, const bool flat = true, TRandom* r=gRandom);
67
68 /// Fills the beam with particles
69 void createBeamParticles(const unsigned int Number_of_particles, const double p_mass=MP, const double p_charge=QP, TRandom* r=gRandom);
70 /// Fills the beam with particles with incremental offset and no initial transverse momentum
71 //@{
72 void createXScanningBeamParticles(const unsigned int, const float);
73 void createYScanningBeamParticles(const unsigned int, const float);
74 //@}
75 /// Fills the beam with particles with incremental initial angle and no offset.
76 //@{
77 void createTXScanningBeamParticles(const unsigned int, const float);
78 void createTYScanningBeamParticles(const unsigned int, const float);
79 //@}
80 /// Returns the ith particle
81 //@{
82 const H_BeamParticle * getBeamParticle(const unsigned int ) const;
83 H_BeamParticle * getBeamParticle(const unsigned int );
84 //@}
85 /// Adds a particle to the beam
86 void add(const H_BeamParticle&);
87 /// Compute the position of each particle in the beamline
88 void computePath(const H_AbstractBeamLine * beamline, const bool NonLinear=false);
89 // Photon emission by the particle
90 void emitGamma(const double gee, const double gq2, const double phimin=0, const double phimax=2*TMath::Pi());
91 // Propagates the beam until a given s
92 void propagate(const float );
93 /// Returns the \f$ \beta \f$ function of the beam
94 //@{
95 float getBetaX(const float , float& );
96 float getBetaY(const float , float& );
97 //@}
98 /// Draws the \f$ \beta \f$ function of the beam
99 //@{
100 TGraphErrors * getBetaX(const float, const unsigned int);
101 TGraphErrors * getBetaY(const float, const unsigned int);
102 //@}
103 /// Returns the position of the beam
104 //@{
105 float getX(const float , float& );
106 float getY(const float , float& );
107 //@}
108 /// Returns the emittance \f$ \epsilon \f$ of the beam in x and y
109 //@{
110 inline const float getEmittanceX() const {
111 if(!x_disp*tx_disp) cout<<"Warning : Degenerate Beam : x-emittance = 0"<<endl;
112 return x_disp * tan(tx_disp/URAD)/URAD;
113 }
114 inline const float getEmittanceY() const {
115 if(!y_disp*ty_disp) cout<<"Warning : Degenerate Beam : y-emittance = 0"<<endl;
116 return y_disp * tan(ty_disp/URAD)/URAD;
117 }
118 //@}
119 /// Sets the initial parameters of the beam \f$ s , x , y , \theta_x , \theta_y \f$
120 //@{
121 void setS(const float fs) { fs_ini = fs;}
122 void setX(const float fx) { fx_ini = fx;}
123 void setY(const float fy) { fy_ini = fy;}
124 void setTX(const float tx) { tx_ini = tx;}
125 void setTY(const float ty) { ty_ini = ty;}
126 void setE(const float fe) {fe_ini = fe;}
127 void setPosition(const float fx, const float fy, const float tx, const float ty, const float fs) {setS(fs); setX(fx); setY(fy); setTX(tx); setTY(ty);}
128 //@}
129 /// Sets the initial divergence of the beam \f$ \sigma_x , \sigma_y , \sigma (\theta_x), \sigma(\theta_y) \f$
130 //@{
131 void setDS(const float ds) {s_disp = ds;}
132 void setDX(const float dx) {x_disp = dx;}
133 void setDY(const float dy) {y_disp = dy;}
134 void setDTX(const float dtx) {tx_disp = dtx;}
135 void setDTY(const float dty) {ty_disp = dty;}
136 void setDE(const float de) {e_disp = de;}
137 void setDispersion(const float dx, const float dy, const float dtx, const float dty, const float ds) {setDS(ds); setDX(dx); setDY(dy); setDTX(dtx); setDTY(dty);}
138 //@}
139 /// Returns the number of particles which have been stopped
140 unsigned int getStoppedNumber(const H_AbstractBeamLine *);
141 /// Returns the list of the stopping elements in the beamline
142 vector<TVectorD> getStoppingElements(const H_AbstractBeamLine *, vector<H_OpticalElement>&, vector<int>&);
143 /// Prints the initial parameters
144 void printInitialState() const;
145 /// Prints the properties for each particle
146 void printProperties() const {cout << *this; return;}
147 /// Prints the list of the stopping elements in the beamline
148 void printStoppingElements(const vector<H_OpticalElement>&, const vector<int>&) const;
149 /// Returns the number of particle is this beam
150 const int getNumberOfBeamParticles() const {return Nparticles;}
151 /// Draws the beam profile at a given s
152 TH2F * drawProfile(const float);
153 TH2F * drawAngleProfile(const float);
154 /// Draws the beam width and height
155 //@{
156 TMultiGraph * drawBeamX(const int) const;
157 TMultiGraph * drawBeamY(const int) const ;
158 //@}
159
160 protected:
161 /// List of particles
162 vector<H_BeamParticle> beamParticles;
163 /// IP position
164 //@{
165 float fx_ini, fy_ini, fs_ini;
166 //@}
167 /// IP \f$ \theta \f$ (angular) position
168 //@{
169 float tx_ini, ty_ini;
170 //@}
171 /// nominal beam energy
172 float fe_ini;
173 /// dispersion in position
174 //@{
175 float x_disp, y_disp, s_disp;
176 //@}
177 /// dispersion in \f$ \theta \f$
178 //@{
179 float tx_disp, ty_disp;
180 //@}
181 /// dispersion in energy
182 float e_disp;
183 /// Number of particles in this beam
184 unsigned int Nparticles;
185
186 friend std::ostream& operator<< (std::ostream& os, const H_Beam& be);
187};
188
189#endif
Note: See TracBrowser for help on using the repository browser.