QMCPACK
RadialOrbitalSetBuilder< COT > Class Template Reference

Build a set of radial orbitals at the origin. More...

+ Inheritance diagram for RadialOrbitalSetBuilder< COT >:
+ Collaboration diagram for RadialOrbitalSetBuilder< COT >:

Public Types

using RealType = typename COT::RealType
 
using RadialOrbitalType = typename COT::RadialOrbital_t
 
using GridType = typename COT::GridType
 
- Public Types inherited from MPIObjectBase
using mpi_comm_type = Communicate::mpi_comm_type
 

Public Member Functions

 RadialOrbitalSetBuilder (Communicate *comm, COT &aos, int radial_grid_size=1001)
 constructor More...
 
bool addGrid (xmlNodePtr cur, const std::string &rad_type)
 add a grid More...
 
bool addGridH5 (hdf_archive &hin)
 
bool addRadialOrbital (xmlNodePtr cur, const std::string &m_infunctype, const QuantumNumberType &nlms)
 add a radial functor More...
 
bool addRadialOrbitalH5 (hdf_archive &hin, const std::string &radtype, const QuantumNumberType &nlms)
 
void finalize ()
 This is when the radial orbitals are actually created. More...
 
- 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)
 

Public Attributes

bool Normalized
 true, if the RadialOrbitalType is normalized More...
 
COT & m_orbitals
 the atomic orbitals More...
 
std::unique_ptr< GridTypeinput_grid
 input grid in case transform is needed More...
 
QuantumNumberType m_nlms
 the quantum number of this node More...
 

Private Member Functions

void addGaussian (xmlNodePtr cur)
 
void addGaussianH5 (hdf_archive &hin)
 
void addSlater (xmlNodePtr cur)
 
template<typename Fin , typename T >
find_cutoff (Fin &in, T rmax)
 compute the safe cutoff radius of a radial functor More...
 

Private Attributes

hdf_archive hin
 hdf file only for numerical basis h5 file generated by SQD More...
 
int radial_grid_size_
 Number of grid points. More...
 
RealType m_rcut
 maximum cutoff for an orbital, could come from XML always set to -1 by constructor More...
 
RealType m_rcut_safe
 safe common cutoff radius More...
 
std::vector< std::unique_ptr< TransformerBase< RealType > > > radTemp
 radial functors to be finalized More...
 
std::tuple< int, double, double > grid_param_in
 

Additional Inherited Members

- Protected Attributes inherited from MPIObjectBase
CommunicatemyComm
 pointer to Communicate More...
 
std::string ClassName
 class Name More...
 
std::string myName
 name of the object More...
 

Detailed Description

template<typename COT>
class qmcplusplus::RadialOrbitalSetBuilder< COT >

Build a set of radial orbitals at the origin.

For a center,

  • only one grid is used
  • any number of radial orbitals
  • using MultiQuinticSpline1D
  • Gaussian and Slater mixing allowed
  • Slater's are assumed to fit in a 100 rmax grid

Definition at line 187 of file MultiFunctorAdapter.h.

Member Typedef Documentation

◆ GridType

using GridType = typename COT::GridType

Definition at line 95 of file RadialOrbitalSetBuilder.h.

◆ RadialOrbitalType

using RadialOrbitalType = typename COT::RadialOrbital_t

Definition at line 94 of file RadialOrbitalSetBuilder.h.

◆ RealType

using RealType = typename COT::RealType

Definition at line 93 of file RadialOrbitalSetBuilder.h.

Constructor & Destructor Documentation

◆ RadialOrbitalSetBuilder()

RadialOrbitalSetBuilder ( Communicate comm,
COT &  aos,
int  radial_grid_size = 1001 
)

constructor

Definition at line 155 of file RadialOrbitalSetBuilder.h.

