qharv.seed package

Submodules

qharv.seed.lcao_h5 module

qharv.seed.lcao_h5.abs_grid(fp, iabs)[source]

Extract <grid> for some <atomicBasisSet>

Parameters

fp (h5py.File) – LCAO h5

Returns

<grid>

Return type

lxml.etree.Element

qharv.seed.lcao_h5.basis_group(fp, iabs, ibg)[source]

Extract <basisGroup>

Parameters

fp (h5py.File) – LCAO h5

Returns

<basisGroup>

Return type

lxml.etree.Element

qharv.seed.lcao_h5.basisset(fp)[source]

Extract <basisset>

Parameters

fp (h5py.File) – LCAO h5

Returns

<basisset>

Return type

lxml.etree.Element

qharv.seed.lcao_h5.bg_radfunc(fp, iabs, ibg, irf)[source]

Extract <radfunc> for some <basisGroup>

Parameters

fp (h5py.File) – LCAO h5

Returns

<radfunc>

Return type

lxml.etree.Element

qharv.seed.lcao_h5.sposet(fp, bsname, cname, nstate=None, ik=0, ispin=0)[source]

qharv.seed.qmcpack_in module

qharv.seed.qmcpack_in.all_electron_hamiltonian(elec_name='e', ion_name='ion0')[source]
qharv.seed.qmcpack_in.assemble_project(nodel, name='qmc', series=0)[source]

assemble QMCPACK input using a list of xml nodes

usually nodel=[qmcsystem, qmc]

Parameters
  • nodel (list) – a list of xml node (lxml.Element)

  • name (str, optional) – project name, default ‘qmc’

qharv.seed.qmcpack_in.bspline_qmcsystem(fh5, tmat=None)[source]

Create Slater-Jastrow system input from pw2qmcpack.x h5 file

Parameters
  • fh5 (str) – path to wf h5 file

  • tmat (np.array) – tile matrix

Returns

<qmcsystem>

Return type

qsys (etree.Element)

qharv.seed.qmcpack_in.bundle_twists(calc_dir, fregex='*twistnum_*.in.xml')[source]

bundle all twist inputs

quick and dirty: cd $calc_dir; ls > prefix.in, then edit prefix.in

Parameters
  • calc_dir (str) – calculation directory

  • fregex (str, optional) – regular expression for all twists

Returns

bundled input text

Return type

str

qharv.seed.qmcpack_in.disperse(ginp_loc, calc_dir, execute=False, overwrite=False)[source]

disperse inputs bundled up in a grouped input

Parameters
  • ginp_loc (str) – location of grouped input e.g. ../runs/dmc/qmc.in

  • calc_dir (str) – folder to output new inputs e.g. dmc1

  • execute (bool,optional) – perform file I/O, default is False i.e. a dry run

  • overwrite (bool,optional) – overwrite existing files, default is False

Returns

a list of new inputs

Return type

list

qharv.seed.qmcpack_in.expand_twists(example_in_xml, twist_list, calc_dir, force=False)[source]

expand example input xml to all twists in twist_list

Naming convention of new inputs:

[prefix].g001.twistnum_[itwist].in.xml

Parameters
  • example_in_xml (str) – example QMCPACK input xml file

  • twist_list (list) – a list of twist indices

  • calc_dir (str) – folder to output new inputs

Returns

None

Examples

>>> expand_twists('./vmc.in.xml',range(64),'.')
>>> expand_twists('./ref/vmc.in.xml',[0,15,17],'./new')
qharv.seed.qmcpack_in.particle_group_from_pos(pos, name, charge, **kwargs)[source]

construct a <group> in the <particleset> xml element

Parameters
  • pos (np.array) – positions, shape (nptcl, ndim)

  • name (str) – name of particle group

  • charge (float) – the amount of charge of this particle species

Returns

<group> including <attrib name=”positions”>

Return type

etree.Element

qharv.seed.qmcpack_in.set_gc_occ(norbl, calc_dir, fregex_fmt='*twistnum_{itwist:d}.in.xml')[source]

edit twist inputs in calc_dir according to occupation vector norbl

Parameters
  • norbl (list) – number of occupied orbitals at each twist

  • calc_dir (str) – location of twist inputs

qharv.seed.qmcpack_in.set_norb(doc, norb)[source]

occupy the lowest norb Kohn-Sham orbitals

