1
0
mirror of https://github.com/TREX-CoE/qmckl.git synced 2024-06-01 02:45:43 +02:00
qmckl/org/qmckl_electron.org
Aurélien Delval d0eb207404
Integration of Verificarlo CI tests (#1)
* comment

* Update distance test code

The distance test has been updated to the latest version, with a first attempt at using vfc_probes inside it

* Functional implementation of vfc_probes in the distance tests

This commit has the first functional vfc_ci tests. Verificarlo tests
should be written over the existing tests, and they can be enabled with
the following configure command:

QMCKL_DEVEL=1 ./configure --prefix=$PWD/_install --enable-maintainer-mode --enable-vfc_ci CC="verificarlo-f -Mpreprocess -D VFC_CI" FC="verificarlo-f -Mpreprocess -D VFC_CI" --host=x86_64

The --enable-vfc_ci flag will trigger the linking of the vfc_ci
library. Moreover, as of now, the "-Mpreprocess" and "-D VFC_CI" flags
have to be specified directly here. There is probably an appropriate
macro to place those flags into but I couldn't find it yet, and could
only manage to build the tests this way.

When the VFC_CI preprocessor is defined, somme additional code to
register and export the test probes can be executed (see
qmckl_distance.org).

As of now, the tests are built as normal, even though they are expected
to fail :

make all
make check

From there, the test_qmckl_distance (and potentially the others)
executable can be called at will. This will typically be done
automatically by vfc_ci, but one could manually execute the executable
by defining the following env variables :

VFC_PROBES_OUTPUT="test.csv" VFC_BACKENDS="libinterflop_ieee.so"

depending on the export file and the Verificarlo backend to be used.

The next steps will be to define more tests such as this one, and to
integrate them into a Verificarlo CI workflow (by writing a
vfc_tests_config.json file and using the automatic CI setup
command).

* Error in FOrtran interface fixed

* Added missing Fortran interfaces

* Modify distance test and install process integration

All probes are now ignored using only the preprocessor (instead
of checking for a facultative argument) in the distance test.
Moreover,preprocessing can now be enabled correctly using FCFLAGS
(the issue seemed to come from the order of the arguments passed
to gfortran/verificarlo-f with the preprocessor arg having to come
first).

* Add vfc_probes to AO tests

vfc_probes have been added to qmckl_ao.org in the same way as
qmckl_distance.org, which means that it can be enabled or disabled at
compile time using the --enable-vfc_ci option.

qmckl_distance.org has been slightly modified with a better indentation,
and configure.ac now adds the "-D VFC_CI" flag to CFLAGS when vfc_ci is
enabled.

Co-authored-by: Anthony Scemama <scemama@irsamc.ups-tlse.fr>
2021-07-07 13:42:42 +02:00

81 KiB

Electrons

In conventional QMC simulations, up-spin and down-spin electrons are different. The electron data structure contains the number of up-spin and down-spin electrons, and the electron coordinates.

Context

The following data stored in the context:

uninitialized int32_t Keeps bit set for uninitialized data
num int64_t Total number of electrons
up_num int64_t Number of up-spin electrons
down_num int64_t Number of down-spin electrons
walk_num int64_t Number of walkers
rescale_factor_kappa_ee double The distance scaling factor
rescale_factor_kappa_en double The distance scaling factor
provided bool If true, electron is valid
coord_new double[walk_num][3][num] New set of electron coordinates
coord_old double[walk_num][3][num] Old set of electron coordinates
coord_new_date uint64_t Last modification date of the coordinates
ee_distance double[walk_num][num][num] Electron-electron distances
ee_distance_date uint64_t Last modification date of the electron-electron distances
en_distance double[walk_num][nucl_num][num] Electron-nucleus distances
en_distance_date uint64_t Last modification date of the electron-electron distances
ee_distance_rescaled double[walk_num][num][num] Electron-electron rescaled distances
ee_distance_rescaled_date uint64_t Last modification date of the electron-electron distances
ee_distance_rescaled_deriv_e double[walk_num][4][num][num] Electron-electron rescaled distances derivatives
ee_distance_rescaled_deriv_e_date uint64_t Last modification date of the electron-electron distance derivatives
en_distance_rescaled double[walk_num][nucl_num][num] Electron-nucleus distances
en_distance_rescaled_date uint64_t Last modification date of the electron-electron distances
en_distance_rescaled_deriv_e double[walk_num][4][num][num] Electron-electron rescaled distances derivatives
en_distance_rescaled_deriv_e_date uint64_t Last modification date of the electron-electron distance derivatives

Data structure

typedef struct qmckl_electron_struct {
int64_t   num;
int64_t   up_num;
int64_t   down_num;
int64_t   walk_num;
double    rescale_factor_kappa_ee;
double    rescale_factor_kappa_en;
int64_t   coord_new_date;
int64_t   ee_distance_date;
int64_t   en_distance_date;
int64_t   ee_distance_rescaled_date;
int64_t   ee_distance_rescaled_deriv_e_date;
int64_t   en_distance_rescaled_date;
int64_t   en_distance_rescaled_deriv_e_date;
double*   coord_new;
double*   coord_old;
double*   ee_distance;
double*   en_distance;
double*   ee_distance_rescaled;
double*   ee_distance_rescaled_deriv_e;
double*   en_distance_rescaled;
double*   en_distance_rescaled_deriv_e;
int32_t   uninitialized;
bool      provided;
} qmckl_electron_struct;

The uninitialized integer contains one bit set to one for each initialization function which has not been called. It becomes equal to zero after all initialization functions have been called. The struct is then initialized and provided == true. Some values are initialized by default, and are not concerned by this mechanism.

qmckl_exit_code qmckl_init_electron(qmckl_context context);
qmckl_exit_code qmckl_init_electron(qmckl_context context) {

if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
 return false;
}

qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
assert (ctx != NULL);

ctx->electron.uninitialized = (1 << 2) - 1;

/* Default values */
ctx->electron.rescale_factor_kappa_ee = 1.0;
ctx->electron.rescale_factor_kappa_en = 1.0;

return QMCKL_SUCCESS;
}
bool qmckl_electron_provided (const qmckl_context context);

