Version: SMASH-1.8
smash::RectangularLattice< T > Class Template Reference

#include <lattice.h>

template<typename T>
class smash::RectangularLattice< T >

A container class to hold all the arrays on the lattice and access them.

Template Parameters
TThe type of the contained values.

Definition at line 47 of file lattice.h.

Collaboration diagram for smash::RectangularLattice< T >:
[legend]

Public Types

using iterator = typename std::vector< T >::iterator
 Iterator of lattice. More...
 
using const_iterator = typename std::vector< T >::const_iterator
 Const interator of lattice. More...
 

Public Member Functions

 RectangularLattice (const std::array< double, 3 > &l, const std::array< int, 3 > &n, const std::array< double, 3 > &orig, bool per, const LatticeUpdate upd)
 Rectangular lattice constructor. More...
 
void reset ()
 Sets all values on lattice to zeros. More...
 
bool out_of_bounds (int ix, int iy, int iz) const
 Checks if 3D index is out of lattice bounds. More...
 
ThreeVector cell_center (int ix, int iy, int iz) const
 Find the coordinates of a given cell. More...
 
ThreeVector cell_center (int index) const
 Find the coordinate of cell center given the 1d index of the cell. More...
 
const std::array< double, 3 > & lattice_sizes () const
 
const std::array< int, 3 > & dimensions () const
 
const std::array< double, 3 > & cell_sizes () const
 
const std::array< double, 3 > & origin () const
 
bool periodic () const
 
LatticeUpdate when_update () const
 
iterator begin ()
 
const_iterator begin () const
 
iterator end ()
 
const_iterator end () const
 
T & operator[] (std::size_t i)
 
const T & operator[] (std::size_t i) const
 
std::size_t size () const
 
T & node (int ix, int iy, int iz)
 Take the value of a cell given its 3-D indices. More...
 
bool value_at (const ThreeVector &r, T &value)
 Interpolates lattice quantity to coordinate r. More...
 
template<typename F >
void iterate_sublattice (const std::array< int, 3 > &lower_bounds, const std::array< int, 3 > &upper_bounds, F &&func)
 A sub-lattice iterator, which iterates in a 3D-structured manner and calls a function on every cell. More...
 
template<typename F >
void iterate_in_radius (const ThreeVector &point, const double r_cut, F &&func)
 Iterates only nodes, whose cell centers lie not further than r_cut in x, y, z directions from the given point and applies a function to each node. More...
 
template<typename L >
bool identical_to_lattice (const L *lat) const
 Checks if lattices of possibly different types have identical structure. More...
 

Protected Attributes

std::vector< T > lattice_
 The lattice itself, array containing physical quantities. More...
 
const std::array< double, 3 > lattice_sizes_
 Lattice sizes in x, y, z directions. More...
 
const std::array< int, 3 > n_cells_
 Number of cells in x,y,z directions. More...
 
const std::array< double, 3 > cell_sizes_
 Cell sizes in x, y, z directions. More...
 
const std::array< double, 3 > origin_
 Coordinates of the left down nearer corner. More...
 
const bool periodic_
 Whether the lattice is periodic. More...
 
const LatticeUpdate when_update_
 When the lattice should be recalculated. More...
 

Private Member Functions

int positive_modulo (int i, int n) const
 Returns division modulo, which is always between 0 and n-1 in is not suitable, because it returns results from -(n-1) to n-1. More...
 

Member Typedef Documentation

◆ iterator

template<typename T>
using smash::RectangularLattice< T >::iterator = typename std::vector<T>::iterator

Iterator of lattice.

Definition at line 160 of file lattice.h.

◆ const_iterator

template<typename T>
using smash::RectangularLattice< T >::const_iterator = typename std::vector<T>::const_iterator

Const interator of lattice.

Definition at line 162 of file lattice.h.

Constructor & Destructor Documentation

◆ RectangularLattice()

template<typename T>
smash::RectangularLattice< T >::RectangularLattice ( const std::array< double, 3 > &  l,
const std::array< int, 3 > &  n,
const std::array< double, 3 > &  orig,
bool  per,
const LatticeUpdate  upd 
)
inline