Parameters
  • doc (etree.Element) – must contain <particleset>, <sposet>, and <det…set>

  • norb (int) – number of orbitals to occupy

Returns

None

Effect:

doc is modified

qharv.seed.qmcpack_in.set_nwalker(doc, nwalker)[source]

set the number of walkers to use in DMC

Parameters
  • doc (lxml.Element) – xml node containing <qmc>

  • nwalker (int) – number of walkers

qharv.seed.qmcpack_in.simulationcell_from_axes(axes, bconds='p p p', rckc=15.0)[source]

construct the <simulationcell> xml element from axes

Parameters
  • axes (np.array) – lattice vectors

  • bconds (str, optional) – boundary conditions in x,y,z directions. p for periodic, n for non-periodic, default to ‘p p p’

  • rckc – long-range cutoff paramter rc*kc, default to 15

Returns

representing <simulationcell>

Return type

etree.Element

qharv.seed.qmcpack_in.ud_electrons(nup, ndown)[source]

construct the <particleset name=”e”> xml element for electrons

Parameters
  • nup (int) – number of up electrons

  • ndown (int) – number of down electrons

Returns

representing <particleset name=”e”>

Return type

etree.Element

qharv.seed.wf_h5 module

qharv.seed.wf_h5.axes_elem_charges_pos(fp)[source]

extract lattice vectors, atomic positions, and element names The main difficulty is constructing the element names of each atomic species. If elem is not needed, use get(fp,’axes’) and get(fp,’pos’) to get the simulation cell and ion positions directly.

Parameters

fp (h5py.File) – hdf5 file object

Returns

(axes, elem, vchg_map, pos)

Return type

(np.array, np.array, dict, np.array)

qharv.seed.wf_h5.get(fp, name)[source]

retrieve data from a known location in pwscf.h5

Parameters
  • fp (h5py.File) – hdf5 file object

  • name (str) – a known name in locations

Returns

whatever fp[loc][()] returns

Return type

array_like

qharv.seed.wf_h5.get_bands(fp, ispin=0)[source]

return the list of available Kohn-Sham eigenvalues

Parameters
  • fp (h5py.File) – wf h5 file

  • ispin (int, optional) – spin index, default 0

Returns

tvecs, twist vectors in reciprocal lattice units (nk, nbnd)

Return type

np.array

qharv.seed.wf_h5.get_cmat(fp, ikpt, ispin, norb=None, npw=None)[source]

get Kohn-Sham orbital coefficients on a list of plane waves (PWs)

Parameters
  • fp (h5py.File) – wf h5 file

  • ikpt (int) – kpoint index

  • ispin (int) – spin index

  • norb (int, optional) – number of orbitals at this kpoint

  • npw (int, optional) – number of PW for each orbital

Returns

cmat orbital coefficient matrix

Return type

np.array

qharv.seed.wf_h5.get_orb_in_pw(fp, ikpt, ispin, istate)[source]

get the plane wave coefficients of a single Kohn-Sham orbital

Parameters
  • fp (h5py.File) – wf h5 file

  • ikpt (int) – kpoint index

  • ispin (int) – spin index

  • istate (int) – band index

Returns

(gvecs, cmat), PWs and coefficient matrix

Return type

(np.array, np.array)

qharv.seed.wf_h5.get_orb_on_grid(fp, ikpt, ispin, istate, grid_shape=None)[source]

get a single Kohn-Sham orbital on a real-space grid

Parameters
  • fp (h5py.File) – wf h5 file

  • ikpt (int) – kpoint index

  • ispin (int) – spin index

  • istate (int) – band index

  • grid_shape (np.array) – (3,) grid shape [nx, ny, nz]

Returns

orbital on grid with shape=grid_shape

Return type

np.array

qharv.seed.wf_h5.get_orbs(fp, orbs, truncate=False, tol=1e-08)[source]

return the list of requested Kohn-Sham orbitals

Parameters
  • fp (h5py.File) – wf h5 file

  • orbs (list) – a list of 3-tuples, each tuple species the KS state by (kpoint/twist, spin, band) i.e. (ik, ispin, ib)

  • truncate (bool, optional) – remove PWs with ``small’’ coefficient

  • tol (float, optional) – define ``small’’ as |ck|^2 < tol

qharv.seed.wf_h5.get_sc_cmat(fp, itwist, ispin, noccl, sorted_orbs=False)[source]

