00001
00002 #ifndef HEPMC_HEAVY_ION_H
00003 #define HEPMC_HEAVY_ION_H
00004
00006
00007
00008
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00029
00030
00031
00033
00034 namespace HepMC {
00035
00036
00038
00045 class HeavyIon {
00046
00047 public:
00048
00049
00051 HeavyIon()
00052 : m_Ncoll_hard(0),
00053 m_Npart_proj(0),
00054 m_Npart_targ(0),
00055 m_Ncoll(0),
00056 m_spectator_neutrons(0),
00057 m_spectator_protons(0),
00058 m_N_Nwounded_collisions(0),
00059 m_Nwounded_N_collisions(0),
00060 m_Nwounded_Nwounded_collisions(0),
00061 m_impact_parameter(0),
00062 m_event_plane_angle(0),
00063 m_eccentricity(0),
00064 m_sigma_inel_NN(0)
00065 {}
00066
00068 HeavyIon( int nh, int np, int nt, int nc, int ns, int nsp,
00069 int nnw=0, int nwn=0, int nwnw=0,
00070 float im=0., float pl=0., float ec=0., float s=0. );
00071
00072 ~HeavyIon() {}
00073
00074
00075
00076 HeavyIon( HeavyIon const & orig );
00077 HeavyIon & operator = ( HeavyIon const & rhs );
00078 void swap( HeavyIon & other );
00079
00080
00081
00082 bool operator==( const HeavyIon& ) const;
00083 bool operator!=( const HeavyIon& ) const;
00084
00085
00087 int Ncoll_hard() const { return m_Ncoll_hard; }
00089 int Npart_proj() const { return m_Npart_proj; }
00091 int Npart_targ() const { return m_Npart_targ; }
00093 int Ncoll() const { return m_Ncoll; }
00095 int spectator_neutrons() const { return m_spectator_neutrons; }
00097 int spectator_protons() const { return m_spectator_protons; }
00099 int N_Nwounded_collisions() const { return m_N_Nwounded_collisions; }
00101 int Nwounded_N_collisions() const { return m_Nwounded_N_collisions; }
00103 int Nwounded_Nwounded_collisions() const { return m_Nwounded_Nwounded_collisions; }
00105 float impact_parameter() const { return m_impact_parameter; }
00107 float event_plane_angle() const { return m_event_plane_angle; }
00110 float eccentricity() const { return m_eccentricity; }
00112 float sigma_inel_NN() const { return m_sigma_inel_NN; }
00113
00114
00116 void set_Ncoll_hard(const int &i) { m_Ncoll_hard=i; }
00118 void set_Npart_proj(const int &i) { m_Npart_proj=i; }
00120 void set_Npart_targ(const int &i) { m_Npart_targ=i; }
00122 void set_Ncoll(const int &i) { m_Ncoll=i; }
00124 void set_spectator_neutrons(const int &i) { m_spectator_neutrons=i; }
00126 void set_spectator_protons(const int &i) { m_spectator_protons=i; }
00128 void set_N_Nwounded_collisions(const int &i) { m_N_Nwounded_collisions=i; }
00130 void set_Nwounded_N_collisions(const int &i) { m_Nwounded_N_collisions=i; }
00132 void set_Nwounded_Nwounded_collisions(const int &i)
00133 { m_Nwounded_Nwounded_collisions=i; }
00135 void set_impact_parameter(const float &f) { m_impact_parameter=f; }
00137 void set_event_plane_angle(const float &f) { m_event_plane_angle=f; }
00139 void set_eccentricity(const float &f) { m_eccentricity=f; }
00141 void set_sigma_inel_NN(const float &f) { m_sigma_inel_NN=f; }
00142
00143 private:
00144 int m_Ncoll_hard;
00145 int m_Npart_proj;
00146 int m_Npart_targ;
00147 int m_Ncoll;
00148 int m_spectator_neutrons;
00149 int m_spectator_protons;
00150 int m_N_Nwounded_collisions;
00151 int m_Nwounded_N_collisions;
00152 int m_Nwounded_Nwounded_collisions;
00153 float m_impact_parameter;
00154 float m_event_plane_angle;
00155 float m_eccentricity;
00156 float m_sigma_inel_NN;
00157
00158 };
00159
00160
00168 inline HeavyIon::HeavyIon( int nh, int np, int nt, int nc, int ns, int nsp,
00169 int nnw, int nwn, int nwnw,
00170 float im, float pl, float ec, float s )
00171 : m_Ncoll_hard(nh),
00172 m_Npart_proj(np),
00173 m_Npart_targ(nt),
00174 m_Ncoll(nc),
00175 m_spectator_neutrons(ns),
00176 m_spectator_protons(nsp),
00177 m_N_Nwounded_collisions(nnw),
00178 m_Nwounded_N_collisions(nwn),
00179 m_Nwounded_Nwounded_collisions(nwnw),
00180 m_impact_parameter(im),
00181 m_event_plane_angle(pl),
00182 m_eccentricity(ec),
00183 m_sigma_inel_NN(s)
00184 {}
00185
00186 inline HeavyIon::HeavyIon( HeavyIon const & orig )
00187 : m_Ncoll_hard(orig.m_Ncoll_hard),
00188 m_Npart_proj(orig.m_Npart_proj),
00189 m_Npart_targ(orig.m_Npart_targ),
00190 m_Ncoll(orig.m_Ncoll),
00191 m_spectator_neutrons(orig.m_spectator_neutrons),
00192 m_spectator_protons(orig.m_spectator_protons),
00193 m_N_Nwounded_collisions(orig.m_N_Nwounded_collisions),
00194 m_Nwounded_N_collisions(orig.m_Nwounded_N_collisions),
00195 m_Nwounded_Nwounded_collisions(orig.m_Nwounded_Nwounded_collisions),
00196 m_impact_parameter(orig.m_impact_parameter),
00197 m_event_plane_angle(orig.m_event_plane_angle),
00198 m_eccentricity(orig.m_eccentricity),
00199 m_sigma_inel_NN(orig.m_sigma_inel_NN)
00200 {}
00201
00202 inline HeavyIon & HeavyIon::operator = ( HeavyIon const & rhs )
00203 {
00204 HeavyIon temp( rhs );
00205 swap( temp );
00206 return *this;
00207 }
00208
00209 inline void HeavyIon::swap( HeavyIon & other )
00210 {
00211 std::swap(m_Ncoll_hard, other.m_Ncoll_hard);
00212 std::swap(m_Npart_proj, other.m_Npart_proj);
00213 std::swap(m_Npart_targ, other.m_Npart_targ);
00214 std::swap(m_Ncoll, other.m_Ncoll);
00215 std::swap(m_N_Nwounded_collisions, other.m_N_Nwounded_collisions);
00216 std::swap(m_Nwounded_N_collisions, other.m_Nwounded_N_collisions);
00217 std::swap(m_Nwounded_Nwounded_collisions, other.m_Nwounded_Nwounded_collisions);
00218 std::swap(m_spectator_neutrons, other.m_spectator_neutrons);
00219 std::swap(m_spectator_protons, other.m_spectator_protons);
00220 std::swap(m_impact_parameter, other.m_impact_parameter);
00221 std::swap(m_event_plane_angle, other.m_event_plane_angle);
00222 std::swap(m_eccentricity, other.m_eccentricity);
00223 std::swap(m_sigma_inel_NN, other.m_sigma_inel_NN);
00224 }
00225
00226 inline bool HeavyIon::operator==( const HeavyIon& a ) const
00227 {
00229 return ( a.Ncoll_hard() == this->Ncoll_hard()
00230 && a.Npart_proj() == this->Npart_proj()
00231 && a.Npart_targ() == this->Npart_targ()
00232 && a.Ncoll() == this->Ncoll()
00233 && a.N_Nwounded_collisions() == this->N_Nwounded_collisions()
00234 && a.Nwounded_N_collisions() == this->Nwounded_N_collisions()
00235 && a.Nwounded_Nwounded_collisions() == this->Nwounded_Nwounded_collisions()
00236 && a.spectator_neutrons() == this->spectator_neutrons()
00237 && a.spectator_protons() == this->spectator_protons()
00238 && a.impact_parameter() == this->impact_parameter()
00239 && a.event_plane_angle() == this->event_plane_angle()
00240 && a.eccentricity() == this->eccentricity()
00241 && a.sigma_inel_NN() == this->sigma_inel_NN() );
00242 }
00243
00244 inline bool HeavyIon::operator!=( const HeavyIon& a ) const
00245 {
00247 return !( a == *this );
00248 }
00249
00250 }
00251
00252 #endif // HEPMC_HEAVY_ION_H