QMCPACK
SplineSetReader< SA > Class Template Reference

General SplineSetReader to handle any unitcell. More...

+ Inheritance diagram for SplineSetReader< SA >:
+ Collaboration diagram for SplineSetReader< SA >:

Public Member Functions

 SplineSetReader (EinsplineSetBuilder *e)
 
std::unique_ptr< SPOSetcreate_spline_set (const std::string &my_name, int spin, const BandInfoGroup &bandgroup) override
 create the actual spline sets More...
 
bool createSplineDataSpaceLookforDumpFile (const BandInfoGroup &bandgroup, SA &bspline) const
 create data space in the spline object and try open spline dump files. More...
 
void readOneOrbitalCoefs (const std::string &s, hdf_archive &h5f, Vector< std::complex< double >> &cG) const
 read planewave coefficients from h5 file More...
 
void initialize_spline_pio_gather (const int spin, const BandInfoGroup &bandgroup, SA &bspline) const
 transforming planewave orbitals to 3D B-spline orbitals in real space. More...
 
- Public Member Functions inherited from BsplineReader
 BsplineReader (EinsplineSetBuilder *e)
 
virtual ~BsplineReader ()
 
std::string getSplineDumpFileName (const BandInfoGroup &bandgroup) const
 
template<typename GT , typename BCT >
bool set_grid (const TinyVector< int, 3 > &halfg, GT *xyz_grid, BCT *xyz_bc) const
 read gvectors and set the mesh, and prepare for einspline More...
 
template<typename SPE >
void check_twists (SPE &bspline, const BandInfoGroup &bandgroup) const
 initialize twist-related data for N orbitals More...
 
std::string psi_g_path (int ti, int spin, int ib) const
 return the path name in hdf5 More...
 
std::string psi_r_path (int ti, int spin, int ib) const
 return the path name in hdf5 More...
 
void setCommon (xmlNodePtr cur)
 setting common parameters More...
 
std::unique_ptr< SPOSetcreate_spline_set (int spin, xmlNodePtr cur, SPOSetInputInfo &input_info)
 create the spline after one of the kind is created More...
 
std::unique_ptr< SPOSetcreate_spline_set (int spin, xmlNodePtr cur)
 create the spline set More...
 
void setCheckNorm (bool new_checknorm)
 Set the checkNorm variable. More...
 
void setRotate (bool new_rotate)
 Set the orbital rotation flag. More...
 
void initialize_spo2band (int spin, const std::vector< BandInfo > &bigspace, SPOSetInfo &sposet, std::vector< int > &band2spo)
 build index tables to map a state to band with k-point folidng More...
 

Additional Inherited Members

- Public Attributes inherited from BsplineReader
EinsplineSetBuildermybuilder
 pointer to the EinsplineSetBuilder More...
 
CommunicatemyComm
 communicator More...
 
bool checkNorm
 mesh size check the norm of orbitals More...
 
bool saveSplineCoefs
 save spline coefficients to storage More...
 
bool rotate
 apply orbital rotations More...
 
std::vector< std::vector< int > > spo2band
 map from spo index to band index More...
 

Detailed Description

template<typename SA>
class qmcplusplus::SplineSetReader< SA >

General SplineSetReader to handle any unitcell.

Definition at line 52 of file SplineSetReader.h.

Constructor & Destructor Documentation

◆ SplineSetReader()

Definition at line 34 of file SplineSetReader.cpp.

34  : BsplineReader(e)
35 {}
BsplineReader(EinsplineSetBuilder *e)

Member Function Documentation

◆ create_spline_set()

std::unique_ptr< SPOSet > create_spline_set ( const std::string &  my_name,
int  spin,
const BandInfoGroup bandgroup 
)
overridevirtual

create the actual spline sets

Implements BsplineReader.

Definition at line 38 of file SplineSetReader.cpp.

