Accumulators

class pyqmc.accumulators.EnergyAccumulator(mol, threshold=10, **kwargs)[source]

Returns local energy of each configuration in a dictionary.

class pyqmc.accumulators.LinearTransform(parameters, to_opt=None)[source]

Linearize a dictionary of wf parameters.

Parameters:
  • parameters (dict) – the wave function parameters
  • to_opt (dict) – is a dictionary with the keys to optimize, and its values are boolean arrays indicating which specific elements to optimize
Note:
to_opt[k] can’t be boolean scalar; it has to be an array with the same dimension as parameters[k]. to_opt doesn’t have to have all the keys of parameters, but all keys of to_opt must be keys of parameters.
deserialize(parameters)[source]

Convert serialized parameters to dictionary

serialize_gradients(pgrad)[source]

Convert the dictionary to a linear list of gradients, mask allows for certain fixed parameters

serialize_parameters(parameters)[source]

Convert the dictionary to a linear list of gradients

class pyqmc.accumulators.PGradTransform(enacc, transform, nodal_cutoff=0.001)[source]
avg(configs, wf, weights=None)[source]

Compute (weighted) average

class pyqmc.accumulators.SqAccumulator(qlist=None, Lvecs=None, nq=4)[source]

Accumulates structure factor

\[S(\vec{q}) = \langle \rho_{\vec{q}} \rho_{-\vec{q}} \rangle = \langle \left| \sum_{j=1}^{N_e} e^{i\vec{q}\cdot\vec{r}_j} \right| \rangle\]

Energy

Effective Core Potential (ECP)

pyqmc.eval_ecp.P_l(x, l)[source]

Legendre functions,

Parameters:
  • x ((nconf,) array) – distances x=r_ea(i)
  • l (int) – angular momentum channel
Returns:

legendre function P_l values for channel \(l\).

Return type:

(nconf, naip) array

pyqmc.eval_ecp.ecp(mol, configs, wf, threshold)[source]
Returns:ECP value, summed over all the electrons and atoms.
pyqmc.eval_ecp.ecp_ea(mol, configs, wf, e, atom, threshold)[source]
Returns:the ECP value between electron e and atom at, local+nonlocal.

TODO: update documentation

pyqmc.eval_ecp.ecp_mask(v_l, threshold)[source]
Returns:a mask for configurations sized nconf based on values of v_l. Also returns acceptance probabilities
pyqmc.eval_ecp.generate_ecp_functors(coeffs)[source]
Parameters:coeffsmol._ecp[atom_name][1] (coefficients of the ECP)
Returns:a functor v_l, with keys as the angular momenta: -1 stands for the nonlocal part, 0,1,2,… are the s,p,d channels, etc.
pyqmc.eval_ecp.get_P_l(r_ea, r_ea_vec, l_list)[source]

The factor \((2l+1)\) and the quadrature weights are included.

Parameters:
  • r_ea ((nconf,)) – distances of electron e and atom a
  • r_ea_vec ((nconf, 3)) – displacements of electron e and atom a
  • l_list (list) – [-1,0,1,…] list of given angular momenta
Returns:

legendre function P_l values for each \(l\) channel.

Return type:

(nconf, naip, nl) array

pyqmc.eval_ecp.get_rot(nconf, naip)[source]
Parameters:
  • nconf (int) – number of configurations
  • naip (int) – number of auxiliary integration points
Returns:

the integration weights, and the positions of the rotated electron e

Return type:

((naip,) array, (nconf, naip, 3) array)

pyqmc.eval_ecp.get_v_l(mol, at_name, r_ea)[source]
Returns:list of the \(l\)’s, and a nconf x nl array, v_l values for each \(l\): l= 0,1,2,…,-1
class pyqmc.eval_ecp.rnExp(n, e, c)[source]

v_l object.

\(cr^{n-2}\cdot\exp(-er^2)\)

Ewald summation

class pyqmc.ewald.Ewald(cell, ewald_gmax=200, nlatvec=1)[source]

The Ewald summation computes the Coulomb energy of a periodic arrangement of charges. The sum (of pair interactions) is separated into real-space (short range) and reciprocal-space (long range) sums, each of which converges quickly. The separation is determined by the parameter \(\alpha\)