get Kohn-Sham orbital coefficients at a supercell twist

Parameters
  • fp (h5py.File) – wf h5 file

  • itwist (int) – twist index

  • ispin (int) – spin index

  • noccl (np.array) – number of occupied orbitals at each kpoint

Returns

cmat orbital coefficient matrix

Return type

np.array

qharv.seed.wf_h5.get_twist(fp, itwist)[source]
qharv.seed.wf_h5.get_twists(fp, ndim=3)[source]

return the list of available twist vectors

Parameters
  • fp (h5py.File) – wf h5 file

  • ndim (int, optional) – number of spatial dimensions, default 3

Returns

tvecs, twist vectors in reciprocal lattice units (nk, ndim)

Return type

np.array

qharv.seed.wf_h5.kpoint_path(ikpt)[source]

construct path to kpoint

e.g. electrons/kpoint_0/spin_0/state_0

Parameters

ikpt (int) – kpoint index

Returns

path in hdf5 file

Return type

str

qharv.seed.wf_h5.ls(handle, r=False, level=0, indent=' ')[source]

List directory structure

Similar to the Linux ls command, but for an hdf5 file

Parameters
  • handle (h5py.Group) – or h5py.File or h5py.Dataset

  • r (bool) – recursive list

  • level (int) – level of indentation, only used if r=True

  • indent (str) – indent string, only used if r=True

Returns

mystr, a string representation of the directory structure

Return type

str

qharv.seed.wf_h5.normalize_cmat(cmat)[source]

normalize PW orbital coefficients

Parameters

cmat (np.array) – coefficient matrix shape (norb, npw)

Effect:

each row of cmat will be normalized to |ci|^2=1

qharv.seed.wf_h5.read(fname, **kwargs)[source]

read h5 file and return a h5py File object

Parameters
  • fname (str) – hdf5 file

  • kwargs (dict) – keyword arguments to pass on to h5py.File, default is {‘mode’: ‘r’}

Returns

h5py File object

Return type

h5py.File

qharv.seed.wf_h5.spin_path(ikpt, ispin)[source]
qharv.seed.wf_h5.state_path(ikpt, ispin, istate)[source]
qharv.seed.wf_h5.write_atoms(fp, elem, pos, pseudized_charge)[source]

create and fill the /atoms group !!!! QMCPACK does NOT check valence_charge, pseudized_charge matters?

Parameters
  • fp (h5py.File) – hdf5 file object

  • elem (np.array) – array of atom names

  • pos (np.array) – array of atomic coordinates in Bohr

  • pseudized_charge (dict) – number of pseudized electrons for each species

Returns

species_id, a list of atomic numbers for each atom

Return type

list

Effect:

fill /atoms group in ‘fp’

Example

>>> fp = h5py.File('pwscf.pwscf.h5', 'a')
>>> wf_h5.write_atoms(fp, ['Li'],
>>>                   np.array([[0, 0, 0]]),
>>>                   pseudized_charge={'Li': 2})
>>> fp.close()
qharv.seed.wf_h5.write_gvecs(fp, gvecs, kpath='/electrons/kpoint_0')[source]

fill the electrons/kpoint_0/gvectors group in wf h5 file

Parameters
  • fp (h5py.File) – hdf5 file object

  • gvecs (np.array) – PW basis as integer vectors

  • kpath (str, optional) – kpoint group to contain gvecs, default is ‘/electrons/kpoint_0’

Example

>>> fp = h5py.File('pwscf.pwscf.h5', 'w')
>>> write_gvecs(fp, gvecs)
>>> fp.close()
qharv.seed.wf_h5.write_kpoint(kgrp, ikpt, utvec, evals, cmats)[source]

fill the electrons/kpoint_$ikpt group in wf h5 file

Parameters
  • kgrp (h5py.Group) – kpoint group

  • ikpt (int) – twist index

  • utvec (np.array) – twist vector in reduced units

  • evals (list) – list of Kohn-Sham eigenvalues to sort orbitals; one real np.array of shape (norb) for each spin

  • cmats (list) – list of Kohn-Sham orbitals in PW basis; one complex np.array of shape (norb, npw) for each spin

Example

