Version: SMASH-1.5
smash::DeformedNucleus Class Reference

#include <deformednucleus.h>

DeformedNucleus: Child of nucleus for deformed nuclei.

All options from the nucleus will still apply. The deformed nucleus adds new or updated features which are outlined below.

Definition at line 27 of file deformednucleus.h.

Inheritance diagram for smash::DeformedNucleus:
[legend]
Collaboration diagram for smash::DeformedNucleus:
[legend]

Public Member Functions

 DeformedNucleus (const std::map< PdgCode, int > &particle_list, int nTest)
 Constructor for DeformedNucles which takes a particle list and the number of testparticles. More...
 
 DeformedNucleus (Configuration &config, int nTest)
 Constructor for DeformedNucleus, that needs the configuration parameters from the inputfile and the number of testparticles. More...
 
double deformed_woods_saxon (double r, double cosx) const
 Return the deformed Woods-Saxon probability for the given position. More...
 
ThreeVector distribute_nucleon () const override
 Deformed Woods-Saxon sampling routine. More...
 
void set_deformation_parameters_automatic ()
 Sets the deformation parameters of the radius according to the current mass number. More...
 
void set_deformation_parameters_from_config (Configuration &config)
 Set parameters for spherical deformation of the nucleus from the values specified in the configuration file. More...
 
void rotate () override
 Rotates the nucleus according to members nucleus_polar_angle_ and nucleus_azimuthal_angle_ and updates nucleon positions. More...
 
void generate_fermi_momenta () override
 Does not allow to generate Fermi-momenta for a deformed nucleus. More...
 
double y_l_0 (int l, double cosx) const
 Spherical harmonics Y_2_0 and Y_4_0. More...
 
void set_beta_2 (double b2)
 Set deformation coefficient for Y_2_0. More...
 
void set_beta_4 (double b4)
 Set deformation coefficient for Y_4_0. More...
 
void set_polar_angle (double theta)
 Set the nucleus polar angle. More...
 
void set_azimuthal_angle (double phi)
 Set the nucleus azimuthal angle. More...
 
- Public Member Functions inherited from smash::Nucleus
 Nucleus (Configuration &config, int nTest)
 Constructor for Nucleus, that needs the configuration parameters from the inputfile and the number of testparticles. More...
 
 Nucleus (const std::map< PdgCode, int > &particle_list, int nTest)
 
virtual ~Nucleus ()=default
 
double mass () const
 
double woods_saxon (double x)
 Woods-Saxon distribution. More...
 
void arrange_nucleons ()
 Sets the positions of the nucleons inside a nucleus. More...
 
virtual void set_parameters_automatic ()
 Sets the deformation parameters of the Woods-Saxon distribution according to the current mass number. More...
 
virtual void set_parameters_from_config (Configuration &config)
 Sets the parameters of the Woods-Saxon according to manually added values in the configuration file. More...
 
void boost (double beta_scalar)
 Boosts the nuclei into the computational frame, such that the nucleons have the appropriate momentum and the nuclei are lorentz-contracted. More...
 
void fill_from_list (const std::map< PdgCode, int > &particle_list, int testparticles)
 Adds particles from a map PDG code => Number_of_particles_with_that_PDG_code to the nucleus. More...
 
void shift (double z_offset, double x_offset, double simulation_time)
 Shifts the nucleus to correct impact parameter and z displacement. More...
 
void copy_particles (Particles *particles)
 Copies the particles from this nucleus into the particle list. More...
 
size_t size () const
 Number of numerical (=test-)particles in the nucleus: More...
 
size_t number_of_particles () const
 Number of physical particles in the nucleus: More...
 
FourVector center () const
 Calculate geometrical center of the nucleus. More...
 
void align_center ()
 Shifts the nucleus so that its center is at (0,0,0) More...
 
std::vector< ParticleData >::iterator begin ()
 For iterators over the particle list: More...
 
std::vector< ParticleData >::iterator end ()
 For iterators over the particle list: More...
 
std::vector< ParticleData >::const_iterator cbegin () const
 For const iterators over the particle list: More...
 
std::vector< ParticleData >::const_iterator cend () const
 For const iterators over the particle list: More...
 
void set_diffusiveness (double diffuse)
 Sets the diffusiveness of the nucleus. More...
 
double get_diffusiveness () const
 