156  : MPIObjectBase(comm), Normalized(true), m_orbitals(aos), radial_grid_size_(radial_grid_size), m_rcut(-1.0)
157 {}
RealType m_rcut
maximum cutoff for an orbital, could come from XML always set to -1 by constructor ...
MPIObjectBase(Communicate *c)
constructor with communicator
bool Normalized
true, if the RadialOrbitalType is normalized
int radial_grid_size_
Number of grid points.

Member Function Documentation

◆ addGaussian()

void addGaussian ( xmlNodePtr  cur)
private

Definition at line 284 of file RadialOrbitalSetBuilder.h.

References GaussianCombo< T >::putBasisGroup().

285 {
286  int L = m_nlms[1];
287  using gto_type = GaussianCombo<OHMMS_PRECISION_FULL>;
288  auto gset = std::make_unique<gto_type>(L, Normalized);
289  gset->putBasisGroup(cur);
290  //Warning::Magic Number for max rmax of gaussians
291  RealType r0 = find_cutoff(*gset, 100.);
292  m_rcut_safe = std::max(m_rcut_safe, r0);
293  radTemp.push_back(std::make_unique<A2NTransformer<RealType, gto_type>>(std::move(gset)));
294  m_orbitals.RnlID.push_back(m_nlms);
295 }
QuantumNumberType m_nlms
the quantum number of this node
T find_cutoff(Fin &in, T rmax)
compute the safe cutoff radius of a radial functor
std::vector< std::unique_ptr< TransformerBase< RealType > > > radTemp
radial functors to be finalized
RealType m_rcut_safe
safe common cutoff radius
QMCTraits::RealType RealType
bool Normalized
true, if the RadialOrbitalType is normalized

◆ addGaussianH5()

void addGaussianH5 ( hdf_archive hin)
private

Definition at line 299 of file RadialOrbitalSetBuilder.h.

References GaussianCombo< T >::putBasisGroupH5().

300 {
301  int L = m_nlms[1];
302  using gto_type = GaussianCombo<OHMMS_PRECISION_FULL>;
303  auto gset = std::make_unique<gto_type>(L, Normalized);
304  gset->putBasisGroupH5(hin, *myComm);
305  //at least gamess derived xml seems to provide the max its grid goes to
306  //So in priniciple this 100 should be coming in from input
307  //m_rcut seems like it once served this purpose but is somehow
308  //a class global variable even though it should apply here and
309  //similar locations on a function by function basis.
310  RealType r0 = find_cutoff(*gset, 100.);
311  m_rcut_safe = 6 * std::max(m_rcut_safe, r0);
312  radTemp.push_back(std::make_unique<A2NTransformer<RealType, gto_type>>(std::move(gset)));
313  m_orbitals.RnlID.push_back(m_nlms);
314 }
QuantumNumberType m_nlms
the quantum number of this node
T find_cutoff(Fin &in, T rmax)
compute the safe cutoff radius of a radial functor
std::vector< std::unique_ptr< TransformerBase< RealType > > > radTemp
radial functors to be finalized
Communicate * myComm
pointer to Communicate
Definition: MPIObjectBase.h:62
RealType m_rcut_safe
safe common cutoff radius
hdf_archive hin
hdf file only for numerical basis h5 file generated by SQD
QMCTraits::RealType RealType
bool Normalized
true, if the RadialOrbitalType is normalized

◆ addGrid()

bool addGrid ( xmlNodePtr  cur,
const std::string &  rad_type 
)

add a grid

Definition at line 160 of file RadialOrbitalSetBuilder.h.

References OneDimGridFactory::createGrid().

Referenced by AOBasisBuilder< COT >::createAOSet().

161 {
162  if (rad_type == "Numerical")
163  {
164  hin.push("grid");
165  addGridH5(hin);
166  hin.pop();
167  }
168  else
170 
171  //set zero to use std::max
172  m_rcut_safe = 0;
173 
174  return true;
175 }
static std::unique_ptr< GridType > createGrid(xmlNodePtr cur)
return a GridType*
RealType m_rcut_safe
safe common cutoff radius
void push(const std::string &gname, bool createit=true)
push a group to the group stack
hdf_archive hin
hdf file only for numerical basis h5 file generated by SQD
std::unique_ptr< GridType > input_grid
input grid in case transform is needed