>>> fp = h5py.File('pwscf.pwscf.h5', 'w')
>>> kgrp = fp.require_group('/electrons/kpoint_0')
>>> evals = [ np.array([0]) ]  # 1 spin, 1 state
>>> cmats = [ np.array([[0]], dtype=complex) ]
>>> write_kpoint(kgrp, 0, [0, 0, 0], evals, cmats)
>>> fp.close()
qharv.seed.wf_h5.write_misc(fp, nup_ndn, nkpt)[source]

fill /electrons/number_of_* and /format

Parameters
  • fp (h5py.File) – hdf5 file object

  • nup_ndn (tuple/list) – number of up and down electrons (nup, ndn)

  • nkpt (int) – number of kpoints/twists

qharv.seed.wf_h5.write_supercell(fp, axes)[source]

create and fill the /supercell group

Parameters
  • fp (h5py.File) – hdf5 file object

  • axes (np.array) – lattice vectors in Bohr

qharv.seed.wf_h5.write_wf(egrp, utvecs, gvecs, evalsl, cmatsl)[source]

fill the wf portion of the electrons group in wf h5 file !!!! WARNING: this function may require too much memory; if so, use write_kpoint directly

Parameters
  • egrp (h5py.Group) – electrons group

  • utvecs (np.array) – all twist vectors in reduced units

  • evalsl (list) – list of Kohn-Sham eigenvalues to sort orbitals; one real np.array of shape (norb) for each spin and twist

  • cmatsl (list) – list of Kohn-Sham orbitals in PW basis; one complex np.array of shape (norb, npw) for each spin and twist

qharv.seed.xml module

qharv.seed.xml.add_backflow(wf_node, bf_node)[source]
qharv.seed.xml.append(root, nodes, copy=True)[source]

append one or more nodes to a root node

Parameters
  • root (lxml.etree._Element) – xml node

  • nodes (list or lxml.etree._Element) – xml node(s)

  • copy (bool, optional) – make copy of nodes to append

qharv.seed.xml.arr2text(arr)[source]

format a numpy array into a text string

qharv.seed.xml.build_coeff(knots, **attribs)[source]

construct an <coefficients/>

Example

build_coeff([1,2]):

<coefficients id=”new” type=”Array”> 1 2 </coefficients>

Parameters

knots (list) – a list of numbers

Returns

<coefficients/>

Return type

lxml.etree._Element

qharv.seed.xml.build_jk2_iso(coeffs, kc)[source]

construct isotropic e-e reciprocal space Jastrow node

Example

build_jk2([1,2], 0.4):

Parameters
  • coefs (list) – a list of numbers at the

  • kc (float) – k space cutoff in a.u.

Returns

<jastrow/>

Return type

lxml.etree._Element

qharv.seed.xml.build_jr2(uuc, udc)[source]
qharv.seed.xml.dset2spo(wf_node, det_map)[source]

change <wavefunction> from old style, <basis> in <determinantset>, to new style, <basis> in <sposet_builder>

Parameters
  • wf_node (etree.Element) – <wavefunction> node

  • det_map (dict) – determinant name -> particle group name e.g. {‘updet’:’u’,’downdet’:’d’}

Returns

None

qharv.seed.xml.get_axes(doc)[source]
qharv.seed.xml.get_nelec(doc)[source]
qharv.seed.xml.get_param(node, pname)[source]
retrieve the str representation of a parameter from:

<parameter name=”pname”> str_rep </parameter>

Parameters
  • node (lxml.etree._Element) – xml node with <parameter>.

  • pname (str) – name of parameter

Returns

string representation of the parameter value

Return type

str

qharv.seed.xml.get_pos(doc, pset='ion0', all_pos=True, group=None)[source]
qharv.seed.xml.ls(node, r=False, level=0, indent=' ')[source]

List directory structure

Similar to the Linux ls command, but for an xml node

Parameters
  • node (lxml.etree._Element) – xml node

  • r (bool) – recursive

  • level (int) – level of indentation, used only if r=True

  • indent (str) – indent string, used only if r=True

Returns

mystr, a string representation of the directory structure

Return type

str

qharv.seed.xml.make_node(tag, attribs=None, text=None, pad=' ')[source]
create etree.Element

<tag **attribs> text </tag>

Args:

tag (str): tag node attribs (dict, optional): attributes, default None text (str, optional): text content, default None pad (str, optional): padding for text, default ‘ ‘

Return:

etree.Element: node

