74 KiB
Slater Determinant
The slater deteminant is required for the calculation of the wavefunction, gradient, and derivatives. These quantities will be used to calculate the local Energy (\[E_L\]).
ψ(x) = det|ϕ₁(x₁)…ϕᵢ(yᵢ)…ϕₙ(xₙ)|
The above slater-matrix is also required and is denoted by Dᵢⱼ(x) such that:
ψ(x) = det|Dᵢⱼ(x)|
We also require the inverse of the slater-matrix which is denoted by D⁻¹ᵢⱼ(x). Using this notation, the acceptance probability which is proportional to ψ(y)/ψ(x) can be calculated as follows:
ψ(yᵢ)/ψ(xᵢ) = ∑ⱼDᵢⱼ(y)D⁻¹ⱼᵢ(x)
Concerning the gradient and laplacian, in fact what is actually calculated is the ratio of the gradient/laplacian and the determinant of the slater matrix:
∇ψ(x)/ψ(x)
and
∇²ψ(x)/ψ(x)
This avoids the unnecessary multiplication and division of by the determinant ψ(x).
Context
The following arrays are stored in the context:
type |
char |
α ('A' ) or β ('B' ) determinant |
det_num_alpha |
int64_t |
Number of determinants per walker |
det_num_beta |
int64_t |
Number of determinants per walker |
mo_index_alpha |
mo_index[det_num_alpha][walker.num][alpha_num] |
Index of MOs for each walker |
mo_index_beta |
mo_index[det_num_beta][walker.num][beta_num] |
Index of MOs for each walker |
Computed data:
up_num |
int64_t |
Number of number of α electrons |
donwn_num |
int64_t |
Number of number of β electrons |
det_value_alpha |
[det_num_alpha][walker.num] |
The α slater matrix for each determinant of each walker. |
det_value_alpha_date |
uint64_t |
Date of The α slater matrix for each determinant of each walker. |
det_value_beta |
[det_num_beta][walker.num] |
The β slater matrix for each determinant of each walker. |
det_value_beta_date |
uint64_t |
Date of The β slater matrix for each determinant of each walker. |
det_adj_matrix_alpha |
[det_num_alpha][walker.num][alpha_num][alpha_num] |
Adjoint of the α slater matrix for each determinant of each walker. |
det_adj_matrix_alpha_date |
uint64_t |
Date of the Adjoint of the α slater matrix for each determinant of each walker. |
det_adj_matrix_beta |
[det_num_beta][walker.num][beta_num][beta_num] |
Adjoint of the β slater matrix for each determinant of each walker. |
det_adj_matrix_beta_date |
uint64_t |
Date of the Adjoint of the β slater matrix for each determinant of each walker. |
det_vgl_alpha |
[5][det_num_alpha][walker.num][alpha_num][alpha_num] |
Value, gradients, Laplacian of Dᵅᵢⱼ(x) at electron positions |
det_vgl_alpha_date |
uint64_t |
Late modification date of Value, gradients, Laplacian of the MOs at electron positions |
det_vgl_beta |
[5][det_num_beta][walker.num][beta_num][beta_num] |
Value, gradients, Laplacian of Dᵝᵢⱼ(x) at electron positions |
det_vgl_beta_date |
uint64_t |
Late modification date of Value, gradients, Laplacian of the MOs at electron positions |
det_inv_matrix_alpha |
[det_num_alpha][walker.num][alpha_num][alpha_num] |
Inverse of the α electron slater matrix for each determinant of each walker. |
det_inv_matrix_alpha_date |
uint64_t |
Date for the Inverse of the α electron slater matrix for each determinant of each walker. |
det_inv_matrix_beta |
[det_num_beta][walker.num][beta_num][beta_num] |
Inverse of the β electron slater matrix for each determinant of each walker. |
det_inv_matrix_beta_date |
uint64_t |
Date for the Inverse of the β electron slater matrix for each determinant of each walker. |
Data structure
typedef struct qmckl_determinant_struct {
char type;
int64_t det_num_alpha;
int64_t det_num_beta ;
int64_t up_num;
int64_t down_num;
int64_t* mo_index_alpha;
int64_t* mo_index_beta;
double * det_value_alpha;
double * det_value_beta;
double * det_vgl_alpha;
double * det_adj_matrix_alpha;
double * det_inv_matrix_alpha;
double * det_vgl_beta;
double * det_adj_matrix_beta;
double * det_inv_matrix_beta;
uint64_t det_value_alpha_date;
uint64_t det_vgl_alpha_date;
uint64_t det_adj_matrix_alpha_date;
uint64_t det_inv_matrix_alpha_date;
uint64_t det_value_beta_date;
uint64_t det_vgl_beta_date;
uint64_t det_adj_matrix_beta_date;
uint64_t det_inv_matrix_beta_date;
int32_t uninitialized;
bool provided;
} qmckl_determinant_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_determinant(qmckl_context context);
qmckl_exit_code qmckl_init_determinant(qmckl_context context) {
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
return false;
}
qmckl_context_struct* const ctx = (qmckl_context_struct*) context;
assert (ctx != NULL);
ctx->det.uninitialized = (1 << 5) - 1;
return QMCKL_SUCCESS;
}
Access functions
When all the data for the slater determinants have been provided, the following
function returns true
.
bool qmckl_determinant_provided (const qmckl_context context);
#+NAME:post
Initialization functions
To set the basis set, all the following functions need to be called.
qmckl_exit_code qmckl_set_determinant_type (const qmckl_context context, const char t);
qmckl_exit_code qmckl_set_determinant_det_num_alpha (const qmckl_context context, const int64_t det_num_alpha);
qmckl_exit_code qmckl_set_determinant_det_num_beta (const qmckl_context context, const int64_t det_num_beta);
qmckl_exit_code qmckl_set_determinant_mo_index_alpha (const qmckl_context context, const int64_t* mo_index_alpha, const int64_t size_max);
qmckl_exit_code qmckl_set_determinant_mo_index_beta (const qmckl_context context, const int64_t* mo_index_beta, const int64_t size_max);
#+NAME:pre2
#+NAME:post2
When the basis set is completely entered, other data structures are computed to accelerate the calculations.
Fortran Interfaces
Test
Computation
Determinant matrix
Get
qmckl_exit_code qmckl_get_det_vgl_alpha(qmckl_context context, double* const det_vgl_alpha);
qmckl_exit_code qmckl_get_det_vgl_beta(qmckl_context context, double* const det_vgl_beta);
Provide
Compute alpha
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
det_num_alpha |
int64_t |
in | Number of determinants |
walk_num |
int64_t |
in | Number of walkers |
alpha_num |
int64_t |
in | Number of electrons |
beta_num |
int64_t |
in | Number of electrons |
elec_num |
int64_t |
in | Number of electrons |
mo_index_alpha |
int64_t[det_num_alpha][walk_num][alpha_num] |
in | MO indices for electrons |
mo_num |
int64_t |
in | Number of MOs |
mo_vgl |
double[5][elec_num][mo_num] |
in | Value, gradients and Laplacian of the MOs |
det_vgl_alpha |
double[det_num_alpha][walk_num][5][alpha_num][alpha_num] |
out | Value, gradients and Laplacian of the Det |
integer function qmckl_compute_det_vgl_alpha_f(context, &
det_num_alpha, walk_num, alpha_num, beta_num, elec_num, &
mo_index_alpha, mo_num, mo_vgl, det_vgl_alpha) &
result(info)
use qmckl
implicit none
integer(qmckl_context) , intent(in) :: context
integer*8, intent(in) :: det_num_alpha
integer*8, intent(in) :: walk_num
integer*8, intent(in) :: alpha_num
integer*8, intent(in) :: beta_num
integer*8, intent(in) :: elec_num
integer*8, intent(in) :: mo_num
integer*8, intent(in) :: mo_index_alpha(alpha_num, walk_num, det_num_alpha)
double precision, intent(in) :: mo_vgl(mo_num, elec_num, 5)
double precision, intent(inout) :: det_vgl_alpha(alpha_num, alpha_num, 5, walk_num, det_num_alpha)
integer*8 :: idet, iwalk, ielec, mo_id, imo
info = QMCKL_SUCCESS
if (context == QMCKL_NULL_CONTEXT) then
info = QMCKL_INVALID_CONTEXT
return
endif
if (walk_num <= 0) then
info = QMCKL_INVALID_ARG_2
return
endif
if (alpha_num <= 0) then
info = QMCKL_INVALID_ARG_3
return
endif
do idet = 1, det_num_alpha
do iwalk = 1, walk_num
do ielec = 1, alpha_num
do imo = 1, alpha_num
mo_id = mo_index_alpha(imo,iwalk,idet)
! Value
det_vgl_alpha(imo, ielec, 1, iwalk, idet) = mo_vgl(mo_id, ielec, 1)
! Grad_x
det_vgl_alpha(imo, ielec, 2, iwalk, idet) = mo_vgl(mo_id, ielec, 2)
! Grad_y
det_vgl_alpha(imo, ielec, 3, iwalk, idet) = mo_vgl(mo_id, ielec, 3)
! Grad_z
det_vgl_alpha(imo, ielec, 4, iwalk, idet) = mo_vgl(mo_id, ielec, 4)
! Lap
det_vgl_alpha(imo, ielec, 5, iwalk, idet) = mo_vgl(mo_id, ielec, 5)
end do
end do
end do
end do
end function qmckl_compute_det_vgl_alpha_f
qmckl_exit_code qmckl_compute_det_vgl_alpha (
const qmckl_context context,
const int64_t det_num_alpha,
const int64_t walk_num,
const int64_t alpha_num,
const int64_t beta_num,
const int64_t elec_num,
const int64_t* mo_index_alpha,
const int64_t mo_num,
const double* mo_vgl,
double* const det_vgl_alpha );
Compute beta
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
det_num_beta |
int64_t |
in | Number of determinants |
walk_num |
int64_t |
in | Number of walkers |
alpha_num |
int64_t |
in | Number of electrons |
beta_num |
int64_t |
in | Number of electrons |
elec_num |
int64_t |
in | Number of electrons |
mo_index_beta |
int64_t[det_num_beta][walk_num][beta_num] |
in | Number of electrons |
mo_num |
int64_t |
in | Number of MOs |
mo_vgl |
double[5][elec_num][mo_num] |
in | Value, gradients and Laplacian of the MOs |
det_vgl_beta |
double[det_num_beta][walk_num][5][beta_num][beta_num] |
out | Value, gradients and Laplacian of the Det |
integer function qmckl_compute_det_vgl_beta_f(context, &
det_num_beta, walk_num, alpha_num, beta_num, elec_num, &
mo_index_beta, mo_num, mo_vgl, det_vgl_beta) &
result(info)
use qmckl
implicit none
integer(qmckl_context) , intent(in) :: context
integer*8, intent(in) :: det_num_beta
integer*8, intent(in) :: walk_num
integer*8, intent(in) :: alpha_num
integer*8, intent(in) :: beta_num
integer*8, intent(in) :: elec_num
integer*8, intent(in) :: mo_num
integer*8, intent(in) :: mo_index_beta(beta_num,walk_num,det_num_beta)
double precision, intent(in) :: mo_vgl(mo_num, elec_num, 5)
double precision, intent(inout) :: det_vgl_beta(beta_num, beta_num, 5, walk_num, det_num_beta)
integer*8 :: idet, iwalk, ielec, mo_id, imo
info = QMCKL_SUCCESS
if (context == QMCKL_NULL_CONTEXT) then
info = QMCKL_INVALID_CONTEXT
return
endif
if (walk_num <= 0) then
info = QMCKL_INVALID_ARG_2
return
endif
if (beta_num <= 0) then
info = QMCKL_INVALID_ARG_3
return
endif
do idet = 1, det_num_beta
do iwalk = 1, walk_num
do ielec = 1, beta_num
do imo = 1, beta_num
mo_id = mo_index_beta(imo, iwalk, idet)
! Value
det_vgl_beta(imo, ielec, 1, iwalk, idet) = mo_vgl(mo_id, alpha_num + ielec, 1)
! Grad_x
det_vgl_beta(imo, ielec, 2, iwalk, idet) = mo_vgl(mo_id, alpha_num + ielec, 2)
! Grad_y
det_vgl_beta(imo, ielec, 3, iwalk, idet) = mo_vgl(mo_id, alpha_num + ielec, 3)
! Grad_z
det_vgl_beta(imo, ielec, 4, iwalk, idet) = mo_vgl(mo_id, alpha_num + ielec, 4)
! Lap
det_vgl_beta(imo, ielec, 5, iwalk, idet) = mo_vgl(mo_id, alpha_num + ielec, 5)
end do
end do
end do
end do
end function qmckl_compute_det_vgl_beta_f
qmckl_exit_code qmckl_compute_det_vgl_beta (
const qmckl_context context,
const int64_t det_num_beta,
const int64_t walk_num,
const int64_t alpha_num,
const int64_t beta_num,
const int64_t elec_num,
const int64_t* mo_index_beta,
const int64_t mo_num,
const double* mo_vgl,
double* const det_vgl_beta );
Test
Inverse of Determinant matrix
Get
qmckl_exit_code qmckl_get_det_inv_matrix_alpha(qmckl_context context, double* const det_inv_matrix_alpha);
qmckl_exit_code qmckl_get_det_inv_matrix_beta(qmckl_context context, double* const det_inv_matrix_beta);
qmckl_exit_code qmckl_get_det_adj_matrix_alpha(qmckl_context context, double* const det_adj_matrix_alpha);
qmckl_exit_code qmckl_get_det_adj_matrix_beta(qmckl_context context, double* const det_adj_matrix_beta);
qmckl_exit_code qmckl_get_det_alpha(qmckl_context context, double* const det_adj_matrix_alpha);
qmckl_exit_code qmckl_get_det_beta(qmckl_context context, double* const det_adj_matrix_beta);
Provide
Compute alpha
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
det_num_alpha |
int64_t |
in | Number of determinants |
walk_num |
int64_t |
in | Number of walkers |
alpha_num |
int64_t |
in | Number of electrons |
det_vgl_alpha |
double[det_num_alpha][walk_num][5][alpha_num][alpha_num] |
in | determinant matrix Value, gradients and Laplacian of the MOs |
det_value_alpha |
double[det_num_alpha][walk_num] |
out | value of determinant matrix |
det_adj_matrix_alpha |
double[det_num_alpha][walk_num][alpha_num][alpha_num] |
out | adjoint of determinant matrix |
det_inv_matrix_alpha |
double[det_num_alpha][walk_num][alpha_num][alpha_num] |
out | inverse of determinant matrix |
integer function qmckl_compute_det_inv_matrix_alpha_f(context, &
det_num_alpha, walk_num, alpha_num, det_vgl_alpha, det_value_alpha, det_adj_matrix_alpha, det_inv_matrix_alpha) &
result(info)
use qmckl
implicit none
integer(qmckl_context) , intent(in) :: context
integer*8, intent(in) :: det_num_alpha
integer*8, intent(in) :: walk_num
integer*8, intent(in) :: alpha_num
double precision, intent(in) :: det_vgl_alpha(alpha_num, alpha_num, 5, walk_num, det_num_alpha)
double precision, intent(inout) :: det_value_alpha(walk_num, det_num_alpha)
double precision, intent(inout) :: det_adj_matrix_alpha(alpha_num, alpha_num, walk_num, det_num_alpha)
double precision, intent(inout) :: det_inv_matrix_alpha(alpha_num, alpha_num, walk_num, det_num_alpha)
double precision,dimension(:,:),allocatable :: matA
double precision :: det_l
integer*8 :: idet, iwalk, ielec, mo_id, imo, LDA, res, i, j
allocate(matA(alpha_num, alpha_num))
info = QMCKL_SUCCESS
if (context == QMCKL_NULL_CONTEXT) then
info = QMCKL_INVALID_CONTEXT
return
endif
if (det_num_alpha <= 0) then
info = QMCKL_INVALID_ARG_2
return
endif
if (walk_num <= 0) then
info = QMCKL_INVALID_ARG_3
return
endif
if (alpha_num <= 0) then
info = QMCKL_INVALID_ARG_4
return
endif
LDA = alpha_num
do idet = 1, det_num_alpha
do iwalk = 1, walk_num
! Value
matA(1:alpha_num,1:alpha_num) = &
det_vgl_alpha(1:alpha_num, 1:alpha_num, 1, iwalk, idet)
res = qmckl_adjugate(context, &
alpha_num, matA, LDA, &
det_adj_matrix_alpha(1:alpha_num, 1:alpha_num, iwalk, idet), &
int(size(det_adj_matrix_alpha,1),8), &
det_l)
det_inv_matrix_alpha(1:alpha_num, 1:alpha_num, iwalk, idet) = &
(1.d0/det_l) * &
det_adj_matrix_alpha(1:alpha_num, 1:alpha_num, iwalk, idet)
det_value_alpha(iwalk, idet) = det_l
end do
end do
deallocate(matA)
end function qmckl_compute_det_inv_matrix_alpha_f
qmckl_exit_code qmckl_compute_det_inv_matrix_alpha (
const qmckl_context context,
const int64_t det_num_alpha,
const int64_t walk_num,
const int64_t alpha_num,
const double* det_vgl_alpha,
double* const det_value_alpha,
double* const det_adj_matrix_alpha,
double* const det_inv_matrix_alpha );
Compute beta
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
det_num_beta |
int64_t |
in | Number of determinants |
walk_num |
int64_t |
in | Number of walkers |
beta_num |
int64_t |
in | Number of electrons |
det_vgl_beta |
double[det_num_beta][walk_num][5][beta_num][beta_num] |
in | determinant matrix Value, gradients and Laplacian of the MOs |
det_value_beta |
double[det_num_beta][walk_num] |
out | value of determinant matrix |
det_adj_matrix_beta |
double[det_num_beta][walk_num][beta_num][beta_num] |
out | adjoint of determinant matrix |
det_inv_matrix_beta |
double[det_num_beta][walk_num][beta_num][beta_num] |
out | inverse of determinant matrix |
integer function qmckl_compute_det_inv_matrix_beta_f(context, &
det_num_beta, walk_num, beta_num, det_vgl_beta, det_value_beta, det_adj_matrix_beta, det_inv_matrix_beta) &
result(info)
use qmckl
implicit none
integer(qmckl_context) , intent(in) :: context
integer*8, intent(in) :: det_num_beta
integer*8, intent(in) :: walk_num
integer*8, intent(in) :: beta_num
double precision, intent(in) :: det_vgl_beta(beta_num, beta_num, 5, walk_num, det_num_beta)
double precision, intent(inout) :: det_value_beta(walk_num, det_num_beta)
double precision, intent(inout) :: det_adj_matrix_beta(beta_num, beta_num, walk_num, det_num_beta)
double precision, intent(inout) :: det_inv_matrix_beta(beta_num, beta_num, walk_num, det_num_beta)
double precision,dimension(:,:),allocatable :: matA
double precision :: det_l
integer*8 :: idet, iwalk, ielec, mo_id, imo, LDA, res
allocate(matA(beta_num, beta_num))
info = QMCKL_SUCCESS
if (context == QMCKL_NULL_CONTEXT) then
info = QMCKL_INVALID_CONTEXT
return
endif
if (det_num_beta <= 0) then
info = QMCKL_INVALID_ARG_2
return
endif
if (walk_num <= 0) then
info = QMCKL_INVALID_ARG_3
return
endif
if (beta_num <= 0) then
info = QMCKL_INVALID_ARG_4
return
endif
LDA = beta_num
do idet = 1, det_num_beta
do iwalk = 1, walk_num
! Value
matA(1:beta_num,1:beta_num) = &
det_vgl_beta(1:beta_num, 1:beta_num, 1, iwalk, idet)
res = qmckl_adjugate(context, &
beta_num, matA, LDA, &
det_adj_matrix_beta(1, 1, iwalk, idet), &
int(size(det_adj_matrix_beta,1),8), &
det_l)
det_inv_matrix_beta(1:beta_num, 1:beta_num, iwalk, idet) = &
(1.d0/det_l) * &
det_adj_matrix_beta(1:beta_num, 1:beta_num, iwalk, idet)
det_value_beta(iwalk, idet) = det_l
end do
end do
deallocate(matA)
end function qmckl_compute_det_inv_matrix_beta_f
qmckl_exit_code qmckl_compute_det_inv_matrix_beta (
const qmckl_context context,
const int64_t det_num_beta,
const int64_t walk_num,
const int64_t beta_num,
const double* det_vgl_beta,
double* const det_value_beta,
double* const det_adj_matrix_beta,
double* const det_inv_matrix_beta );