◆ addGridH5()

bool addGridH5 ( hdf_archive hin)

Definition at line 178 of file RadialOrbitalSetBuilder.h.

References qmcplusplus::app_log(), and hdf_archive::read().

Referenced by AOBasisBuilder< COT >::createAOSetH5().

179 {
180  app_log() << " Grid is created by the input parameters in h5" << std::endl;
181 
182  std::string gridtype;
183  if (myComm->rank() == 0)
184  {
185  hin.read(gridtype, "grid_type");
186  }
187  myComm->bcast(gridtype);
188 
189  int npts = 0;
190  RealType ri = 0.0, rf = 10.0, rmax_safe = 10;
191  if (myComm->rank() == 0)
192  {
193  double tt = 0;
194  hin.read(tt, "grid_ri");
195  ri = tt;
196  hin.read(tt, "grid_rf");
197  rf = tt;
198  // Ye TODO: grid handling will all moved to XML.
199  //hin.read(tt, "rmax_safe");
200  //rmax_safe = tt;
201  hin.read(npts, "grid_npts");
202  }
203  myComm->bcast(ri);
204  myComm->bcast(rf);
205  myComm->bcast(rmax_safe);
206  myComm->bcast(npts);
207 
208  if (gridtype.empty())
209  myComm->barrier_and_abort("Grid type is not specified.");
210 
211  if (gridtype == "log")
212  {
213  app_log() << " Using log grid ri = " << ri << " rf = " << rf << " npts = " << npts << std::endl;
214  input_grid = std::make_unique<LogGrid<RealType>>();
215  }
216  else if (gridtype == "linear")
217  {
218  app_log() << " Using linear grid ri = " << ri << " rf = " << rf << " npts = " << npts << std::endl;
219  input_grid = std::make_unique<LinearGrid<RealType>>();
220  }
221 
222  input_grid->set(ri, rf, npts);
223 
224  //set zero to use std::max
225  m_rcut_safe = 0;
226 
227  return true;
228 }
int rank() const
return the rank
Definition: Communicate.h:116
std::ostream & app_log()
Definition: OutputManager.h:65
Communicate * myComm
pointer to Communicate
Definition: MPIObjectBase.h:62
RealType m_rcut_safe
safe common cutoff radius
hdf_archive hin
hdf file only for numerical basis h5 file generated by SQD
QMCTraits::RealType RealType
void read(T &data, const std::string &aname)
read the data from the group aname and check status runtime error is issued on I/O error ...
Definition: hdf_archive.h:306
void bcast(T &)
void barrier_and_abort(const std::string &msg) const
std::unique_ptr< GridType > input_grid
input grid in case transform is needed

◆ addRadialOrbital()

bool addRadialOrbital ( xmlNodePtr  cur,
const std::string &  m_infunctype,
const QuantumNumberType nlms 
)

add a radial functor

Add a new Slater Type Orbital with quantum numbers $(n,l,m,s)$.

Parameters
curxml element
nlmsquantum number
curthe current xmlNode to be processed
nlmsa vector containing the quantum numbers $(n,l,m,s)$
Returns
true is succeeds

Definition at line 237 of file RadialOrbitalSetBuilder.h.

References OhmmsAttributeSet::add(), and OhmmsAttributeSet::put().

Referenced by AOBasisBuilder< COT >::createAOSet().