41 {
42  auto bspline = std::make_unique<SA>(my_name);
43  app_log() << " ClassName = " << bspline->getClassName() << std::endl;
44  bool foundspline = createSplineDataSpaceLookforDumpFile(bandgroup, *bspline);
45  if (foundspline)
46  {
47  Timer now;
48  hdf_archive h5f(myComm);
49  const auto splinefile = getSplineDumpFileName(bandgroup);
50  h5f.open(splinefile, H5F_ACC_RDONLY);
51  foundspline = bspline->read_splines(h5f);
52  if (foundspline)
53  app_log() << " Successfully restored 3D B-spline coefficients from " << splinefile << ". The reading time is "
54  << now.elapsed() << " sec." << std::endl;
55  }
56 
57  if (!foundspline)
58  {
59  bspline->flush_zero();
60 
61  Timer now;
62  initialize_spline_pio_gather(spin, bandgroup, *bspline);
63  app_log() << " SplineSetReader initialize_spline_pio " << now.elapsed() << " sec" << std::endl;
64 
65  if (saveSplineCoefs && myComm->rank() == 0)
66  {
67  Timer now;
68  const std::string splinefile(getSplineDumpFileName(bandgroup));
69  hdf_archive h5f;
70  h5f.create(splinefile);
71  std::string classname = bspline->getClassName();
72  h5f.write(classname, "class_name");
73  int sizeD = sizeof(typename SA::DataType);
74  h5f.write(sizeD, "sizeof");
75  bspline->write_splines(h5f);
76  h5f.close();
77  app_log() << " Stored spline coefficients in " << splinefile << " for potential reuse. The writing time is "
78  << now.elapsed() << " sec." << std::endl;
79  }
80  }
81 
82  {
83  Timer now;
84  bspline->bcast_tables(myComm);
85  app_log() << " Time to bcast the table = " << now.elapsed() << std::endl;
86  }
87 
88  return bspline;
89 }
double bspline(double x, const std::vector< double > &t, const std::vector< double > &c, int k)
int rank() const
return the rank
Definition: Communicate.h:116
Communicate * myComm
communicator
Definition: BsplineReader.h:44
std::ostream & app_log()
Definition: OutputManager.h:65
bool createSplineDataSpaceLookforDumpFile(const BandInfoGroup &bandgroup, SA &bspline) const
create data space in the spline object and try open spline dump files.
void initialize_spline_pio_gather(const int spin, const BandInfoGroup &bandgroup, SA &bspline) const
transforming planewave orbitals to 3D B-spline orbitals in real space.
std::string getSplineDumpFileName(const BandInfoGroup &bandgroup) const
Definition: BsplineReader.h:59
bool saveSplineCoefs
save spline coefficients to storage
Definition: BsplineReader.h:49

◆ createSplineDataSpaceLookforDumpFile()

bool createSplineDataSpaceLookforDumpFile ( const BandInfoGroup bandgroup,
SA &  bspline 
) const

create data space in the spline object and try open spline dump files.

Parameters
bandgroupband info
bsplinethe spline object being worked on
Returns
true if dumpfile pass class name and data type size check

Definition at line 92 of file SplineSetReader.cpp.

93 {
94  if (bspline.isComplex())
95  app_log() << " Using complex einspline table" << std::endl;
96  else
97  app_log() << " Using real einspline table" << std::endl;
98 
99  //baseclass handles twists
100  check_twists(bspline, bandgroup);
101 
102  Ugrid xyz_grid[3];
103 
104  typename SA::BCType xyz_bc[3];
105  bool havePsig = set_grid(bspline.HalfG, xyz_grid, xyz_bc);
106  if (!havePsig)
107  myComm->barrier_and_abort("SplineSetReader needs psi_g. Set precision=\"double\".");
108  bspline.create_spline(xyz_grid, xyz_bc);
109 
110  int foundspline = 0;
111  Timer now;
112  if (myComm->rank() == 0)
113  {
114  now.restart();
115  hdf_archive h5f(myComm);
116  foundspline = h5f.open(getSplineDumpFileName(bandgroup), H5F_ACC_RDONLY);
117  if (foundspline)
118  {
119  std::string aname("none");
120  foundspline = h5f.readEntry(aname, "class_name");
121  foundspline = (aname.find(bspline.getKeyword()) != std::string::npos);
122  }
123  if (foundspline)
124  {
125  int sizeD = 0;
126  foundspline = h5f.readEntry(sizeD, "sizeof");
127  foundspline = (sizeD == sizeof(typename SA::DataType));
128  }
129  h5f.close();
130  }
131  myComm->bcast(foundspline);
132  return foundspline;
133 }
double bspline(double x, const std::vector< double > &t, const std::vector< double > &c, int k)
int rank() const
return the rank
Definition: Communicate.h:116
Communicate * myComm
communicator
Definition: BsplineReader.h:44
std::ostream & app_log()
Definition: OutputManager.h:65
void check_twists(SPE &bspline, const BandInfoGroup &bandgroup) const
initialize twist-related data for N orbitals
std::string getSplineDumpFileName(const BandInfoGroup &bandgroup) const
Definition: BsplineReader.h:59
bool set_grid(const TinyVector< int, 3 > &halfg, GT *xyz_grid, BCT *xyz_bc) const
read gvectors and set the mesh, and prepare for einspline
Definition: BsplineReader.h:70
void bcast(T &)
void barrier_and_abort(const std::string &msg) const

◆ initialize_spline_pio_gather()

void initialize_spline_pio_gather ( const int  spin,
const BandInfoGroup bandgroup,
SA &  bspline 
) const

transforming planewave orbitals to 3D B-spline orbitals in real space.

Parameters
spinorbital dataset spin index
bandgroupband info
bsplinethe spline object being worked on

