QMCPACK
OneSplineOrbData.hpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////////////
2 // This file is distributed under the University of Illinois/NCSA Open Source License.
3 // See LICENSE file in top directory for details.
4 //
5 // Copyright (c) 2023 QMCPACK developers.
6 //
7 // File developed by: Jeremy McMinnis, jmcminis@gmail.com, University of Illinois at Urbana-Champaign
8 // Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
9 // Ye Luo, yeluo@anl.gov, Argonne National Laboratory
10 // Mark A. Berrill, berrillma@ornl.gov, Oak Ridge National Laboratory
11 // Jeongnim Kim, jeongnim.kim@inte.com, Intel Corp.
12 //
13 // File created by: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
14 //////////////////////////////////////////////////////////////////////////////////////
15 
16 
17 #ifndef QMCPLUSPLUS_ONESPLINEORBDATA_H
18 #define QMCPLUSPLUS_ONESPLINEORBDATA_H
19 
20 #include <complex>
21 #include <fftw3.h>
22 #include <einspline/bspline_structs.h>
23 #include <OhmmsPETE/TinyVector.h>
24 #include <OhmmsPETE/OhmmsArray.h>
25 
26 namespace qmcplusplus
27 {
29 {
33  UBspline_3d_d* spline_r = nullptr;
34  UBspline_3d_d* spline_i = nullptr;
35  fftw_plan FFTplan = nullptr;
36 
38  const bool isComplex_;
39 
40  void create(const TinyVector<int, 3>& halfG);
41 
42 public:
43  OneSplineOrbData(const TinyVector<int, 3>& mesh_size, const TinyVector<int, 3>& halfG, const bool isComplex);
44 
46 
47  auto getRotatePhase() const { return std::complex<double>(rotate_phase_r, rotate_phase_i); }
48  auto& get_spline_r() { return *spline_r; }
49  auto& get_spline_i() { return *spline_i; }
50 
51  void clear();
52 
53  /** fft and spline cG
54  * @param cG psi_g to be processed
55  * @param ti twist index
56  * @param iorb orbital index
57  *
58  * Perform FFT and spline to spline_r and spline_i
59  */
60  void fft_spline(const Vector<std::complex<double>>& cG,
61  const std::vector<TinyVector<int, 3>>& gvecs,
62  const TinyVector<double, 3>& primcell_kpoint,
63  const bool rotate);
64 };
65 } // namespace qmcplusplus
66 #endif
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
void create(const TinyVector< int, 3 > &halfG)
Array< std::complex< double >, 3 > FFTbox
OneSplineOrbData(const TinyVector< int, 3 > &mesh_size, const TinyVector< int, 3 > &halfG, const bool isComplex)
A D-dimensional Array class based on PETE.
Definition: OhmmsArray.h:25
const TinyVector< int, 3 > & mesh_size_
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