mirror of
https://github.com/TREX-CoE/qmckl.git
synced 2024-11-19 20:42:50 +01:00
HPC version of qmckl_compute_jastrow_champ_factor_een_deriv_e
This commit is contained in:
parent
2153cfccf6
commit
5093b2c35c
@ -543,6 +543,7 @@ double qmckl_mat (qmckl_matrix m, int64_t i, int64_t j) // m[j][i]
|
|||||||
double qmckl_ten3(qmckl_tensor t, int64_t i, int64_t j, int64_t k); // t[k][j][i]
|
double qmckl_ten3(qmckl_tensor t, int64_t i, int64_t j, int64_t k); // t[k][j][i]
|
||||||
double qmckl_ten4(qmckl_tensor t, int64_t i, int64_t j, int64_t k, int64_t l);
|
double qmckl_ten4(qmckl_tensor t, int64_t i, int64_t j, int64_t k, int64_t l);
|
||||||
double qmckl_ten5(qmckl_tensor t, int64_t i, int64_t j, int64_t k, int64_t l, int64_t m);
|
double qmckl_ten5(qmckl_tensor t, int64_t i, int64_t j, int64_t k, int64_t l, int64_t m);
|
||||||
|
double qmckl_ten6(qmckl_tensor t, int64_t i, int64_t j, int64_t k, int64_t l, int64_t m, int64_t n);
|
||||||
...
|
...
|
||||||
#+end_src
|
#+end_src
|
||||||
|
|
||||||
@ -553,6 +554,7 @@ double qmckl_ten5(qmckl_tensor t, int64_t i, int64_t j, int64_t k, int64_t l, in
|
|||||||
#define qmckl_ten3(t, i, j, k) t.data[(i) + t.size[0]*((j) + t.size[1]*(k))]
|
#define qmckl_ten3(t, i, j, k) t.data[(i) + t.size[0]*((j) + t.size[1]*(k))]
|
||||||
#define qmckl_ten4(t, i, j, k, l) t.data[(i) + t.size[0]*((j) + t.size[1]*((k) + t.size[2]*(l)))]
|
#define qmckl_ten4(t, i, j, k, l) t.data[(i) + t.size[0]*((j) + t.size[1]*((k) + t.size[2]*(l)))]
|
||||||
#define qmckl_ten5(t, i, j, k, l, m) t.data[(i) + t.size[0]*((j) + t.size[1]*((k) + t.size[2]*((l) + t.size[3]*(m))))]
|
#define qmckl_ten5(t, i, j, k, l, m) t.data[(i) + t.size[0]*((j) + t.size[1]*((k) + t.size[2]*((l) + t.size[3]*(m))))]
|
||||||
|
#define qmckl_ten6(t, i, j, k, l, m, n) t.data[(i) + t.size[0]*((j) + t.size[1]*((k) + t.size[2]*((l) + t.size[3]*((m) + t.size[4]*(n)))))]
|
||||||
#+end_src
|
#+end_src
|
||||||
|
|
||||||
For example:
|
For example:
|
||||||
@ -640,9 +642,7 @@ qmckl_double_of_vector(const qmckl_context context,
|
|||||||
assert (target != NULL);
|
assert (target != NULL);
|
||||||
assert (size_max > (int64_t) 0);
|
assert (size_max > (int64_t) 0);
|
||||||
assert (size_max >= vector.size);
|
assert (size_max >= vector.size);
|
||||||
for (int64_t i=0 ; i<vector.size ; ++i) {
|
memcpy(target, vector.data, vector.size*sizeof(double));
|
||||||
target[i] = vector.data[i];
|
|
||||||
}
|
|
||||||
return QMCKL_SUCCESS;
|
return QMCKL_SUCCESS;
|
||||||
|
|
||||||
}
|
}
|
||||||
|
@ -147,50 +147,50 @@ int main() {
|
|||||||
|
|
||||||
Computed data:
|
Computed data:
|
||||||
|
|
||||||
| Variable | Type | In/Out | |
|
| Variable | Type | In/Out |
|
||||||
|-------------------------------------+-----------------------------------------------------------------+---------------------------------------------------------------------------------------------------------+----------------------------------------------|
|
|-------------------------------------+-----------------------------------------------------------------+---------------------------------------------------------------------------------------------------------|
|
||||||
| ~dim_c_vector~ | ~int64_t~ | Number of unique C coefficients | |
|
| ~dim_c_vector~ | ~int64_t~ | Number of unique C coefficients |
|
||||||
| ~dim_c_vector_date~ | ~uint64_t~ | Number of unique C coefficients | |
|
| ~dim_c_vector_date~ | ~uint64_t~ | Number of unique C coefficients |
|
||||||
| ~asymp_jasa~ | ~double[type_nucl_num]~ | Asymptotic component | |
|
| ~asymp_jasa~ | ~double[type_nucl_num]~ | Asymptotic component |
|
||||||
| ~asymp_jasa_date~ | ~uint64_t~ | Ladt modification of the asymptotic component | |
|
| ~asymp_jasa_date~ | ~uint64_t~ | Ladt modification of the asymptotic component |
|
||||||
| ~asymp_jasb~ | ~double[2]~ | Asymptotic component (up- or down-spin) | |
|
| ~asymp_jasb~ | ~double[2]~ | Asymptotic component (up- or down-spin) |
|
||||||
| ~asymp_jasb_date~ | ~uint64_t~ | Ladt modification of the asymptotic component | |
|
| ~asymp_jasb_date~ | ~uint64_t~ | Ladt modification of the asymptotic component |
|
||||||
| ~c_vector_full~ | ~double[dim_c_vector][nucl_num]~ | vector of non-zero coefficients | |
|
| ~c_vector_full~ | ~double[dim_c_vector][nucl_num]~ | vector of non-zero coefficients |
|
||||||
| ~c_vector_full_date~ | ~uint64_t~ | Keep track of changes here | |
|
| ~c_vector_full_date~ | ~uint64_t~ | Keep track of changes here |
|
||||||
| ~lkpm_combined_index~ | ~int64_t[4][dim_c_vector]~ | Transform l,k,p, and m into consecutive indices | |
|
| ~lkpm_combined_index~ | ~int64_t[4][dim_c_vector]~ | Transform l,k,p, and m into consecutive indices |
|
||||||
| ~lkpm_combined_index_date~ | ~uint64_t~ | Transform l,k,p, and m into consecutive indices | |
|
| ~lkpm_combined_index_date~ | ~uint64_t~ | Transform l,k,p, and m into consecutive indices |
|
||||||
| ~tmp_c~ | ~double[walk_num][cord_num][cord_num+1][nucl_num][elec_num]~ | vector of non-zero coefficients | |
|
| ~tmp_c~ | ~double[walk_num][cord_num][cord_num+1][nucl_num][elec_num]~ | vector of non-zero coefficients |
|
||||||
| ~dtmp_c~ | ~double[walk_num][elec_num][4][nucl_num][cord_num+1][cord_num]~ | vector of non-zero coefficients | |
|
| ~dtmp_c~ | ~double[walk_num][elec_num][4][nucl_num][cord_num+1][cord_num]~ | vector of non-zero coefficients |
|
||||||
| ~ee_distance_rescaled~ | ~double[walk_num][num][num]~ | Electron-electron rescaled 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_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~ | ~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_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~ | ~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_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~ | ~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 | |
|
| ~en_distance_rescaled_deriv_e_date~ | ~uint64_t~ | Last modification date of the electron-electron distance derivatives |
|
||||||
| ~een_rescaled_n~ | ~double[walk_num][cord_num+1][nucl_num][elec_num]~ | The electron-electron rescaled distances raised to the powers defined by cord | |
|
| ~een_rescaled_n~ | ~double[walk_num][cord_num+1][nucl_num][elec_num]~ | The electron-electron rescaled distances raised to the powers defined by cord |
|
||||||
| ~een_rescaled_n_date~ | ~uint64_t~ | Keep track of the date of creation | |
|
| ~een_rescaled_n_date~ | ~uint64_t~ | Keep track of the date of creation |
|
||||||
| ~een_rescaled_e_deriv_e~ | ~double[walk_num][cord_num+1][elec_num][4][elec_num]~ | The electron-electron rescaled distances raised to the powers defined by cord derivatives wrt electrons | |
|
| ~een_rescaled_e_deriv_e~ | ~double[walk_num][cord_num+1][elec_num][4][elec_num]~ | The electron-electron rescaled distances raised to the powers defined by cord derivatives wrt electrons |
|
||||||
| ~een_rescaled_e_deriv_e_date~ | ~uint64_t~ | Keep track of the date of creation | |
|
| ~een_rescaled_e_deriv_e_date~ | ~uint64_t~ | Keep track of the date of creation |
|
||||||
| ~een_rescaled_n_deriv_e~ | ~double[walk_num][cord_num+1][nucl_num][4][elec_num]~ | The electron-electron rescaled distances raised to the powers defined by cord derivatives wrt electrons | |
|
| ~een_rescaled_n_deriv_e~ | ~double[walk_num][cord_num+1][nucl_num][4][elec_num]~ | The electron-electron rescaled distances raised to the powers defined by cord derivatives wrt electrons |
|
||||||
| ~een_rescaled_n_deriv_e_date~ | ~uint64_t~ | Keep track of the date of creation | |
|
| ~een_rescaled_n_deriv_e_date~ | ~uint64_t~ | Keep track of the date of creation |
|
||||||
| ~factor_ee~ | ~double[walk_num]~ | Jastrow factor: electron-electron part | |
|
| ~factor_ee~ | ~double[walk_num]~ | Jastrow factor: electron-electron part |
|
||||||
| ~factor_ee_date~ | ~uint64_t~ | Jastrow factor: electron-electron part | |
|
| ~factor_ee_date~ | ~uint64_t~ | Jastrow factor: electron-electron part |
|
||||||
| ~factor_en~ | ~double[walk_num]~ | Jastrow factor: electron-nucleus part | |
|
| ~factor_en~ | ~double[walk_num]~ | Jastrow factor: electron-nucleus part |
|
||||||
| ~factor_en_date~ | ~uint64_t~ | Jastrow factor: electron-nucleus part | |
|
| ~factor_en_date~ | ~uint64_t~ | Jastrow factor: electron-nucleus part |
|
||||||
| ~factor_een~ | ~double[walk_num]~ | Jastrow factor: electron-electron-nucleus part | |
|
| ~factor_een~ | ~double[walk_num]~ | Jastrow factor: electron-electron-nucleus part |
|
||||||
| ~factor_een_date~ | ~uint64_t~ | Jastrow factor: electron-electron-nucleus part | |
|
| ~factor_een_date~ | ~uint64_t~ | Jastrow factor: electron-electron-nucleus part |
|
||||||
| ~factor_ee_deriv_e~ | ~double[walk_num][4][elec_num]~ | Derivative of the Jastrow factor: electron-electron-nucleus part | |
|
| ~factor_ee_deriv_e~ | ~double[walk_num][4][elec_num]~ | Derivative of the Jastrow factor: electron-electron-nucleus part |
|
||||||
| ~factor_ee_deriv_e_date~ | ~uint64_t~ | Keep track of the date for the derivative | |
|
| ~factor_ee_deriv_e_date~ | ~uint64_t~ | Keep track of the date for the derivative |
|
||||||
| ~factor_en_deriv_e~ | ~double[walk_num][4][elec_num]~ | Derivative of the Jastrow factor: electron-electron-nucleus part | |
|
| ~factor_en_deriv_e~ | ~double[walk_num][4][elec_num]~ | Derivative of the Jastrow factor: electron-electron-nucleus part |
|
||||||
| ~factor_en_deriv_e_date~ | ~uint64_t~ | Keep track of the date for the en derivative | |
|
| ~factor_en_deriv_e_date~ | ~uint64_t~ | Keep track of the date for the en derivative |
|
||||||
| ~factor_een_deriv_e~ | ~double[walk_num][4][elec_num]~ | Derivative of the Jastrow factor: electron-electron-nucleus part | |
|
| ~factor_een_deriv_e~ | ~double[walk_num][4][elec_num]~ | Derivative of the Jastrow factor: electron-electron-nucleus part |
|
||||||
| ~factor_een_deriv_e_date~ | ~uint64_t~ | Keep track of the date for the een derivative | |
|
| ~factor_een_deriv_e_date~ | ~uint64_t~ | Keep track of the date for the een derivative |
|
||||||
| ~value~ | ~double[walk_num]~ | out | Value of the Jastrow factor |
|
| ~value~ | ~double[walk_num]~ | Value of the Jastrow factor |
|
||||||
| ~value_date~ | ~uint64_t~ | out | Keep track of the date |
|
| ~value_date~ | ~uint64_t~ | Keep track of the date |
|
||||||
| ~gl~ | ~double[walk_num][4][elec_num]~ | out | Gradient and Laplacian of the Jastrow factor |
|
| ~gl~ | ~double[walk_num][4][elec_num]~ | Gradient and Laplacian of the Jastrow factor |
|
||||||
| ~value_date~ | ~uint64_t~ | out | Keep track of the date |
|
| ~value_date~ | ~uint64_t~ | Keep track of the date |
|
||||||
|
|
||||||
#+NAME: jastrow_data
|
#+NAME: jastrow_data
|
||||||
#+BEGIN_SRC python :results none :exports none
|
#+BEGIN_SRC python :results none :exports none
|
||||||
@ -9372,7 +9372,7 @@ end function qmckl_compute_jastrow_champ_factor_een_deriv_e_naive_f
|
|||||||
|
|
||||||
|
|
||||||
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
|
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
|
||||||
integer function qmckl_compute_jastrow_champ_factor_een_deriv_e_f( &
|
integer function qmckl_compute_jastrow_champ_factor_een_deriv_e_doc_f( &
|
||||||
context, walk_num, elec_num, nucl_num, &
|
context, walk_num, elec_num, nucl_num, &
|
||||||
cord_num, dim_c_vector, c_vector_full, lkpm_combined_index, &
|
cord_num, dim_c_vector, c_vector_full, lkpm_combined_index, &
|
||||||
tmp_c, dtmp_c, een_rescaled_n, een_rescaled_n_deriv_e, factor_een_deriv_e)&
|
tmp_c, dtmp_c, een_rescaled_n, een_rescaled_n_deriv_e, factor_een_deriv_e)&
|
||||||
@ -9389,7 +9389,7 @@ integer function qmckl_compute_jastrow_champ_factor_een_deriv_e_f( &
|
|||||||
double precision , intent(in) :: een_rescaled_n_deriv_e(elec_num, 4, nucl_num, 0:cord_num, walk_num)
|
double precision , intent(in) :: een_rescaled_n_deriv_e(elec_num, 4, nucl_num, 0:cord_num, walk_num)
|
||||||
double precision , intent(out) :: factor_een_deriv_e(elec_num,4,walk_num)
|
double precision , intent(out) :: factor_een_deriv_e(elec_num,4,walk_num)
|
||||||
|
|
||||||
integer*8 :: i, a, j, l, k, p, m, n, nw, ii
|
integer*8 :: i, a, j, l, k, m, n, nw, ii
|
||||||
double precision :: accu, accu2, cn
|
double precision :: accu, accu2, cn
|
||||||
|
|
||||||
info = QMCKL_SUCCESS
|
info = QMCKL_SUCCESS
|
||||||
@ -9422,121 +9422,294 @@ integer function qmckl_compute_jastrow_champ_factor_een_deriv_e_f( &
|
|||||||
factor_een_deriv_e = 0.0d0
|
factor_een_deriv_e = 0.0d0
|
||||||
|
|
||||||
do nw =1, walk_num
|
do nw =1, walk_num
|
||||||
do n = 1, dim_c_vector
|
do n = 1, dim_c_vector
|
||||||
l = lkpm_combined_index(n, 1)
|
l = lkpm_combined_index(n, 1)
|
||||||
k = lkpm_combined_index(n, 2)
|
k = lkpm_combined_index(n, 2)
|
||||||
p = lkpm_combined_index(n, 3)
|
m = lkpm_combined_index(n, 4)
|
||||||
m = lkpm_combined_index(n, 4)
|
|
||||||
|
|
||||||
do a = 1, nucl_num
|
do a = 1, nucl_num
|
||||||
cn = c_vector_full(a, n)
|
cn = c_vector_full(a, n)
|
||||||
if(cn == 0.d0) cycle
|
if(cn == 0.d0) cycle
|
||||||
|
|
||||||
do ii = 1, 4
|
do ii = 1, 4
|
||||||
do j = 1, elec_num
|
do j = 1, elec_num
|
||||||
factor_een_deriv_e(j,ii,nw) = factor_een_deriv_e(j,ii,nw) + (&
|
factor_een_deriv_e(j,ii,nw) = factor_een_deriv_e(j,ii,nw) + ( &
|
||||||
tmp_c(j,a,m,k,nw) * een_rescaled_n_deriv_e(j,ii,a,m+l,nw) + &
|
tmp_c(j,a,m,k,nw) * een_rescaled_n_deriv_e(j,ii,a,m+l,nw) + &
|
||||||
(dtmp_c(j,ii,a,m,k,nw)) * een_rescaled_n(j,a,m+l,nw) + &
|
(dtmp_c(j,ii,a,m,k,nw)) * een_rescaled_n(j,a,m+l,nw) + &
|
||||||
(dtmp_c(j,ii,a,m+l,k,nw)) * een_rescaled_n(j,a,m ,nw) + &
|
(dtmp_c(j,ii,a,m+l,k,nw)) * een_rescaled_n(j,a,m ,nw) + &
|
||||||
tmp_c(j,a,m+l,k,nw) * een_rescaled_n_deriv_e(j,ii,a,m,nw) &
|
tmp_c(j,a,m+l,k,nw) * een_rescaled_n_deriv_e(j,ii,a,m,nw) &
|
||||||
) * cn
|
) * cn
|
||||||
|
end do
|
||||||
|
end do
|
||||||
|
|
||||||
|
cn = cn + cn
|
||||||
|
do j = 1, elec_num
|
||||||
|
factor_een_deriv_e(j,4,nw) = factor_een_deriv_e(j,4,nw) + ( &
|
||||||
|
(dtmp_c(j,1,a,m ,k,nw)) * een_rescaled_n_deriv_e(j,1,a,m+l,nw) + &
|
||||||
|
(dtmp_c(j,2,a,m ,k,nw)) * een_rescaled_n_deriv_e(j,2,a,m+l,nw) + &
|
||||||
|
(dtmp_c(j,3,a,m ,k,nw)) * een_rescaled_n_deriv_e(j,3,a,m+l,nw) + &
|
||||||
|
(dtmp_c(j,1,a,m+l,k,nw)) * een_rescaled_n_deriv_e(j,1,a,m ,nw) + &
|
||||||
|
(dtmp_c(j,2,a,m+l,k,nw)) * een_rescaled_n_deriv_e(j,2,a,m ,nw) + &
|
||||||
|
(dtmp_c(j,3,a,m+l,k,nw)) * een_rescaled_n_deriv_e(j,3,a,m ,nw) &
|
||||||
|
) * cn
|
||||||
|
end do
|
||||||
end do
|
end do
|
||||||
end do
|
end do
|
||||||
|
|
||||||
cn = cn + cn
|
|
||||||
do j = 1, elec_num
|
|
||||||
factor_een_deriv_e(j,4,nw) = factor_een_deriv_e(j,4,nw) + (&
|
|
||||||
(dtmp_c(j,1,a,m ,k,nw)) * een_rescaled_n_deriv_e(j,1,a,m+l,nw) + &
|
|
||||||
(dtmp_c(j,2,a,m ,k,nw)) * een_rescaled_n_deriv_e(j,2,a,m+l,nw) + &
|
|
||||||
(dtmp_c(j,3,a,m ,k,nw)) * een_rescaled_n_deriv_e(j,3,a,m+l,nw) + &
|
|
||||||
(dtmp_c(j,1,a,m+l,k,nw)) * een_rescaled_n_deriv_e(j,1,a,m ,nw) + &
|
|
||||||
(dtmp_c(j,2,a,m+l,k,nw)) * een_rescaled_n_deriv_e(j,2,a,m ,nw) + &
|
|
||||||
(dtmp_c(j,3,a,m+l,k,nw)) * een_rescaled_n_deriv_e(j,3,a,m ,nw) &
|
|
||||||
) * cn
|
|
||||||
end do
|
|
||||||
end do
|
|
||||||
end do
|
|
||||||
end do
|
end do
|
||||||
|
|
||||||
end function qmckl_compute_jastrow_champ_factor_een_deriv_e_f
|
end function qmckl_compute_jastrow_champ_factor_een_deriv_e_doc_f
|
||||||
#+end_src
|
#+end_src
|
||||||
|
|
||||||
# #+CALL: generate_c_header(table=qmckl_factor_een_deriv_e_args,rettyp=get_value("CRetType"),fname=get_value("Name"))
|
#+CALL: generate_private_c_header(table=qmckl_factor_een_deriv_e_args,rettyp=get_value("CRetType"),fname="qmckl_compute_jastrow_champ_factor_een_deriv_e_doc" )
|
||||||
|
|
||||||
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
|
#+RESULTS:
|
||||||
qmckl_exit_code qmckl_compute_jastrow_champ_factor_een_deriv_e (
|
#+begin_src c :tangle (eval h_private_func) :comments org
|
||||||
const qmckl_context context,
|
qmckl_exit_code qmckl_compute_jastrow_champ_factor_een_deriv_e_doc (
|
||||||
const int64_t walk_num,
|
const qmckl_context context,
|
||||||
const int64_t elec_num,
|
const int64_t walk_num,
|
||||||
const int64_t nucl_num,
|
const int64_t elec_num,
|
||||||
const int64_t cord_num,
|
const int64_t nucl_num,
|
||||||
const int64_t dim_c_vector,
|
const int64_t cord_num,
|
||||||
const double* c_vector_full,
|
const int64_t dim_c_vector,
|
||||||
const int64_t* lkpm_combined_index,
|
const double* c_vector_full,
|
||||||
const double* tmp_c,
|
const int64_t* lkpm_combined_index,
|
||||||
const double* dtmp_c,
|
const double* tmp_c,
|
||||||
const double* een_rescaled_n,
|
const double* dtmp_c,
|
||||||
const double* een_rescaled_n_deriv_e,
|
const double* een_rescaled_n,
|
||||||
double* const factor_een_deriv_e );
|
const double* een_rescaled_n_deriv_e,
|
||||||
#+end_src
|
double* const factor_een_deriv_e );
|
||||||
|
#+end_src
|
||||||
|
|
||||||
#+CALL: generate_c_interface(table=qmckl_factor_een_deriv_e_args,rettyp=get_value("CRetType"),fname=get_value("Name"))
|
#+CALL: generate_c_interface(table=qmckl_factor_een_deriv_e_args,rettyp=get_value("CRetType"),fname="qmckl_compute_jastrow_champ_factor_een_deriv_e_doc"))
|
||||||
|
|
||||||
#+RESULTS:
|
#+RESULTS:
|
||||||
#+begin_src f90 :tangle (eval f) :comments org :exports none
|
#+begin_src f90 :tangle (eval f) :comments org :exports none
|
||||||
integer(c_int32_t) function qmckl_compute_jastrow_champ_factor_een_deriv_e &
|
integer(c_int32_t) function qmckl_compute_jastrow_champ_factor_een_deriv_e_doc &
|
||||||
(context, &
|
|
||||||
walk_num, &
|
|
||||||
elec_num, &
|
|
||||||
nucl_num, &
|
|
||||||
cord_num, &
|
|
||||||
dim_c_vector, &
|
|
||||||
c_vector_full, &
|
|
||||||
lkpm_combined_index, &
|
|
||||||
tmp_c, &
|
|
||||||
dtmp_c, &
|
|
||||||
een_rescaled_n, &
|
|
||||||
een_rescaled_n_deriv_e, &
|
|
||||||
factor_een_deriv_e) &
|
|
||||||
bind(C) result(info)
|
|
||||||
|
|
||||||
use, intrinsic :: iso_c_binding
|
|
||||||
implicit none
|
|
||||||
|
|
||||||
integer (c_int64_t) , intent(in) , value :: context
|
|
||||||
integer (c_int64_t) , intent(in) , value :: walk_num
|
|
||||||
integer (c_int64_t) , intent(in) , value :: elec_num
|
|
||||||
integer (c_int64_t) , intent(in) , value :: nucl_num
|
|
||||||
integer (c_int64_t) , intent(in) , value :: cord_num
|
|
||||||
integer (c_int64_t) , intent(in) , value :: dim_c_vector
|
|
||||||
real (c_double ) , intent(in) :: c_vector_full(nucl_num,dim_c_vector)
|
|
||||||
integer (c_int64_t) , intent(in) :: lkpm_combined_index(dim_c_vector,4)
|
|
||||||
real (c_double ) , intent(in) :: tmp_c(elec_num,nucl_num,0:cord_num,0:cord_num-1,walk_num)
|
|
||||||
real (c_double ) , intent(in) :: dtmp_c(elec_num,4,nucl_num,0:cord_num,0:cord_num-1,walk_num)
|
|
||||||
real (c_double ) , intent(in) :: een_rescaled_n(elec_num,nucl_num,0:cord_num,walk_num)
|
|
||||||
real (c_double ) , intent(in) :: een_rescaled_n_deriv_e(elec_num,4,nucl_num,0:cord_num,walk_num)
|
|
||||||
real (c_double ) , intent(out) :: factor_een_deriv_e(elec_num,4,walk_num)
|
|
||||||
|
|
||||||
integer(c_int32_t), external :: qmckl_compute_jastrow_champ_factor_een_deriv_e_f
|
|
||||||
info = qmckl_compute_jastrow_champ_factor_een_deriv_e_f &
|
|
||||||
(context, &
|
(context, &
|
||||||
walk_num, &
|
walk_num, &
|
||||||
elec_num, &
|
elec_num, &
|
||||||
nucl_num, &
|
nucl_num, &
|
||||||
cord_num, &
|
cord_num, &
|
||||||
dim_c_vector, &
|
dim_c_vector, &
|
||||||
c_vector_full, &
|
c_vector_full, &
|
||||||
lkpm_combined_index, &
|
lkpm_combined_index, &
|
||||||
tmp_c, &
|
tmp_c, &
|
||||||
dtmp_c, &
|
dtmp_c, &
|
||||||
een_rescaled_n, &
|
een_rescaled_n, &
|
||||||
een_rescaled_n_deriv_e, &
|
een_rescaled_n_deriv_e, &
|
||||||
factor_een_deriv_e)
|
factor_een_deriv_e) &
|
||||||
|
bind(C) result(info)
|
||||||
|
|
||||||
end function qmckl_compute_jastrow_champ_factor_een_deriv_e
|
use, intrinsic :: iso_c_binding
|
||||||
|
implicit none
|
||||||
|
|
||||||
|
integer (c_int64_t) , intent(in) , value :: context
|
||||||
|
integer (c_int64_t) , intent(in) , value :: walk_num
|
||||||
|
integer (c_int64_t) , intent(in) , value :: elec_num
|
||||||
|
integer (c_int64_t) , intent(in) , value :: nucl_num
|
||||||
|
integer (c_int64_t) , intent(in) , value :: cord_num
|
||||||
|
integer (c_int64_t) , intent(in) , value :: dim_c_vector
|
||||||
|
real (c_double ) , intent(in) :: c_vector_full(nucl_num,dim_c_vector)
|
||||||
|
integer (c_int64_t) , intent(in) :: lkpm_combined_index(dim_c_vector,4)
|
||||||
|
real (c_double ) , intent(in) :: tmp_c(elec_num,nucl_num,0:cord_num,0:cord_num-1,walk_num)
|
||||||
|
real (c_double ) , intent(in) :: dtmp_c(elec_num,4,nucl_num,0:cord_num,0:cord_num-1,walk_num)
|
||||||
|
real (c_double ) , intent(in) :: een_rescaled_n(elec_num,nucl_num,0:cord_num,walk_num)
|
||||||
|
real (c_double ) , intent(in) :: een_rescaled_n_deriv_e(elec_num,4,nucl_num,0:cord_num,walk_num)
|
||||||
|
real (c_double ) , intent(out) :: factor_een_deriv_e(elec_num,4,walk_num)
|
||||||
|
|
||||||
|
integer(c_int32_t), external :: qmckl_compute_jastrow_champ_factor_een_deriv_e_doc_f
|
||||||
|
info = qmckl_compute_jastrow_champ_factor_een_deriv_e_doc_f &
|
||||||
|
(context, &
|
||||||
|
walk_num, &
|
||||||
|
elec_num, &
|
||||||
|
nucl_num, &
|
||||||
|
cord_num, &
|
||||||
|
dim_c_vector, &
|
||||||
|
c_vector_full, &
|
||||||
|
lkpm_combined_index, &
|
||||||
|
tmp_c, &
|
||||||
|
dtmp_c, &
|
||||||
|
een_rescaled_n, &
|
||||||
|
een_rescaled_n_deriv_e, &
|
||||||
|
factor_een_deriv_e)
|
||||||
|
|
||||||
|
end function qmckl_compute_jastrow_champ_factor_een_deriv_e_doc
|
||||||
#+end_src
|
#+end_src
|
||||||
|
|
||||||
|
#+CALL: generate_private_c_header(table=qmckl_factor_een_deriv_e_args,rettyp=get_value("CRetType"),fname="qmckl_compute_jastrow_champ_factor_een_deriv_e" )
|
||||||
|
|
||||||
|
#+RESULTS:
|
||||||
|
#+begin_src c :tangle (eval h_private_func) :comments org
|
||||||
|
qmckl_exit_code qmckl_compute_jastrow_champ_factor_een_deriv_e (
|
||||||
|
const qmckl_context context,
|
||||||
|
const int64_t walk_num,
|
||||||
|
const int64_t elec_num,
|
||||||
|
const int64_t nucl_num,
|
||||||
|
const int64_t cord_num,
|
||||||
|
const int64_t dim_c_vector,
|
||||||
|
const double* c_vector_full,
|
||||||
|
const int64_t* lkpm_combined_index,
|
||||||
|
const double* tmp_c,
|
||||||
|
const double* dtmp_c,
|
||||||
|
const double* een_rescaled_n,
|
||||||
|
const double* een_rescaled_n_deriv_e,
|
||||||
|
double* const factor_een_deriv_e );
|
||||||
|
#+end_src
|
||||||
|
|
||||||
|
#+begin_src c :tangle (eval c) :comments org
|
||||||
|
qmckl_exit_code
|
||||||
|
qmckl_compute_jastrow_champ_factor_een_deriv_e(const qmckl_context context,
|
||||||
|
const int64_t walk_num,
|
||||||
|
const int64_t elec_num,
|
||||||
|
const int64_t nucl_num,
|
||||||
|
const int64_t cord_num,
|
||||||
|
const int64_t dim_c_vector,
|
||||||
|
const double *c_vector_full,
|
||||||
|
const int64_t *lkpm_combined_index,
|
||||||
|
const double *tmp_c,
|
||||||
|
const double *dtmp_c,
|
||||||
|
const double *een_rescaled_n,
|
||||||
|
const double *een_rescaled_n_deriv_e,
|
||||||
|
double *factor_een_deriv_e)
|
||||||
|
{
|
||||||
|
#ifdef HAVE_HPC
|
||||||
|
return qmckl_compute_jastrow_champ_factor_een_deriv_e_hpc(context, walk_num, elec_num, nucl_num,
|
||||||
|
cord_num, dim_c_vector, c_vector_full,
|
||||||
|
lkpm_combined_index, tmp_c, dtmp_c,
|
||||||
|
een_rescaled_n, een_rescaled_n_deriv_e,
|
||||||
|
factor_een_deriv_e);
|
||||||
|
#else
|
||||||
|
return qmckl_compute_jastrow_champ_factor_een_deriv_e_doc(context, walk_num, elec_num, nucl_num,
|
||||||
|
cord_num, dim_c_vector, c_vector_full,
|
||||||
|
lkpm_combined_index, tmp_c, dtmp_c,
|
||||||
|
een_rescaled_n, een_rescaled_n_deriv_e,
|
||||||
|
factor_een_deriv_e);
|
||||||
|
#endif
|
||||||
|
}
|
||||||
|
#+end_src
|
||||||
|
***** HPC implementation :noexport:
|
||||||
|
#+CALL: generate_private_c_header(table=qmckl_factor_een_deriv_e_args,rettyp=get_value("CRetType"),fname="qmckl_compute_jastrow_champ_factor_een_deriv_e_hpc" )
|
||||||
|
|
||||||
|
#+RESULTS:
|
||||||
|
#+begin_src c :tangle (eval h_private_func) :comments org
|
||||||
|
qmckl_exit_code
|
||||||
|
qmckl_compute_jastrow_champ_factor_een_deriv_e_hpc (
|
||||||
|
const qmckl_context context,
|
||||||
|
const int64_t walk_num,
|
||||||
|
const int64_t elec_num,
|
||||||
|
const int64_t nucl_num,
|
||||||
|
const int64_t cord_num,
|
||||||
|
const int64_t dim_c_vector,
|
||||||
|
const double* c_vector_full,
|
||||||
|
const int64_t* lkpm_combined_index,
|
||||||
|
const double* tmp_c,
|
||||||
|
const double* dtmp_c,
|
||||||
|
const double* een_rescaled_n,
|
||||||
|
const double* een_rescaled_n_deriv_e,
|
||||||
|
double* const factor_een_deriv_e );
|
||||||
|
#+end_src
|
||||||
|
|
||||||
|
#+begin_src c :tangle (eval c) :comments org
|
||||||
|
qmckl_exit_code
|
||||||
|
qmckl_compute_jastrow_champ_factor_een_deriv_e_hpc(const qmckl_context context,
|
||||||
|
const int64_t walk_num,
|
||||||
|
const int64_t elec_num,
|
||||||
|
const int64_t nucl_num,
|
||||||
|
const int64_t cord_num,
|
||||||
|
const int64_t dim_c_vector,
|
||||||
|
const double *c_vector_full,
|
||||||
|
const int64_t *lkpm_combined_index,
|
||||||
|
const double *tmp_c,
|
||||||
|
const double *dtmp_c,
|
||||||
|
const double *een_rescaled_n,
|
||||||
|
const double *een_rescaled_n_deriv_e,
|
||||||
|
double *factor_een_deriv_e)
|
||||||
|
{
|
||||||
|
|
||||||
|
int64_t info = QMCKL_SUCCESS;
|
||||||
|
qmckl_exit_code rc = QMCKL_SUCCESS;
|
||||||
|
|
||||||
|
if (context == QMCKL_NULL_CONTEXT) return QMCKL_INVALID_CONTEXT;
|
||||||
|
if (walk_num <= 0) return QMCKL_INVALID_ARG_2;
|
||||||
|
if (elec_num <= 0) return QMCKL_INVALID_ARG_3;
|
||||||
|
if (nucl_num <= 0) return QMCKL_INVALID_ARG_4;
|
||||||
|
if (cord_num < 0) return QMCKL_INVALID_ARG_5;
|
||||||
|
|
||||||
|
qmckl_matrix c_vector_full_ = qmckl_matrix_alloc(context, nucl_num, dim_c_vector);
|
||||||
|
rc = qmckl_matrix_of_double(context, c_vector_full, nucl_num*dim_c_vector, &c_vector_full_);
|
||||||
|
if (rc != QMCKL_SUCCESS) return rc;
|
||||||
|
|
||||||
|
const int64_t size_tmp_c[5] = {elec_num, nucl_num, cord_num+1, cord_num, walk_num};
|
||||||
|
qmckl_tensor tmp_c_ = qmckl_tensor_alloc(context, 5, &(size_tmp_c[0]));
|
||||||
|
rc = qmckl_tensor_of_double(context, tmp_c, (elec_num*nucl_num*(cord_num+1)*cord_num*walk_num), &tmp_c_);
|
||||||
|
if (rc != QMCKL_SUCCESS) return rc;
|
||||||
|
|
||||||
|
const int64_t size_dtmp_c[6] = {elec_num, 4, nucl_num, cord_num+1, cord_num, walk_num};
|
||||||
|
qmckl_tensor dtmp_c_ = qmckl_tensor_alloc(context, 6, &(size_dtmp_c[0]));
|
||||||
|
rc = qmckl_tensor_of_double(context, dtmp_c, (elec_num*4*nucl_num*(cord_num+1)*cord_num*walk_num), &dtmp_c_);
|
||||||
|
if (rc != QMCKL_SUCCESS) return rc;
|
||||||
|
|
||||||
|
const int64_t size_een_rescaled_n[4] = {elec_num, nucl_num, cord_num, walk_num};
|
||||||
|
qmckl_tensor een_rescaled_n_ = qmckl_tensor_alloc(context, 4, &(size_een_rescaled_n[0]));
|
||||||
|
rc = qmckl_tensor_of_double(context, een_rescaled_n, (elec_num*nucl_num*cord_num*walk_num), &een_rescaled_n_);
|
||||||
|
if (rc != QMCKL_SUCCESS) return rc;
|
||||||
|
|
||||||
|
const int64_t size_een_rescaled_n_deriv_e[5] = {elec_num, 4, nucl_num, cord_num, walk_num};
|
||||||
|
qmckl_tensor een_rescaled_n_deriv_e_ = qmckl_tensor_alloc(context, 5, &(size_een_rescaled_n_deriv_e[0]));
|
||||||
|
rc = qmckl_tensor_of_double(context, een_rescaled_n_deriv_e, (elec_num*4*nucl_num*cord_num*walk_num), &een_rescaled_n_deriv_e_);
|
||||||
|
if (rc != QMCKL_SUCCESS) return rc;
|
||||||
|
|
||||||
|
const int64_t size_factor_een_deriv_e[3] = {elec_num, 4, walk_num};
|
||||||
|
qmckl_tensor factor_een_deriv_e_ = qmckl_tensor_alloc(context, 3, &(size_factor_een_deriv_e[0]));
|
||||||
|
factor_een_deriv_e_ = qmckl_tensor_set(factor_een_deriv_e_,0.);
|
||||||
|
|
||||||
|
|
||||||
|
#ifdef HAVE_OPENMP
|
||||||
|
#pragma omp parallel for
|
||||||
|
#endif
|
||||||
|
for (int64_t nw = 0; nw < walk_num; ++nw) {
|
||||||
|
for (int64_t n = 0; n < dim_c_vector; ++n) {
|
||||||
|
int64_t l = lkpm_combined_index[n];
|
||||||
|
int64_t k = lkpm_combined_index[n+ dim_c_vector];
|
||||||
|
int64_t m = lkpm_combined_index[n+3*dim_c_vector];
|
||||||
|
|
||||||
|
for (int64_t a = 0; a < nucl_num; a++) {
|
||||||
|
double cn = qmckl_mat(c_vector_full_, a, n);
|
||||||
|
if (cn == 0.0) continue;
|
||||||
|
|
||||||
|
for (int64_t ii = 0; ii < 4; ++ii) {
|
||||||
|
for (int64_t j = 0; j < elec_num; ++j) {
|
||||||
|
qmckl_ten3(factor_een_deriv_e_,j,ii,nw) += cn *
|
||||||
|
(qmckl_ten5(tmp_c_,j,a,m,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,ii,a,m+l,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,ii,a,m ,k,nw) * qmckl_ten4(een_rescaled_n_,j,a,m+l,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,ii,a,m+l,k,nw) * qmckl_ten4(een_rescaled_n_,j,a,m ,nw) +
|
||||||
|
qmckl_ten5(tmp_c_,j,a,m+l,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,ii,a,m,nw));
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
cn += cn;
|
||||||
|
|
||||||
|
for (int64_t j = 0; j < elec_num; j++) {
|
||||||
|
qmckl_ten3(factor_een_deriv_e_,j,3,nw) += cn *
|
||||||
|
(qmckl_ten6(dtmp_c_,j,0,a,m ,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,0,a,m+l,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,1,a,m ,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,1,a,m+l,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,2,a,m ,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,2,a,m+l,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,0,a,m+l,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,0,a,m ,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,1,a,m+l,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,1,a,m ,nw) +
|
||||||
|
qmckl_ten6(dtmp_c_,j,2,a,m+l,k,nw) * qmckl_ten5(een_rescaled_n_deriv_e_,j,2,a,m ,nw));
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
info = qmckl_double_of_tensor(context, factor_een_deriv_e_, factor_een_deriv_e, (4*elec_num*walk_num));
|
||||||
|
qmckl_matrix_free(context, &c_vector_full_);
|
||||||
|
qmckl_tensor_free(context, &tmp_c_);
|
||||||
|
qmckl_tensor_free(context, &dtmp_c_);
|
||||||
|
qmckl_tensor_free(context, &een_rescaled_n_);
|
||||||
|
qmckl_tensor_free(context, &een_rescaled_n_deriv_e_);
|
||||||
|
return info;
|
||||||
|
}
|
||||||
|
#+end_src
|
||||||
**** Test
|
**** Test
|
||||||
#+begin_src python :results output :exports none :noweb yes
|
#+begin_src python :results output :exports none :noweb yes
|
||||||
import numpy as np
|
import numpy as np
|
||||||
@ -9558,7 +9731,6 @@ print(een_rescaled_e_deriv_e_t.shape)
|
|||||||
for n in range(0, dim_c_vector):
|
for n in range(0, dim_c_vector):
|
||||||
l = lkpm_of_cindex[0,n]
|
l = lkpm_of_cindex[0,n]
|
||||||
k = lkpm_of_cindex[1,n]
|
k = lkpm_of_cindex[1,n]
|
||||||
p = lkpm_of_cindex[2,n]
|
|
||||||
m = lkpm_of_cindex[3,n]
|
m = lkpm_of_cindex[3,n]
|
||||||
|
|
||||||
for a in range(0, nucl_num):
|
for a in range(0, nucl_num):
|
||||||
@ -9596,7 +9768,17 @@ assert(qmckl_jastrow_champ_provided(context));
|
|||||||
double factor_een_deriv_e[4][walk_num][elec_num];
|
double factor_een_deriv_e[4][walk_num][elec_num];
|
||||||
rc = qmckl_get_jastrow_champ_factor_een_deriv_e(context, &(factor_een_deriv_e[0][0][0]),4*walk_num*elec_num);
|
rc = qmckl_get_jastrow_champ_factor_een_deriv_e(context, &(factor_een_deriv_e[0][0][0]),4*walk_num*elec_num);
|
||||||
|
|
||||||
assert(fabs(factor_een_deriv_e[0][0][0] + 0.0005481671107226865) < 1e-12);
|
printf("%20.15e\n", factor_een_deriv_e[0][0][0]);
|
||||||
|
assert(fabs(factor_een_deriv_e[0][0][0] - (-5.481671107220383e-04)) < 1e-12);
|
||||||
|
|
||||||
|
printf("%20.15e\n", factor_een_deriv_e[1][0][1]);
|
||||||
|
assert(fabs(factor_een_deriv_e[1][0][1] - (-5.402107832095666e-02)) < 1e-12);
|
||||||
|
|
||||||
|
printf("%20.15e\n", factor_een_deriv_e[2][0][2]);
|
||||||
|
assert(fabs(factor_een_deriv_e[2][0][2] - (-1.648945927082279e-01)) < 1e-12);
|
||||||
|
|
||||||
|
printf("%20.15e\n", factor_een_deriv_e[3][0][3]);
|
||||||
|
assert(fabs(factor_een_deriv_e[3][0][3] - (-1.269746119491287e+00)) < 1e-12);
|
||||||
#+end_src
|
#+end_src
|
||||||
|
|
||||||
** Total Jastrow
|
** Total Jastrow
|
||||||
|
Loading…
Reference in New Issue
Block a user