Electrons
Table of Contents
- 1. Context
- 2. Computation
- 2.1. Electron-electron distances
- 2.2. Electron-electron rescaled distances
- 2.3. Electron-electron rescaled distance gradients and laplacian with respect to electron coords
- 2.4. Electron-electron potential
- 2.5. Electron-nucleus distances
- 2.6. Electron-nucleus rescaled distances
- 2.7. Electron-nucleus rescaled distance gradients and laplacian with respect to electron coords
- 2.8. Electron-nucleus potential
- 2.9. Generate initial coordinates
1 Context
The following data stored in the context:
Variable | Type | Description |
---|---|---|
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 |
qmckl_matrix |
Current set of electron coordinates. Pointer to ctx->points |
coord_old |
qmckl_matrix |
Old set of electron coordinates |
coord_new_date |
uint64_t |
Last modification date of the coordinates |
Computed data:
Variable | Type | Description |
---|---|---|
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 |
ee_pot |
double[walk_num] |
Electron-electron rescaled distances derivatives |
ee_pot_date |
uint64_t |
Last modification date of the electron-electron distance derivatives |
en_pot |
double[walknum] | Electron-nucleus potential energy |
en_pot_date |
int64t | Date when the electron-nucleus potential energy was computed |
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][nucl_num][num] |
Electron-electron rescaled distances derivatives |
en_distance_rescaled_deriv_e_date |
uint64_t |
Last modification date of the electron-electron distance derivatives |
1.1 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_pot_date; int64_t en_pot_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; qmckl_matrix coord_new; qmckl_matrix coord_old; double* ee_distance; double* en_distance; double* ee_pot; double* en_pot; 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);
1.2 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.
1.2.1 Number of electrons
1.2.2 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.
1.2.3 Scaling factors Kappa
1.2.4 Electron coordinates
Returns the current electron coordinates. The pointer is assumed
to point on a memory block of size size_max
≥ 3 * elec_num * walk_num
.
The order of the indices is:
Normal | Transposed | |
---|---|---|
C | [walk_num*elec_num][3] |
[3][walk_num*elec_num] |
Fortran | (3,walk_num*elec_num) |
(walk_num*elec_num, 3) |
As the coord_new
attribute is a pointer equal to points
,
returning the current electron coordinates is equivalent to
returning the current points.
1.3 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, const int64_t size_max); 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);
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.
size_max
should be equal to elec_num * walk_num * 3
, to be symmetric
with qmckl_get_electron_coord
.
Important: changing the electron coordinates increments the date in the context.
interface integer(c_int32_t) function qmckl_set_electron_coord(context, transp, coord, size_max) 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(*) integer (c_int64_t) , intent(in) , value :: size_max end function end interface
1.4 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, walk_num*elec_num*3); assert(rc == QMCKL_SUCCESS); double elec_coord2[walk_num*3*elec_num]; rc = qmckl_get_electron_coord (context, 'N', elec_coord2, walk_num*3*elec_num); assert(rc == QMCKL_SUCCESS); for (int64_t i=0 ; i<3*elec_num*walk_num ; ++i) { printf("%f %f\n", elec_coord[i], elec_coord2[i]); assert( elec_coord[i] == elec_coord2[i] ); }
2 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.
2.1 Electron-electron distances
2.1.1 Get
qmckl_exit_code qmckl_get_electron_ee_distance(qmckl_context context, double* const distance);
2.1.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
walk_num |
int64_t |
in | Number of walkers |
coord |
double[3][walk_num][elec_num] |
in | Electron coordinates |
ee_distance |
double[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,walk_num,3) double precision , intent(out) :: ee_distance(elec_num,elec_num,walk_num) integer*8 :: k, i, j double precision :: x, y, z 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,k,1), elec_num * walk_num, & coord(1,k,1), elec_num * walk_num, & ee_distance(1,1,k), elec_num) if (info /= QMCKL_SUCCESS) then exit endif end do end function qmckl_compute_ee_distance_f
2.1.3 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);
2.2 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.
2.2.1 Get
qmckl_exit_code qmckl_get_electron_ee_distance_rescaled(qmckl_context context, double* const distance_rescaled);
2.2.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
rescale_factor_kappa_ee |
double |
in | Factor to rescale ee distances |
walk_num |
int64_t |
in | Number of walkers |
coord |
double[walk_num][3][elec_num] |
in | Electron coordinates |
ee_distance |
double[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,k,1), elec_num * walk_num, & coord(1,k,1), elec_num * walk_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
2.2.3 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);
2.3 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\).
2.3.1 Get
qmckl_exit_code qmckl_get_electron_ee_distance_rescaled_deriv_e(qmckl_context context, double* const distance_rescaled_deriv_e);
2.3.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
rescale_factor_kappa_ee |
double |
in | Factor to rescale ee distances |
walk_num |
int64_t |
in | Number of walkers |
coord |
double[walk_num][3][elec_num] |
in | Electron coordinates |
ee_distance_deriv_e |
double[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
2.3.3 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);
2.4 Electron-electron potential
ee_pot
calculates the ee
potential energy.
\[ \mathcal{V}_{ee} = \sum_{i=1}^{N_e}\sum_{j>i}^{N_e}\frac{1}{r_{ij}} \]
where \(\mathcal{V}_{ee}\) is the ee
potential and \[r_{ij}\] the ee
distance.
2.4.1 Get
qmckl_exit_code qmckl_get_electron_ee_potential(qmckl_context context, double* const ee_pot);
2.4.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
walk_num |
int64_t |
in | Number of walkers |
ee_distance |
double[walk_num][elec_num][elec_num] |
in | Electron-electron rescaled distances |
ee_pot |
double[walk_num] |
out | Electron-electron potential |
integer function qmckl_compute_ee_potential_f(context, elec_num, walk_num, & ee_distance, ee_pot) & 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) :: ee_distance(elec_num,elec_num,walk_num) double precision , intent(out) :: ee_pot(walk_num) integer*8 :: nw, i, j 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 ee_pot = 0.0d0 do nw=1,walk_num do j=2,elec_num do i=1,j-1 if (dabs(ee_distance(i,j,nw)) > 1e-5) then ee_pot(nw) = ee_pot(nw) + 1.0d0/(ee_distance(i,j,nw)) endif end do end do end do end function qmckl_compute_ee_potential_f
qmckl_exit_code qmckl_compute_ee_potential ( const qmckl_context context, const int64_t elec_num, const int64_t walk_num, const double* ee_distance, double* const ee_pot );
2.4.3 Test
double ee_pot[walk_num]; rc = qmckl_get_electron_ee_potential(context, &(ee_pot[0])); assert (rc == QMCKL_SUCCESS);
2.5 Electron-nucleus distances
2.5.1 Get
qmckl_exit_code qmckl_get_electron_en_distance(qmckl_context context, double* distance);
2.5.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
nucl_num |
int64_t |
in | Number of nuclei |
walk_num |
int64_t |
in | Number of walkers |
elec_coord |
double[walk_num][3][elec_num] |
in | Electron coordinates |
nucl_coord |
double[3][elec_num] |
in | Nuclear coordinates |
en_distance |
double[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,walk_num,3) 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,k,1), elec_num * walk_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
2.5.3 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, nucl_num); assert (rc == QMCKL_SUCCESS); rc = qmckl_set_nucleus_coord (context, 'T', nucl_coord, 3*nucl_num); 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);
2.6 Electron-nucleus rescaled distances
en_distance_rescaled
stores the matrix of the rescaled distances between
electrons and nuclei.
\[ C_{ij} = \left( 1 - \exp{-\kappa C_{ij}}\right)/\kappa \]
where \(C_{ij}\) is the matrix of electron-nucleus distances.
2.6.1 Get
qmckl_exit_code qmckl_get_electron_en_distance_rescaled(qmckl_context context, double* distance_rescaled);
2.6.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
nucl_num |
int64_t |
in | Number of nuclei |
rescale_factor_kappa_en |
double |
in | The factor for rescaled distances |
walk_num |
int64_t |
in | Number of walkers |
elec_coord |
double[walk_num][3][elec_num] |
in | Electron coordinates |
nucl_coord |
double[3][elec_num] |
in | Nuclear coordinates |
en_distance_rescaled |
double[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,walk_num,3) 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,k,1), elec_num*walk_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
2.6.3 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, nucl_num); assert (rc == QMCKL_SUCCESS); rc = qmckl_set_nucleus_coord (context, 'T', nucl_coord, 3*nucl_num); 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);
2.7 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][nucl_num][elec_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\).
2.7.1 Get
qmckl_exit_code qmckl_get_electron_en_distance_rescaled_deriv_e(qmckl_context context, double* distance_rescaled_deriv_e);
2.7.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
nucl_num |
int64_t |
in | Number of nuclei |
rescale_factor_kappa_en |
double |
in | The factor for rescaled distances |
walk_num |
int64_t |
in | Number of walkers |
elec_coord |
double[walk_num][3][elec_num] |
in | Electron coordinates |
nucl_coord |
double[3][elec_num] |
in | Nuclear coordinates |
en_distance_rescaled_deriv_e |
double[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,walk_num,3) double precision , intent(in) :: nucl_coord(nucl_num,3) double precision , intent(out) :: en_distance_rescaled_deriv_e(4,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,k,1), elec_num*walk_num, & nucl_coord, nucl_num, & en_distance_rescaled_deriv_e(1,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
2.7.3 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, nucl_num); assert (rc == QMCKL_SUCCESS); rc = qmckl_set_nucleus_coord (context, 'T', nucl_coord, 3*nucl_num); 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);
2.8 Electron-nucleus potential
en_potential
stores the en
potential energy
\[ \mathcal{V}_{en} = -\sum_{i=1}^{N_e}\sum_{A=1}^{N_n}\frac{Z_A}{r_{iA}} \]
where \(\mathcal{V}_{en}\) is the en
potential, \[r_{iA}\] the en
distance and \[Z_A\] is the nuclear charge.
2.8.1 Get
qmckl_exit_code qmckl_get_electron_en_potential(qmckl_context context, double* const en_pot);
2.8.2 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
elec_num |
int64_t |
in | Number of electrons |
nucl_num |
int64_t |
in | Number of nuclei |
walk_num |
int64_t |
in | Number of walkers |
charge |
double[nucl_num] |
in | charge of nucleus |
en_distance |
double[walk_num][nucl_num][elec_num] |
in | Electron-electron rescaled distances |
en_pot |
double[walk_num] |
out | Electron-electron potential |
integer function qmckl_compute_en_potential_f(context, elec_num, nucl_num, walk_num, & charge, en_distance, en_pot) & 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) :: charge(nucl_num) double precision , intent(in) :: en_distance(elec_num,nucl_num,walk_num) double precision , intent(out) :: en_pot(walk_num) integer*8 :: nw, i, j 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 en_pot = 0.0d0 do nw=1,walk_num do j=1,nucl_num do i=1,elec_num if (dabs(en_distance(i,j,nw)) > 1e-5) then en_pot(nw) = en_pot(nw) - charge(j)/(en_distance(i,j,nw)) endif end do end do end do end function qmckl_compute_en_potential_f
qmckl_exit_code qmckl_compute_en_potential ( const qmckl_context context, const int64_t elec_num, const int64_t nucl_num, const int64_t walk_num, const double* charge, const double* en_distance, double* const en_pot );
2.8.3 Test
double en_pot[walk_num]; rc = qmckl_get_electron_en_potential(context, &(en_pot[0])); assert (rc == QMCKL_SUCCESS);