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

#include <matching_coefficients.hpp>

Public Member Functions

 Matching_coefficients (Unit_cell const &unit_cell__, fft::Gvec const &gkvec__)
 Constructor. More...
 
template<bool conjugate, typename T , typename = std::enable_if_t<!std::is_scalar<T>::value>>
void generate (Atom const &atom__, sddk::mdarray< T, 2 > &alm__) const
 Generate plane-wave matching coefficients for the radial solutions of a given atom. More...
 
auto const & gkvec () const
 

Private Member Functions

template<int N, bool conjugate, typename T , typename = std::enable_if_t<!std::is_scalar<T>::value>>
void generate (int ngk, std::vector< std::complex< double > > const &phase_factors__, int iat, int l, int lm, int nu, r3::matrix< double > const &A, T *alm) const
 Generate matching coefficients for a specific \( \ell \) and order. More...
 

Private Attributes

Unit_cell const & unit_cell_
 Description of the unit cell. More...
 
fft::Gvec const & gkvec_
 Description of the G+k vectors. More...
 
std::vector< double > gkvec_len_
 
sddk::mdarray< std::complex< double >, 2 > gkvec_ylm_
 Spherical harmonics Ylm(theta, phi) of the G+k vectors. More...
 
sddk::mdarray< std::complex< double >, 4 > alm_b_
 Precomputed values for the linear equations for matching coefficients. More...
 

Detailed Description

The following matching conditions must be fulfilled:

\[ \frac{\partial^j}{\partial r^j} \sum_{L \nu} A_{L \nu}^{\bf k}({\bf G})u_{\ell \nu}(r) Y_{L}(\hat {\bf r}) \bigg|_{R^{MT}} = \frac{\partial^j}{\partial r^j} \frac{4 \pi}{\sqrt \Omega} e^{i{\bf (G+k)\tau}} \sum_{L}i^{\ell} j_{\ell}(|{\bf G+k}|r) Y_{L}^{*}(\widehat {\bf G+k}) Y_{L}(\hat {\bf r}) \bigg|_{R^{MT}} \]

where \( L = \{ \ell, m \} \). Dropping sum over L we arrive to the following system of linear equations:

\[ \sum_{\nu} \frac{\partial^j u_{\ell \nu}(r)}{\partial r^j} \bigg|_{R^{MT}} A_{L \nu}^{\bf k}({\bf G}) = \frac{4 \pi}{\sqrt \Omega} e^{i{\bf (G+k)\tau}} i^{\ell} \frac{\partial^j j_{\ell}(|{\bf G+k}|r)}{\partial r^j} \bigg|_{R^{MT}} Y_{L}^{*}(\widehat {\bf G+k}) \]

The matching coefficients are then equal to:

\[ A_{L \nu}^{\bf k}({\bf G}) = \sum_{j} \bigg[ \frac{\partial^j u_{\ell \nu}(r)}{\partial r^j} \bigg|_{R^{MT}} \bigg]_{\nu j}^{-1} \frac{\partial^j j_{\ell}(|{\bf G+k}|r)}{\partial r^j} \bigg|_{R^{MT}} \frac{4 \pi}{\sqrt \Omega} i^{\ell} e^{i{\bf (G+k)\tau}} Y_{L}^{*}(\widehat {\bf G+k}) \]

Definition at line 53 of file matching_coefficients.hpp.

Constructor & Destructor Documentation

◆ Matching_coefficients()

sirius::Matching_coefficients::Matching_coefficients ( Unit_cell const &  unit_cell__,
fft::Gvec const &  gkvec__ 
)
inline

Constructor.

Definition at line 112 of file matching_coefficients.hpp.

Member Function Documentation

◆ generate() [1/2]

template<int N, bool conjugate, typename T , typename = std::enable_if_t<!std::is_scalar<T>::value>>
void sirius::Matching_coefficients::generate ( int  ngk,
std::vector< std::complex< double > > const &  phase_factors__,
int  iat,
int  l,
int  lm,
int  nu,
r3::matrix< double > const &  A,
T *  alm 
) const
inlineprivate

Generate matching coefficients for a specific \( \ell \) and order.

Parameters
[in]ngkNumber of G+k vectors.
[in]phase_factorsPhase factors of G+k vectors.
[in]iatIndex of atom type.
[in]lOrbital quantum nuber.
[in]lmComposite l,m index.
[in]nuOrder of radial function \( u_{\ell \nu}(r) \) for which coefficients are generated.
[in]AInverse matrix of radial derivatives.
[out]almPointer to alm coefficients.

Definition at line 81 of file matching_coefficients.hpp.

◆ generate() [2/2]

template<bool conjugate, typename T , typename = std::enable_if_t<!std::is_scalar<T>::value>>
void sirius::Matching_coefficients::generate ( Atom const &  atom__,
sddk::mdarray< T, 2 > &  alm__ 
) const
inline

Generate plane-wave matching coefficients for the radial solutions of a given atom.

Parameters
[in]atomAtom, for which matching coefficients are generated.
[out]almArray of matching coefficients with dimension indices \( ({\bf G+k}, \xi) \).

Definition at line 193 of file matching_coefficients.hpp.

◆ gkvec()

auto const & sirius::Matching_coefficients::gkvec ( ) const
inline

Definition at line 296 of file matching_coefficients.hpp.

Member Data Documentation

◆ unit_cell_

Unit_cell const& sirius::Matching_coefficients::unit_cell_
private

Description of the unit cell.

Definition at line 57 of file matching_coefficients.hpp.

◆ gkvec_

fft::Gvec const& sirius::Matching_coefficients::gkvec_
private

Description of the G+k vectors.

Definition at line 60 of file matching_coefficients.hpp.

◆ gkvec_len_

std::vector<double> sirius::Matching_coefficients::gkvec_len_
private

Definition at line 62 of file matching_coefficients.hpp.

◆ gkvec_ylm_

sddk::mdarray<std::complex<double>, 2> sirius::Matching_coefficients::gkvec_ylm_
private

Spherical harmonics Ylm(theta, phi) of the G+k vectors.

Definition at line 65 of file matching_coefficients.hpp.

◆ alm_b_

sddk::mdarray<std::complex<double>, 4> sirius::Matching_coefficients::alm_b_
private

Precomputed values for the linear equations for matching coefficients.

Definition at line 68 of file matching_coefficients.hpp.


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