QMCPACK
HybridRepSetReader.h
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: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
8 //
9 // File created by: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
10 //////////////////////////////////////////////////////////////////////////////////////
11 
12 
13 /** @file
14  *
15  * derived from BsplineReader
16  */
17 
18 #ifndef QMCPLUSPLUS_HYBRIDREP_READER_H
19 #define QMCPLUSPLUS_HYBRIDREP_READER_H
20 
21 #include "BsplineReader.h"
22 #include "SplineSetReader.h"
23 
24 #if !defined(QMC_COMPLEX)
25 #include "HybridRepReal.h"
26 #endif
27 #include "HybridRepCplx.h"
28 
29 #if defined(QMC_COMPLEX)
30 #include "SplineC2C.h"
31 #include "SplineC2COMPTarget.h"
32 #else
33 #include "SplineR2R.h"
34 #include "SplineC2R.h"
35 #include "SplineC2ROMPTarget.h"
36 #endif
37 
38 
39 namespace qmcplusplus
40 {
41 /** General HybridRepSetReader to handle any unitcell
42  */
43 template<typename SA>
44 class HybridRepSetReader : public BsplineReader
45 {
47  using DataType = typename SA::DataType;
49 
50 public:
52 
53  std::unique_ptr<SPOSet> create_spline_set(const std::string& my_name,
54  int spin,
55  const BandInfoGroup& bandgroup) override;
56 
57  /** initialize basic parameters of atomic orbitals */
59 
60  /** initialize construct atomic orbital radial functions from plane waves */
61  inline void create_atomic_centers_Gspace(const Vector<std::complex<double>>& cG,
62  Communicate& band_group_comm,
63  const int iorb,
64  const std::complex<double>& rotate_phase,
65  SA& bspline) const;
66 
67  /** transforming planewave orbitals to 3D B-spline orbitals and 1D B-spline radial orbitals in real space.
68  * @param spin orbital dataset spin index
69  * @param bandgroup band info
70  * @param bspline the spline object being worked on
71  */
72  void initialize_hybrid_pio_gather(const int spin, const BandInfoGroup& bandgroup, SA& bspline) const;
73 };
74 
75 #if defined(QMC_COMPLEX)
78 #if !defined(QMC_MIXED_PRECISION)
81 #endif
82 #else
86 #if !defined(QMC_MIXED_PRECISION)
90 #endif
91 #endif
92 } // namespace qmcplusplus
93 #endif
double bspline(double x, const std::vector< double > &t, const std::vector< double > &c, int k)
hold HybridRepCplx
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
HybridRepSetReader(EinsplineSetBuilder *e)
Wrapping information on parallelism.
Definition: Communicate.h:68
a group of bands
Definition: BandInfo.h:66
class to handle complex splines to real orbitals with splines of arbitrary precision ...
void initialize_hybridrep_atomic_centers(SA &bspline) const
initialize basic parameters of atomic orbitals
base class to read data and manage spline tables
void initialize_hybrid_pio_gather(const int spin, const BandInfoGroup &bandgroup, SA &bspline) const
transforming planewave orbitals to 3D B-spline orbitals and 1D B-spline radial orbitals in real space...
class to handle complex splines to real orbitals with splines of arbitrary precision splines storage ...
void create_atomic_centers_Gspace(const Vector< std::complex< double >> &cG, Communicate &band_group_comm, const int iorb, const std::complex< double > &rotate_phase, SA &bspline) const
initialize construct atomic orbital radial functions from plane waves
The most general reader class for the following classes using the full single grid for the supercell...
class to handle complex splines to complex orbitals with splines of arbitrary precision splines stora...
hold HybridRepReal
General HybridRepSetReader to handle any unitcell.
std::unique_ptr< SPOSet > create_spline_set(const std::string &my_name, int spin, const BandInfoGroup &bandgroup) override
create the actual spline sets
class to handle complex splines to complex orbitals with splines of arbitrary precision ...