qmctorch.wavefunction package

Subpackages

Submodules

Module contents

class qmctorch.wavefunction.WaveFunction(*args: Any, **kwargs: Any)[source]

Bases: Module

forward(x)[source]

Compute the value of the wave function. for a multiple conformation of the electrons

Parameters:
  • parameters – variational param of the wf

  • pos – position of the electrons

Returns: values of psi

electronic_potential(pos)[source]

Computes the electron-electron term

Parameters:

x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the electon-electron energy at each sampling points

Return type:

torch.tensor

nuclear_potential(pos)[source]

Computes the electron-nuclear term

Parameters:

x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the electon-nuclear energy at each sampling points

Return type:

torch.tensor

nuclear_repulsion()[source]

Computes the nuclear-nuclear repulsion term

Returns:

values of the nuclear-nuclear energy at each sampling points

Return type:

torch.tensor

gradients_autograd(pos, pdf=False)[source]

Computes the gradients of the wavefunction (or density) w.r.t the values of the pos.

Parameters:
  • pos (torch.tensor) – positions of the walkers

  • pdf (bool, optional) – if true compute the grads of the density

Returns:

values of the gradients

Return type:

torch.tensor

kinetic_energy_autograd(pos)[source]

Compute the kinetic energy through the 2nd derivative w.r.t the value of the pos.

Parameters:

pos (torch.tensor) – positions of the walkers

Returns:

values of nabla^2 * Psi

local_energy(pos)[source]

Computes the local energy

\[E = K(R) + V_{ee}(R) + V_{en}(R) + V_{nn}\]
Parameters:

pos (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the local enrgies at each sampling points

Return type:

[torch.tensor]

Examples::
>>> mol = Molecule('h2.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrow(mol, configs='cas(2,2)')
>>> pos = torch.rand(500,6)
>>> vals = wf.local_energy(pos)

Note

by default kinetic_energy refers to kinetic_energy_autograd users can overwrite it to poit to any other methods see kinetic_energy_jacobi in wf_orbital

energy(pos)[source]

Total energy for the sampling points.

variance(pos)[source]

Variance of the energy at the sampling points.

sampling_error(eloc)[source]

Compute the statistical uncertainty. Assuming the samples are uncorrelated.

pdf(pos, return_grad=False)[source]

density of the wave function.

get_number_parameters()[source]

Computes the total number of parameters.

load(filename, group='wf_opt', model='best')[source]

Load trained parameters

Parameters:
  • filename (str) – hdf5 filename

  • group (str, optional) – group in the hdf5 file where the model is stored. Defaults to ‘wf_opt’.

  • model (str, optional) – ‘best’ or ‘ last’. Defaults to ‘best’.

class qmctorch.wavefunction.SlaterJastrow(*args: Any, **kwargs: Any)[source]

Bases: SlaterJastrowBase

Slater Jastrow wave function with electron-electron Jastrow factor

\[\Psi(R_{at}, r) = J(r)\sum_n c_n D^\uparrow_n(r^\uparrow)D^\downarrow_n(r^\downarrow)\]

with

\[J(r) = \exp\left( K_{ee}(r) \right)\]

with K, a kernel function depending only on the electron-eletron distances

Parameters:
  • mol (Molecule) – a QMCTorch molecule object

  • configs (str, optional) – defines the CI configurations to be used. Defaults to ‘ground_state’. - ground_state : only the ground state determinant in the wave function - single(n,m) : only single excitation with n electrons and m orbitals - single_double(n,m) : single and double excitation with n electrons and m orbitals - cas(n, m) : all possible configuration using n eletrons and m orbitals

  • kinetic (str, optional) – method to compute the kinetic energy. Defaults to ‘jacobi’. - jacobi : use the Jacobi formula to compute the kinetic energy - auto : use automatic differentiation to compute the kinetic energy

  • jastrow_kernel (JastrowKernelBase, optional) – Class that computes the jastrow kernels

  • jastrow_kernel_kwargs (dict, optional) – keyword arguments for the jastrow kernel contructor

  • cuda (bool, optional) – turns GPU ON/OFF Defaults to Fals e.

  • include_all_mo (bool, optional) – include either all molecular orbitals or only the ones that are popualted in the configs. Defaults to False

Examples::
>>> from qmctorch.scf import Molecule
>>> from qmctorch.wavefunction import SlaterJastrow
>>> mol = Molecule('h2o.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrow(mol, configs='cas(2,2)')
forward(x, ao=None)[source]

computes the value of the wave function for the sampling points

\[\Psi(R) = J(R) \sum_{n} c_n D^{u}_n(r^u) \times D^{d}_n(r^d)\]
Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • ao (torch.tensor, optional) – values of the atomic orbitals (Nbatch, Nelec, Nao)

Returns:

values of the wave functions at each sampling point (Nbatch, 1)

Return type:

torch.tensor

Examples::
>>> mol = Molecule('h2.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrow(mol, configs='cas(2,2)')
>>> pos = torch.rand(500,6)
>>> vals = wf(pos)
ao2mo(ao)[source]

Get the values of the MO from the values of AO.

pos2mo(x, derivative=0)[source]

Get the values of MOs

Parameters:
  • [nbatch (x {torch.tensor} -- positions of the electrons) –

  • nelec*ndim]

Keyword Arguments:

(default (derivative {int} -- order of the derivative) – {0})

Returns:

torch.tensor – MO matrix [nbatch, nelec, nmo]

kinetic_energy_jacobi(x, **kwargs)[source]

Compute the value of the kinetic enery using the Jacobi Formula. C. Filippi, Simple Formalism for Efficient Derivatives .

\[\frac{\Delta \Psi(R)}{\Psi(R)} = \Psi(R)^{-1} \sum_n c_n (\frac{\Delta D_n^u}{D_n^u} + \frac{\Delta D_n^d}{D_n^d}) D_n^u D_n^d\]

We compute the laplacian of the determinants through the Jacobi formula

\[\frac{\Delta \det(A)}{\det(A)} = Tr(A^{-1} \Delta A)\]

Here \(A = J(R) \phi\) and therefore :

\[\Delta A = (\Delta J) D + 2 \nabla J \nabla D + (\Delta D) J\]
Parameters:

x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the kinetic energy at each sampling points

Return type:

torch.tensor

gradients_jacobi(x, sum_grad=False, pdf=False)[source]

Compute the gradients of the wave function (or density) using the Jacobi Formula C. Filippi, Simple Formalism for Efficient Derivatives.

\[\frac{K(R)}{\Psi(R)} = Tr(A^{-1} B_{grad})\]

The gradients of the wave function

\[\Psi(R) = J(R) \sum_n c_n D^{u}_n D^{d}_n = J(R) \Sigma\]

are computed following

\[\nabla \Psi(R) = \left( \nabla J(R) \right) \Sigma + J(R) \left(\nabla \Sigma \right)\]

with

\[\nabla \Sigma = \sum_n c_n (\frac{\nabla D^u_n}{D^u_n} + \frac{\nabla D^d_n}{D^d_n}) D^u_n D^d_n\]

that we compute with the Jacobi formula as:

\[\nabla \Sigma = \sum_n c_n (Tr( (D^u_n)^{-1} \nabla D^u_n) + Tr( (D^d_n)^{-1} \nabla D^d_n)) D^u_n D^d_n\]
Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • pdf (bool, optional) – if true compute the grads of the density

Returns:

values of the gradients wrt the walker pos at each sampling points

Return type:

torch.tensor

get_kinetic_operator(x, ao, dao, d2ao, mo)[source]

Compute the Bkin matrix

Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • mo (torch.tensor, optional) – precomputed values of the MOs

Returns:

matrix of the kinetic operator

Return type:

torch.tensor

class qmctorch.wavefunction.SlaterManyBodyJastrow(*args: Any, **kwargs: Any)[source]

Bases: SlaterJastrow

Slater Jastrow wave function with many body Jastrow factor

\[\Psi(R_{at}, r) = J(r)\sum_n c_n D^\uparrow_n(r^\uparrow)D^\downarrow_n(r^\downarrow)\]

with

\[J(r) = \exp\left( K_{ee}(r) + K_{en}(R_{at},r) + K_{een}(R_{at}, r) \right)\]

with the different kernels representing electron-electron, electron-nuclei and electron-electron-nuclei terms

Parameters:
  • mol (Molecule) – a QMCTorch molecule object

  • configs (str, optional) – defines the CI configurations to be used. Defaults to ‘ground_state’. - ground_state : only the ground state determinant in the wave function - single(n,m) : only single excitation with n electrons and m orbitals - single_double(n,m) : single and double excitation with n electrons and m orbitals - cas(n, m) : all possible configuration using n eletrons and m orbitals

  • kinetic (str, optional) – method to compute the kinetic energy. Defaults to ‘jacobi’. - jacobi : use the Jacobi formula to compute the kinetic energy - auto : use automatic differentiation to compute the kinetic energy

  • jastrow_kernel (dict, optional) – different Jastrow kernels for the different terms. By default only electron-electron and electron-nuclei terms are used

  • jastrow_kernel_kwargs (dict, optional) – keyword arguments for the jastrow kernels contructor

  • cuda (bool, optional) – turns GPU ON/OFF Defaults to False.

  • include_all_mo (bool, optional) – include either all molecular orbitals or only the ones that are popualted in the configs. Defaults to False

Examples::
>>> from qmctorch.scf import Molecule
>>> from qmctorch.wavefunction import SlaterManyBodyJastrow
>>> mol = Molecule('h2o.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterManyBodyJastrow(mol, configs='cas(2,2)')
class qmctorch.wavefunction.SlaterJastrowBackFlow(*args: Any, **kwargs: Any)[source]

Bases: SlaterJastrowBase

Slater Jastrow wave function with electron-electron Jastrow factor and backflow

\[\Psi(R_{at}, r) = J(r)\sum_n c_n D^\uparrow_n(q^\uparrow)D^\downarrow_n(q^\downarrow)\]

with

\[J(r) = \exp\left( K_{ee}(r) \right)\]

with K, a kernel function depending only on the electron-eletron distances, and

\[q(r_i) = r_i + \sum){j\neq i} K_{BF}(r_{ij})(r_i-r_j)\]

is a backflow transformation defined by the kernel K_{BF}. Note that different transformation can be used for different orbital via the orbital_dependent_backflow option.

Args: :param mol: a QMCTorch molecule object :type mol: Molecule :param configs: defines the CI configurations to be used. Defaults to ‘ground_state’.

  • ground_state : only the ground state determinant in the wave function

  • single(n,m) : only single excitation with n electrons and m orbitals

  • single_double(n,m) : single and double excitation with n electrons and m orbitals

  • cas(n, m) : all possible configuration using n eletrons and m orbitals

Parameters:
  • kinetic (str, optional) – method to compute the kinetic energy. Defaults to ‘jacobi’. - jacobi : use the Jacobi formula to compute the kinetic energy - auto : use automatic differentiation to compute the kinetic energy

  • jastrow_kernel (JastrowKernelBase, optional) – Class that computes the jastrow kernels

  • jastrow_kernel_kwargs (dict, optional) – keyword arguments for the jastrow kernel contructor

  • backflow_kernel (BackFlowKernelBase, optional) – kernel function of the backflow transformation. - By default an inverse kernel K(r_{ij}) = w/r_{ij} is used

  • backflow_kernel_kwargs (dict, optional) – keyword arguments for the backflow kernel contructor

  • orbital_dependent_backflow (bool, optional) – every orbital has a different transformation if True. Default to False

  • cuda (bool, optional) – turns GPU ON/OFF Defaults to False.

  • include_all_mo (bool, optional) – include either all molecular orbitals or only the ones that are popualted in the configs. Defaults to False

Examples::
>>> from qmctorch.scf import Molecule
>>> from qmctorch.wavefunction import SlaterJastrowBackFlow
>>> mol = Molecule('h2o.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrowBackFlow(mol, configs='cas(2,2)')
forward(x, ao=None)[source]

computes the value of the wave function for the sampling points

\[J(R) \Psi(R) = J(R) \sum_{n} c_n D^{u}_n(r^u) \times D^{d}_n(r^d)\]
Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • ao (torch.tensor, optional) – values of the atomic orbitals (Nbatch, Nelec, Nao)

Returns:

values of the wave functions at each sampling point (Nbatch, 1)

Return type:

torch.tensor

Examples::
>>> mol = Molecule('h2.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrow(mol, configs='cas(2,2)')
>>> pos = torch.rand(500,6)
>>> vals = wf(pos)
ao2mo(ao)[source]

transforms AO values in to MO values.

pos2mo(x, derivative=0, sum_grad=True)[source]

Compute the MO vals from the pos

Parameters:
  • x ([type]) – [description]

  • derivative (int, optional) – [description]. Defaults to 0.

  • sum_grad (bool, optional) – [description]. Defaults to True.

Returns:

[description]

Return type:

[type]

kinetic_energy_jacobi(x, **kwargs)[source]

Compute the value of the kinetic enery using the Jacobi Formula.

\[\frac{\Delta (J(R) \Psi(R))}{ J(R) \Psi(R)} = \frac{\Delta J(R)}{J(R)} + 2 \frac{\nabla J(R)}{J(R)} \frac{\nabla \Psi(R)}{\Psi(R)} + \frac{\Delta \Psi(R)}{\Psi(R)}\]

The lapacian of the determinental part is computed via

\[\Delta_i \Psi(R) \sum_n c_n ( \frac{\Delta_i D_n^{u}}{D_n^{u}} + \frac{\Delta_i D_n^{d}}{D_n^{d}} + 2 \frac{\nabla_i D_n^{u}}{D_n^{u}} \frac{\nabla_i D_n^{d}}{D_n^{d}} ) D_n^{u} D_n^{d}\]

Since the backflow orbitals are multi-electronic the laplacian of the determinants are obtained

\[\frac{\Delta det(A)}{det(A)} = Tr(A^{-1} \Delta A) + Tr(A^{-1} \nabla A) Tr(A^{-1} \nabla A) + Tr( (A^{-1} \nabla A) (A^{-1} \nabla A ))\]
Parameters:

x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the kinetic energy at each sampling points

Return type:

torch.tensor

gradients_jacobi(x, sum_grad=True)[source]

Computes the gradients of the wf using Jacobi’s Formula

Parameters:

x ([type]) – [description]

class qmctorch.wavefunction.SlaterOrbitalDependentJastrow(*args: Any, **kwargs: Any)[source]

Bases: SlaterJastrowBase

Slater Jastrow Wave function with an orbital dependent Electron-Electron Jastrow Factor

\[\Psi(R_{at}, r) = \sum_n c_n D^\uparrow_n(r^\uparrow)D^\downarrow_n(r^\downarrow)\]

where each molecular orbital of the determinants is multiplied with a different electron-electron Jastrow

\[\phi_i(r) \rightarrow J_i(r) \phi_i(r)\]
Parameters:
  • mol (Molecule) – a QMCTorch molecule object

  • configs (str, optional) – defines the CI configurations to be used. Defaults to ‘ground_state’. - ground_state : only the ground state determinant in the wave function - single(n,m) : only single excitation with n electrons and m orbitals - single_double(n,m) : single and double excitation with n electrons and m orbitals - cas(n, m) : all possible configuration using n eletrons and m orbitals

  • kinetic (str, optional) – method to compute the kinetic energy. Defaults to ‘jacobi’. - jacobi : use the Jacobi formula to compute the kinetic energy - auto : use automatic differentiation to compute the kinetic energy

  • jastrow_kernel (JastrowKernelBase, optional) – Class that computes the jastrow kernels

  • jastrow_kernel_kwargs (dict, optional) – keyword arguments for the jastrow kernel contructor

  • cuda (bool, optional) – turns GPU ON/OFF Defaults to False.

  • include_all_mo (bool, optional) – include either all molecular orbitals or only the ones that are popualted in the configs. Defaults to False

Examples::
>>> from qmctorch.scf import Molecule
>>> from qmctorch.wavefunction import SlaterOrbitalDependentJastrow
>>> mol = Molecule('h2o.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterOrbitalDependentJastrow(mol, configs='cas(2,2)')
ordered_jastrow(pos, derivative=0, sum_grad=True)[source]

Returns the value of the jastrow with the correct dimensions

Parameters:
  • pos (torch.tensor) – Positions of the electrons Size : Nbatch, Nelec x Ndim

  • derivative (int, optional) – order of the derivative (0,1,2,). Defaults to 0.

  • sum_grad (bool, optional) – Return the sum_grad (i.e. the sum of the derivatives) or the individual terms. Defaults to True. False only for derivative=1

Returns:

value of the jastrow parameter for all confs

Nbatch, Nelec, Nmo (sum_grad = True) Nbatch, Nelec, Nmo, Ndim (sum_grad = False)

Return type:

torch.tensor

forward(x, ao=None)[source]

computes the value of the wave function for the sampling points

\[\Psi(R) = \sum_{n} c_n D^{u}_n(r^u) \times D^{d}_n(r^d)\]
Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • ao (torch.tensor, optional) – values of the atomic orbitals (Nbatch, Nelec, Nao)

Returns:

values of the wave functions at each sampling point (Nbatch, 1)

Return type:

torch.tensor

Examples::
>>> mol = Molecule('h2.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrow(mol, configs='cas(2,2)')
>>> pos = torch.rand(500,6)
>>> vals = wf(pos)
ao2mo(ao)[source]

Get the values of the MO from the values of AO.

ao2cmo(ao, jastrow)[source]
pos2mo(x, derivative=0, sum_grad=True)[source]

Compute the uncorrelated MOs from the positions.

pos2cmo(x, derivative=0, sum_grad=True)[source]

Get the values of correlated MOs

Parameters:
  • [nbatch (x {torch.tensor} -- positions of the electrons) –

  • nelec*ndim]

Returns:

torch.tensor – MO matrix [nbatch, nelec, nmo]

kinetic_energy_jacobi(x, **kwargs)[source]

Compute the value of the kinetic enery using the Jacobi Formula. C. Filippi, Simple Formalism for Efficient Derivatives .

\[\frac{K(R)}{\Psi(R)} = Tr(A^{-1} B_{kin})\]
Parameters:

x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the kinetic energy at each sampling points

Return type:

torch.tensor

gradients_jacobi(x, sum_grad=True, pdf=False)[source]

Computes the gradients of the wf using Jacobi’s Formula

Parameters:

x ([type]) – [description]

get_hessian_operator(x)[source]

Compute the Bkin matrix

Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • mo (torch.tensor, optional) – precomputed values of the MOs

Returns:

matrix of the kinetic operator

Return type:

torch.tensor

get_gradient_operator(x)[source]

Compute the gradient operator

Parameters:
  • x ([type]) – [description]

  • ao ([type]) – [description]

  • dao ([type]) – [description]

class qmctorch.wavefunction.SlaterManyBodyJastrowBackflow(*args: Any, **kwargs: Any)[source]

Bases: SlaterJastrow

Slater Jastrow wave function with many-body Jastrow factor and backflow

\[\Psi(R_{at}, r) = J(R_{at}, r)\sum_n c_n D^\uparrow_n(q^\uparrow)D^\downarrow_n(q^\downarrow)\]

with

\[J(r) = \exp\left( K_{ee}(r) + K_{en}(R_{at},r) + K_{een}(R_{at}, r) \right)\]

with the different kernels representing electron-electron, electron-nuclei and electron-electron-nuclei terms and

\[q(r_i) = r_i + \sum){j\neq i} K_{BF}(r_{ij})(r_i-r_j)\]

is a backflow transformation defined by the kernel K_{BF}. Note that different transformation can be used for different orbital via the orbital_dependent_backflow option.

Args: :param mol: a QMCTorch molecule object :type mol: Molecule :param configs: defines the CI configurations to be used. Defaults to ‘ground_state’.

  • ground_state : only the ground state determinant in the wave function

  • single(n,m) : only single excitation with n electrons and m orbitals

  • single_double(n,m) : single and double excitation with n electrons and m orbitals

  • cas(n, m) : all possible configuration using n eletrons and m orbitals

Parameters:
  • kinetic (str, optional) – method to compute the kinetic energy. Defaults to ‘jacobi’. - jacobi : use the Jacobi formula to compute the kinetic energy - auto : use automatic differentiation to compute the kinetic energy

  • jastrow_kernel (dict, optional) – different Jastrow kernels for the different terms. By default only electron-electron and electron-nuclei terms are used

  • jastrow_kernel_kwargs (dict, optional) – keyword arguments for the jastrow kernels contructor

  • backflow_kernel (BackFlowKernelBase, optional) – kernel function of the backflow transformation. - By default an inverse kernel K(r_{ij}) = w/r_{ij} is used

  • backflow_kernel_kwargs (dict, optional) – keyword arguments for the backflow kernel contructor

  • orbital_dependent_backflow (bool, optional) – every orbital has a different transformation if True. Default to False

  • cuda (bool, optional) – turns GPU ON/OFF Defaults to False.

  • include_all_mo (bool, optional) – include either all molecular orbitals or only the ones that are popualted in the configs. Defaults to False

Examples::
>>> from qmctorch.scf import Molecule
>>> from qmctorch.wavefunction import SlaterManyBodyJastrowBackflow
>>> mol = Molecule('h2o.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterManyBodyJastrowBackflow(mol, configs='cas(2,2)')
forward(x, ao=None)[source]

computes the value of the wave function for the sampling points

\[J(R) \Psi(R) = J(R) \sum_{n} c_n D^{u}_n(r^u) \times D^{d}_n(r^d)\]
Parameters:
  • x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

  • ao (torch.tensor, optional) – values of the atomic orbitals (Nbatch, Nelec, Nao)

Returns:

values of the wave functions at each sampling point (Nbatch, 1)

Return type:

torch.tensor

Examples::
>>> mol = Molecule('h2.xyz', calculator='adf', basis = 'dzp')
>>> wf = SlaterJastrow(mol, configs='cas(2,2)')
>>> pos = torch.rand(500,6)
>>> vals = wf(pos)
ao2mo(ao)[source]

transforms AO values in to MO values.

pos2mo(x, derivative=0, sum_grad=True)[source]

Compute the MO vals from the pos

Parameters:
  • x ([type]) – [description]

  • derivative (int, optional) – [description]. Defaults to 0.

  • sum_grad (bool, optional) – [description]. Defaults to True.

Returns:

[description]

Return type:

[type]

kinetic_energy_jacobi(x, **kwargs)[source]

Compute the value of the kinetic enery using the Jacobi Formula.

\[\begin{split}\\frac{\Delta (J(R) \Psi(R))}{ J(R) \Psi(R)} = \\frac{\\Delta J(R)}{J(R} + 2 \\frac{\\nabla J(R)}{J(R)} \\frac{\\nabla \\Psi(R)}{\\Psi(R)} + \\frac{\\Delta \\Psi(R)}{\\Psi(R)}\end{split}\]

The lapacian of the determinental part is computed via

\[\begin{split}\\Delta_i \\Psi(R) \\sum_n c_n ( \\frac{\\Delta_i D_n^{u}}{D_n^{u}} + \\frac{\\Delta_i D_n^{d}}{D_n^{d}} + 2 \\frac{\\nabla_i D_n^{u}}{D_n^{u}} \\frac{\\nabla_i D_n^{d}}{D_n^{d}} ) D_n^{u} D_n^{d}\end{split}\]

Since the backflow orbitals are multi-electronic the laplacian of the determinants are obtained

\[\begin{split}\\frac{\\Delta det(A)}{det(A)} = Tr(A^{-1} \\Delta A) + Tr(A^{-1} \\nabla A) Tr(A^{-1} \\nabla A) + Tr( (A^{-1} \\nabla A) (A^{-1} \\nabla A ))\end{split}\]
Parameters:

x (torch.tensor) – sampling points (Nbatch, 3*Nelec)

Returns:

values of the kinetic energy at each sampling points

Return type:

torch.tensor

gradients_jacobi(x, sum_grad=True)[source]

Computes the gradients of the wf using Jacobi’s Formula

Parameters:

x ([type]) – [description]