SIRIUS 7.5.0
Electronic structure library and applications
|
Compute gradient of beta-projectors over atomic positions \( d \langle {\bf G+k} | \beta \rangle / d \tau_{\alpha} \). More...
#include <beta_projectors_gradient.hpp>
Inherits sirius::Beta_projectors_base< T >.
Public Member Functions | |
Beta_projectors_gradient (Simulation_context &ctx__, fft::Gvec const &gkvec__, Beta_projectors< T > &beta__) | |
Public Member Functions inherited from sirius::Beta_projectors_base< T > | |
Beta_projectors_base (Simulation_context &ctx__, fft::Gvec const &gkvec__, int N__) | |
Beta_projector_generator< T > | make_generator () const |
Beta_projector_generator< T > | make_generator (sddk::device_t pu) const |
Beta_projector_generator< T > | make_generator (sddk::memory_t mem) const |
auto const & | ctx () const |
auto | num_gkvec_loc () const |
int | num_total_beta () const |
int | num_comp () const |
auto const & | unit_cell () const |
auto & | pw_coeffs_t (int ig__, int n__, int j__) |
auto | pw_coeffs_t (int j__) |
auto const & | pw_coeffs_a () const |
int | num_beta_t () const |
int | num_chunks () const |
int | max_num_beta () const |
int | nrows () const |
const mpi::Communicator & | comm () const |
Private Member Functions | |
void | generate_pw_coefs_t (Beta_projectors< T > &beta__) |
Additional Inherited Members | |
Protected Member Functions inherited from sirius::Beta_projectors_base< T > | |
void | split_in_chunks () |
Split beta-projectors into chunks. More... | |
Protected Attributes inherited from sirius::Beta_projectors_base< T > | |
Simulation_context & | ctx_ |
fft::Gvec const & | gkvec_ |
List of G+k vectors. More... | |
sddk::mdarray< double, 2 > | gkvec_coord_ |
Coordinates of G+k vectors used by GPU kernel. More... | |
int | N_ |
Number of different components: 1 for beta-projectors, 3 for gradient, 9 for strain derivatives. More... | |
sddk::mdarray< std::complex< T >, 3 > | pw_coeffs_t_ |
Phase-factor independent coefficients of |beta> functions for atom types. More... | |
bool | reallocate_pw_coeffs_t_on_gpu_ {true} |
sddk::matrix< std::complex< T > > | pw_coeffs_a_ |
Set of beta PW coefficients for a chunk of atoms. More... | |
sddk::matrix< std::complex< T > > | beta_pw_all_atoms_ |
Set of beta PW coefficients for all atoms. More... | |
std::vector< beta_chunk_t > | beta_chunks_ |
int | max_num_beta_ |
int | num_total_beta_ |
total number of beta-projectors (=number of columns) More... | |
int | num_beta_t_ |
Total number of beta-projectors among atom types. More... | |
Compute gradient of beta-projectors over atomic positions \( d \langle {\bf G+k} | \beta \rangle / d \tau_{\alpha} \).
Definition at line 35 of file beta_projectors_gradient.hpp.
|
inline |
Definition at line 56 of file beta_projectors_gradient.hpp.
|
inlineprivate |
Definition at line 38 of file beta_projectors_gradient.hpp.