ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
ModulePW::PW_Basis Class Reference

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>

Inheritance diagram for ModulePW::PW_Basis:
Collaboration diagram for ModulePW::PW_Basis:

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
 

Detailed Description

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).

Author
qianrui, Sunliang on 2021-10-15

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);

  1. setup FFT grids for PW_Basis pwtest.initgrids(lat0,latvec,gridecut); pwtest.initgrids(lat0,latvec,N1,N2,N3); //double lat0: unit length, (unit: bohr) //ModuleBase::Matrix3 latvec: lattice vector, (unit: lat0), e.g. ModuleBase::Matrix3 latvec(1, 1, 0, 0, 2, 0, 0, 0, 2); //double gridecut: cutoff energy to generate FFT grids, (unit: Ry) //int N1,N2,N3: FFT grids
  2. init parameters pwtest.initparameters(gamma_only,ggecut,dividemthd); //bool gamma_only: if use gamma_only //double ggecut: cutoff kinetic energy for planewaves,(unit in Ry) G^2 < ggecut //int dividemthd: method to divide planewaves to different cores
  3. Setup transforms from real space to reciprocal space or from reciprocal space to real space. pwtest.setuptransform(); pwtest.recip2real(rhog,rhor); //rhog to rhor pwtest.real2recip(rhor,rhog); //rhor to rhog
  4. Generate the wave vector for planewaves pwtest.collect_local_pw(); //then we can use pwtest.gg, pwtest.gdirect, pwtest.gcar, (unit in lat0^-1 or lat0^-2)

Member Typedef Documentation

◆ delmem_int_op

◆ resmem_int_op

◆ syncmem_int_h2d_op

using ModulePW::PW_Basis::syncmem_int_h2d_op = base_device::memory::synchronize_memory_op<int, base_device::DEVICE_GPU, base_device::DEVICE_CPU>

Constructor & Destructor Documentation

◆ PW_Basis() [1/2]

ModulePW::PW_Basis::PW_Basis ( )

◆ PW_Basis() [2/2]

ModulePW::PW_Basis::PW_Basis ( std::string  device_,
std::string  precision_ 
)
Here is the call graph for this function:

◆ ~PW_Basis()

ModulePW::PW_Basis::~PW_Basis ( )
virtual

Member Function Documentation

◆ collect_local_pw()

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

Here is the call graph for this function:

◆ collect_st()

void ModulePW::PW_Basis::collect_st ( int *  st_length2D,
int *  st_bottom2D,
int *  st_i,
int *  st_j,
int *  st_length 
)
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.

Parameters
intot_npw, this->nstot, st_length2D, st_bottom2D
outst_i, st_j, st_length
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ collect_uniqgg()

void ModulePW::PW_Basis::collect_uniqgg ( )

Collect modulus of planewaves on current cores known: ig2isz, is2fftixy output: ig2igg, gg_uniq, ngg

Here is the call graph for this function:
Here is the caller graph for this function:

◆ count_pw_st()

void ModulePW::PW_Basis::count_pw_st ( int *  st_length2D,
int *  st_bottom2D 
)
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).

Parameters
infftnx, fftny, nz, ggecut, GGT
outtot_npw, this->nstot, st_length2D, st_bottom2D, this->riy, this->liy
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp, and /home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method2.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ create_maps()

void ModulePW::PW_Basis::create_maps ( int *  st_length2D)
protected

(3) Create the maps from ixy to ip, istot, and from istot to ixy, and construt npw_per simultaneously.

Parameters
inst_length2D
outthis->fftixy2ip, this->istot2ixy, npw_per
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method2.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ distribute_g()

void ModulePW::PW_Basis::distribute_g ( )
protected

distribute plane waves to different cores

Parameters
inG, GT, GGT, fftnx, fftny, nz, poolnproc, poolrank, ggecut
outig2isz[ig], istot2ixy[is], is2fftixy[is], fftixy2ip[ixy], gg[ig], gcar[ig], gdirect[ig], nst, nstot
Here is the call graph for this function:
Here is the caller graph for this function:

◆ distribute_r()

void ModulePW::PW_Basis::distribute_r ( )
protectedvirtual

distribute real-space grids to different processors

Parameters
innx, ny, nz, poolnproc, poolrank
outnrxx, startz, numz

Reimplemented in ModulePW::PW_Basis_Big, and ModulePW::PW_Basis_K_Big.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ distribution_method1()

void ModulePW::PW_Basis::distribution_method1 ( )
protected

(1) Count the total number of planewaves (tot_npw) and sticks (this->nstot).

Note
the funcion here is defined in pw_distributeg.cpp Actually we will scan [(2 * ibox[0] + 1) * (2 * ibox[1] + 1)] points on x-y plane, but we define st_length2D with (fftny * fftnx) points here, because the diameter of the sphere should be shorter than the sides of the cube. calculate this->nstot and this->npwtot, liy, riy
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ distribution_method2()

void ModulePW::PW_Basis::distribution_method2 ( )
protected

(1) Count the total number of planewaves (tot_npw) and sticks (this->nstot).

Note
the funcion here is defined in pw_distributeg.cpp Actually we will scan [(2 * ibox[0] + 1) * (2 * ibox[1] + 1)] points on x-y plane, but we define st_length2D with (fftny * fftnx) points here, because the diameter of the sphere should be shorter than the sides of the cube. calculate this->nstot and this->npwtot, liy, riy
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method2.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ divide_sticks_1()

void ModulePW::PW_Basis::divide_sticks_1 ( int *  st_i,
int *  st_j,
int *  st_length 
)
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.

Parameters
intot_npw, this->nstot, st_i, st_j, st_length
outnpw_per, nst_per, this->fftixy2ip, this->startnsz_per
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ divide_sticks_2()

void ModulePW::PW_Basis::divide_sticks_2 ( )
protected

(2) Devide the sticks to each core according to the number of sticks Sticks are in the order of ixy increasing.

Parameters
inthis->nstot, this->poolnproc
outnst_per, this->startnsz_per
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method2.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gatherp_scatters()

template<typename T >
void ModulePW::PW_Basis::gatherp_scatters ( std::complex< T > *  in,
std::complex< T > *  out 
) const
protected

gather planes and scatter sticks

Parameters
in(nplane,fftny,fftnx)
out(nz,nst)
Note
in and out should be in different places
in[] will be changed
Here is the caller graph for this function:

◆ gathers_scatterp()

template<typename T >
void ModulePW::PW_Basis::gathers_scatterp ( std::complex< T > *  in,
std::complex< T > *  out 
) const
protected

gather sticks and scatter planes

Parameters
in(nz,nst)
out(nplane,fftny,fftnx)
Note
in and out should be in different places
in[] will be changed
Here is the caller graph for this function:

◆ get_device()

std::string ModulePW::PW_Basis::get_device ( ) const
inline
Here is the caller graph for this function:

◆ get_ig2isz_is2fftixy()

void ModulePW::PW_Basis::get_ig2isz_is2fftixy ( int *  st_bottom2D,
int *  st_length2D 
)
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.

Parameters
inthis->nstot, st_bottom2D, st_length2D
outig2isz, is2fftixy
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp, and /home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method2.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_istot2ixy()

void ModulePW::PW_Basis::get_istot2ixy ( int *  st_i,
int *  st_j 
)
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

Parameters
inthis->nstot, st_i, st_j, this->startnsz_per
outistot2ixy
Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_precision()

std::string ModulePW::PW_Basis::get_precision ( ) const
inline
Here is the caller graph for this function:

◆ getfftixy2is()

void ModulePW::PW_Basis::getfftixy2is ( int *  fftixy2is) const
Here is the caller graph for this function:

◆ getstartgr()

void ModulePW::PW_Basis::getstartgr ( )
protected
Here is the caller graph for this function:

◆ initgrids() [1/2]

void ModulePW::PW_Basis::initgrids ( const double  lat0_in,
const ModuleBase::Matrix3  latvec_in,
const double  gridecut 
)
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.

Here is the call graph for this function:

◆ initgrids() [2/2]

void ModulePW::PW_Basis::initgrids ( const double  lat0_in,
const ModuleBase::Matrix3  latvec_in,
const int  nx_in,
int  ny_in,
int  nz_in 
)
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.

Here is the call graph for this function:

◆ initmpi()

void ModulePW::PW_Basis::initmpi ( const int  poolnproc_in,
const int  poolrank_in,
MPI_Comm  pool_world_in 
)

◆ initparameters()

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 
)

◆ real2recip() [1/2]

template<typename FPTYPE >
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)

Parameters
in(nplane,ny,nx), double data
out(nz, ns), std::complex<double> data
Here is the call graph for this function:

◆ real2recip() [2/2]

template<typename FPTYPE >
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)

Parameters
in(nplane,ny,nx), std::complex<double> data
out(nz, ns), std::complex<double> data
Here is the call graph for this function:

◆ real2recip_gpu() [1/2]

template<typename FPTYPE >
void ModulePW::PW_Basis::real2recip_gpu ( const FPTYPE *  in,
std::complex< FPTYPE > *  out,
const bool  add = false,
const FPTYPE  factor = 1.0 
) const
Here is the caller graph for this function:

◆ real2recip_gpu() [2/2]

template<typename FPTYPE >
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

◆ real_to_recip() [1/2]

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 ModulePW::PW_Basis::real_to_recip ( TR *  in,
TK *  out,
const bool  add = false,
const typename GetTypeReal< TK >::type  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.

Template Parameters
FPTYPEThe underlying fundamental type of the input data (e.g., float, double). SFINAE constraint ensures that FPTYPE is a fundamental type, not a complex type.
DeviceThe device type, which must be base_device::DEVICE_CPU.
std::enable_if<std::is_same<FPTYPE,typenameGetTypeReal<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>::typeSFINAE constraint to ensure that Device is base_device::DEVICE_CPU.
Parameters
inPointer to the input data array in real space, of type FPTYPE*.
outPointer to the output data array in reciprocal space, of type std::complex<FPTYPE>*.
ikOptional parameter, default value 0, representing some index or identifier.
addOptional boolean flag, default value false, indicating whether to add the result to the existing values in the output array.
factorOptional scaling factor, default value 1.0, applied to the output values.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ real_to_recip() [2/2]

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 ModulePW::PW_Basis::real_to_recip ( TR *  in,
TK *  out,
const bool  add = false,
const typename GetTypeReal< TK >::type  factor = 1.0 
) const
inline
Here is the call graph for this function:

◆ recip2real() [1/2]

template<typename FPTYPE >
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)

Parameters
in(nz,ns), std::complex<double>
out(nplane, ny, nx), double
Here is the call graph for this function:

◆ recip2real() [2/2]

template<typename FPTYPE >
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)

Parameters
in(nz,ns), std::complex<double>
out(nplane, ny, nx), std::complex<double>
Here is the call graph for this function:

◆ recip2real_gpu() [1/2]

template<typename FPTYPE >
void ModulePW::PW_Basis::recip2real_gpu ( const std::complex< FPTYPE > *  in,
FPTYPE *  out,
const bool  add = false,
const FPTYPE  factor = 1.0 
) const
Here is the caller graph for this function:

◆ recip2real_gpu() [2/2]

template<typename FPTYPE >
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

◆ recip_to_real() [1/2]

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 ModulePW::PW_Basis::recip_to_real ( TK *  in,
TR *  out,
const bool  add = false,
const typename GetTypeReal< TK >::type  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.

Template Parameters
FPTYPEThe type of the input data, which can only be a compelx type (e.g., std::complex<float>, std::complex<double>)
DeviceThe device type, must be base_device::DEVICE_CPU.
std::enable_if<!std::is_same<FPTYPE,typenameGetTypeReal<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>::typeSFINAE constraint to ensure that Device is base_device::DEVICE_CPU.
Parameters
inPointer to the input data array in reciprocal space.
outPointer 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.
addBoolean flag indicating whether to add the result to the existing values in the output array.
factorA scaling factor applied to the output values.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ recip_to_real() [2/2]

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 ModulePW::PW_Basis::recip_to_real ( TK *  in,
TR *  out,
const bool  add = false,
const typename GetTypeReal< TK >::type  factor = 1.0 
) const
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.

Template Parameters
FPTYPEThe underlying fundamental type of the input data (e.g., float, double). Note that the actual type passed should be std::complex<FPTYPE>.
DeviceThe device type, which can be any supported device type (e.g., base_device::DEVICE_CPU).
Parameters
inPointer to the input data array in reciprocal space, of type std::complex<FPTYPE>*.
outPointer to the output data array in real space, of type FPTYPE*.
addOptional boolean flag, default value false, indicating whether to add the result to the existing values in the output array.
factorOptional scaling factor, default value 1.0, applied to the output values.
Here is the call graph for this function:

◆ set_device()

void ModulePW::PW_Basis::set_device ( std::string  device_)
Here is the caller graph for this function:

◆ set_precision()

void ModulePW::PW_Basis::set_precision ( std::string  precision_)
Here is the caller graph for this function:

◆ setfullpw()

void ModulePW::PW_Basis::setfullpw ( const bool  inpt_full_pw = false,
const int  inpt_full_pw_dim = 0 
)
Here is the caller graph for this function:

◆ setuptransform()

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

Here is the call graph for this function:

Member Data Documentation

◆ classname

std::string ModulePW::PW_Basis::classname

◆ d_is2fftixy

int * ModulePW::PW_Basis::d_is2fftixy = nullptr

◆ device

std::string ModulePW::PW_Basis::device = "cpu"
protected

cpu or gpu

◆ distribution_type

int ModulePW::PW_Basis::distribution_type = 1

distribution method

◆ double_data_

bool ModulePW::PW_Basis::double_data_ = true
protected

if has double data

◆ fft_bundle

FFT_Bundle ModulePW::PW_Basis::fft_bundle

◆ fftixy2ip

int* ModulePW::PW_Basis::fftixy2ip =nullptr

◆ fftnx

int ModulePW::PW_Basis::fftnx =0

◆ fftnxy

int ModulePW::PW_Basis::fftnxy =0

◆ fftnxyz

int ModulePW::PW_Basis::fftnxyz =0

◆ fftny

int ModulePW::PW_Basis::fftny =0

◆ fftnz

int ModulePW::PW_Basis::fftnz =0

◆ float_data_

bool ModulePW::PW_Basis::float_data_ = false
protected

if has float data

◆ full_pw

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.

Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp.

◆ full_pw_dim

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.

◆ G

ModuleBase::Matrix3 ModulePW::PW_Basis::G

reciprocal lattice vector, unit in 1/lat0

◆ gamma_only

bool ModulePW::PW_Basis::gamma_only = false

◆ gcar

ModuleBase::Vector3<double>* ModulePW::PW_Basis::gcar =nullptr

◆ gdirect

ModuleBase::Vector3<double>* ModulePW::PW_Basis::gdirect =nullptr

◆ gg

double* ModulePW::PW_Basis::gg =nullptr

◆ gg_uniq

double* ModulePW::PW_Basis::gg_uniq =nullptr

◆ ggecut

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.

Examples
/home/runner/work/abacus-develop/abacus-develop/source/source_basis/module_pw/pw_distributeg_method1.cpp.

◆ GGT

ModuleBase::Matrix3 ModulePW::PW_Basis::GGT

◆ gridecut_lat

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.

◆ GT

ModuleBase::Matrix3 ModulePW::PW_Basis::GT

traspose of G

◆ ig2igg

int* ModulePW::PW_Basis::ig2igg =nullptr

◆ ig2isz

int* ModulePW::PW_Basis::ig2isz =nullptr

◆ ig2ixyz_gpu

int* ModulePW::PW_Basis::ig2ixyz_gpu = nullptr

◆ ig_gge0

int ModulePW::PW_Basis::ig_gge0 =-1

◆ is2fftixy

int* ModulePW::PW_Basis::is2fftixy =nullptr

◆ istot2ixy

int* ModulePW::PW_Basis::istot2ixy =nullptr

◆ lat0

double ModulePW::PW_Basis::lat0 = 1

unit length for lattice, unit in bohr

◆ latvec

ModuleBase::Matrix3 ModulePW::PW_Basis::latvec

Unitcell lattice vectors, unit in lat0.

◆ lix

int ModulePW::PW_Basis::lix =0

◆ liy

int ModulePW::PW_Basis::liy =0

◆ ng_xeq0

int ModulePW::PW_Basis::ng_xeq0 = 0

◆ ngg

int ModulePW::PW_Basis::ngg =0

◆ nmaxgr

int ModulePW::PW_Basis::nmaxgr = 0

◆ nplane

int ModulePW::PW_Basis::nplane =0

◆ npw

int ModulePW::PW_Basis::npw =0

◆ npw_per

int* ModulePW::PW_Basis::npw_per =nullptr

◆ npwtot

int ModulePW::PW_Basis::npwtot =0

◆ nrxx

int ModulePW::PW_Basis::nrxx =0

◆ nst

int ModulePW::PW_Basis::nst =0

◆ nst_per

int* ModulePW::PW_Basis::nst_per =nullptr

◆ nstnz

int ModulePW::PW_Basis::nstnz =0

◆ nstot

int ModulePW::PW_Basis::nstot =0

◆ numg

int* ModulePW::PW_Basis::numg =nullptr

◆ numr

int* ModulePW::PW_Basis::numr =nullptr

◆ numz

int* ModulePW::PW_Basis::numz =nullptr

◆ nx

int ModulePW::PW_Basis::nx =0

◆ nxy

int ModulePW::PW_Basis::nxy =0

◆ nxyz

int ModulePW::PW_Basis::nxyz =0

◆ ny

int ModulePW::PW_Basis::ny =0

◆ nz

int ModulePW::PW_Basis::nz =0

◆ omega

double ModulePW::PW_Basis::omega = 1.0

volume of the cell

◆ pool_world

MPI_Comm ModulePW::PW_Basis::pool_world =MPI_COMM_NULL

◆ poolnproc

int ModulePW::PW_Basis::poolnproc = 1

◆ poolrank

int ModulePW::PW_Basis::poolrank = 0

◆ precision

std::string ModulePW::PW_Basis::precision = "double"
protected

single, double, mixing

◆ rix

int ModulePW::PW_Basis::rix =0

◆ riy

int ModulePW::PW_Basis::riy =0

◆ startg

int* ModulePW::PW_Basis::startg =nullptr

◆ startnsz_per

int* ModulePW::PW_Basis::startnsz_per =nullptr
protected

◆ startr

int* ModulePW::PW_Basis::startr =nullptr

◆ startz

int* ModulePW::PW_Basis::startz =nullptr

◆ startz_current

int ModulePW::PW_Basis::startz_current =0

◆ tpiba

double ModulePW::PW_Basis::tpiba = 1

2pi/lat0

◆ tpiba2

double ModulePW::PW_Basis::tpiba2 = 1

4pi^2/lat0^2

◆ xprime

bool ModulePW::PW_Basis::xprime = true

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