Definition at line 160 of file SplineSetReader.cpp.

163 {
164  //distribute bands over processor groups
165  int Nbands = bandgroup.getNumDistinctOrbitals();
166  const int Nprocs = myComm->size();
167  const int Nbandgroups = std::min(Nbands, Nprocs);
168  Communicate band_group_comm(*myComm, Nbandgroups);
169  std::vector<int> band_groups(Nbandgroups + 1, 0);
170  FairDivideLow(Nbands, Nbandgroups, band_groups);
171  int iorb_first = band_groups[band_group_comm.getGroupID()];
172  int iorb_last = band_groups[band_group_comm.getGroupID() + 1];
173 
174  app_log() << "Start transforming plane waves to 3D B-Splines." << std::endl;
175  OneSplineOrbData oneband(mybuilder->MeshSize, bspline.HalfG, bspline.isComplex());
176  hdf_archive h5f(&band_group_comm, false);
177  Vector<std::complex<double>> cG(mybuilder->Gvecs[0].size());
178  const std::vector<BandInfo>& cur_bands = bandgroup.myBands;
179  if (band_group_comm.isGroupLeader())
180  {
181  h5f.open(mybuilder->H5FileName, H5F_ACC_RDONLY);
182  for (int iorb = iorb_first; iorb < iorb_last; iorb++)
183  {
184  const auto& cur_band = cur_bands[bspline.BandIndexMap[iorb]];
185  const int ti = cur_band.TwistIndex;
186  readOneOrbitalCoefs(psi_g_path(ti, spin, cur_band.BandIndex), h5f, cG);
187  oneband.fft_spline(cG, mybuilder->Gvecs[0], mybuilder->primcell_kpoints[ti], rotate);
188  bspline.set_spline(&oneband.get_spline_r(), &oneband.get_spline_i(), cur_band.TwistIndex, iorb, 0);
189  }
190 
191  {
192  band_group_comm.getGroupLeaderComm()->barrier();
193  Timer now;
194  bspline.gather_tables(band_group_comm.getGroupLeaderComm());
195  app_log() << " Time to gather the table = " << now.elapsed() << std::endl;
196  }
197  }
198 }
double bspline(double x, const std::vector< double > &t, const std::vector< double > &c, int k)
Communicate * myComm
communicator
Definition: BsplineReader.h:44
std::ostream & app_log()
Definition: OutputManager.h:65
int size() const
return the number of tasks
Definition: Communicate.h:118
T min(T a, T b)
void readOneOrbitalCoefs(const std::string &s, hdf_archive &h5f, Vector< std::complex< double >> &cG) const
read planewave coefficients from h5 file
void FairDivideLow(int ntot, int npart, IV &adist)
partition ntot elements among npart
Definition: FairDivide.h:114
Wrapping information on parallelism.
Definition: Communicate.h:68
std::string psi_g_path(int ti, int spin, int ib) const
return the path name in hdf5
std::vector< std::vector< TinyVector< int, 3 > > > Gvecs
std::vector< TinyVector< double, OHMMS_DIM > > primcell_kpoints
EinsplineSetBuilder * mybuilder
pointer to the EinsplineSetBuilder
Definition: BsplineReader.h:42
bool rotate
apply orbital rotations
Definition: BsplineReader.h:51

◆ readOneOrbitalCoefs()

void readOneOrbitalCoefs ( const std::string &  s,
hdf_archive h5f,
Vector< std::complex< double >> &  cG 
) const

read planewave coefficients from h5 file

Parameters
sdata set full path in h5
h5fhdf5 file handle
cGvector to store coefficients

Definition at line 136 of file SplineSetReader.cpp.

139 {
140  if (!h5f.readEntry(cG, s))
141  {
142  std::ostringstream msg;
143  msg << "SplineSetReader Failed to read band(s) from h5 file. " << "Attempted dataset " << s << " with " << cG.size()
144  << " complex numbers." << std::endl;
145  throw std::runtime_error(msg.str());
146  }
147  double total_norm = compute_norm(cG);
148  if ((checkNorm) && (std::abs(total_norm - 1.0) > PW_COEFF_NORM_TOLERANCE))
149  {
150  std::ostringstream msg;
151  msg << "SplineSetReader The orbital dataset " << s << " has a wrong norm " << total_norm
152  << ", computed from plane wave coefficients!" << std::endl
153  << "This may indicate a problem with the HDF5 library versions used "
154  << "during wavefunction conversion or read." << std::endl;
155  throw std::runtime_error(msg.str());
156  }
157 }
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
T compute_norm(const Vector< std::complex< T >> &cG)
Compute the norm of an orbital.
#define PW_COEFF_NORM_TOLERANCE
bool checkNorm
mesh size check the norm of orbitals
Definition: BsplineReader.h:47

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