Rectangular lattice constructor.

Parameters
[in]l3-dimensional array (lx,ly,lz) indicates the size of the lattice in x, y, z directions respectively [fm].
[in]n3-dimensional array (nx,ny,nz) indicates the number of cells of the lattice in x, y, z directions respectively. Each cell in the lattice is labeled by three integers i, j, k where \(i\in[0, nx-1]\), \(j\in[0, ny-1]\), \(k\in[0, nz-1]\). The sizes of each cell are given by lx/nx, ly/ny, lz/nz in x,y,z directions respectively.
[in]origA 3-dimensional array indicating the coordinates of the origin [fm].
[in]perBoolean indicating whether a periodic boundary condition is applied.
[in]updEnum indicating how frequently the lattice is updated.

Definition at line 66 of file lattice.h.

70  : lattice_sizes_(l),
71  n_cells_(n),
72  cell_sizes_{l[0] / n[0], l[1] / n[1], l[2] / n[2]},
73  origin_(orig),
74  periodic_(per),
75  when_update_(upd) {
76  lattice_.resize(n_cells_[0] * n_cells_[1] * n_cells_[2]);
77  logg[LLattice].debug(
78  "Rectangular lattice created: sizes[fm] = (", lattice_sizes_[0], ",",
79  lattice_sizes_[1], ",", lattice_sizes_[2], "), dims = (", n_cells_[0],
80  ",", n_cells_[1], ",", n_cells_[2], "), origin = (", origin_[0], ",",
81  origin_[1], ",", origin_[2], "), periodic: ", periodic_);
82  if (n_cells_[0] < 1 || n_cells_[1] < 1 || n_cells_[2] < 1 ||
83  lattice_sizes_[0] < 0.0 || lattice_sizes_[1] < 0.0 ||
84  lattice_sizes_[2] < 0.0) {
85  throw std::invalid_argument(
86  "Lattice sizes should be positive, "
87  "lattice dimensions should be > 0.");
88  }
89  }

Member Function Documentation

◆ reset()

template<typename T>
void smash::RectangularLattice< T >::reset ( )
inline

Sets all values on lattice to zeros.

Definition at line 92 of file lattice.h.

92 { std::fill(lattice_.begin(), lattice_.end(), T()); }
Here is the caller graph for this function:

◆ out_of_bounds()

template<typename T>
bool smash::RectangularLattice< T >::out_of_bounds ( int  ix,
int  iy,
int  iz 
) const
inline

Checks if 3D index is out of lattice bounds.

Parameters
[in]ixThe index of the cell in x direction.
[in]iyThe index of the cell in y direction.
[in]izThe index of the cell in z direction.
Returns
Whether the cell is out of the lattice.

Definition at line 102 of file lattice.h.

102  {
103  // clang-format off
104  return !periodic_ &&
105  (ix < 0 || ix >= n_cells_[0] ||
106  iy < 0 || iy >= n_cells_[1] ||
107  iz < 0 || iz >= n_cells_[2]);
108  // clang-format on
109  }
Here is the caller graph for this function:

◆ cell_center() [1/2]

template<typename T>
ThreeVector smash::RectangularLattice< T >::cell_center ( int  ix,
int  iy,
int  iz 
) const
inline

Find the coordinates of a given cell.

Parameters
[in]ixThe index of the cell in x direction.
[in]iyThe index of the cell in y direction.
[in]izThe index of the cell in z direction.
Returns
Coordinates of the center of the given cell [fm].

Definition at line 119 of file lattice.h.

119  {
120  return ThreeVector(origin_[0] + cell_sizes_[0] * (ix + 0.5),
121  origin_[1] + cell_sizes_[1] * (iy + 0.5),
122  origin_[2] + cell_sizes_[2] * (iz + 0.5));
123  }
Here is the caller graph for this function:

◆ cell_center() [2/2]

template<typename T>
ThreeVector smash::RectangularLattice< T >::cell_center ( int  index) const
inline

