ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
|
A class which can convert a function of "r" to the corresponding linear superposition of plane waves (real space to reciprocal space) or convert a linear superposition of plane waves to the function of "r" (reciprocal to real). More...
#include <pw_basis.h>
Public Types | |
using | resmem_int_op = base_device::memory::resize_memory_op< int, base_device::DEVICE_GPU > |
using | delmem_int_op = base_device::memory::delete_memory_op< int, base_device::DEVICE_GPU > |
using | syncmem_int_h2d_op = base_device::memory::synchronize_memory_op< int, base_device::DEVICE_GPU, base_device::DEVICE_CPU > |
Public Member Functions | |
PW_Basis () | |
PW_Basis (std::string device_, std::string precision_) | |
virtual | ~PW_Basis () |
void | initmpi (const int poolnproc_in, const int poolrank_in, MPI_Comm pool_world_in) |
virtual void | initgrids (const double lat0_in, const ModuleBase::Matrix3 latvec_in, const double gridecut) |
virtual void | initgrids (const double lat0_in, const ModuleBase::Matrix3 latvec_in, const int nx_in, int ny_in, int nz_in) |
void | initparameters (const bool gamma_only_in, const double pwecut_in, const int distribution_type_in=1, const bool xprime_in=true) |
void | setfullpw (const bool inpt_full_pw=false, const int inpt_full_pw_dim=0) |
void | setuptransform () |
void | collect_local_pw () |
void | collect_uniqgg () |
template<typename FPTYPE > | |
void | real2recip (const FPTYPE *in, std::complex< FPTYPE > *out, const bool add=false, const FPTYPE factor=1.0) const |
transform real space to reciprocal space | |
template<typename FPTYPE > | |
void | real2recip (const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const bool add=false, const FPTYPE factor=1.0) const |
transform real space to reciprocal space | |
template<typename FPTYPE > | |
void | recip2real (const std::complex< FPTYPE > *in, FPTYPE *out, const bool add=false, const FPTYPE factor=1.0) const |
transform reciprocal space to real space | |
template<typename FPTYPE > | |
void | recip2real (const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const bool add=false, const FPTYPE factor=1.0) const |
transform reciprocal space to real space | |
template<typename FPTYPE > | |
void | real2recip_gpu (const FPTYPE *in, std::complex< FPTYPE > *out, const bool add=false, const FPTYPE factor=1.0) const |
template<typename FPTYPE > | |
void | real2recip_gpu (const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const bool add=false, const FPTYPE factor=1.0) const |
template<typename FPTYPE > | |
void | recip2real_gpu (const std::complex< FPTYPE > *in, FPTYPE *out, const bool add=false, const FPTYPE factor=1.0) const |
template<typename FPTYPE > | |
void | recip2real_gpu (const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const bool add=false, const FPTYPE factor=1.0) const |
template<typename TK , typename TR , typename Device , typename std::enable_if<!std::is_same< TK, typename GetTypeReal< TK >::type >::value &&(std::is_same< TR, typename GetTypeReal< TK >::type >::value||std::is_same< TR, TK >::value) &&std::is_same< Device, base_device::DEVICE_CPU >::value, int >::type = 0> | |
void | recip_to_real (TK *in, TR *out, const bool add=false, const typename GetTypeReal< TK >::type factor=1.0) const |
Converts data from reciprocal space to real space on Cpu. | |
template<typename TK , typename TR , typename Device , typename std::enable_if<!std::is_same< TK, typename GetTypeReal< TK >::type >::value &&(std::is_same< TR, typename GetTypeReal< TK >::type >::value||std::is_same< TR, TK >::value) &&std::is_same< Device, base_device::DEVICE_GPU >::value, int >::type = 0> | |
void | recip_to_real (TK *in, TR *out, const bool add=false, const typename GetTypeReal< TK >::type factor=1.0) const |
Converts data from reciprocal space (Fourier space) to real space. | |
template<typename TR , typename TK , typename Device , typename std::enable_if<!std::is_same< TK, typename GetTypeReal< TK >::type >::value &&(std::is_same< TR, typename GetTypeReal< TK >::type >::value||std::is_same< TR, TK >::value) &&std::is_same< Device, base_device::DEVICE_CPU >::value, int >::type = 0> | |
void | real_to_recip (TR *in, TK *out, const bool add=false, const typename GetTypeReal< TK >::type factor=1.0) const |
Converts data from real space to reciprocal space (Fourier space). | |
template<typename TR , typename TK , typename Device , typename std::enable_if<!std::is_same< TK, typename GetTypeReal< TK >::type >::value &&(std::is_same< TR, typename GetTypeReal< TK >::type >::value||std::is_same< TR, TK >::value) &&std::is_same< Device, base_device::DEVICE_GPU >::value, int >::type = 0> | |
void | real_to_recip (TR *in, TK *out, const bool add=false, const typename GetTypeReal< TK >::type factor=1.0) const |
void | getfftixy2is (int *fftixy2is) const |
void | set_device (std::string device_) |
void | set_precision (std::string precision_) |
std::string | get_device () const |
std::string | get_precision () const |
Public Attributes | |
std::string | classname |
MPI_Comm | pool_world =MPI_COMM_NULL |
int * | ig2isz =nullptr |
int * | ig2ixyz_gpu = nullptr |
int * | istot2ixy =nullptr |
int * | is2fftixy =nullptr |
int * | d_is2fftixy = nullptr |
int * | fftixy2ip =nullptr |
int | nst =0 |
int * | nst_per =nullptr |
int | nstnz =0 |
int | nstot =0 |
int | npw =0 |
int * | npw_per =nullptr |
int | npwtot =0 |
int | nrxx =0 |
int * | startz =nullptr |
int * | numz =nullptr |
int * | numg =nullptr |
int * | numr =nullptr |
int * | startg =nullptr |
int * | startr =nullptr |
int | startz_current =0 |
int | nplane =0 |
ModuleBase::Vector3< double > * | gdirect =nullptr |
ModuleBase::Vector3< double > * | gcar =nullptr |
double * | gg =nullptr |
int | ig_gge0 =-1 |
int | ngg =0 |
int * | ig2igg =nullptr |
double * | gg_uniq =nullptr |
bool | gamma_only = false |
only half g are used. | |
bool | full_pw = false |
double | ggecut = 0.0 |
Energy cut off for g^2/2 = ecutwfc(Ry)*lat0^2/4pi^2, unit in 1/lat0^2. | |
double | gridecut_lat = 0.0 |
Energy cut off for all fft grids = ecutrho(Ry)*lat0^2/4pi^2, unit in 1/lat0^2. | |
double | lat0 = 1 |
unit length for lattice, unit in bohr | |
double | tpiba = 1 |
2pi/lat0 | |
double | tpiba2 = 1 |
4pi^2/lat0^2 | |
ModuleBase::Matrix3 | latvec |
Unitcell lattice vectors, unit in lat0. | |
ModuleBase::Matrix3 | G |
reciprocal lattice vector, unit in 1/lat0 | |
ModuleBase::Matrix3 | GT |
traspose of G | |
ModuleBase::Matrix3 | GGT |
GGT = G*GT. | |
double | omega = 1.0 |
volume of the cell | |
int | distribution_type = 1 |
distribution method | |
int | full_pw_dim = 0 |
int | poolnproc = 1 |
int | poolrank = 0 |
int | fftnx =0 |
int | fftny =0 |
int | fftnz =0 |
int | fftnxyz =0 |
int | fftnxy =0 |
int | nx =0 |
int | ny =0 |
int | nz =0 |
int | nxyz =0 |
int | nxy =0 |
int | liy =0 |
int | riy =0 |
int | lix =0 |
int | rix =0 |
bool | xprime = true |
int | ng_xeq0 = 0 |
int | nmaxgr = 0 |
FFT_Bundle | fft_bundle |
Protected Member Functions | |
void | distribute_g () |
distribute plane waves to different cores | |
virtual void | distribute_r () |
distribute real-space grids to different processors | |
void | getstartgr () |
void | distribution_method1 () |
void | divide_sticks_1 (int *st_i, int *st_j, int *st_length) |
(3-1) Distribute sticks to cores according to the number of plane waves. | |
void | distribution_method2 () |
void | divide_sticks_2 () |
(2) Devide the sticks to each core according to the number of sticks Sticks are in the order of ixy increasing. | |
void | count_pw_st (int *st_length2D, int *st_bottom2D) |
(1) We count the total number of planewaves (tot_npw) and sticks (this->nstot) here. | |
void | get_ig2isz_is2fftixy (int *st_bottom, int *st_length) |
(5) Construct ig2isz, and is2fftixy. | |
void | collect_st (int *st_length2D, int *st_bottom2D, int *st_i, int *st_j, int *st_length) |
(2) Collect the x, y indexs, length of the sticks. | |
void | get_istot2ixy (int *st_i, int *st_j) |
(3-2) Rearrange sticks in the order of the ip of core increasing, in each core, sticks are sorted in the order of ixy increasing. | |
void | create_maps (int *st_length2D) |
(3) Create the maps from ixy to ip, istot, and from istot to ixy, and construt npw_per simultaneously. | |
template<typename T > | |
void | gatherp_scatters (std::complex< T > *in, std::complex< T > *out) const |
gather planes and scatter sticks | |
template<typename T > | |
void | gathers_scatterp (std::complex< T > *in, std::complex< T > *out) const |
gather sticks and scatter planes | |
Protected Attributes | |
int * | startnsz_per =nullptr |
std::string | device = "cpu" |
cpu or gpu | |
std::string | precision = "double" |
single, double, mixing | |
bool | double_data_ = true |
if has double data | |
bool | float_data_ = false |
if has float data | |
A class which can convert a function of "r" to the corresponding linear superposition of plane waves (real space to reciprocal space) or convert a linear superposition of plane waves to the function of "r" (reciprocal to real).
Math: plane waves: <r|g>=1/sqrt(V) * exp(igr) f(r) = 1/sqrt(V) * \sum_g{c(g)*exp(igr)} c(g) = \int f(r)*exp(-igr) dr USAGE: ModulePW::PW_Basis pwtest; 0. init mpi for PW_Basis pwtest.inimpi(nproc_in_pool,rank_in_pool,POOL_WORLD);
using ModulePW::PW_Basis::delmem_int_op = base_device::memory::delete_memory_op<int, base_device::DEVICE_GPU> |
using ModulePW::PW_Basis::resmem_int_op = base_device::memory::resize_memory_op<int, base_device::DEVICE_GPU> |
using ModulePW::PW_Basis::syncmem_int_h2d_op = base_device::memory::synchronize_memory_op<int, base_device::DEVICE_GPU, base_device::DEVICE_CPU> |
ModulePW::PW_Basis::PW_Basis | ( | ) |
ModulePW::PW_Basis::PW_Basis | ( | std::string | device_, |
std::string | precision_ | ||
) |
|
virtual |
void ModulePW::PW_Basis::collect_local_pw | ( | ) |
Collect planewaves on current core, and construct gg, gdirect, gcar according to ig2isz and is2fftixy. known: ig2isz, is2fftixy output: gg, gdirect, gcar
|
protected |
(2) Collect the x, y indexs, length of the sticks.
Firstly, we scan the area and construct temp_st_*. Then, as we will distribute the longest sticks preferentially in Step(3), we will sort temp_st_length from largest to smallest, and reaarange st_* to the same order.
in | tot_npw, this->nstot, st_length2D, st_bottom2D |
out | st_i, st_j, st_length |
void ModulePW::PW_Basis::collect_uniqgg | ( | ) |
Collect modulus of planewaves on current cores known: ig2isz, is2fftixy output: ig2igg, gg_uniq, ngg
|
protected |
(1) We count the total number of planewaves (tot_npw) and sticks (this->nstot) here.
Meanwhile, we record the number of planewaves on (x, y) in st_length2D, and store the smallest z-coordinate of each stick in st_bottom2D, so that we can scan a much smaller area in step(2).
in | fftnx, fftny, nz, ggecut, GGT |
out | tot_npw, this->nstot, st_length2D, st_bottom2D, this->riy, this->liy |
|
protected |
(3) Create the maps from ixy to ip, istot, and from istot to ixy, and construt npw_per simultaneously.
in | st_length2D |
out | this->fftixy2ip, this->istot2ixy, npw_per |
|
protected |
distribute plane waves to different cores
in | G, GT, GGT, fftnx, fftny, nz, poolnproc, poolrank, ggecut |
out | ig2isz[ig], istot2ixy[is], is2fftixy[is], fftixy2ip[ixy], gg[ig], gcar[ig], gdirect[ig], nst, nstot |
|
protectedvirtual |
distribute real-space grids to different processors
in | nx, ny, nz, poolnproc, poolrank |
out | nrxx, startz, numz |
Reimplemented in ModulePW::PW_Basis_Big, and ModulePW::PW_Basis_K_Big.
|
protected |
(1) Count the total number of planewaves (tot_npw) and sticks (this->nstot).
|
protected |
(1) Count the total number of planewaves (tot_npw) and sticks (this->nstot).
|
protected |
(3-1) Distribute sticks to cores according to the number of plane waves.
We have rearranged sticks in the order of length decreasing, so that we will distribute the longest stick preferentially here. For each stick, we find the core that contains the least planewaves firstly, and distribute the stick to it, then update npw_per, this->fftixy2ip, and this->startnsz_per.
in | tot_npw, this->nstot, st_i, st_j, st_length |
out | npw_per, nst_per, this->fftixy2ip, this->startnsz_per |
|
protected |
(2) Devide the sticks to each core according to the number of sticks Sticks are in the order of ixy increasing.
in | this->nstot, this->poolnproc |
out | nst_per, this->startnsz_per |
|
protected |
gather planes and scatter sticks
in | (nplane,fftny,fftnx) |
out | (nz,nst) |
|
protected |
gather sticks and scatter planes
in | (nz,nst) |
out | (nplane,fftny,fftnx) |
|
inline |
|
protected |
(5) Construct ig2isz, and is2fftixy.
is2fftixy contains the x-coordinate and y-coordinate of sticks on current core. ig2isz contains the z-coordinate of planewaves on current core. We will scan all the sticks and find the planewaves on them, then store the information into ig2isz and is2fftixy.
in | this->nstot, st_bottom2D, st_length2D |
out | ig2isz, is2fftixy |
|
protected |
(3-2) Rearrange sticks in the order of the ip of core increasing, in each core, sticks are sorted in the order of ixy increasing.
(st_start + st_move) is the new index of sticks. Then get istot2ixy (istot2ixy[is]: iy + ix * fftny of is^th stick among all sticks) on the first core
in | this->nstot, st_i, st_j, this->startnsz_per |
out | istot2ixy |
|
inline |
void ModulePW::PW_Basis::getfftixy2is | ( | int * | fftixy2is | ) | const |
|
protected |
|
virtual |
Init the grids for FFT Input: lattice vectors of the cell, Energy cut off for G^2/2 Output: fftnx, fftny, fftnz, fftnxyz, latvec, G, GT, GGT
Reimplemented in ModulePW::PW_Basis_Big.
|
virtual |
Init the grids for FFT Input: lattice vectors of the cell, nx, ny, nz Output: nx, ny, nz, nxyz, latvec, G, GT, GGT
Reimplemented in ModulePW::PW_Basis_Big.
void ModulePW::PW_Basis::initmpi | ( | const int | poolnproc_in, |
const int | poolrank_in, | ||
MPI_Comm | pool_world_in | ||
) |
void ModulePW::PW_Basis::initparameters | ( | const bool | gamma_only_in, |
const double | pwecut_in, | ||
const int | distribution_type_in = 1 , |
||
const bool | xprime_in = true |
||
) |
void ModulePW::PW_Basis::real2recip | ( | const FPTYPE * | in, |
std::complex< FPTYPE > * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
transform real space to reciprocal space
c(g)=\int dr*f(r)*exp(-ig*r) Here we calculate c(g)
in | (nplane,ny,nx), double data |
out | (nz, ns), std::complex<double> data |
void ModulePW::PW_Basis::real2recip | ( | const std::complex< FPTYPE > * | in, |
std::complex< FPTYPE > * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
transform real space to reciprocal space
c(g)=\int dr*f(r)*exp(-ig*r) Here we calculate c(g)
in | (nplane,ny,nx), std::complex<double> data |
out | (nz, ns), std::complex<double> data |
void ModulePW::PW_Basis::real2recip_gpu | ( | const FPTYPE * | in, |
std::complex< FPTYPE > * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
void ModulePW::PW_Basis::real2recip_gpu | ( | const std::complex< FPTYPE > * | in, |
std::complex< FPTYPE > * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
|
inline |
Converts data from real space to reciprocal space (Fourier space).
This function handles the conversion of data from real space to reciprocal space (typically after performing a Fourier transform), supporting scenarios where the input is of a fundamental type (e.g., float, double) and the output is of a complex type.
FPTYPE | The underlying fundamental type of the input data (e.g., float, double). SFINAE constraint ensures that FPTYPE is a fundamental type, not a complex type. |
Device | The device type, which must be base_device::DEVICE_CPU. |
std::enable_if<std::is_same<FPTYPE,typename | GetTypeReal<FPTYPE>::type>::value, int>::type SFINAE constraint to ensure that FPTYPE is a fundamental type. |
std::enable_if<std::is_same<Device,base_device::DEVICE_CPU>::value,int>::type | SFINAE constraint to ensure that Device is base_device::DEVICE_CPU. |
in | Pointer to the input data array in real space, of type FPTYPE*. |
out | Pointer to the output data array in reciprocal space, of type std::complex<FPTYPE>*. |
ik | Optional parameter, default value 0, representing some index or identifier. |
add | Optional boolean flag, default value false, indicating whether to add the result to the existing values in the output array. |
factor | Optional scaling factor, default value 1.0, applied to the output values. |
|
inline |
void ModulePW::PW_Basis::recip2real | ( | const std::complex< FPTYPE > * | in, |
FPTYPE * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
transform reciprocal space to real space
f(r)=1/V * \sum_{g} c(g)*exp(ig*r) Here we calculate f(r)
in | (nz,ns), std::complex<double> |
out | (nplane, ny, nx), double |
void ModulePW::PW_Basis::recip2real | ( | const std::complex< FPTYPE > * | in, |
std::complex< FPTYPE > * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
transform reciprocal space to real space
f(r)=1/V * \sum_{g} c(g)*exp(ig*r) Here we calculate f(r)
in | (nz,ns), std::complex<double> |
out | (nplane, ny, nx), std::complex<double> |
void ModulePW::PW_Basis::recip2real_gpu | ( | const std::complex< FPTYPE > * | in, |
FPTYPE * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
void ModulePW::PW_Basis::recip2real_gpu | ( | const std::complex< FPTYPE > * | in, |
std::complex< FPTYPE > * | out, | ||
const bool | add = false , |
||
const FPTYPE | factor = 1.0 |
||
) | const |
|
inline |
Converts data from reciprocal space to real space on Cpu.
This function handles the conversion of data from reciprocal space (Fourier space) to real space. It supports complex types as input. The output can be either the same fundamental type or the underlying real type of a complex type.
FPTYPE | The type of the input data, which can only be a compelx type (e.g., std::complex<float>, std::complex<double>) |
Device | The device type, must be base_device::DEVICE_CPU. |
std::enable_if<!std::is_same<FPTYPE,typename | GetTypeReal<FPTYPE>::type>::value, int>::type SFINAE constraint to ensure that FPTYPE is a complex type. |
std::enable_if<std::is_same<Device,base_device::DEVICE_CPU>::value,int>::type | SFINAE constraint to ensure that Device is base_device::DEVICE_CPU. |
in | Pointer to the input data array in reciprocal space. |
out | Pointer to the output data array in real space. If FPTYPE is a complex type, this should point to an array of the underlying real type. |
add | Boolean flag indicating whether to add the result to the existing values in the output array. |
factor | A scaling factor applied to the output values. |
|
inline |
Converts data from reciprocal space (Fourier space) to real space.
This function handles the conversion of data from reciprocal space (typically after a Fourier transform) to real space, supporting scenarios where the input is of a complex type and the output is of the underlying fundamental type.
FPTYPE | The underlying fundamental type of the input data (e.g., float, double). Note that the actual type passed should be std::complex<FPTYPE>. |
Device | The device type, which can be any supported device type (e.g., base_device::DEVICE_CPU). |
in | Pointer to the input data array in reciprocal space, of type std::complex<FPTYPE>*. |
out | Pointer to the output data array in real space, of type FPTYPE*. |
add | Optional boolean flag, default value false, indicating whether to add the result to the existing values in the output array. |
factor | Optional scaling factor, default value 1.0, applied to the output values. |
void ModulePW::PW_Basis::set_device | ( | std::string | device_ | ) |
void ModulePW::PW_Basis::set_precision | ( | std::string | precision_ | ) |
void ModulePW::PW_Basis::setfullpw | ( | const bool | inpt_full_pw = false , |
const int | inpt_full_pw_dim = 0 |
||
) |
void ModulePW::PW_Basis::setuptransform | ( | ) |
distribute plane wave basis and real-space grids to different processors set up maps for fft and create arrays for MPI_Alltoall set up ffts
std::string ModulePW::PW_Basis::classname |
int * ModulePW::PW_Basis::d_is2fftixy = nullptr |
|
protected |
cpu or gpu
int ModulePW::PW_Basis::distribution_type = 1 |
distribution method
|
protected |
if has double data
FFT_Bundle ModulePW::PW_Basis::fft_bundle |
int* ModulePW::PW_Basis::fftixy2ip =nullptr |
int ModulePW::PW_Basis::fftnx =0 |
int ModulePW::PW_Basis::fftnxy =0 |
int ModulePW::PW_Basis::fftnxyz =0 |
int ModulePW::PW_Basis::fftny =0 |
int ModulePW::PW_Basis::fftnz =0 |
|
protected |
if has float data
bool ModulePW::PW_Basis::full_pw = false |
If set to 1, ecut will be ignored while collecting planewaves, so that all planewaves will be used. !! Note this parameter is not used in PW_BASIS_K !! sunliang added 2022-08-30.
int ModulePW::PW_Basis::full_pw_dim = 0 |
If full_pw = 1, the dimention of FFT will be testricted to be (0) either odd or even; (1) odd only; (2) even only. sunliang added 2022-08-30.
ModuleBase::Matrix3 ModulePW::PW_Basis::G |
reciprocal lattice vector, unit in 1/lat0
bool ModulePW::PW_Basis::gamma_only = false |
only half g are used.
ModuleBase::Vector3<double>* ModulePW::PW_Basis::gcar =nullptr |
ModuleBase::Vector3<double>* ModulePW::PW_Basis::gdirect =nullptr |
double* ModulePW::PW_Basis::gg =nullptr |
double* ModulePW::PW_Basis::gg_uniq =nullptr |
double ModulePW::PW_Basis::ggecut = 0.0 |
Energy cut off for g^2/2 = ecutwfc(Ry)*lat0^2/4pi^2, unit in 1/lat0^2.
ModuleBase::Matrix3 ModulePW::PW_Basis::GGT |
double ModulePW::PW_Basis::gridecut_lat = 0.0 |
Energy cut off for all fft grids = ecutrho(Ry)*lat0^2/4pi^2, unit in 1/lat0^2.
ModuleBase::Matrix3 ModulePW::PW_Basis::GT |
traspose of G
int* ModulePW::PW_Basis::ig2igg =nullptr |
int* ModulePW::PW_Basis::ig2isz =nullptr |
int* ModulePW::PW_Basis::ig2ixyz_gpu = nullptr |
int ModulePW::PW_Basis::ig_gge0 =-1 |
int* ModulePW::PW_Basis::is2fftixy =nullptr |
int* ModulePW::PW_Basis::istot2ixy =nullptr |
double ModulePW::PW_Basis::lat0 = 1 |
unit length for lattice, unit in bohr
ModuleBase::Matrix3 ModulePW::PW_Basis::latvec |
Unitcell lattice vectors, unit in lat0.
int ModulePW::PW_Basis::lix =0 |
int ModulePW::PW_Basis::liy =0 |
int ModulePW::PW_Basis::ng_xeq0 = 0 |
int ModulePW::PW_Basis::ngg =0 |
int ModulePW::PW_Basis::nmaxgr = 0 |
int ModulePW::PW_Basis::nplane =0 |
int ModulePW::PW_Basis::npw =0 |
int* ModulePW::PW_Basis::npw_per =nullptr |
int ModulePW::PW_Basis::npwtot =0 |
int ModulePW::PW_Basis::nrxx =0 |
int ModulePW::PW_Basis::nst =0 |
int* ModulePW::PW_Basis::nst_per =nullptr |
int ModulePW::PW_Basis::nstnz =0 |
int ModulePW::PW_Basis::nstot =0 |
int* ModulePW::PW_Basis::numg =nullptr |
int* ModulePW::PW_Basis::numr =nullptr |
int* ModulePW::PW_Basis::numz =nullptr |
int ModulePW::PW_Basis::nx =0 |
int ModulePW::PW_Basis::nxy =0 |
int ModulePW::PW_Basis::nxyz =0 |
int ModulePW::PW_Basis::ny =0 |
int ModulePW::PW_Basis::nz =0 |
double ModulePW::PW_Basis::omega = 1.0 |
volume of the cell
MPI_Comm ModulePW::PW_Basis::pool_world =MPI_COMM_NULL |
int ModulePW::PW_Basis::poolnproc = 1 |
int ModulePW::PW_Basis::poolrank = 0 |
|
protected |
single, double, mixing
int ModulePW::PW_Basis::rix =0 |
int ModulePW::PW_Basis::riy =0 |
int* ModulePW::PW_Basis::startg =nullptr |
|
protected |
int* ModulePW::PW_Basis::startr =nullptr |
int* ModulePW::PW_Basis::startz =nullptr |
int ModulePW::PW_Basis::startz_current =0 |
double ModulePW::PW_Basis::tpiba = 1 |
2pi/lat0
double ModulePW::PW_Basis::tpiba2 = 1 |
4pi^2/lat0^2
bool ModulePW::PW_Basis::xprime = true |