Access functions

Access functions return QMCKL_SUCCESS when the data has been successfully retrieved. It returnes QMCKL_INVALID_CONTEXT when the context is not a valid context, and QMCKL_NOT_PROVIDED when the data has not been provided. If the function returns successfully, the variable pointed by the pointer given in argument contains the requested data. Otherwise, this variable is untouched.

#+NAME:post

Number of electrons

Number of walkers

A walker is a set of electron coordinates that are arguments of the wave function. walk_num is the number of walkers.

Scaling factors Kappa

Electron coordinates

Returns the current electron coordinates. The pointer is assumed to point on a memory block of size 3 * elec_num * walk_num. The order of the indices is:

Normal Transposed
C [walk_num][elec_num][3] [walk_num][3][elec_num]
Fortran (3,elec_num,walk_num) (elec_num,3,walk_num)

Initialization functions

To set the data relative to the electrons in the context, the following functions need to be called. When the data structure is initialized, the internal coord_new and coord_old arrays are both allocated.

qmckl_exit_code qmckl_set_electron_num      (qmckl_context context, const int64_t up_num, const int64_t down_num);
qmckl_exit_code qmckl_set_electron_walk_num (qmckl_context context, const int64_t walk_num);
qmckl_exit_code qmckl_set_electron_coord    (qmckl_context context, const char transp, const double* coord);

qmckl_exit_code qmckl_set_electron_rescale_factor_ee (qmckl_context context, const double kappa_ee);
qmckl_exit_code qmckl_set_electron_rescale_factor_en (qmckl_context context, const double kappa_en);

#+NAME:pre2

#+NAME:post2

To set the number of electrons, we give the number of up-spin and down-spin electrons to the context and we set the number of walkers.

The following function sets the number of walkers.

Next we set the rescale parameter for the rescaled distance metric.

interface
integer(c_int32_t) function qmckl_set_electron_num(context, alpha, beta) bind(C)
 use, intrinsic :: iso_c_binding
 import
 implicit none

 integer (c_int64_t) , intent(in)  , value :: context                                  
 integer (c_int64_t) , intent(in)  , value :: alpha
 integer (c_int64_t) , intent(in)  , value :: beta 