The Ewald separation is:

\[E_{\rm Coulomb} = E_{\rm real\ space} + E_{\rm reciprocal\ space} + E_{\rm self} + E_{\rm charged}\]
\[E_{\rm real\ space} = \frac{1}{2} {\sum_{\vec{n}}}^\dagger \sum_{i=1}^N \sum_{j=1}^N q_i q_j \frac{{\rm erfc}(\alpha |\vec{x}_{ij}+\vec{n}|)}{|\vec{x}_{ij}+\vec{n}|}\]

\(\qquad\qquad\) (The \({\sum}^\dagger\) means to exclude the self-terms, i.e. \(\vec{n}=0, \, i=j\))

\[E_{\rm reciprocal\ space} = \frac{4\pi}{V} \frac{1}{2} \sum_{k \ne 0} \frac{1}{k^2} e^{-\frac{k^2}{4\alpha^2}} \left| \sum_{i=1}^N q_i e^{-i\vec{k}\cdot\vec{x}_i} \right|^2\]
\[E_{\rm self} = -\frac{\alpha}{\sqrt{\pi}} \sum_{i=1}^N q_i^2 \qquad E_{\rm charged} = -\frac{\pi}{2V\alpha^2} \left| \sum_{i=1}^N q_i \right|^2\]

The self energy corrects for a self-interaction included in the reciprocal-space term, and the charged-system correction is only necessary for systems with nonzero net charge.

In our implementation, the parts are further split into electron-electron, electron-ion, and ion-ion contributions. We use lower-case summation indices for electrons, and upper case for ions.

For ease of notation (and reading the code), let pair distances be denoted by

\[r_{ijn} = |\vec{x}_{ij}+\vec{n}|\]
\[r_{iIn} = |\vec{x}_{iI}+\vec{n}|\]
\[r_{IJn} = |\vec{x}_{IJ}+\vec{n}|\]

Real space terms, arranged to sum over each pair only once:

\[E_{\rm real\ space}^{\text{ion-ion}} = \sum_{\vec{n}} \sum_{I<J}^{N_{ion}} Z_I Z_J \frac{{\rm erfc}(\alpha r_{IJn})}{r_{IJn}} + \frac{1}{2} \sum_{I=1}^{N_{ion}} Z_I^2 C_{\rm self\ image}\]
\[E_{\rm real\ space}^{ee} = \sum_{\vec{n}} \sum_{i<j}^{N_e} \frac{{\rm erfc}(\alpha r_{ijn})}{r_{ijn}} + \frac{N_e}{2} C_{\rm self\ image}\]
\[E_{\rm real\ space}^{e\text{-ion}} = {\sum_{\vec{n}}} \sum_{i=1}^{N_e} \sum_{I=1}^{N_{ion}} -Z_I \frac{{\rm erfc}(\alpha r_{iIn})}{r_{iIn}}\]

where the interactions between a particle and its own image in other cells is represented by the sum

\[C_{\rm self\ image} = \frac{1}{2} \sum_{\vec{n} \ne 0} \frac{{\rm erfc}(\alpha |\vec{n}|)}{|\vec{n}|}\]

Reciprocal space terms, summing over \(\vec{G}>0\) – reciprocal lattice vectors in positive octants. In other words, only one of \(\vec{G}\) and \(-\vec{G}\) is included in the sum, and \(\vec{G}=0\) is not included.

\[E_{\rm reciprocal\ space}^{\text{ion-ion}} = \sum_{\vec{G} > 0 } W_G \left| \sum_{I=1}^{N_{ion}} Z_I e^{-i\vec{G}\cdot\vec{x}_I} \right|^2\]
\[E_{\rm reciprocal\ space}^{ee} = \sum_{\vec{G}>0} W_G \left| \sum_{i=1}^{N_e} e^{-i\vec{k}\cdot\vec{x}_i} \right|^2\]
\[E_{\rm reciprocal\ space}^{e\text{-ion}} = \sum_{\vec{G}>0} W_G {\rm Re} \left[ 2 \sum_{i=1}^{N_e} \sum_{I=1}^{N_{ion}} -Z_I e^{-i\vec{k}\cdot\vec{x}_i} e^{i\vec{k}\cdot\vec{x}_I} \right]\]

where gweight is a factor that doesn’t depend on the coordinates:

\[W_G = \frac{4\pi}{V |\vec{G}|^2} e^{- \frac{|\vec{G}|^2}{ 4\alpha^2}}\]

Self energy:

\[E_{\rm self}^{e} = - \frac{\alpha}{\sqrt{\pi}} N_e \qquad E_{\rm self}^{\rm ion} = - \frac{\alpha}{\sqrt{\pi}} \sum_{I=1}^{N_{ion}} Z_I^2\]

Charged-system energy:

\[E_{\rm charged}^{ee} = - \frac{\pi}{2V\alpha^2} N_e^2 \qquad E_{\rm charged}^{e\text{-ion}} = \frac{\pi}{2V\alpha^2} 2 N_e \sum_{I=1}^{N_{ion}} Z_I\]
\[E_{\rm charged}^{\text{ion-ion}} = - \frac{\pi}{2V\alpha^2} \left[ \sum_{I=1}^{N_{ion}} Z_I^2 + 2 \sum_{I<J}^{N_{ion}} Z_I Z_J \right]\]
energy(configs)[source]

Compute Coulomb energy for a set of configs.

\[\begin{split}E_{\rm Coulomb} &= E_{\rm real+reciprocal}^{ee} + E_{\rm self+charged}^{ee} \\&+ E_{\rm real+reciprocal}^{e\text{-ion}} + E_{\rm self+charged}^{e\text{-ion}} \\&+ E_{\rm real+reciprocal}^{\text{ion-ion}} + E_{\rm self+charged}^{\text{ion-ion}}\end{split}\]
Parameters:configs ((nconf, nelec, 3) PeriodicConfigs object) – electron positions (walkers)
Returns:
  • ee: electron-electron part
  • ei: electron-ion part
  • ii: ion-ion part
Return type:float, float, float
energy_separated(configs)[source]

Compute Coulomb energy contribution from each electron (does not include ion-ion energy).

NOTE: energy() needs to be called first to update the separated energy values

Parameters:configs ((nconf, nelec, 3) PeriodicConfigs object) – electron positions (walkers)
Returns:energies
Return type:(nelec,) array
ewald_electron(configs)[source]

Compute the Ewald sum for e-e and e-ion

Note: We ignore the constant term \(\frac{N_e}{2} C_{\rm self\ image}\) in the real-space e-e sum corresponding to the interaction of an electron with its own image in other cells.

Real space e-e:

\[E_{\rm real\ space}^{ee} = \sum_{\vec{n}} \sum_{i<j}^{N_e} \frac{{\rm erfc}(\alpha r_{ijn})}{r_{ijn}}\]

Real space e-i:

\[E_{\rm real\ space}^{e\text{-ion}} = {\sum_{\vec{n}}} \sum_{i=1}^{N_e} \sum_{I=1}^{N_{ion}} -Z_I \frac{{\rm erfc}(\alpha r_{iIn})}{r_{iIn}}\]

Reciprocal space e-e:

\[E_{\rm reciprocal\ space}^{ee} = \sum_{\vec{G}>0} W_G \left| \sum_{i=1}^{N_e} e^{-i\vec{k}\cdot\vec{x}_i} \right|^2\]

Reciprocal space e-i:

\[E_{\rm reciprocal\ space}^{e\text{-ion}} = \sum_{\vec{G}>0} W_G {\rm Re} \left[ 2 \sum_{i=1}^{N_e} \sum_{I=1}^{N_{ion}} -Z_I e^{-i\vec{k}\cdot\vec{x}_i} e^{i\vec{k}\cdot\vec{x}_I} \right]\]
Parameters:configs ((nconf, nelec, 3) PeriodicConfigs object) – electron positions (walkers)
Returns:
  • ee: electron-electron part
  • ei: electron-ion part
Return type:float, float
ewald_ion()[source]

Compute ion contribution to Ewald sums. Since the ions don’t move in our calculations, the ion-ion term only needs to be computed once.

Note: We ignore the constant term \(\frac{1}{2} \sum_{I} Z_I^2 C_{\rm self\ image}\) in the real-space ion-ion sum corresponding to the interaction of an ion with its own image in other cells.

The real-space part:

\[E_{\rm real\ space}^{\text{ion-ion}} = \sum_{\vec{n}} \sum_{I<J}^{N_{ion}} Z_I Z_J \frac{{\rm erfc}(\alpha |\vec{x}_{IJ}+\vec{n}|)}{|\vec{x}_{IJ}+\vec{n}|}\]

The reciprocal-space part:

\[E_{\rm reciprocal\ space}^{\text{ion-ion}} = \sum_{\vec{G} > 0 } W_G \left| \sum_{I=1}^{N_{ion}} Z_I e^{-i\vec{G}\cdot\vec{x}_I} \right|^2\]
Returns:ion-ion component of Ewald sum
Return type:float
set_ewald_constants(cellvolume)[source]

Compute Ewald constants (independent of particle positions): self energy and charged-system energy. Here we compute the combined terms. These terms are independent of the convergence parameters gmax and nlatvec, but do depend on the partitioning parameter \(\alpha\).

We define two constants, squareconst, the coefficient of the squared charges, and ijconst, the coefficient of the pairs:

\[C_{ij} = - \frac{\pi}{V\alpha^2}\]
\[C_{\rm square} = - \frac{\alpha}{\sqrt{\pi}} - \frac{\pi}{2V\alpha^2} = - \frac{\alpha}{\sqrt{\pi}} - \frac{C_{ij}}{2}\]

The Ewald object doesn’t retain information about the configurations, including number of electrons, so the electron constants are defined as functions of \(N_e\).

Self plus charged-system energy:

\[E_{\rm self+charged}^{ee} = N_e C_{\rm square} + \frac{N_e(N_e-1)}{2} C_{ij}\]
\[E_{\rm self+charged}^{e\text{-ion}} = - N_e \sum_{I=1}^{N_{ion}} Z_I C_{ij}\]
\[E_{\rm self+charged}^{\text{ion-ion}} = \sum_{I=1}^{N_{ion}} Z_I^2 C_{\rm square} + \sum_{I<J}^{N_{ion}} Z_I Z_J C_{ij}\]

We also compute contributions from a single electron, to separate the Ewald sum by electron.

\[E_{\rm self+charged}^{\rm single} = C_{\rm square} + \frac{N_e-1}{2} C_{ij} - \sum_{I=1}^{N_{ion}} Z_I C_{ij}\]
\[E_{\rm self+charged}^{\text{single-test}} = C_{\rm square} - \sum_{I=1}^{N_{ion}} Z_I C_{ij}\]
set_lattice_displacements(nlatvec)[source]

Generates list of lattice-vector displacements to add together for real-space sum

Parameters:nlatvec (int) – sum goes from -nlatvec to nlatvec in each lattice direction.
set_up_reciprocal_ewald_sum(ewald_gmax)[source]

Determine parameters for Ewald sums.

\(\alpha\) determines the partitioning of the real and reciprocal-space parts.

We define a weight gweight for the part of the reciprocal-space sums that doesn’t depend on the coordinates:

\[W_G = \frac{4\pi}{V |\vec{G}|^2} e^{- \frac{|\vec{G}|^2}{ 4\alpha^2}}\]
Parameters:ewald_gmax (int) – max number of reciprocal lattice vectors to check away from 0

1-RDM

Evaluate the OBDM for a wave function object.

class pyqmc.obdm.OBDMAccumulator(mol, orb_coeff, nsweeps=5, tstep=0.5, warmup=100, naux=None, spin=None, electrons=None, kpts=None)[source]

Return the obdm as an array with indices rho[spin][i][j] = <c_{spin,i}c^+_{spin,j}>

\[\rho^\sigma_{ij} = \langle c_{\sigma, i} c^\dagger_{\sigma, j} \rangle\]

We are measuring the amplitude of moving one electron (e.g. the first one) from orbital \(\phi_i\) to orbital \(\phi_j\) (Eq (7) of DOI:10.1063/1.4793531)

