QMCPACK
OneSplineOrbData Class Reference
+ Collaboration diagram for OneSplineOrbData:

Public Member Functions

 OneSplineOrbData (const TinyVector< int, 3 > &mesh_size, const TinyVector< int, 3 > &halfG, const bool isComplex)
 
 ~OneSplineOrbData ()
 
auto getRotatePhase () const
 
auto & get_spline_r ()
 
auto & get_spline_i ()
 
void clear ()
 
void fft_spline (const Vector< std::complex< double >> &cG, const std::vector< TinyVector< int, 3 >> &gvecs, const TinyVector< double, 3 > &primcell_kpoint, const bool rotate)
 fft and spline cG More...
 

Private Member Functions

void create (const TinyVector< int, 3 > &halfG)
 

Private Attributes

Array< std::complex< double >, 3 > FFTbox
 
Array< double, 3 > splineData_r
 
Array< double, 3 > splineData_i
 
double rotate_phase_r
 
double rotate_phase_i
 
UBspline_3d_d * spline_r = nullptr
 
UBspline_3d_d * spline_i = nullptr
 
fftw_plan FFTplan = nullptr
 
const TinyVector< int, 3 > & mesh_size_
 
const bool isComplex_
 

Detailed Description

Definition at line 28 of file OneSplineOrbData.hpp.

Constructor & Destructor Documentation

◆ OneSplineOrbData()

OneSplineOrbData ( const TinyVector< int, 3 > &  mesh_size,
const TinyVector< int, 3 > &  halfG,
const bool  isComplex 
)

Definition at line 43 of file OneSplineOrbData.cpp.

References OneSplineOrbData::create().

46  : mesh_size_(mesh_size), isComplex_(isComplex)
47 {
48  create(halfG);
49 }
void create(const TinyVector< int, 3 > &halfG)
const TinyVector< int, 3 > & mesh_size_

◆ ~OneSplineOrbData()

Definition at line 51 of file OneSplineOrbData.cpp.

References OneSplineOrbData::clear().

Member Function Documentation

◆ clear()

void clear ( )

Definition at line 53 of file OneSplineOrbData.cpp.

References OneSplineOrbData::FFTplan, OneSplineOrbData::spline_i, and OneSplineOrbData::spline_r.

Referenced by OneSplineOrbData::~OneSplineOrbData().

54 {
55  einspline::destroy(spline_r);
56  einspline::destroy(spline_i);
57  if (FFTplan != nullptr)
58  fftw_destroy_plan(FFTplan);
59  FFTplan = nullptr;
60 }

◆ create()

void create ( const TinyVector< int, 3 > &  halfG)
private

Definition at line 23 of file OneSplineOrbData.cpp.

References qmcplusplus::create(), Array< T, D, ALLOC >::data(), OneSplineOrbData::FFTbox, OneSplineOrbData::FFTplan, OneSplineOrbData::isComplex_, OneSplineOrbData::mesh_size_, Array< T, D, ALLOC >::resize(), OneSplineOrbData::spline_i, OneSplineOrbData::spline_r, OneSplineOrbData::splineData_i, and OneSplineOrbData::splineData_r.

Referenced by OneSplineOrbData::OneSplineOrbData().

24 {
25  const int nx = mesh_size_[0];
26  const int ny = mesh_size_[1];
27  const int nz = mesh_size_[2];
28  //perform FFT using FFTW
29  FFTbox.resize(nx, ny, nz);
30  FFTplan = fftw_plan_dft_3d(nx, ny, nz, reinterpret_cast<fftw_complex*>(FFTbox.data()),
31  reinterpret_cast<fftw_complex*>(FFTbox.data()), +1, FFTW_ESTIMATE);
32  splineData_r.resize(nx, ny, nz);
33  if (isComplex_)
34  splineData_i.resize(nx, ny, nz);
35 
36  TinyVector<double, 3> start(0.0);
37  TinyVector<double, 3> end(1.0);
38  spline_r = einspline::create(spline_r, start, end, mesh_size_, halfG);
39  if (isComplex_)
40  spline_i = einspline::create(spline_i, start, end, mesh_size_, halfG);
41 }
Type_t * data()
Definition: OhmmsArray.h:87
void resize(const std::array< SIZET, D > &dims)
Resize the container.
Definition: OhmmsArray.h:65
Array< std::complex< double >, 3 > FFTbox
hFile create(filename)
const TinyVector< int, 3 > & mesh_size_

◆ fft_spline()

