![]() |
QMCPACK
|
Each SplineC2X needs a reader derived from BsplineReader. More...
Public Member Functions | |
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... | |
virtual std::unique_ptr< SPOSet > | create_spline_set (const std::string &my_name, int spin, const BandInfoGroup &bandgroup)=0 |
create the actual spline sets More... | |
void | setCommon (xmlNodePtr cur) |
setting common parameters More... | |
std::unique_ptr< SPOSet > | create_spline_set (int spin, xmlNodePtr cur, SPOSetInputInfo &input_info) |
create the spline after one of the kind is created More... | |
std::unique_ptr< SPOSet > | create_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... | |
Public Attributes | |
EinsplineSetBuilder * | mybuilder |
pointer to the EinsplineSetBuilder More... | |
Communicate * | myComm |
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... | |
Each SplineC2X needs a reader derived from BsplineReader.
This base class handles common chores
Definition at line 39 of file BsplineReader.h.
BsplineReader | ( | EinsplineSetBuilder * | e | ) |
Definition at line 30 of file BsplineReader.cpp.
References MPIObjectBase::getCommunicator(), BsplineReader::mybuilder, and BsplineReader::myComm.
|
virtualdefault |
|
inline |
initialize twist-related data for N orbitals
Definition at line 101 of file BsplineReader.h.
References qmcplusplus::abs(), qmcplusplus::app_log(), bspline(), qmcplusplus::dot(), qmcplusplus::Units::charge::e, BandInfoGroup::getFirstSPO(), BandInfoGroup::getLastSPO(), ParticleSet::getLattice(), BandInfoGroup::getNumDistinctOrbitals(), BandInfoGroup::getNumSPOs(), CrystalLattice< T, D >::k_cart(), BandInfoGroup::myBands, BsplineReader::mybuilder, qmcplusplus::Units::force::N, EinsplineSetBuilder::PrimCell, EinsplineSetBuilder::primcell_kpoints, EinsplineSetBuilder::TargetPtcl, qmcplusplus::transpose(), and BandInfoGroup::TwistIndex.
|
pure virtual |
create the actual spline sets
Implemented in SplineSetReader< SA >, SplineSetReader< typename SA::SplineBase >, and HybridRepSetReader< SA >.
Referenced by BsplineReader::create_spline_set().
std::unique_ptr< SPOSet > create_spline_set | ( | int | spin, |
xmlNodePtr | cur, | ||
SPOSetInputInfo & | input_info | ||
) |
create the spline after one of the kind is created
Definition at line 123 of file BsplineReader.cpp.
References OhmmsAttributeSet::add(), SPOSetBuilder::clear_states(), BsplineReader::create_spline_set(), EinsplineSetBuilder::FullBands, MPIObjectBase::getName(), BandInfoGroup::GroupID, BsplineReader::initialize_spo2band(), qmcplusplus::make_bandgroup_name(), SPOSetInputInfo::max_index(), SPOSetInputInfo::min_index(), BsplineReader::mybuilder, BandInfoGroup::myName, OhmmsAttributeSet::put(), BandInfoGroup::selectBands(), BsplineReader::spo2band, SPOSetBuilder::states, EinsplineSetBuilder::TileMatrix, EinsplineSetBuilder::twist_num_, and BandInfoGroup::TwistIndex.
std::unique_ptr< SPOSet > create_spline_set | ( | int | spin, |
xmlNodePtr | cur | ||
) |
create the spline set
Definition at line 87 of file BsplineReader.cpp.
References OhmmsAttributeSet::add(), APP_ABORT_TRACE, SPOSetBuilder::clear_states(), BsplineReader::create_spline_set(), EinsplineSetBuilder::FullBands, MPIObjectBase::getName(), BandInfoGroup::GroupID, BsplineReader::initialize_spo2band(), qmcplusplus::make_bandgroup_name(), BsplineReader::mybuilder, BandInfoGroup::myName, qmcplusplus::Units::time::ns, OhmmsAttributeSet::put(), BandInfoGroup::selectBands(), BsplineReader::spo2band, SPOSetBuilder::states, EinsplineSetBuilder::TileMatrix, EinsplineSetBuilder::twist_num_, and BandInfoGroup::TwistIndex.
|
inline |
Definition at line 59 of file BsplineReader.h.
References EinsplineSetBuilder::MeshSize, BsplineReader::mybuilder, and BandInfoGroup::myName.
void initialize_spo2band | ( | int | spin, |
const std::vector< BandInfo > & | bigspace, | ||
SPOSetInfo & | sposet, | ||
std::vector< int > & | spo2band | ||
) |
build index tables to map a state to band with k-point folidng
bigspace | full BandInfo constructed by EinsplineSetBuilder |
sposet | SPOSetInfo owned by someone, most likely EinsplinseSetBuilder |
spo2band | spo2band[i] is the index in bigspace |
At gamma or arbitrary kpoints with complex wavefunctions, spo2band[i]==i
Definition at line 163 of file BsplineReader.cpp.
References SPOSetInfo::add(), qmcplusplus::comm, qmcplusplus::Units::charge::e, Communicate::getGroupID(), MPIObjectBase::getName(), CrystalLattice< T, D >::k_cart(), qmcplusplus::make_bandinfo_filename(), BsplineReader::mybuilder, BsplineReader::myComm, qmcplusplus::Units::time::ns, EinsplineSetBuilder::PrimCell, EinsplineSetBuilder::primcell_kpoints, Communicate::rank(), qmcplusplus::Units::time::s, BsplineReader::spo2band, EinsplineSetBuilder::TileMatrix, and EinsplineSetBuilder::twist_num_.
Referenced by BsplineReader::create_spline_set().
|
inline |
return the path name in hdf5
ti | twist index |
spin | spin index |
ib | band index |
Definition at line 150 of file BsplineReader.h.
|
inline |
return the path name in hdf5
ti | twist index |
spin | spin index |
ib | band index |
Definition at line 162 of file BsplineReader.h.
|
inline |
read gvectors and set the mesh, and prepare for einspline
Definition at line 70 of file BsplineReader.h.
References EinsplineSetBuilder::MeshSize, BsplineReader::mybuilder, and EinsplineSetBuilder::ReadGvectors_ESHDF().
|
inline |
Set the checkNorm variable.
Definition at line 186 of file BsplineReader.h.
References BsplineReader::checkNorm.
void setCommon | ( | xmlNodePtr | cur | ) |
setting common parameters
Definition at line 68 of file BsplineReader.cpp.
References OhmmsAttributeSet::add(), qmcplusplus::app_log(), BsplineReader::checkNorm, OhmmsAttributeSet::put(), and BsplineReader::saveSplineCoefs.
|
inline |
Set the orbital rotation flag.
Rotations are applied to balance the real/imaginary components.
Definition at line 189 of file BsplineReader.h.
References BsplineReader::rotate.
bool checkNorm |
mesh size check the norm of orbitals
Definition at line 47 of file BsplineReader.h.
Referenced by BsplineReader::setCheckNorm(), and BsplineReader::setCommon().
EinsplineSetBuilder* mybuilder |
pointer to the EinsplineSetBuilder
Definition at line 42 of file BsplineReader.h.
Referenced by BsplineReader::BsplineReader(), BsplineReader::check_twists(), BsplineReader::create_spline_set(), BsplineReader::getSplineDumpFileName(), BsplineReader::initialize_spo2band(), and BsplineReader::set_grid().
Communicate* myComm |
communicator
Definition at line 44 of file BsplineReader.h.
Referenced by BsplineReader::BsplineReader(), and BsplineReader::initialize_spo2band().
bool rotate |
apply orbital rotations
Definition at line 51 of file BsplineReader.h.
Referenced by BsplineReader::setRotate().
bool saveSplineCoefs |
save spline coefficients to storage
Definition at line 49 of file BsplineReader.h.
Referenced by BsplineReader::setCommon().
std::vector<std::vector<int> > spo2band |
map from spo index to band index
Definition at line 53 of file BsplineReader.h.
Referenced by BsplineReader::create_spline_set(), and BsplineReader::initialize_spo2band().