\[\rho_{i,k} = \int dR dr' \Psi^*(R') \phi_k(r') \phi_i^*(r) \Psi(R)\]

(The complex conjugate is on the wrong orbital in Eq (7) in the paper.) Sampling \(R\) from \(|\Psi(R)^2|\) and \(r'\) from \(f(r) = \sum_i |\phi(r)|^2\)

\[\rho_{i,k} = \int dR dr' \frac{\Psi^*(R')}{\Psi^*(R)} \left[\Psi^*(R) \Psi(R)\right] \frac{\phi_k(r') \phi_i^*(r)}{f(r)} \left[f(r)\right]\]

The distributions (in square brackets) are accounted for by the Monte Carlo integration

\[\rho_{i,k} = \left\langle \frac{\Psi^*(R')}{\Psi^*(R)} \frac{\phi_k(r') \phi_i^*(r)}{f(r)} \right\rangle\]

Eq (9) in the paper is the complex conjugate of this

\[\rho_{i,k}^* = \left\langle \frac{\Psi(R')}{\Psi(R)} \frac{\phi_k^*(r') \phi_i(r)}{f(r)} \right\rangle\]

Args:

mol (Mole): PySCF Mole object.

configs (array): electron positions.

wf (pyqmc wave function object): wave function to evaluate on.

orb_coeff (array): coefficients with size (nbasis,norb) relating mol basis to basis
of 1-RDM desired.
tstep (float): width of the Gaussian to update a walker position for the
extra coordinate.

spin: 0 or 1 for up or down. Defaults to all electrons.

pyqmc.obdm.sample_onebody(configs, orbitals, spin, nsamples=1, tstep=0.5)[source]

For a set of orbitals defined by orb_coeff, return samples from \(f(r) = \sum_i \phi_i(r)^2\).

2-RDM

Evaluate the TBDM for a wave function object.

class pyqmc.tbdm.TBDMAccumulator(mol, orb_coeff, spin, nsweeps=4, tstep=0.5, warmup=200, naux=None, ijkl=None, kpts=None)[source]

Returns one spin sector of the tbdm[s1,s2] as an array (norb_s1,norb_s1,norb_s2,norb_s2) with indices (using pySCF’s convention): tbdm[s1,s2][i,j,k,l] = < c^+_{s1,i} c^+_{s2,k} c_{s2,l} c_{s1,j} > = phi*_{s1,j} phi*_{s2,l} phi_{s2,k} phi_{s1,i}.

We use pySCF’s index convention (while Eq. 10 in DOI:10.1063/1.4793531 uses QWalk’s) QWalk -> tbdm[s1,s2,i,j,k,l] = < c^+_{s1,i} c^+_{s2,j} c_{s2,l} c_{s1,k} > = phi*_{s1,k} phi*_{s2,l} phi_{s2,j} phi_{s1,i} pySCF -> tbdm[s1,s2,i,j,k,l] = < c^+_{s1,i} c^+_{s2,k} c_{s2,l} c_{s1,j} > = phi*_{s1,j} phi*_{s2,l} phi_{s2,k} phi_{s1,i}

Args:

mol (Mole): PySCF Mole object.

orb_coeff (array): coefficients with size (2,nbasis,norb) relating mol basis to basis
of the 2-RDM.

nsweeps (int): number of sweeps over all the electron pairs.

tstep (float): width of the Gaussian to update a walker position for each extra coordinate.

warmup (int): number of warmup steps for single-particle local orbital sampling.

naux (int): number of auxiliary configurations for extra moves sampling the local
orbitals.
spin (tuple): size 2 indicates spin sector to be computed, either (0,0), (0,1), (1,0) or (1,1)
for up-up, up-down, down-up or down-down.

ijkl (array): contains M tbdm matrix elements to calculate with dim (M,4).

get_configurations(nconf)[source]

Obtain a sequence of auxilliary configurations. This function returns one auxilliary configuration for each nconf. Changes internal state: self._aux_configs is updated to the last sampled location.

This will resample the auxilliary configurations to match the number of walkers.

returns a dictionary with the following elements, separated by spin:
assignments: [nsweeps, nconf]: assignment of configurations for each sweep to an auxilliary walker. orbs: [nsweeps, conf, norb]: orbital values configs: [nsweeps] Configuration object with nconf configurations of 1 electron acceptance: [nsweeps, naux] acceptance probability for each auxilliary walker
pyqmc.tbdm.normalize_tbdm(tbdm, norm_a, norm_b)[source]

Returns tbdm by taking the ratio of the averages in Eq. (10) of DOI:10.1063/1.4793531.