Find the coordinate of cell center given the 1d index of the cell.

Parameters
[in]index1-dimensional index of the given cell. It can be related to a 3-dimensional one by index = ix + nx (iy + iz * ny).
Returns
Coordinates of the center of the given cell [fm].

Definition at line 133 of file lattice.h.

133  {
134  const int ix = index % n_cells_[0];
135  index = index / n_cells_[0];
136  const int iy = index % n_cells_[1];
137  const int iz = index / n_cells_[1];
138  return cell_center(ix, iy, iz);
139  }

◆ lattice_sizes()

template<typename T>
const std::array<double, 3>& smash::RectangularLattice< T >::lattice_sizes ( ) const
inline
Returns
Lengths of the lattice in x, y, z directions.

Definition at line 142 of file lattice.h.

142 { return lattice_sizes_; }

◆ dimensions()

template<typename T>
const std::array<int, 3>& smash::RectangularLattice< T >::dimensions ( ) const
inline
Returns
Number of cells in x, y, z directions.

Definition at line 145 of file lattice.h.

145 { return n_cells_; }
Here is the caller graph for this function:

◆ cell_sizes()

template<typename T>
const std::array<double, 3>& smash::RectangularLattice< T >::cell_sizes ( ) const
inline
Returns
Lengths of one cell in x, y, z directions.

Definition at line 148 of file lattice.h.

148 { return cell_sizes_; }
Here is the caller graph for this function:

◆ origin()

template<typename T>
const std::array<double, 3>& smash::RectangularLattice< T >::origin ( ) const
inline
Returns
Lattice origin: left, down, near corner coordinates.

Definition at line 151 of file lattice.h.

151 { return origin_; }
Here is the caller graph for this function:

◆ periodic()

template<typename T>
bool smash::RectangularLattice< T >::periodic ( ) const
inline
Returns
If lattice is periodic or not.

Definition at line 154 of file lattice.h.

154 { return periodic_; }

◆ when_update()

template<typename T>
LatticeUpdate smash::RectangularLattice< T >::when_update ( ) const
inline
Returns
The enum, which tells at which time lattice needs to be updated.

Definition at line 157 of file lattice.h.

157 { return when_update_; }
Here is the caller graph for this function:

◆ begin() [1/2]

template<typename T>
iterator smash::RectangularLattice< T >::begin ( )
inline
Returns
First element of lattice.

Definition at line 164 of file lattice.h.

164 { return lattice_.begin(); }

◆ begin() [2/2]

template<typename T>
const_iterator smash::RectangularLattice< T >::begin ( ) const
inline
Returns
First element of lattice (const).

Definition at line 166 of file lattice.h.

166 { return lattice_.begin(); }

◆ end() [1/2]

template<typename T>
iterator smash::RectangularLattice< T >::end ( )
inline
Returns
Last element of lattice.

Definition at line 168 of file lattice.h.

168 { return lattice_.end(); }

◆ end() [2/2]

template<typename T>
const_iterator smash::RectangularLattice< T >::end ( ) const
inline
Returns
Last element of lattice (const).

Definition at line 170 of file lattice.h.

170 { return lattice_.end(); }

◆ operator[]() [1/2]

template<typename T>
T& smash::RectangularLattice< T >::operator[] ( std::size_t  i)
inline
Returns
ith element of lattice.

Definition at line 172 of file lattice.h.

172 { return lattice_[i]; }

◆ operator[]() [2/2]

template<typename T>
const T& smash::RectangularLattice< T >::operator[] ( std::size_t  i) const
inline
Returns
ith element of lattice (const).

Definition at line 174 of file lattice.h.

174 { return lattice_[i]; }

◆ size()

template<typename T>
std::size_t smash::RectangularLattice< T >::size ( ) const
inline
Returns
Size of lattice.

Definition at line 176 of file lattice.h.

176 { return lattice_.size(); }
Here is the caller graph for this function:

◆ node()

template<typename T>
T& smash::RectangularLattice< T >::node ( int  ix,
int  iy,
int  iz 
)
inline