void set_saturation_density (double density)
 Sets the saturation density of the nucleus. More...
 
double get_saturation_density () const
 
double default_nuclear_radius ()
 
void set_nuclear_radius (double rad)
 Sets the nuclear radius. More...
 
double get_nuclear_radius () const
 

Private Attributes

double beta2_ = 0.0
 Deformation parameter for angular momentum l=2. More...
 
double beta4_ = 0.0
 Deformation parameter for angular momentum l=4. More...
 
Angles nuclear_orientation_
 Nucleus orientation (initial profile in xz plane) in terms of a pair of angles (theta, phi) More...
 

Constructor & Destructor Documentation

◆ DeformedNucleus() [1/2]

smash::DeformedNucleus::DeformedNucleus ( const std::map< PdgCode, int > &  particle_list,
int  nTest 
)

Constructor for DeformedNucles which takes a particle list and the number of testparticles.

This constructor is only used for testing purposes.

Parameters
[in]particle_listMap with PDGCode and number of particles which make up the nucleus
[in]nTestnumber of testparticles

Definition at line 74 of file deformednucleus.cc.

76  : Nucleus(particle_list, nTest) {}
Nucleus(Configuration &config, int nTest)
Constructor for Nucleus, that needs the configuration parameters from the inputfile and the number of...
Definition: nucleus.cc:32

◆ DeformedNucleus() [2/2]

smash::DeformedNucleus::DeformedNucleus ( Configuration config,
int  nTest 
)

Constructor for DeformedNucleus, that needs the configuration parameters from the inputfile and the number of testparticles.

Parameters
[in]configcontains the parameters from the inputfile on the numbers of particles with a certain PDG code
[in]nTestnumber of testparticles

Definition at line 78 of file deformednucleus.cc.

79  : Nucleus(config, nTest) {
80  if (config.has_value({"Deformed", "Automatic"}) &&
81  config.take({"Deformed", "Automatic"})) {
83  } else {
85  }
86 }
void set_deformation_parameters_from_config(Configuration &config)
Set parameters for spherical deformation of the nucleus from the values specified in the configuratio...
Nucleus(Configuration &config, int nTest)
Constructor for Nucleus, that needs the configuration parameters from the inputfile and the number of...
Definition: nucleus.cc:32
void set_deformation_parameters_automatic()
Sets the deformation parameters of the radius according to the current mass number.
Here is the call graph for this function:

Member Function Documentation

◆ deformed_woods_saxon()

double smash::DeformedNucleus::deformed_woods_saxon ( double  r,
double  cosx 
) const

Return the deformed Woods-Saxon probability for the given position.

Parameters
[in]rThe radius at which to sample
[in]cosxThe cosine of the polar angle at which to sample
Returns
The Woods-Saxon probability

Definition at line 88 of file deformednucleus.cc.