240 {
241  std::string radtype(m_infunctype);
242  std::string dsname("0");
243  OhmmsAttributeSet aAttrib;
244  aAttrib.add(radtype, "type");
245  aAttrib.add(m_rcut, "rmax");
246  aAttrib.add(dsname, "ds");
247  aAttrib.put(cur);
248  m_nlms = nlms;
249  if (radtype == "Gaussian" || radtype == "GTO")
250  addGaussian(cur);
251  else if (radtype == "Slater" || radtype == "STO")
252  addSlater(cur);
253  else
254  myComm->barrier_and_abort("Purely numerical atomic orbitals are not supported any longer.");
255  return true;
256 }
QuantumNumberType m_nlms
the quantum number of this node
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
RealType m_rcut
maximum cutoff for an orbital, could come from XML always set to -1 by constructor ...
Communicate * myComm
pointer to Communicate
Definition: MPIObjectBase.h:62
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
void barrier_and_abort(const std::string &msg) const
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

◆ addRadialOrbitalH5()

bool addRadialOrbitalH5 ( hdf_archive hin,
const std::string &  radtype,
const QuantumNumberType nlms 
)

Definition at line 259 of file RadialOrbitalSetBuilder.h.

References hdf_archive::read().

Referenced by AOBasisBuilder< COT >::createAOSetH5().

262 {
263  std::string dsname("0");
264  std::string radtype(radtype_atomicBasisSet);
265  if (myComm->rank() == 0)
266  hin.read(radtype, "type");
267  myComm->bcast(radtype);
268 
269  m_nlms = nlms;
270  if (radtype == "Gaussian" || radtype == "GTO")
272  else if (radtype == "Slater" || radtype == "STO")
273  // addSlaterH5(hin);
275  " RadType: Slater. Any type other than Gaussian not implemented in H5 format. Please contact developers.");
276  else
278  " RadType: Numerical. Any type other than Gaussian not implemented in H5 format. Please contact developers.");
279  return true;
280 }
QuantumNumberType m_nlms
the quantum number of this node
int rank() const
return the rank
Definition: Communicate.h:116
Communicate * myComm
pointer to Communicate
Definition: MPIObjectBase.h:62
hdf_archive hin
hdf file only for numerical basis h5 file generated by SQD
void read(T &data, const std::string &aname)
read the data from the group aname and check status runtime error is issued on I/O error ...
Definition: hdf_archive.h:306
void bcast(T &)
void barrier_and_abort(const std::string &msg) const

◆ addSlater()

void addSlater ( xmlNodePtr  cur)
private

Definition at line 352 of file RadialOrbitalSetBuilder.h.

References SlaterCombo< T >::putBasisGroup().

353 {
354  using sto_type = SlaterCombo<OHMMS_PRECISION_FULL>;
355  auto gset = std::make_unique<sto_type>(m_nlms[1], Normalized);
356 
357  gset->putBasisGroup(cur);
358 
359  //need a find_cutoff for STO's, but this was previously in finalize and wiping out GTO's m_rcut_safe
360  m_rcut_safe = std::max(m_rcut_safe, static_cast<RealType>(100));
361  radTemp.push_back(std::make_unique<A2NTransformer<RealType, sto_type>>(std::move(gset)));
362  m_orbitals.RnlID.push_back(m_nlms);
363 }
QuantumNumberType m_nlms
the quantum number of this node
std::vector< std::unique_ptr< TransformerBase< RealType > > > radTemp
radial functors to be finalized
RealType m_rcut_safe
safe common cutoff radius
bool Normalized
true, if the RadialOrbitalType is normalized

◆ finalize()

void finalize ( )

This is when the radial orbitals are actually created.

Definition at line 328 of file RadialOrbitalSetBuilder.h.

References qmcplusplus::app_log(), qmcplusplus::convert(), and qmcplusplus::Units::charge::e.

Referenced by AOBasisBuilder< COT >::createAOSet(), and AOBasisBuilder< COT >::createAOSetH5().

329 {
330  // This is a temporary grid used in conversion, at the full precision of the calculation
331  // to reduce error. It doesn't need to be on the heap but this is the result of a
332  // series of design decisions requiring a base class pointer here.
333  std::unique_ptr<OneDimGridBase<OHMMS_PRECISION_FULL>> grid_prec;
334  grid_prec = std::make_unique<LogGrid<OHMMS_PRECISION_FULL>>();
335  // FIXME: should not hard-coded, probably should be input grid
336  grid_prec->set(1.e-6, m_rcut_safe, 1001);
337 
338  auto& multiset = m_orbitals.MultiRnl;
339  const int norbs = radTemp.size();
340  multiset.initialize(*grid_prec, norbs);
341 
342  for (int ib = 0; ib < norbs; ++ib)
343  radTemp[ib]->convert(*grid_prec, multiset, ib, 5);
344 
345  multiset.finalize();
346 
347  app_log() << " Setting cutoff radius " << m_rcut_safe << std::endl << std::endl;
348  m_orbitals.setRmax(static_cast<RealType>(m_rcut_safe));
349 }
void convert(const PL &lat, const PV &pin, PV &pout)
std::ostream & app_log()
Definition: OutputManager.h:65
std::vector< std::unique_ptr< TransformerBase< RealType > > > radTemp
radial functors to be finalized
RealType m_rcut_safe
safe common cutoff radius

◆ find_cutoff()

T find_cutoff ( Fin &  in,
rmax 
)
private

compute the safe cutoff radius of a radial functor

temporary function to compute the cutoff without constructing NGFunctor

Definition at line 371 of file RadialOrbitalSetBuilder.h.

References qmcplusplus::abs(), qmcplusplus::Units::charge::e, and LogGridLight< T >::set().

372 {
373  LogGridLight<OHMMS_PRECISION_FULL> agrid;
374  //WARNING Magic number, should come from input or be set somewhere more cnetral.
375  const OHMMS_PRECISION_FULL eps = 1e-6;
376  bool too_small = true;
378  int i = radial_grid_size_ - 1;
379  T r = rmax;
380  while (too_small && i > 0)
381  {
382  r = agrid(i--);
383  T x = in.f(r);
384  too_small = (std::abs(x) < eps);
385  }
386  return static_cast<OHMMS_PRECISION_FULL>(r);
387 }
MakeReturn< UnaryNode< FnFabs, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t abs(const Vector< T1, C1 > &l)
int radial_grid_size_
Number of grid points.

Member Data Documentation

◆ grid_param_in

std::tuple<int, double, double> grid_param_in
private

Definition at line 151 of file RadialOrbitalSetBuilder.h.

◆ hin

hdf_archive hin
private

hdf file only for numerical basis h5 file generated by SQD

Definition at line 136 of file RadialOrbitalSetBuilder.h.

Referenced by RadialOrbitalSetBuilder< SoaAtomicBasisSet< MultiFunctorAdapter< FN >, SH > >::addRadialOrbitalH5().

◆ input_grid

std::unique_ptr<GridType> input_grid

input grid in case transform is needed

Definition at line 103 of file RadialOrbitalSetBuilder.h.

◆ m_nlms

the quantum number of this node

Definition at line 105 of file RadialOrbitalSetBuilder.h.

◆ m_orbitals

◆ m_rcut

RealType m_rcut
private

maximum cutoff for an orbital, could come from XML always set to -1 by constructor

Definition at line 142 of file RadialOrbitalSetBuilder.h.

◆ m_rcut_safe

RealType m_rcut_safe
private

safe common cutoff radius

Definition at line 145 of file RadialOrbitalSetBuilder.h.

◆ Normalized

◆ radial_grid_size_

int radial_grid_size_
private

Number of grid points.

Definition at line 139 of file RadialOrbitalSetBuilder.h.

◆ radTemp

std::vector<std::unique_ptr<TransformerBase<RealType> > > radTemp
private

radial functors to be finalized

Definition at line 149 of file RadialOrbitalSetBuilder.h.


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