Take the value of a cell given its 3-D indices.

Parameters
[in]ixThe index of the cell in x direction.
[in]iyThe index of the cell in y direction.
[in]izThe index of the cell in z direction.
Returns
Physical quantity evaluated at the cell center.

Definition at line 186 of file lattice.h.

186  {
187  return periodic_
188  ? lattice_[positive_modulo(ix, n_cells_[0]) +
189  n_cells_[0] *
190  (positive_modulo(iy, n_cells_[1]) +
191  n_cells_[1] * positive_modulo(iz, n_cells_[2]))]
192  : lattice_[ix + n_cells_[0] * (iy + n_cells_[1] * iz)];
193  }
Here is the caller graph for this function:

◆ value_at()

template<typename T>
bool smash::RectangularLattice< T >::value_at ( const ThreeVector r,
T &  value 
)
inline

Interpolates lattice quantity to coordinate r.

Result is stored in the value variable. Returns true if coordinate r is on the lattice, false if out of the lattice. In the latter case, the value is set to the default value (usually 0).

Parameters
[in]rPosition where the physical quantity would be evaluated.
[out]valuePhysical quantity evaluated at the nearest cell to the given position.
Returns
Boolean indicates whether the position r is located inside the lattice.
Todo:
(oliiny): maybe 1-order interpolation instead of 0-order?

Definition at line 209 of file lattice.h.

209  {
210  const int ix = std::floor((r.x1() - origin_[0]) / cell_sizes_[0]);
211  const int iy = std::floor((r.x2() - origin_[1]) / cell_sizes_[1]);
212  const int iz = std::floor((r.x3() - origin_[2]) / cell_sizes_[2]);
213  if (out_of_bounds(ix, iy, iz)) {
214  value = T();
215  return false;
216  } else {
217  value = node(ix, iy, iz);
218  return true;
219  }
220  }

◆ iterate_sublattice()

template<typename T>
template<typename F >
void smash::RectangularLattice< T >::iterate_sublattice ( const std::array< int, 3 > &  lower_bounds,
const std::array< int, 3 > &  upper_bounds,
F &&  func 
)
inline

A sub-lattice iterator, which iterates in a 3D-structured manner and calls a function on every cell.

Template Parameters
FType of the function. Arguments are the current node and the 3 integer indices of the cell.
Parameters
[in]lower_boundsStarting numbers for iterating ix, iy, iz.
[in]upper_boundsEnding numbers for iterating ix, iy, iz.
[in]funcFunction acting on the cells (such as taking value).

Definition at line 233 of file lattice.h.

234  {
235  logg[LLattice].debug(
236  "Iterating sublattice with lower bound index (", lower_bounds[0], ",",
237  lower_bounds[1], ",", lower_bounds[2], "), upper bound index (",
238  upper_bounds[0], ",", upper_bounds[1], ",", upper_bounds[2], ")");
239 
240  if (periodic_) {
241  for (int iz = lower_bounds[2]; iz < upper_bounds[2]; iz++) {
242  const int z_offset = positive_modulo(iz, n_cells_[2]) * n_cells_[1];
243  for (int iy = lower_bounds[1]; iy < upper_bounds[1]; iy++) {
244  const int y_offset =
245  n_cells_[0] * (positive_modulo(iy, n_cells_[1]) + z_offset);
246  for (int ix = lower_bounds[0]; ix < upper_bounds[0]; ix++) {
247  const int index = positive_modulo(ix, n_cells_[0]) + y_offset;
248  func(lattice_[index], ix, iy, iz);
249  }
250  }
251  }
252  } else {
253  for (int iz = lower_bounds[2]; iz < upper_bounds[2]; iz++) {
254  const int z_offset = iz * n_cells_[1];
255  for (int iy = lower_bounds[1]; iy < upper_bounds[1]; iy++) {
256  const int y_offset = n_cells_[0] * (iy + z_offset);
257  for (int ix = lower_bounds[0]; ix < upper_bounds[0]; ix++) {
258  func(lattice_[ix + y_offset], ix, iy, iz);
259  }
260  }
261  }
262  }
263  }
Here is the caller graph for this function:

◆ iterate_in_radius()

template<typename T>
template<typename F >
void smash::RectangularLattice< T >::iterate_in_radius ( const ThreeVector point,
const double  r_cut,
F &&  func 
)
inline

Iterates only nodes, whose cell centers lie not further than r_cut in x, y, z directions from the given point and applies a function to each node.

Useful for adding quantities from one particle to the lattice.

Template Parameters
FType of the function. Arguments are the current node and the 3 integer indices of the cell.
Parameters
[in]pointPosition, usually the position of particle [fm].
[in]r_cutMaximum distance from the cell center to the given position. [fm]
[in]funcFunction acting on the cells (such as taking value).

Definition at line 278 of file lattice.h.

279  {
280  std::array<int, 3> l_bounds, u_bounds;
281 
282  /* Array holds value at the cell center: r_center = r_0 + (i+0.5)cell_size,
283  * where i is index in any direction. Therefore we want cells with condition
284  * (r-r_cut)*csize - 0.5 < i < (r+r_cut)*csize - 0.5, r = r_center - r_0 */
285  for (int i = 0; i < 3; i++) {
286  l_bounds[i] =
287  std::ceil((point[i] - origin_[i] - r_cut) / cell_sizes_[i] - 0.5);
288  u_bounds[i] =
289  std::ceil((point[i] - origin_[i] + r_cut) / cell_sizes_[i] - 0.5);
290  }
291 
292  if (!periodic_) {
293  for (int i = 0; i < 3; i++) {
294  if (l_bounds[i] < 0) {
295  l_bounds[i] = 0;
296  }
297  if (u_bounds[i] > n_cells_[i]) {
298  u_bounds[i] = n_cells_[i];
299  }
300  if (l_bounds[i] > n_cells_[i] || u_bounds[i] < 0) {
301  return;
302  }
303  }
304  }
305  iterate_sublattice(l_bounds, u_bounds, std::forward<F>(func));
306  }
Here is the caller graph for this function:

◆ identical_to_lattice()

template<typename T>
template<typename L >
bool smash::RectangularLattice< T >::identical_to_lattice ( const L *  lat) const
inline

Checks if lattices of possibly different types have identical structure.

Template Parameters
LType of the other lattice.
Parameters
[in]latThe other lattice being compared with the current one
Returns
Whether the two lattices have the same sizes, cell numbers, origins, and boundary conditions.

Definition at line 317 of file lattice.h.

317  {
318  return n_cells_[0] == lat->dimensions()[0] &&
319  n_cells_[1] == lat->dimensions()[1] &&
320  n_cells_[2] == lat->dimensions()[2] &&
321  std::abs(lattice_sizes_[0] - lat->lattice_sizes()[0]) <
322  really_small &&
323  std::abs(lattice_sizes_[1] - lat->lattice_sizes()[1]) <
324  really_small &&
325  std::abs(lattice_sizes_[2] - lat->lattice_sizes()[2]) <
326  really_small &&
327  std::abs(origin_[0] - lat->origin()[0]) < really_small &&
328  std::abs(origin_[1] - lat->origin()[1]) < really_small &&
329  std::abs(origin_[2] - lat->origin()[2]) < really_small &&
330  periodic_ == lat->periodic();
331  }

◆ positive_modulo()

template<typename T>
int smash::RectangularLattice< T >::positive_modulo ( int  i,
int  n 
) const
inlineprivate

Returns division modulo, which is always between 0 and n-1 in is not suitable, because it returns results from -(n-1) to n-1.

Parameters
[in]iDividend.
[in]nDivisor.
Returns
Positive remainder.

Definition at line 358 of file lattice.h.

358  {
359  /* (i % n + n) % n would be correct, but slow.
360  * Instead I rely on the fact that i should never go too far
361  * in negative region and replace i%n + n by i + 256 * n = i + (n << 8) */
362  // FIXME: This should use asserts, also checking for under- or overflows.
363  return (i + (n << 8)) % n;
364  }
Here is the caller graph for this function:

Member Data Documentation

◆ lattice_

template<typename T>
std::vector<T> smash::RectangularLattice< T >::lattice_
protected

The lattice itself, array containing physical quantities.

Definition at line 335 of file lattice.h.

◆ lattice_sizes_

template<typename T>
const std::array<double, 3> smash::RectangularLattice< T >::lattice_sizes_
protected

Lattice sizes in x, y, z directions.

Definition at line 337 of file lattice.h.

◆ n_cells_

template<typename T>
const std::array<int, 3> smash::RectangularLattice< T >::n_cells_
protected

Number of cells in x,y,z directions.

Definition at line 339 of file lattice.h.

◆ cell_sizes_

template<typename T>
const std::array<double, 3> smash::RectangularLattice< T >::cell_sizes_
protected

Cell sizes in x, y, z directions.

Definition at line 341 of file lattice.h.

◆ origin_

template<typename T>
const std::array<double, 3> smash::RectangularLattice< T >::origin_
protected

Coordinates of the left down nearer corner.

Definition at line 343 of file lattice.h.

◆ periodic_

template<typename T>
const bool smash::RectangularLattice< T >::periodic_
protected

Whether the lattice is periodic.

Definition at line 345 of file lattice.h.

◆ when_update_

template<typename T>
const LatticeUpdate smash::RectangularLattice< T >::when_update_
protected

When the lattice should be recalculated.

Definition at line 347 of file lattice.h.


The documentation for this class was generated from the following file:
smash::RectangularLattice::out_of_bounds
bool out_of_bounds(int ix, int iy, int iz) const
Checks if 3D index is out of lattice bounds.
Definition: lattice.h:102
smash::RectangularLattice::n_cells_
const std::array< int, 3 > n_cells_
Number of cells in x,y,z directions.
Definition: lattice.h:339
smash::RectangularLattice::positive_modulo
int positive_modulo(int i, int n) const
Returns division modulo, which is always between 0 and n-1 in is not suitable, because it returns res...
Definition: lattice.h:358
smash::RectangularLattice::lattice_sizes_
const std::array< double, 3 > lattice_sizes_
Lattice sizes in x, y, z directions.
Definition: lattice.h:337
smash::logg
std::array< einhard::Logger<>, std::tuple_size< LogArea::AreaTuple >::value > logg
An array that stores all pre-configured Logger objects.
Definition: logging.cc:39
smash::really_small
constexpr double really_small
Numerical error tolerance.
Definition: constants.h:37
smash::RectangularLattice::cell_center
ThreeVector cell_center(int ix, int iy, int iz) const
Find the coordinates of a given cell.
Definition: lattice.h:119
smash::RectangularLattice::cell_sizes_
const std::array< double, 3 > cell_sizes_
Cell sizes in x, y, z directions.
Definition: lattice.h:341
smash::RectangularLattice::when_update_
const LatticeUpdate when_update_
When the lattice should be recalculated.
Definition: lattice.h:347
smash::RectangularLattice::lattice_
std::vector< T > lattice_
The lattice itself, array containing physical quantities.
Definition: lattice.h:335
smash::RectangularLattice::node
T & node(int ix, int iy, int iz)
Take the value of a cell given its 3-D indices.
Definition: lattice.h:186
smash::pdg::n
constexpr int n
Neutron.
Definition: pdgcode_constants.h:30
smash::RectangularLattice::periodic_
const bool periodic_
Whether the lattice is periodic.
Definition: lattice.h:345
smash::LLattice
static constexpr int LLattice
Definition: lattice.h:25
smash::RectangularLattice::origin_
const std::array< double, 3 > origin_
Coordinates of the left down nearer corner.
Definition: lattice.h:343
smash::RectangularLattice::iterate_sublattice
void iterate_sublattice(const std::array< int, 3 > &lower_bounds, const std::array< int, 3 > &upper_bounds, F &&func)
A sub-lattice iterator, which iterates in a 3D-structured manner and calls a function on every cell.
Definition: lattice.h:233