AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
alma::Gamma_grid Class Reference

Objects of this class represent a regular grid with the Gamma point in one corner. More...

#include <qpoint_grid.hpp>

Public Member Functions

 Gamma_grid (const Crystal_structure &poscar, const Symmetry_operations &symms, const Harmonic_ifcs &force_constants, int _na, int _nb, int _nc)
 Constructor: initialize all internal variables and compute the spectrum at each q point. More...
 
 Gamma_grid (const Crystal_structure &poscar, const Symmetry_operations &symms, const Harmonic_ifcs &force_constants, const Dielectric_parameters &born, int _na, int _nb, int _nc)
 Constructor: initialize all internal variables and compute the spectrum at each q point. More...
 
void enforce_asr ()
 Set the three lowest frequencies at Gamma to zero.
 
std::size_t three_to_one (const std::array< int, 3 > &indices) const
 Return the index of a q point identified by its position along the three axes. More...
 
std::array< int, 3 > one_to_three (std::size_t iq) const
 Return the coordinates of a q point identified by its index. More...
 
std::size_t get_nequivalences () const
 
const Spectrum_at_pointget_spectrum_at_q (int iq) const
 Access the harmonic properties at a point in the grid. More...
 
std::size_t get_cardinal (std::size_t ic) const
 Return the cardinal of an equivalence class. More...
 
std::size_t get_representative (std::size_t ic) const
 Return a representative of an equivalence class. More...
 
std::vector< std::size_t > get_equivalence (std::size_t ic) const
 Return the elements in an equivalence class. More...
 
Eigen::VectorXd get_q (std::size_t iq) const
 Return the Cartesian coordinates of a q point. More...
 
double base_sigma (const Eigen::VectorXd &v) const
 Return the base broadening (without any prefactor) for a mode. More...
 
 Gamma_grid (const Crystal_structure &poscar, const Symmetry_operations &symms, int _na, int _nb, int _nc)
 Very basic constructor that builds a stub object. More...
 
std::size_t polar_opposite (std::size_t q)
 Find the index of the polar opposite q point. More...
 
std::vector< size_t > equivalent_qpoints (std::size_t original) const
 Return the images of a q point through all the symmetry operations, including inversions. More...
 
std::vector< std::array< std::size_t, 2 > > equivalent_qpairs (const std::array< std::size_t, 2 > &original) const
 Find all q-point pairs equivalent to the input. More...
 
std::vector< std::array< std::size_t, 3 > > equivalent_qtriplets (const std::array< std::size_t, 3 > &original) const
 Find all q-point triplets equivalent to the input. More...
 
std::vector< Tetrahedronget_tetrahedra (std::size_t q) const
 Decompose the q-th microcell in five tetrahedra. More...
 
std::vector< Triangleget_triangles (std::size_t ia) const
 
std::size_t getParentIdx (std::size_t iq) const
 
std::size_t getSymIdxToParent (std::size_t iq) const
 
std::vector< std::vector< std::size_t > > getSymmetryMap ()
 
template<typename T >
auto copy_symmetry (std::size_t iq, const Symmetry_operations &symms, const Eigen::MatrixBase< T > &x) const -> Eigen::Matrix< typename T::Scalar, Eigen::Dynamic, Eigen::Dynamic >
 Remove the component of a Cartesian vector which does transform as a q-point in the grid. More...
 

Public Attributes

const int na
 Size of the grid along the first reciprocal axis.
 
const int nb
 Size of the grid along the second reciprocal axis.
 
const int nc
 Size of the grid along the third reciprocal axis.
 
const std::size_t nqpoints
 Total number of q points in the grid.
 
const Eigen::MatrixXd rlattvec
 Reciprocal lattice basis vectors.
 
const Eigen::MatrixXd dq
 Side vectors of each element in reciprocal space.
 

Friends

void save_bulk_hdf5 (const char *filename, const std::string &description, const Crystal_structure &cell, const Symmetry_operations &symmetries, const Gamma_grid &grid, const std::vector< Threeph_process > &processes, const boost::mpi::communicator &comm)
 Save all relevant information about a bulk material to an HDF5 file. More...
 
std::tuple< std::string, std::unique_ptr< Crystal_structure >, std::unique_ptr< Symmetry_operations >, std::unique_ptr< Gamma_grid >, std::unique_ptr< std::vector< Threeph_process > > > load_bulk_hdf5 (const char *filename, const boost::mpi::communicator &comm)
 Reconstruct all data structures from a file saved with save_bulk_hdf5(). More...
 

Detailed Description

Objects of this class represent a regular grid with the Gamma point in one corner.

Constructor & Destructor Documentation

◆ Gamma_grid() [1/3]

alma::Gamma_grid::Gamma_grid ( const Crystal_structure poscar,
const Symmetry_operations symms,
const Harmonic_ifcs force_constants,
int  _na,
int  _nb,
int  _nc 
)

Constructor: initialize all internal variables and compute the spectrum at each q point.

Do not correct the dynamical matrix for the effect of long-range interactions.

◆ Gamma_grid() [2/3]

alma::Gamma_grid::Gamma_grid ( const Crystal_structure poscar,
const Symmetry_operations symms,
const Harmonic_ifcs force_constants,
const Dielectric_parameters born,
int  _na,
int  _nb,
int  _nc 
)

Constructor: initialize all internal variables and compute the spectrum at each q point.

Correct the dynamical matrix for the effect of long-range interactions.

◆ Gamma_grid() [3/3]

alma::Gamma_grid::Gamma_grid ( const Crystal_structure poscar,
const Symmetry_operations symms,
int  _na,
int  _nb,
int  _nc 
)

Very basic constructor that builds a stub object.

Useful for deserialization or for obtaining an equivalences vector without computing the spectrum.

Member Function Documentation

◆ base_sigma()

double alma::Gamma_grid::base_sigma ( const Eigen::VectorXd &  v) const
inline

Return the base broadening (without any prefactor) for a mode.

Parameters
[in]v- group velocity, or difference of group velocities in the case of three-phonon processes
Returns
- the standard deviation of a Gaussian

◆ copy_symmetry()

template<typename T >
auto alma::Gamma_grid::copy_symmetry ( std::size_t  iq,
const Symmetry_operations symms,
const Eigen::MatrixBase< T > &  x 
) const -> Eigen::Matrix<typename T::Scalar, Eigen::Dynamic, Eigen::Dynamic>
inline

Remove the component of a Cartesian vector which does transform as a q-point in the grid.

Parameters
[in]iq- index of the q point
[in]symms- set of symmetry operations of the crystal
[in]x- vector in Cartesian coordinates. Several vectors can be provided by making v a matrix and each vector a column
Returns
the symmetrized version of x

◆ equivalent_qpairs()

std::vector< std::array< std::size_t, 2 > > alma::Gamma_grid::equivalent_qpairs ( const std::array< std::size_t, 2 > &  original) const

Find all q-point pairs equivalent to the input.

Given a pair of indices, obtain all equivalent pairs after looking the up in the symmetry_map.

Parameters
[in]paira pair of q point indices
Returns
a vector of pairs, including the input

◆ equivalent_qpoints()

std::vector<size_t> alma::Gamma_grid::equivalent_qpoints ( std::size_t  original) const
inline

Return the images of a q point through all the symmetry operations, including inversions.

Parameters
[in]q- a q point
Returns
a vector of q point indices. Elements with even indices correspond to operations in the space group, while elements with odd indices compound them with time reversal.

◆ equivalent_qtriplets()

std::vector< std::array< std::size_t, 3 > > alma::Gamma_grid::equivalent_qtriplets ( const std::array< std::size_t, 3 > &  original) const

Find all q-point triplets equivalent to the input.

Given a triplet of indices, obtain all equivalent triplets after looking the up in the symmetry_map.

Parameters
[in]atriplet of q point indices
Returns
a vector of triplets, including the input

◆ get_cardinal()

std::size_t alma::Gamma_grid::get_cardinal ( std::size_t  ic) const
inline

Return the cardinal of an equivalence class.

Parameters
[in]ic- index of the equivalence class.
Returns
the number of points in an equivalence class.

◆ get_equivalence()

std::vector<std::size_t> alma::Gamma_grid::get_equivalence ( std::size_t  ic) const
inline

Return the elements in an equivalence class.

Parameters
[in]ic- index of the equivalence class.
Returns
a vector with the elements in the equivalence class.

◆ get_nequivalences()

std::size_t alma::Gamma_grid::get_nequivalences ( ) const
inline
Returns
the number of q-point equivalence classes in the grid.

◆ get_q()

Eigen::VectorXd alma::Gamma_grid::get_q ( std::size_t  iq) const
inline

Return the Cartesian coordinates of a q point.

Parameters
[in]iq- a q point index
Returns
three Cartesian coordinates

◆ get_representative()

std::size_t alma::Gamma_grid::get_representative ( std::size_t  ic) const
inline

Return a representative of an equivalence class.

The method is guaranteed to always return the same point.

Parameters
[in]ic- index of the equivalence class.
Returns
a representative of the equivalence class.

◆ get_spectrum_at_q()

const Spectrum_at_point& alma::Gamma_grid::get_spectrum_at_q ( int  iq) const
inline

Access the harmonic properties at a point in the grid.

The index is interpreted using modular arithmetic, so even negative values are accepted.

Parameters
[in]iq- index of the q point
Returns
the harmonic properties at the q point.

◆ get_tetrahedra()

std::vector<Tetrahedron> alma::Gamma_grid::get_tetrahedra ( std::size_t  q) const
inline

Decompose the q-th microcell in five tetrahedra.

Parameters
[in]q- the index of the q point
Returns
a vector of five Tetrahedron objects

◆ get_triangles()

std::vector< Triangle > alma::Gamma_grid::get_triangles ( std::size_t  ia) const
Returns
a vector of triangle objects with each containing three indices

◆ getParentIdx()

std::size_t alma::Gamma_grid::getParentIdx ( std::size_t  iq) const
Returns
the index of the representative of the equivalence class of the provided q-point

◆ getSymIdxToParent()

std::size_t alma::Gamma_grid::getSymIdxToParent ( std::size_t  iq) const
Returns
the index of the symmetry operation that maps the provided q-point to the representative of the equivalence class to which it belongs

◆ getSymmetryMap()

std::vector<std::vector<std::size_t> > alma::Gamma_grid::getSymmetryMap ( )
inline
Returns
the grid's symmetry map

◆ one_to_three()

std::array<int, 3> alma::Gamma_grid::one_to_three ( std::size_t  iq) const
inline

Return the coordinates of a q point identified by its index.

Parameters
[in]index- index of the q point
Returns
- an array with the positions of the q point along the three axes.

◆ polar_opposite()

std::size_t alma::Gamma_grid::polar_opposite ( std::size_t  q)
inline

Find the index of the polar opposite q point.

Parameters
[in]q- a q point
Returns
the index of the polar opposite of q

◆ three_to_one()

std::size_t alma::Gamma_grid::three_to_one ( const std::array< int, 3 > &  indices) const
inline

Return the index of a q point identified by its position along the three axes.

Indices are interpreted modulo-{na, nb, nc}, so even negative values are accepted.

Parameters
[in]indices- positions of the q point along the three axes
Returns
the index of the q point in this grid

Friends And Related Function Documentation

◆ load_bulk_hdf5

std::tuple<std::string, std::unique_ptr<Crystal_structure>, std::unique_ptr<Symmetry_operations>, std::unique_ptr<Gamma_grid>, std::unique_ptr<std::vector<Threeph_process> > > load_bulk_hdf5 ( const char *  filename,
const boost::mpi::communicator &  comm 
)
friend

Reconstruct all data structures from a file saved with save_bulk_hdf5().

Parameters
[in]filename- path to the HDF5 file
[in]comm- MPI communicator used to determine which part of the data to load in this process.
Returns
a tuple with all the information in the file, namely: description - a free-format description cell - description of the unit cell symmetries - symmetry operations of the unit cell grid - phonon spectrum on a refular grid processes - allowed three-phonon processes

◆ save_bulk_hdf5

void save_bulk_hdf5 ( const char *  filename,
const std::string &  description,
const Crystal_structure cell,
const Symmetry_operations symmetries,
const Gamma_grid grid,
const std::vector< Threeph_process > &  processes,
const boost::mpi::communicator &  comm 
)
friend

Save all relevant information about a bulk material to an HDF5 file.

Parameters
[in]filename- path to the HDF5 file. It can be NULL for all processes except the one with rank == 0.
[in]description- a free-format description
[in]cell- description of the unit cell
[in]symmetries- symmetry operations of the unit cell
[in]grid- phonon spectrum on a refular grid
[in]processes- allowed three-phonon processes
[in]comm- MPI communicator used to coordinate with all other processes

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