QMCPACK
EinsplineSpinorSetBuilder Class Reference
+ Inheritance diagram for EinsplineSpinorSetBuilder:
+ Collaboration diagram for EinsplineSpinorSetBuilder:

Public Member Functions

 EinsplineSpinorSetBuilder (ParticleSet &p, const PSetMap &psets, Communicate *comm, xmlNodePtr cur)
 constructor More...
 
 ~EinsplineSpinorSetBuilder () override
 destructor More...
 
std::unique_ptr< SPOSetcreateSPOSetFromXML (xmlNodePtr cur) override
 initialize the Antisymmetric wave function for electrons More...
 
- Public Member Functions inherited from EinsplineSetBuilder
 EinsplineSetBuilder (ParticleSet &p, const PSetMap &psets, Communicate *comm, xmlNodePtr cur)
 constructor More...
 
 ~EinsplineSetBuilder () override
 destructor More...
 
std::unique_ptr< SPOSetcreateSPOSet (xmlNodePtr cur, SPOSetInputInfo &input_info) override
 initialize with the existing SPOSet More...
 
bool ReadOrbitalInfo (bool skipChecks=false)
 
bool ReadOrbitalInfo_ESHDF (bool skipChecks=false)
 
void BroadcastOrbitalInfo ()
 
bool CheckLattice ()
 
bool ReadGvectors_ESHDF ()
 read gvectors for each twist More...
 
bool TwistPair (PosType a, PosType b) const
 
void TileIons ()
 
void OccupyBands (int spin, int sortBands, int numOrbs, bool skipChecks=false)
 
void OccupyBands_ESHDF (int spin, int sortBands, int numOrbs)
 
std::string OrbitalPath (int ti, int bi)
 
- Public Member Functions inherited from SPOSetBuilder
 SPOSetBuilder (const std::string &type_name, Communicate *comm)
 
virtual ~SPOSetBuilder ()
 
void reserve_states (int nsets=1)
 reserve space for states (usually only one set, multiple for e.g. spin dependent einspline) More...
 
void modify_states (int index=0)
 allow modification of state information More...
 
void clear_states (int index=0)
 clear state information More...
 
std::unique_ptr< SPOSetcreateSPOSet (xmlNodePtr cur)
 create an sposet from xml and save the resulting SPOSet More...
 
std::unique_ptr< SPOSetcreateRotatedSPOSet (xmlNodePtr cur)
 create orbital rotation transformation from xml and save the resulting SPOSet More...
 
const std::string & getTypeName () const
 
- Public Member Functions inherited from MPIObjectBase
 MPIObjectBase (Communicate *c)
 constructor with communicator More...
 
int rank () const
 return the rank of the communicator More...
 
int getGroupID () const
 return the group id of the communicator More...
 
CommunicategetCommunicator () const
 return myComm More...
 
CommunicategetCommRef () const
 return a TEMPORARY reference to Communicate More...
 
mpi_comm_type getMPI () const
 return MPI communicator if one wants to use MPI directly More...
 
bool is_manager () const
 return true if the rank == 0 More...
 
const std::string & getName () const
 return the name More...
 
void setName (const std::string &aname)
 

Private Types

using PSetMap = std::map< std::string, const std::unique_ptr< ParticleSet > >
 

Additional Inherited Members

- Public Types inherited from EinsplineSetBuilder
enum  FormatType { QMCPACK, ESHDF }
 
using PSetMap = std::map< std::string, const std::unique_ptr< ParticleSet > >
 
using UnitCellType = CrystalLattice< ParticleSet::Scalar_t, DIM >
 
- Public Types inherited from SPOSetBuilder
using indices_t = std::vector< int >
 
using energies_t = std::vector< RealType >
 
- Public Types inherited from QMCTraits
enum  { DIM = OHMMS_DIM, DIM_VGL = OHMMS_DIM + 2 }
 
using QTBase = QMCTypes< OHMMS_PRECISION, DIM >
 
using QTFull = QMCTypes< OHMMS_PRECISION_FULL, DIM >
 
using RealType = QTBase::RealType
 
using ComplexType = QTBase::ComplexType
 
using ValueType = QTBase::ValueType
 
using PosType = QTBase::PosType
 
using GradType = QTBase::GradType
 
using TensorType = QTBase::TensorType
 
using IndexType = OHMMS_INDEXTYPE
 define other types More...
 
using FullPrecRealType = QTFull::RealType
 
using FullPrecValueType = QTFull::ValueType
 
using PropertySetType = RecordNamedProperty< FullPrecRealType >
 define PropertyList_t More...
 
using PtclGrpIndexes = std::vector< std::pair< int, int > >
 
- Public Types inherited from MPIObjectBase
using mpi_comm_type = Communicate::mpi_comm_type
 
- Public Attributes inherited from EinsplineSetBuilder
const PSetMapParticleSets
 reference to the particleset pool More...
 
ParticleSetTargetPtcl
 quantum particle set More...
 
ParticleSetSourcePtcl
 ionic system More...
 
std::vector< std::unique_ptr< std::vector< BandInfo > > > FullBands
 Helper vector for sorting bands. More...
 
std::unique_ptr< BsplineReaderMixedSplineReader
 reader to use BsplineReader More...
 
bool HaveOrbDerivs
 This is true if we have the orbital derivatives w.r.t. the ion positions. More...
 
xmlNodePtr XMLRoot
 root XML node with href, sort, tilematrix, twistnum, source, precision,truncate,version More...
 
std::map< H5OrbSet, SPOSet *, H5OrbSetSPOSetMap
 
hdf_archive H5File
 
std::filesystem::path H5FileName
 
FormatType Format
 
TinyVector< int, 3 > Version
 
std::string parameterGroup
 
std::string ionsGroup
 
std::string eigenstatesGroup
 
std::vector< int > Occ
 
Tensor< double, OHMMS_DIMLattice
 
Tensor< double, OHMMS_DIMRecipLattice
 
Tensor< double, OHMMS_DIMLatticeInv
 
Tensor< double, OHMMS_DIMSuperLattice
 
Tensor< double, OHMMS_DIMGGt
 
UnitCellType SuperCell
 
UnitCellType PrimCell
 
UnitCellType PrimCellInv
 
int NumBands
 
int NumElectrons
 
int NumSpins
 
int NumTwists
 
int MaxNumGvecs
 
double MeshFactor
 
RealType MatchingTol
 
TinyVector< int, 3 > MeshSize
 
std::vector< std::vector< TinyVector< int, 3 > > > Gvecs
 
Vector< int > IonTypes
 
Vector< TinyVector< double, OHMMS_DIM > > IonPos
 
std::vector< int > Super2Prim
 
int twist_num_
 
std::vector< TinyVector< double, OHMMS_DIM > > primcell_kpoints
 
Tensor< int, OHMMS_DIMTileMatrix
 
std::vector< TinyVector< int, OHMMS_DIM > > UseTwists
 
std::vector< int > IncludeTwists
 
std::vector< int > DistinctTwists
 
bool use_real_splines_
 if false, splines are conceptually complex valued More...
 
int NumDistinctOrbitals
 
std::vector< bool > MakeTwoCopies
 
std::map< TinyVector< int, OHMMS_DIM >, int, Int3lessTwistMap
 
struct qmcplusplus::EinsplineSetBuilder::CenterInfo AtomicCentersInfo
 
int LastSpinSet
 
int NumOrbitalsRead
 
std::string occ_format
 
int particle_hole_pairs
 
bool makeRotations
 
- Public Attributes inherited from SPOSetBuilder
bool legacy
 whether implementation conforms only to legacy standard More...
 
std::vector< std::unique_ptr< SPOSetInfo > > states
 state info of all possible states available in the basis More...
 
- Protected Member Functions inherited from EinsplineSetBuilder
void bcastSortBands (int splin, int N, bool root)
 broadcast SortBands More...
 
void set_metadata (int numOrbs, int twist_num_inp, const TinyVector< double, OHMMS_DIM > &twist_inp, bool skipChecks=false)
 a specific but clean code path in createSPOSetFromXML, for PBC, double, ESHDF More...
 
void AnalyzeTwists2 (const int twist_num_inp, const TinyVector< double, OHMMS_DIM > &twist_inp)
 analyze twists of orbitals in h5 and determinine twist_num_ More...
 
- Protected Attributes inherited from SPOSetBuilder
const std::string type_name_
 type name of the SPO objects built by this builder. More...
 
- Protected Attributes inherited from MPIObjectBase
CommunicatemyComm
 pointer to Communicate More...
 
std::string ClassName
 class Name More...
 
std::string myName
 name of the object More...
 
- Static Protected Attributes inherited from EinsplineSetBuilder
static constexpr int TWISTNUM_NO_INPUT = -9999
 twistnum_inp == -9999 to indicate no given input after parsing XML More...
 
static constexpr double TWIST_NO_INPUT = -9999
 twist_inp[i] <= -9999 to indicate no given input after parsing XML More...
 

Detailed Description

Definition at line 29 of file EinsplineSpinorSetBuilder.h.

Member Typedef Documentation

◆ PSetMap

using PSetMap = std::map<std::string, const std::unique_ptr<ParticleSet> >
private

Definition at line 31 of file EinsplineSpinorSetBuilder.h.

Constructor & Destructor Documentation

◆ EinsplineSpinorSetBuilder()

EinsplineSpinorSetBuilder ( ParticleSet p,
const PSetMap psets,
Communicate comm,
xmlNodePtr  cur 
)
inline

constructor

Definition at line 35 of file EinsplineSpinorSetBuilder.h.

36  : EinsplineSetBuilder(p, psets, comm, cur){};
EinsplineSetBuilder(ParticleSet &p, const PSetMap &psets, Communicate *comm, xmlNodePtr cur)
constructor

◆ ~EinsplineSpinorSetBuilder()

~EinsplineSpinorSetBuilder ( )
inlineoverride

destructor

Definition at line 39 of file EinsplineSpinorSetBuilder.h.

39 {};

Member Function Documentation

◆ createSPOSetFromXML()

std::unique_ptr< SPOSet > createSPOSetFromXML ( xmlNodePtr  cur)
overridevirtual

initialize the Antisymmetric wave function for electrons

Parameters
curthe current xml node

Reimplemented from EinsplineSetBuilder.

Definition at line 31 of file EinsplineSpinorSetBuilder.cpp.

References OhmmsAttributeSet::add(), qmcplusplus::app_warning(), Communicate::barrier_and_abort(), EinsplineSetBuilder::bcastSortBands(), qmcplusplus::createBsplineComplex(), qmcplusplus::createBsplineReal(), qmcplusplus::createGlobalTimer(), DELETED, EinsplineSetBuilder::FullBands, ParticleSet::get(), EinsplineSetBuilder::H5FileName, EinsplineSetBuilder::MeshFactor, EinsplineSetBuilder::MixedSplineReader, MPIObjectBase::myComm, MPIObjectBase::myName, EinsplineSetBuilder::NumDistinctOrbitals, EinsplineSetBuilder::Occ, EinsplineSetBuilder::occ_format, EinsplineSetBuilder::OccupyBands(), EinsplineSetBuilder::particle_hole_pairs, EinsplineSetBuilder::ParticleSets, OhmmsAttributeSet::put(), putContent(), Communicate::rank(), Timer::restart(), EinsplineSetBuilder::set_metadata(), EinsplineSetBuilder::SourcePtcl, EinsplineSetBuilder::SPOSetMap, SPOSetBuilder::states, EinsplineSetBuilder::TileIons(), EinsplineSetBuilder::TileMatrix, qmcplusplus::timer_level_medium, EinsplineSetBuilder::TWIST_NO_INPUT, EinsplineSetBuilder::TWISTNUM_NO_INPUT, EinsplineSetBuilder::use_real_splines_, and EinsplineSetBuilder::XMLRoot.

32 {
33  int numOrbs = 0;
34  int sortBands(1);
35  int spinSet = 0;
36  int spinSet2 = 1;
37  int twist_num_inp = TWISTNUM_NO_INPUT;
39 
40  //There have to be two "spin states"... one for the up channel and one for the down channel.
41  // We force this for spinors and manually resize states and FullBands.
42  states.clear();
43  states.resize(2);
44 
45  FullBands.resize(2);
46 
47  SPOSet* UpOrbitalSet;
48  std::string sourceName;
49  std::string spo_prec("double");
50  std::string truncate("no");
51  std::string hybrid_rep("no");
52  std::string spo_object_name;
53 
54  ScopedTimer spo_timer_scope(createGlobalTimer("einspline::CreateSpinorSetFromXML", timer_level_medium));
55 
56  {
58  TinyVector<int, OHMMS_DIM> TileFactor_do_not_use;
59  a.add(H5FileName, "href");
60  a.add(TileFactor_do_not_use, "tile", {}, TagStatus::DELETED);
61  a.add(sortBands, "sort");
62  a.add(TileMatrix, "tilematrix");
63  a.add(twist_num_inp, "twistnum");
64  a.add(twist_inp, "twist");
65  a.add(sourceName, "source");
66  a.add(MeshFactor, "meshfactor");
67  a.add(hybrid_rep, "hybridrep");
68  a.add(spo_prec, "precision");
69  a.add(truncate, "truncate");
70  a.add(myName, "tag");
71 
72  a.put(XMLRoot);
73  a.add(numOrbs, "size");
74  a.add(numOrbs, "norbs");
75  a.add(spinSet, "spindataset");
76  a.add(spinSet, "group");
77  a.put(cur);
78 
79  if (myName.empty())
80  myName = "einspline.spinor";
81  }
82 
83  auto pit(ParticleSets.find(sourceName));
84  if (pit == ParticleSets.end())
85  myComm->barrier_and_abort("Einspline needs the source particleset");
86  else
87  SourcePtcl = pit->second.get();
88 
89  ///////////////////////////////////////////////
90  // Read occupation information from XML file //
91  ///////////////////////////////////////////////
92  const std::vector<int> last_occ(Occ);
93  Occ.resize(0, 0); // correspond to ground
94  bool NewOcc(false);
95 
96  {
97  OhmmsAttributeSet oAttrib;
98  oAttrib.add(spinSet, "spindataset");
99  oAttrib.add(spo_object_name, "name");
100  oAttrib.add(spo_object_name, "id");
101  oAttrib.put(cur);
102  }
103 
104  xmlNodePtr spo_cur = cur;
105  cur = cur->children;
106  while (cur != NULL)
107  {
108  std::string cname((const char*)(cur->name));
109  if (cname == "occupation")
110  {
111  std::string occ_mode("ground");
112  occ_format = "energy";
114  OhmmsAttributeSet oAttrib;
115  oAttrib.add(occ_mode, "mode");
116  oAttrib.add(spinSet, "spindataset");
117  oAttrib.add(occ_format, "format");
118  oAttrib.add(particle_hole_pairs, "pairs");
119  oAttrib.put(cur);
120  if (occ_mode == "excited")
121  putContent(Occ, cur);
122  else if (occ_mode != "ground")
123  myComm->barrier_and_abort("EinsplineSetBuilder::createSPOSet Only ground state occupation currently "
124  "supported in EinsplineSetBuilder.");
125  }
126  cur = cur->next;
127  }
128 
129  if (Occ != last_occ)
130  {
131  NewOcc = true;
132  }
133  else
134  NewOcc = false;
135 
136  H5OrbSet aset(H5FileName, spinSet, numOrbs);
137  const auto iter = SPOSetMap.find(aset);
138  if ((iter != SPOSetMap.end()) && (!NewOcc))
139  app_warning() << "!!!!!!! Identical SPOSets are detected by EinsplineSpinorSetBuilder! "
140  "Implicit sharing one SPOSet for spin-up and spin-down electrons has been removed. "
141  "Each determinant creates its own SPOSet with dedicated memory for spline coefficients. "
142  "To avoid increasing the memory footprint of spline coefficients, "
143  "create a single SPOset outside the determinantset using 'sposet_collection' "
144  "and reference it by name on the determinant line."
145  << std::endl;
146 
147  if (FullBands[spinSet] == 0)
148  FullBands[spinSet] = std::make_unique<std::vector<BandInfo>>();
149 
150  if (FullBands[spinSet2] == 0)
151  FullBands[spinSet2] = std::make_unique<std::vector<BandInfo>>();
152 
153  //This is to skip checks on ion-ID's, spin types, etc. If we've made it here, we assume we know better
154  //than Einspline on what the data means...
155  bool skipChecks = true;
156 
157  set_metadata(numOrbs, twist_num_inp, twist_inp, skipChecks);
158 
159  //////////////////////////////////
160  // Create the OrbitalSet object
161  //////////////////////////////////
162  Timer mytimer;
163  mytimer.restart();
164  OccupyBands(spinSet, sortBands, numOrbs, skipChecks);
165  if (spinSet == 0)
166  TileIons();
167 
168  bool use_single = (spo_prec == "single" || spo_prec == "float");
169 
170  // safeguard for a removed feature
171  if (truncate == "yes")
173  "The 'truncate' feature of spline SPO has been removed. Please use hybrid orbital representation.");
174 
175  std::string useGPU("no");
176 #if !defined(QMC_COMPLEX)
177  if (use_real_splines_)
178  {
179  if (MixedSplineReader == 0)
180  MixedSplineReader = createBsplineReal(this, use_single, hybrid_rep == "yes", useGPU);
181  }
182  else
183 #endif
184  {
185  if (MixedSplineReader == 0)
186  MixedSplineReader = createBsplineComplex(this, use_single, hybrid_rep == "yes", useGPU);
187  }
188 
189  MixedSplineReader->setCommon(XMLRoot);
190  //Norm for spinor wavefunctions is different from SPO's by a factor of sqrt(2). Disable the unit norm check.
191  MixedSplineReader->setCheckNorm(false);
192  //Set no rotation to the orbitals
193  MixedSplineReader->setRotate(false);
194 
195  //Make the up spin set.
196  bcastSortBands(spinSet, NumDistinctOrbitals, myComm->rank() == 0);
197  auto bspline_zd_u = MixedSplineReader->create_spline_set(spinSet, spo_cur);
198  bspline_zd_u->finalizeConstruction();
199 
200  //Make the down spin set.
201  OccupyBands(spinSet2, sortBands, numOrbs, skipChecks);
202  bcastSortBands(spinSet2, NumDistinctOrbitals, myComm->rank() == 0);
203  auto bspline_zd_d = MixedSplineReader->create_spline_set(spinSet2, spo_cur);
204  bspline_zd_d->finalizeConstruction();
205 
206  //register with spin set and we're off to the races.
207  auto spinor_set = std::make_unique<SpinorSet>(spo_object_name);
208  spinor_set->set_spos(std::move(bspline_zd_u), std::move(bspline_zd_d));
209  return spinor_set;
210 };
const PSetMap & ParticleSets
reference to the particleset pool
Tensor< int, OHMMS_DIM > TileMatrix
std::ostream & app_warning()
Definition: OutputManager.h:69
int rank() const
return the rank
Definition: Communicate.h:116
std::map< H5OrbSet, SPOSet *, H5OrbSet > SPOSetMap
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
std::string myName
name of the object
Definition: MPIObjectBase.h:67
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
ParticleSet * SourcePtcl
ionic system
std::vector< std::unique_ptr< SPOSetInfo > > states
state info of all possible states available in the basis
Definition: SPOSetBuilder.h:57
std::vector< std::unique_ptr< std::vector< BandInfo > > > FullBands
Helper vector for sorting bands.
static constexpr double TWIST_NO_INPUT
twist_inp[i] <= -9999 to indicate no given input after parsing XML
void bcastSortBands(int splin, int N, bool root)
broadcast SortBands
Communicate * myComm
pointer to Communicate
Definition: MPIObjectBase.h:62
std::unique_ptr< BsplineReader > createBsplineReal(EinsplineSetBuilder *e, bool use_single, bool hybrid_rep, const std::string &useGPU)
create a reader which handles real splines, R2R case spline storage and computation precision is doub...
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
NewTimer & createGlobalTimer(const std::string &myname, timer_levels mylevel)
std::unique_ptr< BsplineReader > createBsplineComplex(EinsplineSetBuilder *e, bool hybrid_rep, const std::string &useGPU)
create a reader which handles complex (double size real) splines, C2R or C2C case spline storage and ...
void OccupyBands(int spin, int sortBands, int numOrbs, bool skipChecks=false)
std::unique_ptr< BsplineReader > MixedSplineReader
reader to use BsplineReader
static constexpr int TWISTNUM_NO_INPUT
twistnum_inp == -9999 to indicate no given input after parsing XML
bool putContent(T &a, xmlNodePtr cur)
replaces a&#39;s value with the first "element" in the "string" returned by XMLNodeString{cur}.
Definition: libxmldefs.h:88
void set_metadata(int numOrbs, int twist_num_inp, const TinyVector< double, OHMMS_DIM > &twist_inp, bool skipChecks=false)
a specific but clean code path in createSPOSetFromXML, for PBC, double, ESHDF
xmlNodePtr XMLRoot
root XML node with href, sort, tilematrix, twistnum, source, precision,truncate,version ...
void barrier_and_abort(const std::string &msg) const
bool get(std::ostream &os) const override
dummy. For satisfying OhmmsElementBase.
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42
bool use_real_splines_
if false, splines are conceptually complex valued

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