end function
end interface
interface
integer(c_int32_t) function qmckl_set_electron_walk_num(context, walk_num) bind(C)
 use, intrinsic :: iso_c_binding
 import
 implicit none

 integer (c_int64_t) , intent(in)  , value :: context                                  
 integer (c_int64_t) , intent(in)  , value :: walk_num
end function
end interface

The following function sets the electron coordinates of all the walkers. When this is done, the pointers to the old and new sets of coordinates are swapped, and the new coordinates are overwritten. This can be done only when the data relative to electrons have been set.

Important: changing the electron coordinates increments the date in the context.

interface
integer(c_int32_t) function qmckl_set_electron_coord(context, transp, coord) bind(C)
 use, intrinsic :: iso_c_binding
 import
 implicit none

 integer (c_int64_t) , intent(in)  , value :: context                                  
 character           , intent(in)  , value :: transp 
 double precision    , intent(in)          :: coord(*)
end function
end interface

Test

/* Reference input data */
int64_t walk_num      = chbrclf_walk_num;
int64_t elec_num      = chbrclf_elec_num;
int64_t elec_up_num   = chbrclf_elec_up_num;
int64_t elec_dn_num   = chbrclf_elec_dn_num;
double  rescale_factor_kappa_ee   = 1.0;
double  rescale_factor_kappa_en   = 1.0;
double  nucl_rescale_factor_kappa = 1.0;
double* elec_coord    = &(chbrclf_elec_coord[0][0][0]);

int64_t  nucl_num      = chbrclf_nucl_num;
double*  charge        = chbrclf_charge;
double*  nucl_coord    = &(chbrclf_nucl_coord[0][0]);

/* --- */

qmckl_exit_code rc;

assert(!qmckl_electron_provided(context));

int64_t n;
rc = qmckl_get_electron_num (context, &n);
assert(rc == QMCKL_NOT_PROVIDED);

rc = qmckl_get_electron_up_num (context, &n);
assert(rc == QMCKL_NOT_PROVIDED);

rc = qmckl_get_electron_down_num (context, &n);
assert(rc == QMCKL_NOT_PROVIDED);


rc = qmckl_set_electron_num (context, elec_up_num, elec_dn_num);
assert(rc == QMCKL_SUCCESS);
assert(!qmckl_electron_provided(context));

rc = qmckl_get_electron_up_num (context, &n);
assert(rc == QMCKL_SUCCESS);
assert(n == elec_up_num);

rc = qmckl_get_electron_down_num (context, &n);
assert(rc == QMCKL_SUCCESS);
assert(n == elec_dn_num);

rc = qmckl_get_electron_num (context, &n);
assert(rc == QMCKL_SUCCESS);
assert(n == elec_num);

double k_ee = 0.;
double k_en = 0.;
rc = qmckl_get_electron_rescale_factor_ee (context, &k_ee);
assert(rc == QMCKL_SUCCESS);
assert(k_ee == 1.0);

rc = qmckl_get_electron_rescale_factor_en (context, &k_en);
assert(rc == QMCKL_SUCCESS);
assert(k_en == 1.0);

rc = qmckl_set_electron_rescale_factor_en(context, rescale_factor_kappa_en);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_set_electron_rescale_factor_ee(context, rescale_factor_kappa_ee);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_get_electron_rescale_factor_ee (context, &k_ee);
assert(rc == QMCKL_SUCCESS);
assert(k_ee == rescale_factor_kappa_ee);

rc = qmckl_get_electron_rescale_factor_en (context, &k_en);
assert(rc == QMCKL_SUCCESS);
assert(k_en == rescale_factor_kappa_en);


int64_t w;
rc = qmckl_get_electron_walk_num (context, &w);
assert(rc == QMCKL_NOT_PROVIDED);


rc = qmckl_set_electron_walk_num (context, walk_num);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_get_electron_walk_num (context, &w);
assert(rc == QMCKL_SUCCESS);
assert(w == walk_num);

assert(qmckl_electron_provided(context));

rc = qmckl_set_electron_coord (context, 'N', elec_coord);
assert(rc == QMCKL_SUCCESS);

double elec_coord2[walk_num*3*elec_num];

rc = qmckl_get_electron_coord (context, 'N', elec_coord2);
assert(rc == QMCKL_SUCCESS);
for (int64_t i=0 ; i<3*elec_num ; ++i) {
assert( elec_coord[i] == elec_coord2[i] );
}

Computation

The computed data is stored in the context so that it can be reused by different kernels. To ensure that the data is valid, for each computed data the date of the context is stored when it is computed. To know if some data needs to be recomputed, we check if the date of the dependencies are more recent than the date of the data to compute. If it is the case, then the data is recomputed and the current date is stored.

Electron-electron distances

Get

qmckl_exit_code qmckl_get_electron_ee_distance(qmckl_context context, double* const distance);

Compute

qmckl_context context in Global state
int64_t elec_num in Number of electrons
int64_t walk_num in Number of walkers
double coord[walk_num][3][elec_num] in Electron coordinates
double ee_distance[walk_num][elec_num][elec_num] out Electron-electron distances
integer function qmckl_compute_ee_distance_f(context, elec_num, walk_num, coord, ee_distance) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: elec_num
integer*8             , intent(in)  :: walk_num
double precision      , intent(in)  :: coord(elec_num,3,walk_num)
double precision      , intent(out) :: ee_distance(elec_num,elec_num,walk_num)

integer*8 :: k

info = QMCKL_SUCCESS

if (context == QMCKL_NULL_CONTEXT) then
 info = QMCKL_INVALID_CONTEXT
 return
endif

if (elec_num <= 0) then
 info = QMCKL_INVALID_ARG_2
 return
endif

if (walk_num <= 0) then
 info = QMCKL_INVALID_ARG_3
 return
endif

do k=1,walk_num
 info = qmckl_distance(context, 'T', 'T', elec_num, elec_num, &
      coord(1,1,k), elec_num, &
      coord(1,1,k), elec_num, &
      ee_distance(1,1,k), elec_num)
 if (info /= QMCKL_SUCCESS) then
    exit
 endif
end do

end function qmckl_compute_ee_distance_f

Test

assert(qmckl_electron_provided(context));


double ee_distance[walk_num * elec_num * elec_num];
rc = qmckl_get_electron_ee_distance(context, ee_distance);

// (e1,e2,w)
// (0,0,0) == 0.
assert(ee_distance[0] == 0.);

// (1,0,0) == (0,1,0)
assert(ee_distance[1] == ee_distance[elec_num]);

// value of (1,0,0)
assert(fabs(ee_distance[1]-7.152322512964209) < 1.e-12);

// (0,0,1) == 0.
assert(ee_distance[elec_num*elec_num] == 0.);

// (1,0,1) == (0,1,1)
assert(ee_distance[elec_num*elec_num+1] == ee_distance[elec_num*elec_num+elec_num]);

// value of (1,0,1)
assert(fabs(ee_distance[elec_num*elec_num+1]-6.5517646321055665) < 1.e-12);

Electron-electron rescaled distances

ee_distance_rescaled stores the matrix of the rescaled distances between all pairs of electrons:

\[ C_{ij} = \left( 1 - \exp{-\kappa C_{ij}}\right)/\kappa \]

where \(C_{ij}\) is the matrix of electron-electron distances.

Get

qmckl_exit_code qmckl_get_electron_ee_distance_rescaled(qmckl_context context, double* const distance_rescaled);

Compute

qmckl_context context in Global state
int64_t elec_num in Number of electrons
double rescale_factor_kappa_ee in Factor to rescale ee distances
int64_t walk_num in Number of walkers
double coord[walk_num][3][elec_num] in Electron coordinates
double ee_distance[walk_num][elec_num][elec_num] out Electron-electron rescaled distances
integer function qmckl_compute_ee_distance_rescaled_f(context, elec_num, rescale_factor_kappa_ee, walk_num, &
 coord, ee_distance_rescaled) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: elec_num
double precision      , intent(in)  :: rescale_factor_kappa_ee
integer*8             , intent(in)  :: walk_num
double precision      , intent(in)  :: coord(elec_num,3,walk_num)
double precision      , intent(out) :: ee_distance_rescaled(elec_num,elec_num,walk_num)

integer*8 :: k

info = QMCKL_SUCCESS

if (context == QMCKL_NULL_CONTEXT) then
 info = QMCKL_INVALID_CONTEXT
 return
endif

if (elec_num <= 0) then
 info = QMCKL_INVALID_ARG_2
 return
endif

if (walk_num <= 0) then
 info = QMCKL_INVALID_ARG_3
 return
endif

do k=1,walk_num
 info = qmckl_distance_rescaled(context, 'T', 'T', elec_num, elec_num, &
      coord(1,1,k), elec_num, &
      coord(1,1,k), elec_num, &
      ee_distance_rescaled(1,1,k), elec_num, rescale_factor_kappa_ee)
 if (info /= QMCKL_SUCCESS) then
    exit
 endif
end do

end function qmckl_compute_ee_distance_rescaled_f

Test

assert(qmckl_electron_provided(context));


double ee_distance_rescaled[walk_num * elec_num * elec_num];
rc = qmckl_get_electron_ee_distance_rescaled(context, ee_distance_rescaled);

// (e1,e2,w)
// (0,0,0) == 0.
assert(ee_distance_rescaled[0] == 0.);

// (1,0,0) == (0,1,0)
assert(ee_distance_rescaled[1] == ee_distance_rescaled[elec_num]);

// value of (1,0,0)
assert(fabs(ee_distance_rescaled[1]-0.9992169566605263) < 1.e-12);

// (0,0,1) == 0.
assert(ee_distance_rescaled[elec_num*elec_num] == 0.);

// (1,0,1) == (0,1,1)
assert(ee_distance_rescaled[elec_num*elec_num+1] == ee_distance_rescaled[elec_num*elec_num+elec_num]);

// value of (1,0,1)
assert(fabs(ee_distance_rescaled[elec_num*elec_num+1]-0.9985724058042633) < 1.e-12);

Electron-electron rescaled distance gradients and laplacian with respect to electron coords

The rescaled distances which is given as $R = (1 - \exp{-\kappa r})/\kappa$ needs to be perturbed with respect to the electorn coordinates. This data is stored in the ee_distance_rescaled_deriv_e tensor. The The first three elements of this three index tensor [4][num][num] gives the derivatives in the x, y, and z directions $dx, dy, dz$ and the last index gives the Laplacian $\partial x^2 + \partial y^2 + \partial z^2$.

Get

qmckl_exit_code qmckl_get_electron_ee_distance_rescaled_deriv_e(qmckl_context context, double* const distance_rescaled_deriv_e);

Compute

qmckl_context context in Global state
int64_t elec_num in Number of electrons
double rescale_factor_kappa_ee in Factor to rescale ee distances
int64_t walk_num in Number of walkers
double coord[walk_num][3][elec_num] in Electron coordinates
double ee_distance_deriv_e[walk_num][4][elec_num][elec_num] out Electron-electron rescaled distance derivatives
integer function qmckl_compute_ee_distance_rescaled_deriv_e_f(context, elec_num, rescale_factor_kappa_ee, walk_num, &
 coord, ee_distance_rescaled_deriv_e) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: elec_num
double precision      , intent(in)  :: rescale_factor_kappa_ee
integer*8             , intent(in)  :: walk_num
double precision      , intent(in)  :: coord(elec_num,3,walk_num)
double precision      , intent(out) :: ee_distance_rescaled_deriv_e(4,elec_num,elec_num,walk_num)

integer*8 :: k

info = QMCKL_SUCCESS

if (context == QMCKL_NULL_CONTEXT) then
 info = QMCKL_INVALID_CONTEXT
 return
endif

if (elec_num <= 0) then
 info = QMCKL_INVALID_ARG_2
 return
endif

if (walk_num <= 0) then
 info = QMCKL_INVALID_ARG_3
 return
endif

do k=1,walk_num
 info = qmckl_distance_rescaled_deriv_e(context, 'T', 'T', elec_num, elec_num, &
      coord(1,1,k), elec_num, &
      coord(1,1,k), elec_num, &
      ee_distance_rescaled_deriv_e(1,1,1,k), elec_num, rescale_factor_kappa_ee)
 if (info /= QMCKL_SUCCESS) then
    exit
 endif
end do

end function qmckl_compute_ee_distance_rescaled_deriv_e_f

Test

assert(qmckl_electron_provided(context));


double ee_distance_rescaled_deriv_e[4 * walk_num * elec_num * elec_num];
rc = qmckl_get_electron_ee_distance_rescaled_deriv_e(context, ee_distance_rescaled_deriv_e);

// TODO: Get exact values
//// (e1,e2,w)
//// (0,0,0) == 0.
//assert(ee_distance[0] == 0.);
//
//// (1,0,0) == (0,1,0)
//assert(ee_distance[1] == ee_distance[elec_num]);
//
//// value of (1,0,0)
//assert(fabs(ee_distance[1]-7.152322512964209) < 1.e-12);
//
//// (0,0,1) == 0.
//assert(ee_distance[elec_num*elec_num] == 0.);
//
//// (1,0,1) == (0,1,1)
//assert(ee_distance[elec_num*elec_num+1] == ee_distance[elec_num*elec_num+elec_num]);
//
//// value of (1,0,1)
//assert(fabs(ee_distance[elec_num*elec_num+1]-6.5517646321055665) < 1.e-12);

Electron-nucleus distances

Get

qmckl_exit_code qmckl_get_electron_en_distance(qmckl_context context, double* distance);

Compute

qmckl_context context in Global state
int64_t elec_num in Number of electrons
int64_t nucl_num in Number of nuclei
int64_t walk_num in Number of walkers
double elec_coord[walk_num][3][elec_num] in Electron coordinates
double nucl_coord[3][elec_num] in Nuclear coordinates
double en_distance[walk_num][nucl_num][elec_num] out Electron-nucleus distances
integer function qmckl_compute_en_distance_f(context, elec_num, nucl_num, walk_num, elec_coord, nucl_coord, en_distance) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: elec_num
integer*8             , intent(in)  :: nucl_num
integer*8             , intent(in)  :: walk_num
double precision      , intent(in)  :: elec_coord(elec_num,3,walk_num)
double precision      , intent(in)  :: nucl_coord(nucl_num,3)
double precision      , intent(out) :: en_distance(elec_num,nucl_num,walk_num)

integer*8 :: k

info = QMCKL_SUCCESS

if (context == QMCKL_NULL_CONTEXT) then
 info = QMCKL_INVALID_CONTEXT
 return
endif

if (elec_num <= 0) then
 info = QMCKL_INVALID_ARG_2
 return
endif

if (nucl_num <= 0) then
 info = QMCKL_INVALID_ARG_3
 return
endif

if (walk_num <= 0) then
 info = QMCKL_INVALID_ARG_4
 return
endif

do k=1,walk_num
 info = qmckl_distance(context, 'T', 'T', elec_num, nucl_num, &
      elec_coord(1,1,k), elec_num, &
      nucl_coord, nucl_num, &
      en_distance(1,1,k), elec_num)
 if (info /= QMCKL_SUCCESS) then
    exit
 endif
end do

end function qmckl_compute_en_distance_f

Test

assert(!qmckl_nucleus_provided(context));
assert(qmckl_electron_provided(context));

rc = qmckl_set_nucleus_num (context, nucl_num);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_charge (context, charge);
assert (rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_coord (context, 'T', nucl_coord);
assert (rc == QMCKL_SUCCESS);

assert(qmckl_nucleus_provided(context));

double en_distance[walk_num][nucl_num][elec_num];

rc = qmckl_get_electron_en_distance(context, &(en_distance[0][0][0]));
assert (rc == QMCKL_SUCCESS);

// (e,n,w) in Fortran notation
// (1,1,1)
assert(fabs(en_distance[0][0][0] - 7.546738741619978) < 1.e-12);

// (1,2,1)
assert(fabs(en_distance[0][1][0] - 8.77102435246984) < 1.e-12);

// (2,1,1)
assert(fabs(en_distance[0][0][1] - 3.698922010513608) < 1.e-12);

// (1,1,2)
assert(fabs(en_distance[1][0][0] - 5.824059436060509) < 1.e-12);

// (1,2,2)
assert(fabs(en_distance[1][1][0] - 7.080482110317645) < 1.e-12);

// (2,1,2)
assert(fabs(en_distance[1][0][1] - 3.1804527583077356) < 1.e-12);

Electron-nucleus rescaled distances

en_distance_rescaled stores the matrix of the rescaled distances between electrons and nucleii.

\[ C_{ij} = \left( 1 - \exp{-\kappa C_{ij}}\right)/\kappa \]

where \(C_{ij}\) is the matrix of electron-nucleus distances.

Get

qmckl_exit_code qmckl_get_electron_en_distance_rescaled(qmckl_context context, double* distance_rescaled);

Compute

qmckl_context context in Global state
int64_t elec_num in Number of electrons
int64_t nucl_num in Number of nuclei
double rescale_factor_kappa_en in The factor for rescaled distances
int64_t walk_num in Number of walkers
double elec_coord[walk_num][3][elec_num] in Electron coordinates
double nucl_coord[3][elec_num] in Nuclear coordinates
double en_distance_rescaled_date[walk_num][nucl_num][elec_num] out Electron-nucleus distances
integer function qmckl_compute_en_distance_rescaled_f(context, elec_num, nucl_num, rescale_factor_kappa_en, walk_num, elec_coord, &
 nucl_coord, en_distance_rescaled) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: elec_num
integer*8             , intent(in)  :: nucl_num
double precision      , intent(in)  :: rescale_factor_kappa_en
integer*8             , intent(in)  :: walk_num
double precision      , intent(in)  :: elec_coord(elec_num,3,walk_num)
double precision      , intent(in)  :: nucl_coord(nucl_num,3)
double precision      , intent(out) :: en_distance_rescaled(elec_num,nucl_num,walk_num)

integer*8 :: k

info = QMCKL_SUCCESS

if (context == QMCKL_NULL_CONTEXT) then
 info = QMCKL_INVALID_CONTEXT
 return
endif

if (elec_num <= 0) then
 info = QMCKL_INVALID_ARG_2
 return
endif

if (nucl_num <= 0) then
 info = QMCKL_INVALID_ARG_3
 return
endif

! TODO: comparison with 0
!if (rescale_factor_kappa_en <= 0) then
!   info = QMCKL_INVALID_ARG_4
!   return
!endif

if (walk_num <= 0) then
 info = QMCKL_INVALID_ARG_5
 return
endif

do k=1,walk_num
 info = qmckl_distance_rescaled(context, 'T', 'T', elec_num, nucl_num, &
      elec_coord(1,1,k), elec_num, &
      nucl_coord, nucl_num, &
      en_distance_rescaled(1,1,k), elec_num, rescale_factor_kappa_en)
 if (info /= QMCKL_SUCCESS) then
    exit
 endif
end do

end function qmckl_compute_en_distance_rescaled_f

Test

assert(qmckl_electron_provided(context));

rc = qmckl_set_nucleus_num (context, nucl_num);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_charge (context, charge);
assert (rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_coord (context, 'T', nucl_coord);
assert (rc == QMCKL_SUCCESS);

assert(qmckl_nucleus_provided(context));

double en_distance_rescaled[walk_num][nucl_num][elec_num];

rc = qmckl_get_electron_en_distance_rescaled(context, &(en_distance_rescaled[0][0][0]));

assert (rc == QMCKL_SUCCESS);

// (e,n,w) in Fortran notation
// (1,1,1)
assert(fabs(en_distance_rescaled[0][0][0] - 0.9994721712909764) < 1.e-12);

// (1,2,1)
assert(fabs(en_distance_rescaled[0][1][0] - 0.9998448354439821) < 1.e-12);

// (2,1,1)
assert(fabs(en_distance_rescaled[0][0][1] - 0.9752498074577688) < 1.e-12);

// (1,1,2)
assert(fabs(en_distance_rescaled[1][0][0] - 0.9970444172399963) < 1.e-12);

// (1,2,2)
assert(fabs(en_distance_rescaled[1][1][0] - 0.9991586325813303) < 1.e-12);

// (2,1,2)
assert(fabs(en_distance_rescaled[1][0][1] - 0.9584331688679852) < 1.e-12);

Electron-nucleus rescaled distance gradients and laplacian with respect to electron coords

The rescaled distances which is given as $R = (1 - \exp{-\kappa r})/\kappa$ needs to be perturbed with respect to the nuclear coordinates. This data is stored in the en_distance_rescaled_deriv_e tensor. The The first three elements of this three index tensor [4][num][num] gives the derivatives in the x, y, and z directions $dx, dy, dz$ and the last index gives the Laplacian $\partial x^2 + \partial y^2 + \partial z^2$.

Get

qmckl_exit_code qmckl_get_electron_en_distance_rescaled_deriv_e(qmckl_context context, double* distance_rescaled_deriv_e);

Compute

qmckl_context context in Global state
int64_t elec_num in Number of electrons
int64_t nucl_num in Number of nuclei
double rescale_factor_kappa_en in The factor for rescaled distances
int64_t walk_num in Number of walkers
double elec_coord[walk_num][3][elec_num] in Electron coordinates
double nucl_coord[3][elec_num] in Nuclear coordinates
double en_distance_rescaled_deriv_e_date[walk_num][4][nucl_num][elec_num] out Electron-nucleus distance derivatives
integer function qmckl_compute_en_distance_rescaled_deriv_e_f(context, elec_num, nucl_num, &
 rescale_factor_kappa_en, walk_num, elec_coord, &
 nucl_coord, en_distance_rescaled_deriv_e) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: elec_num
integer*8             , intent(in)  :: nucl_num
double precision      , intent(in)  :: rescale_factor_kappa_en
integer*8             , intent(in)  :: walk_num
double precision      , intent(in)  :: elec_coord(elec_num,3,walk_num)
double precision      , intent(in)  :: nucl_coord(nucl_num,3)
double precision      , intent(out) :: en_distance_rescaled_deriv_e(elec_num,nucl_num,walk_num)

integer*8 :: k

info = QMCKL_SUCCESS

if (context == QMCKL_NULL_CONTEXT) then
 info = QMCKL_INVALID_CONTEXT
 return
endif

if (elec_num <= 0) then
 info = QMCKL_INVALID_ARG_2
 return
endif

if (nucl_num <= 0) then
 info = QMCKL_INVALID_ARG_3
 return
endif

! TODO: comparison with 0
!if (rescale_factor_kappa_en <= 0) then
!   info = QMCKL_INVALID_ARG_4
!   return
!endif

if (walk_num <= 0) then
 info = QMCKL_INVALID_ARG_5
 return
endif

do k=1,walk_num
 info = qmckl_distance_rescaled_deriv_e(context, 'T', 'T', elec_num, nucl_num, &
      elec_coord(1,1,k), elec_num, &
      nucl_coord, nucl_num, &
      en_distance_rescaled_deriv_e(1,1,k), elec_num, rescale_factor_kappa_en)
 if (info /= QMCKL_SUCCESS) then
    exit
 endif
end do

end function qmckl_compute_en_distance_rescaled_deriv_e_f

Test

assert(qmckl_electron_provided(context));

rc = qmckl_set_nucleus_num (context, nucl_num);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_rescale_factor (context, nucl_rescale_factor_kappa);
assert(rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_charge (context, charge);
assert (rc == QMCKL_SUCCESS);

rc = qmckl_set_nucleus_coord (context, 'T', nucl_coord);
assert (rc == QMCKL_SUCCESS);

assert(qmckl_nucleus_provided(context));

double en_distance_rescaled_deriv_e[walk_num][4][nucl_num][elec_num];

rc = qmckl_get_electron_en_distance_rescaled_deriv_e(context, &(en_distance_rescaled_deriv_e[0][0][0][0]));

assert (rc == QMCKL_SUCCESS);

// TODO: check exact values
//// (e,n,w) in Fortran notation
//// (1,1,1)
//assert(fabs(en_distance_rescaled[0][0][0] - 7.546738741619978) < 1.e-12);
//
//// (1,2,1)
//assert(fabs(en_distance_rescaled[0][1][0] - 8.77102435246984) < 1.e-12);
//
//// (2,1,1)
//assert(fabs(en_distance_rescaled[0][0][1] - 3.698922010513608) < 1.e-12);
//
//// (1,1,2)
//assert(fabs(en_distance_rescaled[1][0][0] - 5.824059436060509) < 1.e-12);
//
//// (1,2,2)
//assert(fabs(en_distance_rescaled[1][1][0] - 7.080482110317645) < 1.e-12);
//
//// (2,1,2)
//assert(fabs(en_distance_rescaled[1][0][1] - 3.1804527583077356) < 1.e-12);