void fft_spline ( const Vector< std::complex< double >> &  cG,
const std::vector< TinyVector< int, 3 >> &  gvecs,
const TinyVector< double, 3 > &  primcell_kpoint,
const bool  rotate 
)

fft and spline cG

Parameters
cGpsi_g to be processed
titwist index
iorborbital index

Perform FFT and spline to spline_r and spline_i

Definition at line 69 of file OneSplineOrbData.cpp.

References Array< T, D, ALLOC >::data(), OneSplineOrbData::FFTbox, OneSplineOrbData::FFTplan, qmcplusplus::fix_phase_rotate_c2c(), qmcplusplus::fix_phase_rotate_c2r(), OneSplineOrbData::isComplex_, OneSplineOrbData::mesh_size_, OneSplineOrbData::rotate_phase_i, OneSplineOrbData::rotate_phase_r, OneSplineOrbData::spline_i, OneSplineOrbData::spline_r, OneSplineOrbData::splineData_i, OneSplineOrbData::splineData_r, qmcplusplus::split_real_components_c2c(), and qmcplusplus::unpack4fftw().

73 {
74  unpack4fftw(cG, gvecs, mesh_size_, FFTbox);
75  fftw_execute(FFTplan);
76  if (isComplex_)
77  {
78  if (rotate)
80  else
81  {
83  rotate_phase_r = 1.0;
84  rotate_phase_i = 0.0;
85  }
86  einspline::set(spline_r, splineData_r.data());
87  einspline::set(spline_i, splineData_i.data());
88  }
89  else
90  {
92  einspline::set(spline_r, splineData_r.data());
93  }
94 }
Type_t * data()
Definition: OhmmsArray.h:87
void fix_phase_rotate_c2c(const Array< std::complex< T >, 3 > &in, Array< std::complex< T1 >, 3 > &out, const TinyVector< T2, 3 > &twist)
Array< std::complex< double >, 3 > FFTbox
void fix_phase_rotate_c2r(Array< std::complex< T >, 3 > &in, Array< T1, 3 > &out, const TinyVector< T2, 3 > &twist, T &phase_r, T &phase_i)
rotate the state after 3dfft
void split_real_components_c2c(const Array< std::complex< T >, 3 > &in, Array< T1, 3 > &out_r, Array< T1, 3 > &out_i)
Split FFTs into real/imaginary components.
const TinyVector< int, 3 > & mesh_size_
void unpack4fftw(const Vector< std::complex< T >> &cG, const std::vector< TinyVector< int, 3 >> &gvecs, const TinyVector< int, 3 > &maxg, Array< std::complex< T >, 3 > &fftbox)
unpack packed cG to fftbox

◆ get_spline_i()

auto& get_spline_i ( )
inline

Definition at line 49 of file OneSplineOrbData.hpp.

References OneSplineOrbData::spline_i.

49 { return *spline_i; }

◆ get_spline_r()

auto& get_spline_r ( )
inline

Definition at line 48 of file OneSplineOrbData.hpp.

References OneSplineOrbData::spline_r.

48 { return *spline_r; }

◆ getRotatePhase()

auto getRotatePhase ( ) const
inline

Member Data Documentation

◆ FFTbox

Array<std::complex<double>, 3> FFTbox
private

Definition at line 30 of file OneSplineOrbData.hpp.

Referenced by OneSplineOrbData::create(), and OneSplineOrbData::fft_spline().

◆ FFTplan

fftw_plan FFTplan = nullptr
private

◆ isComplex_

const bool isComplex_
private

Definition at line 38 of file OneSplineOrbData.hpp.

Referenced by OneSplineOrbData::create(), and OneSplineOrbData::fft_spline().

◆ mesh_size_

const TinyVector<int, 3>& mesh_size_
private

Definition at line 37 of file OneSplineOrbData.hpp.

Referenced by OneSplineOrbData::create(), and OneSplineOrbData::fft_spline().

◆ rotate_phase_i

double rotate_phase_i
private

◆ rotate_phase_r

double rotate_phase_r
private

◆ spline_i

UBspline_3d_d* spline_i = nullptr
private

◆ spline_r

UBspline_3d_d* spline_r = nullptr
private

◆ splineData_i

Array<double, 3> splineData_i
private

Definition at line 31 of file OneSplineOrbData.hpp.

Referenced by OneSplineOrbData::create(), and OneSplineOrbData::fft_spline().

◆ splineData_r

Array<double, 3> splineData_r
private

Definition at line 31 of file OneSplineOrbData.hpp.

Referenced by OneSplineOrbData::create(), and OneSplineOrbData::fft_spline().


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