88  {
90  (1 + std::exp(r - Nucleus::get_nuclear_radius() *
91  (1 + beta2_ * y_l_0(2, cosx) +
92  beta4_ * y_l_0(4, cosx)) /
94 }
double y_l_0(int l, double cosx) const
Spherical harmonics Y_2_0 and Y_4_0.
double beta4_
Deformation parameter for angular momentum l=4.
double beta2_
Deformation parameter for angular momentum l=2.
double get_nuclear_radius() const
Definition: nucleus.h:261
double get_diffusiveness() const
Definition: nucleus.h:231
double get_saturation_density() const
Definition: nucleus.h:243
Here is the call graph for this function:
Here is the caller graph for this function:

◆ distribute_nucleon()

ThreeVector smash::DeformedNucleus::distribute_nucleon ( ) const
overridevirtual

Deformed Woods-Saxon sampling routine.

Returns
Spatial position from uniformly sampling the deformed woods-saxon distribution

Reimplemented from smash::Nucleus.

Definition at line 96 of file deformednucleus.cc.

96  {
97  double a_radius;
98  Angles a_direction;
99  // Set a sensible maximum bound for radial sampling.
100  double radius_max =
103 
104  // Sample the distribution.
105  do {
106  a_direction.distribute_isotropically();
107  // sample r**2 dr
108  a_radius = radius_max * std::cbrt(random::canonical());
109  } while (random::canonical() >
110  deformed_woods_saxon(a_radius, a_direction.costheta()) /
112 
113  // Update (x, y, z) positions.
114  return a_direction.threevec() * a_radius;
115 }
double get_nuclear_radius() const
Definition: nucleus.h:261
T canonical()
Definition: random.h:110
double get_diffusiveness() const
Definition: nucleus.h:231
double get_saturation_density() const
Definition: nucleus.h:243
double deformed_woods_saxon(double r, double cosx) const
Return the deformed Woods-Saxon probability for the given position.
Here is the call graph for this function:

◆ set_deformation_parameters_automatic()

void smash::DeformedNucleus::set_deformation_parameters_automatic ( )

Sets the deformation parameters of the radius according to the current mass number.

The deformation parameters are taken from Moller:1993ed. Corrections to the deformation parameter beta2 in Uranium come from Kuhlman:2005ts. For finite nucleon size corrections to the nuclear density and radius for copper and gold, see Hirano:2009ah, and Hirano:2010jg for uranium.

Definition at line 117 of file deformednucleus.cc.

117  {
118  // Set the deformation parameters extracted from \iref{Moller:1993ed}.
119  switch (Nucleus::number_of_particles()) {
120  case 238: // Uranium
121  set_beta_2(0.28);
122  set_beta_4(0.093);
123  break;
124  case 208: // Lead
125  set_beta_2(0.0);
126  set_beta_4(0.0);
127  break;
128  case 197: // Gold
129  set_beta_2(-0.131);
130  set_beta_4(-0.031);
131  break;
132  case 63: // Copper
133  set_beta_2(0.162);
134  set_beta_4(-0.006);
135  break;
136  default:
137  throw std::domain_error(
138  "Mass number not listed for automatically setting deformation "
139  "parameters. Please specify at least \"Beta_2\" and \"Beta_4\" "
140  "manually and set \"Automatic: False.\" ");
141  }
142 
143  // Set a random nuclear rotation.
145 }
void set_beta_4(double b4)
Set deformation coefficient for Y_4_0.
Angles nuclear_orientation_
Nucleus orientation (initial profile in xz plane) in terms of a pair of angles (theta, phi)
void distribute_isotropically()
Populate the object with a new direction.
Definition: angles.h:188
void set_beta_2(double b2)
Set deformation coefficient for Y_2_0.
size_t number_of_particles() const
Number of physical particles in the nucleus:
Definition: nucleus.h:152
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_deformation_parameters_from_config()

void smash::DeformedNucleus::set_deformation_parameters_from_config ( Configuration config)

Set parameters for spherical deformation of the nucleus from the values specified in the configuration file.

Parameters
configThe configuration for the deformation of this nucleus (projectile or target).

Definition at line 147 of file deformednucleus.cc.

148  {
149  // Deformation parameters.
150  if (config.has_value({"Deformed", "Beta_2"})) {
151  set_beta_2(static_cast<double>(config.take({"Deformed", "Beta_2"})));
152  }
153  if (config.has_value({"Deformed", "Beta_4"})) {
154  set_beta_4(static_cast<double>(config.take({"Deformed", "Beta_4"})));
155  }
156  if (config.has_value({"Deformed", "Theta"})) {
157  set_polar_angle(static_cast<double>(config.take({"Deformed", "Theta"})));
158  }
159  if (config.has_value({"Deformed", "Phi"})) {
160  set_azimuthal_angle(static_cast<double>(config.take({"Deformed", "Phi"})));
161  }
162 }
void set_polar_angle(double theta)
Set the nucleus polar angle.
void set_beta_4(double b4)
Set deformation coefficient for Y_4_0.
void set_azimuthal_angle(double phi)
Set the nucleus azimuthal angle.
void set_beta_2(double b2)
Set deformation coefficient for Y_2_0.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ rotate()

void smash::DeformedNucleus::rotate ( )
overridevirtual

Rotates the nucleus according to members nucleus_polar_angle_ and nucleus_azimuthal_angle_ and updates nucleon positions.

Reimplemented from smash::Nucleus.

Definition at line 164 of file deformednucleus.cc.

164  {
165  for (auto &particle : *this) {
166  /* Rotate every vector by the nuclear azimuth phi and polar angle
167  * theta (the Euler angles). This means applying the matrix for a
168  * rotation of phi about z, followed by the matrix for a rotation
169  * theta about the rotated x axis. The third angle psi is 0 by symmetry.*/
170  ThreeVector three_pos = particle.position().threevec();
171  three_pos.rotate(nuclear_orientation_.phi(), nuclear_orientation_.theta(),
172  0.);
173  particle.set_3position(three_pos);
174  }
175 }
double theta() const
Definition: angles.h:272
Angles nuclear_orientation_
Nucleus orientation (initial profile in xz plane) in terms of a pair of angles (theta, phi)
double phi() const
Definition: angles.h:260
Here is the call graph for this function:

◆ generate_fermi_momenta()

void smash::DeformedNucleus::generate_fermi_momenta ( )
overridevirtual

Does not allow to generate Fermi-momenta for a deformed nucleus.

Exceptions
domain_errorif this function is ever called

Reimplemented from smash::Nucleus.

Definition at line 177 of file deformednucleus.cc.

177  {
178  throw std::domain_error(
179  "Fermi momenta currently not implemented"
180  " for a deformed nucleus.");
181 }

◆ y_l_0()

double smash::DeformedNucleus::y_l_0 ( int  l,
double  cosx 
) const

Spherical harmonics Y_2_0 and Y_4_0.

Parameters
[in]lAngular momentum value (2 and 4 are supported)
[in]cosxCosine of the polar angle
Returns
Value of the corresponding spherical harmonic
Exceptions
domain_errorif unsupported l is encountered

Definition at line 183 of file deformednucleus.cc.

183  {
184  if (l == 2) {
185  return (1. / 4) * std::sqrt(5 / M_PI) * (3. * (cosx * cosx) - 1);
186  } else if (l == 4) {
187  return (3. / 16) * std::sqrt(1 / M_PI) *
188  (35. * (cosx * cosx) * (cosx * cosx) - 30. * (cosx * cosx) + 3);
189  } else {
190  throw std::domain_error(
191  "Not a valid angular momentum quantum number in"
192  "DeformedNucleus::y_l_0.");
193  }
194 }
Here is the caller graph for this function:

◆ set_beta_2()

void smash::DeformedNucleus::set_beta_2 ( double  b2)
inline

Set deformation coefficient for Y_2_0.

Parameters
[in]b2deformation coefficient for l=2

Definition at line 108 of file deformednucleus.h.

108 { beta2_ = b2; }
double beta2_
Deformation parameter for angular momentum l=2.
Here is the caller graph for this function:

◆ set_beta_4()

void smash::DeformedNucleus::set_beta_4 ( double  b4)
inline

Set deformation coefficient for Y_4_0.

Parameters
[in]b4deformation coefficient for l=4

Definition at line 113 of file deformednucleus.h.

113 { beta4_ = b4; }
double beta4_
Deformation parameter for angular momentum l=4.
Here is the caller graph for this function:

◆ set_polar_angle()

void smash::DeformedNucleus::set_polar_angle ( double  theta)
inline

Set the nucleus polar angle.

Parameters
[in]thetaPolar angle of position inside nucleus

Definition at line 118 of file deformednucleus.h.

118  {
120  }
Angles nuclear_orientation_
Nucleus orientation (initial profile in xz plane) in terms of a pair of angles (theta, phi)
void set_theta(const double theta)
Set the polar angle.
Definition: angles.h:218
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_azimuthal_angle()

void smash::DeformedNucleus::set_azimuthal_angle ( double  phi)
inline

Set the nucleus azimuthal angle.

Parameters
[in]phiAzimuthal angle of position inside nucleus

Definition at line 125 of file deformednucleus.h.

125  {
127  }
Angles nuclear_orientation_
Nucleus orientation (initial profile in xz plane) in terms of a pair of angles (theta, phi)
void set_phi(const double phi)
Sets the azimuthal angle.
Definition: angles.h:194
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ beta2_

double smash::DeformedNucleus::beta2_ = 0.0
private

Deformation parameter for angular momentum l=2.

Definition at line 131 of file deformednucleus.h.

◆ beta4_

double smash::DeformedNucleus::beta4_ = 0.0
private

Deformation parameter for angular momentum l=4.

Definition at line 133 of file deformednucleus.h.

◆ nuclear_orientation_

Angles smash::DeformedNucleus::nuclear_orientation_
private

Nucleus orientation (initial profile in xz plane) in terms of a pair of angles (theta, phi)

Definition at line 138 of file deformednucleus.h.


The documentation for this class was generated from the following files: