SIRIUS 7.5.0
Electronic structure library and applications
Public Member Functions | Private Member Functions | Private Attributes | List of all members
sirius::Unit_cell Class Reference

Representation of a unit cell. More...

#include <unit_cell.hpp>

Public Member Functions

 Unit_cell (Simulation_parameters const &parameters__, mpi::Communicator const &comm__)
 
void initialize ()
 Initialize the unit cell data. More...
 
Atom_typeadd_atom_type (const std::string label__, const std::string file_name__="")
 Add new atom type to the list of atom types and read necessary data from the .json file. More...
 
void add_atom (const std::string label, r3::vector< double > position, r3::vector< double > vector_field)
 Add new atom to the list of atom types. More...
 
void add_atom (const std::string label, r3::vector< double > position)
 Add new atom without vector field to the list of atom types. More...
 
void init_paw ()
 Add PAW atoms. More...
 
int num_paw_atoms () const
 Return number of atoms with PAW pseudopotential. More...
 
auto const & paw_atoms () const
 
auto const & spl_num_paw_atoms () const
 Get split index of PAW atoms. More...
 
auto spl_num_paw_atoms (paw_atom_index_t::local idx__) const
 
auto paw_atom_index (paw_atom_index_t::global ipaw__) const
 Return global index of atom by the index in the list of PAW atoms. More...
 
void print_info (std::ostream &out__, int verbosity__) const
 Print basic info. More...
 
void print_geometry_info (std::ostream &out__, int verbosity__) const
 
void print_nearest_neighbours (std::ostream &out__) const
 
unit_cell_parameters_descriptor unit_cell_parameters ()
 
void get_symmetry ()
 Get crystal symmetries and equivalent atoms. More...
 
void write_cif ()
 Write structure to CIF file. More...
 
json serialize (bool cart_pos__=false) const
 Write structure to JSON file. More...
 
void set_lattice_vectors (r3::matrix< double > lattice_vectors__)
 Set matrix of lattice vectors. More...
 
void set_lattice_vectors (r3::vector< double > a0__, r3::vector< double > a1__, r3::vector< double > a2__)
 Set lattice vectors. More...
 
void find_nearest_neighbours (double cluster_radius)
 Find the cluster of nearest neighbours around each atom. More...
 
bool is_point_in_mt (r3::vector< double > vc, int &ja, int &jr, double &dr, double tp[2]) const
 
void generate_radial_functions (std::ostream &out__)
 
void generate_radial_integrals ()
 
std::string chemical_formula ()
 Get a simple simple chemical formula bases on the total unit cell. More...
 
void update ()
 Update the parameters that depend on atomic positions or lattice vectors. More...
 
void import (config_t::unit_cell_t const &inp__)
 Import unit cell description from the input data structure. More...
 
int atom_id_by_position (r3::vector< double > position__)
 Get atom ID (global index) by it's position in fractional coordinates. More...
 
double min_bond_length () const
 Find the minimum bond length. More...
 
std::vector< double > find_mt_radii (int auto_rmt__, bool inflate__)
 Automatically determine new muffin-tin radii as a half distance between neighbor atoms. More...
 
std::pair< int, std::vector< int > > num_ps_atomic_wf () const
 Get the total number of pseudo atomic wave-functions. More...
 
auto const & num_hubbard_wf () const
 Return number of Hubbard wave-functions. More...
 
template<typename T >
auto get_cartesian_coordinates (r3::vector< T > a__) const
 Get Cartesian coordinates of the vector by its fractional coordinates. More...
 
auto get_fractional_coordinates (r3::vector< double > a__) const
 Get fractional coordinates of the vector by its Cartesian coordinates. More...
 
double omega () const
 Unit cell volume. More...
 
int num_atom_types () const
 Number of atom types. More...
 
Atom_typeatom_type (int id__)
 Return atom type instance by id. More...
 
Atom_type const & atom_type (int id__) const
 Return const atom type instance by id. More...
 
auto const & atom_type (std::string const label__) const
 Return const atom type instance by label. More...
 
auto & atom_type (std::string const label__)
 Return atom type instance by label. More...
 
int num_atom_symmetry_classes () const
 Number of atom symmetry classes. More...
 
Atom_symmetry_class const & atom_symmetry_class (int id__) const
 Return const symmetry class instance by class id. More...
 
Atom_symmetry_classatom_symmetry_class (int id__)
 Return symmetry class instance by class id. More...
 
int num_atoms () const
 Number of atoms in the unit cell. More...
 
Atom const & atom (int id__) const
 Return const atom instance by id. More...
 
Atomatom (int id__)
 Return atom instance by id. More...
 
int total_nuclear_charge () const
 Total nuclear charge. More...
 
double num_electrons () const
 Total number of electrons (core + valence). More...
 
double num_valence_electrons () const
 Number of valence electrons. More...
 
double num_core_electrons () const
 Number of core electrons. More...
 
int max_num_mt_points () const
 Maximum number of muffin-tin points among all atom types. More...
 
int mt_aw_basis_size () const
 Total number of the augmented wave basis functions over all atoms. More...
 
int mt_lo_basis_size () const
 Total number of local orbital basis functions over all atoms. More...
 
int max_mt_basis_size () const
 Maximum number of basis functions among all atom types. More...
 
int max_mt_radial_basis_size () const
 Maximum number of radial functions among all atom types. More...
 
double min_mt_radius () const
 Minimum muffin-tin radius. More...
 
double max_mt_radius () const
 Maximum muffin-tin radius. More...
 
int max_mt_aw_basis_size () const
 Maximum number of AW basis functions among all atom types. More...
 
int max_mt_lo_basis_size () const
 Maximum local orbital basis size among all atoms. More...
 
auto max_num_atoms () const
 Maximum number of atoms across all atom types. More...
 
void set_equivalent_atoms (int const *equivalent_atoms__)
 
auto const & spl_num_atoms () const
 
auto const & spl_num_atom_symmetry_classes () const
 
double volume_mt () const
 
double volume_it () const
 
int lmax () const
 Maximum orbital quantum number of radial functions between all atom types. More...
 
int lmax_apw () const
 
int num_nearest_neighbours (int ia) const
 
auto const & nearest_neighbour (int i, int ia) const
 
auto const & symmetry () const
 
auto const & lattice_vectors () const
 
auto const & inverse_lattice_vectors () const
 
auto const & reciprocal_lattice_vectors () const
 
auto lattice_vector (int idx__) const
 Return a single lattice vector. More...
 
auto const & parameters () const
 
auto const & comm () const
 
auto const & atom_coord (int iat__) const
 
bool augment () const
 Return 'True' if at least one atom in the unit cell has an augmentation charge. More...
 

Private Member Functions

bool check_mt_overlap (int &ia__, int &ja__)
 Check if MT spheres overlap. More...
 
int next_atom_type_id (std::string label__)
 

Private Attributes

Simulation_parameters const & parameters_
 Basic parameters of the simulation. More...
 
std::map< std::string, int > atom_type_id_map_
 Mapping between atom type label and an ordered internal id in the range [0, \( N_{types} \)). More...
 
std::vector< std::shared_ptr< Atom_type > > atom_types_
 List of atom types. More...
 
std::vector< std::shared_ptr< Atom_symmetry_class > > atom_symmetry_classes_
 List of atom classes. More...
 
std::vector< std::shared_ptr< Atom > > atoms_
 List of atoms. More...
 
splindex_block< atom_index_tspl_num_atoms_
 Split index of atoms. More...
 
std::vector< int > paw_atom_index_
 Global index of atom by index of PAW atom. More...
 
splindex_block< paw_atom_index_tspl_num_paw_atoms_
 Split index of PAW atoms. More...
 
splindex_block< atom_symmetry_class_index_tspl_num_atom_symmetry_classes_
 Split index of atom symmetry classes. More...
 
r3::matrix< double > lattice_vectors_
 Bravais lattice vectors in column order. More...
 
r3::matrix< double > inverse_lattice_vectors_
 Inverse matrix of Bravais lattice vectors. More...
 
r3::matrix< double > reciprocal_lattice_vectors_
 Reciprocal lattice vectors in column order. More...
 
double omega_ {0}
 Volume \( \Omega \) of the unit cell. Volume of Brillouin zone is then \( (2\Pi)^3 / \Omega \). More...
 
double volume_mt_ {0}
 Total volume of the muffin-tin spheres. More...
 
double volume_it_ {0}
 Volume of the interstitial region. More...
 
std::vector< int > equivalent_atoms_
 List of equivalent atoms, provided externally. More...
 
std::vector< std::vector< nearest_neighbour_descriptor > > nearest_neighbours_
 List of nearest neighbours for each atom. More...
 
std::unique_ptr< Crystal_symmetrysymmetry_
 
std::vector< sddk::mdarray< double, 2 > > atom_coord_
 Atomic coordinates in GPU-friendly ordering packed in arrays for each atom type. More...
 
mpi::Communicator const & comm_
 
std::pair< int, std::vector< int > > num_hubbard_wf_
 

Detailed Description

Representation of a unit cell.

Definition at line 42 of file unit_cell.hpp.

Constructor & Destructor Documentation

◆ Unit_cell()

sirius::Unit_cell::Unit_cell ( Simulation_parameters const &  parameters__,
mpi::Communicator const &  comm__ 
)

Definition at line 32 of file unit_cell.cpp.

Member Function Documentation

◆ check_mt_overlap()

bool sirius::Unit_cell::check_mt_overlap ( int &  ia__,
int &  ja__ 
)
inlineprivate

Check if MT spheres overlap.

Definition at line 144 of file unit_cell.cpp.

◆ next_atom_type_id()

int sirius::Unit_cell::next_atom_type_id ( std::string  label__)
private

Definition at line 944 of file unit_cell.cpp.

◆ initialize()

void sirius::Unit_cell::initialize ( )

Initialize the unit cell data.

Several things must be done during this phase:

  1. Compute number of electrons
  2. Compute MT basis function indices
  3. [if needed] Scale MT radii
  4. Check MT overlap
  5. Create radial grid for each atom type
  6. Find symmetry and assign symmetry class to each atom
  7. Create split indices for atoms and atom classes

Definition at line 669 of file unit_cell.cpp.

◆ add_atom_type()

Atom_type & sirius::Unit_cell::add_atom_type ( const std::string  label__,
const std::string  file_name__ = "" 
)

Add new atom type to the list of atom types and read necessary data from the .json file.

Definition at line 629 of file unit_cell.cpp.

◆ add_atom() [1/2]

void sirius::Unit_cell::add_atom ( const std::string  label,
r3::vector< double >  position,
r3::vector< double >  vector_field 
)

Add new atom to the list of atom types.

Definition at line 637 of file unit_cell.cpp.

◆ add_atom() [2/2]

void sirius::Unit_cell::add_atom ( const std::string  label,
r3::vector< double >  position 
)
inline

Add new atom without vector field to the list of atom types.

Definition at line 152 of file unit_cell.hpp.

◆ init_paw()

void sirius::Unit_cell::init_paw ( )

Add PAW atoms.

Definition at line 958 of file unit_cell.cpp.

◆ num_paw_atoms()

int sirius::Unit_cell::num_paw_atoms ( ) const
inline

Return number of atoms with PAW pseudopotential.

Definition at line 161 of file unit_cell.hpp.

◆ paw_atoms()

auto const & sirius::Unit_cell::paw_atoms ( ) const
inline

Definition at line 166 of file unit_cell.hpp.

◆ spl_num_paw_atoms() [1/2]

auto const & sirius::Unit_cell::spl_num_paw_atoms ( ) const
inline

Get split index of PAW atoms.

Definition at line 172 of file unit_cell.hpp.

◆ spl_num_paw_atoms() [2/2]

auto sirius::Unit_cell::spl_num_paw_atoms ( paw_atom_index_t::local  idx__) const
inline

Definition at line 177 of file unit_cell.hpp.

◆ paw_atom_index()

auto sirius::Unit_cell::paw_atom_index ( paw_atom_index_t::global  ipaw__) const
inline

Return global index of atom by the index in the list of PAW atoms.

Definition at line 183 of file unit_cell.hpp.

◆ print_info()

void sirius::Unit_cell::print_info ( std::ostream &  out__,
int  verbosity__ 
) const

Print basic info.

Definition at line 170 of file unit_cell.cpp.

◆ print_geometry_info()

void sirius::Unit_cell::print_geometry_info ( std::ostream &  out__,
int  verbosity__ 
) const

Definition at line 266 of file unit_cell.cpp.

◆ print_nearest_neighbours()

void sirius::Unit_cell::print_nearest_neighbours ( std::ostream &  out__) const

Definition at line 463 of file unit_cell.cpp.

◆ unit_cell_parameters()

unit_cell_parameters_descriptor sirius::Unit_cell::unit_cell_parameters ( )

Definition at line 322 of file unit_cell.cpp.

◆ get_symmetry()

void sirius::Unit_cell::get_symmetry ( )

Get crystal symmetries and equivalent atoms.

Makes a call to spglib providing the basic unit cell information: lattice vectors and atomic types and positions. Gets back symmetry operations and a table of equivalent atoms. The table of equivalent atoms is then used to make a list of atom symmetry classes and related data.

Definition at line 726 of file unit_cell.cpp.

◆ write_cif()

void sirius::Unit_cell::write_cif ( )

Write structure to CIF file.

Definition at line 342 of file unit_cell.cpp.

◆ serialize()

json sirius::Unit_cell::serialize ( bool  cart_pos__ = false) const

Write structure to JSON file.

Definition at line 374 of file unit_cell.cpp.

◆ set_lattice_vectors() [1/2]

void sirius::Unit_cell::set_lattice_vectors ( r3::matrix< double >  lattice_vectors__)

Set matrix of lattice vectors.

Initializes lattice vectors, inverse lattice vector matrix, reciprocal lattice vectors and the unit cell volume.

Definition at line 911 of file unit_cell.cpp.

◆ set_lattice_vectors() [2/2]

void sirius::Unit_cell::set_lattice_vectors ( r3::vector< double >  a0__,
r3::vector< double >  a1__,
r3::vector< double >  a2__ 
)

Set lattice vectors.

Definition at line 920 of file unit_cell.cpp.

◆ find_nearest_neighbours()

void sirius::Unit_cell::find_nearest_neighbours ( double  cluster_radius)

Find the cluster of nearest neighbours around each atom.

Definition at line 412 of file unit_cell.cpp.

◆ is_point_in_mt()

bool sirius::Unit_cell::is_point_in_mt ( r3::vector< double >  vc,
int &  ja,
int &  jr,
double &  dr,
double  tp[2] 
) const

Definition at line 489 of file unit_cell.cpp.

◆ generate_radial_functions()

void sirius::Unit_cell::generate_radial_functions ( std::ostream &  out__)

Definition at line 538 of file unit_cell.cpp.

◆ generate_radial_integrals()

void sirius::Unit_cell::generate_radial_integrals ( )

Definition at line 570 of file unit_cell.cpp.

◆ chemical_formula()

std::string sirius::Unit_cell::chemical_formula ( )

Get a simple simple chemical formula bases on the total unit cell.

Atoms of each type are counted and packed in a string. For example, O2Ni2 or La2O4Cu

Definition at line 608 of file unit_cell.cpp.

◆ update()

void sirius::Unit_cell::update ( )

Update the parameters that depend on atomic positions or lattice vectors.

Definition at line 839 of file unit_cell.cpp.

◆ import()

void sirius::Unit_cell::import ( config_t::unit_cell_t const &  inp__)

Import unit cell description from the input data structure.

Set lattice vectors, atom types and coordinates of atoms. The "atom_coordinate_units" parameter by default is assumed to be "lattice" which means that the atomic coordinates are provided in lattice (fractional) units. It can also be specified in "A" or "au" which means that the input atomic coordinates are Cartesian and provided in Angstroms or atomic units of length. This is useful in setting up the molecule calculation.

Definition at line 793 of file unit_cell.cpp.

◆ atom_id_by_position()

int sirius::Unit_cell::atom_id_by_position ( r3::vector< double >  position__)

Get atom ID (global index) by it's position in fractional coordinates.

Definition at line 932 of file unit_cell.cpp.

◆ min_bond_length()

double sirius::Unit_cell::min_bond_length ( ) const

Find the minimum bond length.

This is useful to check the sanity of the crystal structure.

Definition at line 656 of file unit_cell.cpp.

◆ find_mt_radii()

std::vector< double > sirius::Unit_cell::find_mt_radii ( int  auto_rmt__,
bool  inflate__ 
)

Automatically determine new muffin-tin radii as a half distance between neighbor atoms.

In order to guarantee a unique solution muffin-tin radii are dermined as a half distance bethween nearest atoms. Initial values of the muffin-tin radii are ignored.

Definition at line 41 of file unit_cell.cpp.

◆ num_ps_atomic_wf()

std::pair< int, std::vector< int > > sirius::Unit_cell::num_ps_atomic_wf ( ) const

Get the total number of pseudo atomic wave-functions.

Definition at line 971 of file unit_cell.cpp.

◆ num_hubbard_wf()

auto const & sirius::Unit_cell::num_hubbard_wf ( ) const
inline

Return number of Hubbard wave-functions.

Definition at line 256 of file unit_cell.hpp.

◆ get_cartesian_coordinates()

template<typename T >
auto sirius::Unit_cell::get_cartesian_coordinates ( r3::vector< T >  a__) const
inline

Get Cartesian coordinates of the vector by its fractional coordinates.

Definition at line 263 of file unit_cell.hpp.

◆ get_fractional_coordinates()

auto sirius::Unit_cell::get_fractional_coordinates ( r3::vector< double >  a__) const
inline

Get fractional coordinates of the vector by its Cartesian coordinates.

Definition at line 269 of file unit_cell.hpp.

◆ omega()

double sirius::Unit_cell::omega ( ) const
inline

Unit cell volume.

Definition at line 275 of file unit_cell.hpp.

◆ num_atom_types()

int sirius::Unit_cell::num_atom_types ( ) const
inline

Number of atom types.

Definition at line 281 of file unit_cell.hpp.

◆ atom_type() [1/4]

Atom_type & sirius::Unit_cell::atom_type ( int  id__)
inline

Return atom type instance by id.

Definition at line 288 of file unit_cell.hpp.

◆ atom_type() [2/4]

Atom_type const & sirius::Unit_cell::atom_type ( int  id__) const
inline

Return const atom type instance by id.

Definition at line 295 of file unit_cell.hpp.

◆ atom_type() [3/4]

auto const & sirius::Unit_cell::atom_type ( std::string const  label__) const
inline

Return const atom type instance by label.

Definition at line 302 of file unit_cell.hpp.

◆ atom_type() [4/4]

auto & sirius::Unit_cell::atom_type ( std::string const  label__)
inline

Return atom type instance by label.

Definition at line 314 of file unit_cell.hpp.

◆ num_atom_symmetry_classes()

int sirius::Unit_cell::num_atom_symmetry_classes ( ) const
inline

Number of atom symmetry classes.

Definition at line 320 of file unit_cell.hpp.

◆ atom_symmetry_class() [1/2]

Atom_symmetry_class const & sirius::Unit_cell::atom_symmetry_class ( int  id__) const
inline

Return const symmetry class instance by class id.

Definition at line 326 of file unit_cell.hpp.

◆ atom_symmetry_class() [2/2]

Atom_symmetry_class & sirius::Unit_cell::atom_symmetry_class ( int  id__)
inline

Return symmetry class instance by class id.

Definition at line 332 of file unit_cell.hpp.

◆ num_atoms()

int sirius::Unit_cell::num_atoms ( ) const
inline

Number of atoms in the unit cell.

Definition at line 338 of file unit_cell.hpp.

◆ atom() [1/2]

Atom const & sirius::Unit_cell::atom ( int  id__) const
inline

Return const atom instance by id.

Definition at line 344 of file unit_cell.hpp.

◆ atom() [2/2]

Atom & sirius::Unit_cell::atom ( int  id__)
inline

Return atom instance by id.

Definition at line 351 of file unit_cell.hpp.

◆ total_nuclear_charge()

int sirius::Unit_cell::total_nuclear_charge ( ) const
inline

Total nuclear charge.

Definition at line 358 of file unit_cell.hpp.

◆ num_electrons()

double sirius::Unit_cell::num_electrons ( ) const
inline

Total number of electrons (core + valence).

Definition at line 369 of file unit_cell.hpp.

◆ num_valence_electrons()

double sirius::Unit_cell::num_valence_electrons ( ) const
inline

Number of valence electrons.

Definition at line 375 of file unit_cell.hpp.

◆ num_core_electrons()

double sirius::Unit_cell::num_core_electrons ( ) const
inline

Number of core electrons.

Definition at line 386 of file unit_cell.hpp.

◆ max_num_mt_points()

int sirius::Unit_cell::max_num_mt_points ( ) const
inline

Maximum number of muffin-tin points among all atom types.

Definition at line 397 of file unit_cell.hpp.

◆ mt_aw_basis_size()

int sirius::Unit_cell::mt_aw_basis_size ( ) const
inline

Total number of the augmented wave basis functions over all atoms.

This is equal to the total number of matching coefficients for each plane-wave.

Definition at line 408 of file unit_cell.hpp.

◆ mt_lo_basis_size()

int sirius::Unit_cell::mt_lo_basis_size ( ) const
inline

Total number of local orbital basis functions over all atoms.

Definition at line 419 of file unit_cell.hpp.

◆ max_mt_basis_size()

int sirius::Unit_cell::max_mt_basis_size ( ) const
inline

Maximum number of basis functions among all atom types.

Definition at line 430 of file unit_cell.hpp.

◆ max_mt_radial_basis_size()

int sirius::Unit_cell::max_mt_radial_basis_size ( ) const
inline

Maximum number of radial functions among all atom types.

Definition at line 440 of file unit_cell.hpp.

◆ min_mt_radius()

double sirius::Unit_cell::min_mt_radius ( ) const
inline

Minimum muffin-tin radius.

Definition at line 450 of file unit_cell.hpp.

◆ max_mt_radius()

double sirius::Unit_cell::max_mt_radius ( ) const
inline

Maximum muffin-tin radius.

Definition at line 460 of file unit_cell.hpp.

◆ max_mt_aw_basis_size()

int sirius::Unit_cell::max_mt_aw_basis_size ( ) const
inline

Maximum number of AW basis functions among all atom types.

Definition at line 470 of file unit_cell.hpp.

◆ max_mt_lo_basis_size()

int sirius::Unit_cell::max_mt_lo_basis_size ( ) const
inline

Maximum local orbital basis size among all atoms.

Definition at line 480 of file unit_cell.hpp.

◆ max_num_atoms()

auto sirius::Unit_cell::max_num_atoms ( ) const
inline

Maximum number of atoms across all atom types.

Definition at line 490 of file unit_cell.hpp.

◆ set_equivalent_atoms()

void sirius::Unit_cell::set_equivalent_atoms ( int const *  equivalent_atoms__)
inline

Definition at line 499 of file unit_cell.hpp.

◆ spl_num_atoms()

auto const & sirius::Unit_cell::spl_num_atoms ( ) const
inline

Definition at line 505 of file unit_cell.hpp.

◆ spl_num_atom_symmetry_classes()

auto const & sirius::Unit_cell::spl_num_atom_symmetry_classes ( ) const
inline

Definition at line 510 of file unit_cell.hpp.

◆ volume_mt()

double sirius::Unit_cell::volume_mt ( ) const
inline

Definition at line 515 of file unit_cell.hpp.

◆ volume_it()

double sirius::Unit_cell::volume_it ( ) const
inline

Definition at line 520 of file unit_cell.hpp.

◆ lmax()

int sirius::Unit_cell::lmax ( ) const
inline

Maximum orbital quantum number of radial functions between all atom types.

Definition at line 526 of file unit_cell.hpp.

◆ lmax_apw()

int sirius::Unit_cell::lmax_apw ( ) const
inline

Definition at line 535 of file unit_cell.hpp.

◆ num_nearest_neighbours()

int sirius::Unit_cell::num_nearest_neighbours ( int  ia) const
inline

Definition at line 544 of file unit_cell.hpp.

◆ nearest_neighbour()

auto const & sirius::Unit_cell::nearest_neighbour ( int  i,
int  ia 
) const
inline

Definition at line 549 of file unit_cell.hpp.

◆ symmetry()

auto const & sirius::Unit_cell::symmetry ( ) const
inline

Definition at line 554 of file unit_cell.hpp.

◆ lattice_vectors()

auto const & sirius::Unit_cell::lattice_vectors ( ) const
inline

Definition at line 560 of file unit_cell.hpp.

◆ inverse_lattice_vectors()

auto const & sirius::Unit_cell::inverse_lattice_vectors ( ) const
inline

Definition at line 565 of file unit_cell.hpp.

◆ reciprocal_lattice_vectors()

auto const & sirius::Unit_cell::reciprocal_lattice_vectors ( ) const
inline

Definition at line 570 of file unit_cell.hpp.

◆ lattice_vector()

auto sirius::Unit_cell::lattice_vector ( int  idx__) const
inline

Return a single lattice vector.

Definition at line 576 of file unit_cell.hpp.

◆ parameters()

auto const & sirius::Unit_cell::parameters ( ) const
inline

Definition at line 581 of file unit_cell.hpp.

◆ comm()

auto const & sirius::Unit_cell::comm ( ) const
inline

Definition at line 586 of file unit_cell.hpp.

◆ atom_coord()

auto const & sirius::Unit_cell::atom_coord ( int  iat__) const
inline

Definition at line 591 of file unit_cell.hpp.

◆ augment()

bool sirius::Unit_cell::augment ( ) const
inline

Return 'True' if at least one atom in the unit cell has an augmentation charge.

Definition at line 597 of file unit_cell.hpp.

Member Data Documentation

◆ parameters_

Simulation_parameters const& sirius::Unit_cell::parameters_
private

Basic parameters of the simulation.

Definition at line 46 of file unit_cell.hpp.

◆ atom_type_id_map_

std::map<std::string, int> sirius::Unit_cell::atom_type_id_map_
private

Mapping between atom type label and an ordered internal id in the range [0, \( N_{types} \)).

Definition at line 49 of file unit_cell.hpp.

◆ atom_types_

std::vector<std::shared_ptr<Atom_type> > sirius::Unit_cell::atom_types_
private

List of atom types.

Definition at line 52 of file unit_cell.hpp.

◆ atom_symmetry_classes_

std::vector<std::shared_ptr<Atom_symmetry_class> > sirius::Unit_cell::atom_symmetry_classes_
private

List of atom classes.

Definition at line 55 of file unit_cell.hpp.

◆ atoms_

std::vector<std::shared_ptr<Atom> > sirius::Unit_cell::atoms_
private

List of atoms.

Definition at line 58 of file unit_cell.hpp.

◆ spl_num_atoms_

splindex_block<atom_index_t> sirius::Unit_cell::spl_num_atoms_
private

Split index of atoms.

Definition at line 61 of file unit_cell.hpp.

◆ paw_atom_index_

std::vector<int> sirius::Unit_cell::paw_atom_index_
private

Global index of atom by index of PAW atom.

Definition at line 64 of file unit_cell.hpp.

◆ spl_num_paw_atoms_

splindex_block<paw_atom_index_t> sirius::Unit_cell::spl_num_paw_atoms_
private

Split index of PAW atoms.

Definition at line 67 of file unit_cell.hpp.

◆ spl_num_atom_symmetry_classes_

splindex_block<atom_symmetry_class_index_t> sirius::Unit_cell::spl_num_atom_symmetry_classes_
private

Split index of atom symmetry classes.

Definition at line 70 of file unit_cell.hpp.

◆ lattice_vectors_

r3::matrix<double> sirius::Unit_cell::lattice_vectors_
private

Bravais lattice vectors in column order.

The following convention is used to transform fractional coordinates to Cartesian:

\[ \vec v_{C} = {\bf L} \vec v_{f} \]

Definition at line 78 of file unit_cell.hpp.

◆ inverse_lattice_vectors_

r3::matrix<double> sirius::Unit_cell::inverse_lattice_vectors_
private

Inverse matrix of Bravais lattice vectors.

This matrix is used to find fractional coordinates by Cartesian coordinates:

\[ \vec v_{f} = {\bf L}^{-1} \vec v_{C} \]

Definition at line 86 of file unit_cell.hpp.

◆ reciprocal_lattice_vectors_

r3::matrix<double> sirius::Unit_cell::reciprocal_lattice_vectors_
private

Reciprocal lattice vectors in column order.

The following convention is used:

\[ \vec a_{i} \vec b_{j} = 2 \pi \delta_{ij} \]

or in matrix notation

\[ {\bf A} {\bf B}^{T} = 2 \pi {\bf I} \]

Definition at line 98 of file unit_cell.hpp.

◆ omega_

double sirius::Unit_cell::omega_ {0}
private

Volume \( \Omega \) of the unit cell. Volume of Brillouin zone is then \( (2\Pi)^3 / \Omega \).

Definition at line 101 of file unit_cell.hpp.

◆ volume_mt_

double sirius::Unit_cell::volume_mt_ {0}
private

Total volume of the muffin-tin spheres.

Definition at line 104 of file unit_cell.hpp.

◆ volume_it_

double sirius::Unit_cell::volume_it_ {0}
private

Volume of the interstitial region.

Definition at line 107 of file unit_cell.hpp.

◆ equivalent_atoms_

std::vector<int> sirius::Unit_cell::equivalent_atoms_
private

List of equivalent atoms, provided externally.

Definition at line 110 of file unit_cell.hpp.

◆ nearest_neighbours_

std::vector<std::vector<nearest_neighbour_descriptor> > sirius::Unit_cell::nearest_neighbours_
private

List of nearest neighbours for each atom.

Definition at line 113 of file unit_cell.hpp.

◆ symmetry_

std::unique_ptr<Crystal_symmetry> sirius::Unit_cell::symmetry_
private

Definition at line 115 of file unit_cell.hpp.

◆ atom_coord_

std::vector<sddk::mdarray<double, 2> > sirius::Unit_cell::atom_coord_
private

Atomic coordinates in GPU-friendly ordering packed in arrays for each atom type.

Definition at line 118 of file unit_cell.hpp.

◆ comm_

mpi::Communicator const& sirius::Unit_cell::comm_
private

Definition at line 120 of file unit_cell.hpp.

◆ num_hubbard_wf_

std::pair<int, std::vector<int> > sirius::Unit_cell::num_hubbard_wf_
private

Definition at line 122 of file unit_cell.hpp.


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