Examples:
>>> sim = make_node('simulation')
>>> epset = make_node('particleset', {'name': 'e'})
>>> lrnode = make_node('parameter', {'name': 'LR_dim_cutoff'}, str(20))
>>> bconds = make_node('parameter', {'name': 'bconds'}, 'p p p', pad='
‘)
>>> seed = make_node('seed', text=str(31415))
qharv.seed.xml.parse(text)[source]

parse the text representation of an xml node delegate to read()

Parameters

text (str) – string representation of an xml node

Returns

root, parsed xml node

Return type

lxml.etree._Element

qharv.seed.xml.read(fname)[source]

read an xml file wrap around lxml.etree.parse

Parameters

fname (str) – filename to read from

Returns

doc, parsed xml document

Return type

lxml.etree._ElementTree

qharv.seed.xml.remove(nodes)[source]

remove nodes from the xml tree

Parameters

nodes (list) – xml nodes

qharv.seed.xml.set_param(node, pname, pval, new=False, pad=' ')[source]

set <parameter> with name ‘pname’ to ‘pval’ if new=True, then <parameter name=”pname”> does not exist. create it

Parameters
  • node (lxml.etree._Element) – xml node with children having tag ‘parameter’

  • pname (str) – name of parameter

  • pval (str) – value of parameter

  • new (bool) – create new <paremter> node, default is false

Effect:

the text of <parameter> with ‘pname’ will be set to ‘pval’

qharv.seed.xml.show(node)[source]
qharv.seed.xml.str_rep(node)[source]

return the string representation of an xml node

Parameters

node (lxml.etree._Element) – xml node

Returns

string representation of node

Return type

str

qharv.seed.xml.swap_node(node0, node1)[source]

replace the node0 with node1 node0 must have a parent

Parameters
  • node0 (etree.Element) – node to be swapped out

  • node1 (etree.Element) – replacement node

Effect:

node0 is replaced by node1 in node0’s owning tree

qharv.seed.xml.text2arr(text, dtype=<class 'float'>, flatten=False)[source]

convert a text string into a numpy array

qharv.seed.xml.text2vec(text, dtype=<class 'float'>)[source]

convert a text string into a 1D numpy array

qharv.seed.xml.turn_off_jas_opt(wf_node)[source]
qharv.seed.xml.write(fname, doc)[source]

write an xml file wrap around lxml.etree._ElementTree.write

Parameters
  • fname (str) – filename to write to

  • doc (lxml.etree._ElementTree) – xml file in memory

Effect:

write fname using contents of doc

qharv.seed.xml_examples module

qharv.seed.xml_examples.bcc54_dynamic_backflow()[source]
qharv.seed.xml_examples.bcc54_static_backflow()[source]
qharv.seed.xml_examples.dynamic_ae_ham()[source]
qharv.seed.xml_examples.ee_ham()[source]
qharv.seed.xml_examples.heg_system(rs, nshell_up, polarized)[source]

construct QMCPACK input xml <qmcsystem> node for the homogeneous electron gas (HEG). Momentum shells must be fully filled. The HEG must either be fully polarized or fully unpolarized.

number of particles is determined by the number of filled shells.

0 -> 1 1 -> 7 2 -> 19 3 -> 27 4 -> 33

e.g. heg_system(1.0, 2, False) for unpolarized HEG with 2 shells (38 electrons) at rs=1.0

  • note: this function aims to return the simplest functional input,

for more flexibility please edit the returned xml element instead of adding to this function.

Parameters
  • rs (float) – Wigner-Seitz radius

  • nshell_up (int) – number k shells filled by up electrons

  • polarized (bool) – True: polarized (ndn=0); False: unpolarized (nup=ndn)

Returns

qsys_node containing the <qmcsystem> xml node

Return type

lxml.etree.Element

qharv.seed.xml_examples.i4_dynamic_jastrow()[source]
qharv.seed.xml_examples.pbyp_dmc()[source]
qharv.seed.xml_examples.pbyp_var_optimize()[source]
qharv.seed.xml_examples.pbyp_vmc()[source]
qharv.seed.xml_examples.static_ae_ham()[source]
qharv.seed.xml_examples.static_coul_ae_ham()[source]
qharv.seed.xml_examples.wbyw_dmc()[source]
qharv.seed.xml_examples.wbyw_optimize()[source]
qharv.seed.xml_examples.wbyw_vmc()[source]

Module contents