1
0
mirror of https://github.com/TREX-CoE/qmckl.git synced 2024-06-26 15:12:24 +02:00

Fixed Jastrow GL for spin-independent Jastrow

This commit is contained in:
Anthony Scemama 2023-09-27 15:56:37 +02:00
parent 5118359099
commit c70b7b246b
2 changed files with 166 additions and 160 deletions

View File

@ -1430,11 +1430,7 @@ function qmckl_distance_rescaled_gl(context, transa, transb, m, n, &
x = A(1,i) - B(1,j)
y = A(2,i) - B(2,j)
z = A(3,i) - B(3,j)
dist = dsqrt(x*x + y*y + z*z)
! Avoid floating-point exception
if (dist == 0.d0) then
dist = 69.d0/rescale_factor_kappa
endif
dist = max(1.d-20, dsqrt(x*x + y*y + z*z))
dist_inv = 1.0d0/dist
rij = dexp(-rescale_factor_kappa * dist)
C(1,i,j) = x * dist_inv * rij
@ -1451,11 +1447,7 @@ function qmckl_distance_rescaled_gl(context, transa, transb, m, n, &
x = A(i,1) - B(1,j)
y = A(i,2) - B(2,j)
z = A(i,3) - B(3,j)
dist = dsqrt(x*x + y*y + z*z)
! Avoid floating-point exception
if (dist == 0.d0) then
dist = 69.d0/rescale_factor_kappa
endif
dist = max(1.d-20, dsqrt(x*x + y*y + z*z))
dist_inv = 1.0d0/dist
rij = dexp(-rescale_factor_kappa * dist)
C(1,i,j) = x * dist_inv * rij
@ -1472,11 +1464,7 @@ function qmckl_distance_rescaled_gl(context, transa, transb, m, n, &
x = A(1,i) - B(j,1)
y = A(2,i) - B(j,2)
z = A(3,i) - B(j,3)
dist = dsqrt(x*x + y*y + z*z)
! Avoid floating-point exception
if (dist == 0.d0) then
dist = 69.d0/rescale_factor_kappa
endif
dist = max(1.d-20, dsqrt(x*x + y*y + z*z))
dist_inv = 1.0d0/dist
rij = dexp(-rescale_factor_kappa * dist)
C(1,i,j) = x * dist_inv * rij
@ -1493,11 +1481,7 @@ function qmckl_distance_rescaled_gl(context, transa, transb, m, n, &
x = A(i,1) - B(j,1)
y = A(i,2) - B(j,2)
z = A(i,3) - B(j,3)
dist = dsqrt(x*x + y*y + z*z)
! Avoid floating-point exception
if (dist == 0.d0) then
dist = 69.d0/rescale_factor_kappa
endif
dist = max(1.d-20, dsqrt(x*x + y*y + z*z))
dist_inv = 1.0d0/dist
rij = dexp(-rescale_factor_kappa * dist)
C(1,i,j) = x * dist_inv * rij

View File

@ -149,7 +149,7 @@ int main() {
| ~b_vector~ | ~double[bord_num + 1]~ | b polynomial coefficients |
| ~c_vector~ | ~double[dim_c_vector][type_nucl_num]~ | c polynomial coefficients |
| ~c_vector~ | ~double[dim_c_vector][type_nucl_num]~ | c polynomial coefficients |
| ~spin_independent~ | ~bool~ | If true, use same parameters for parallel and anti-parallel pairs |
| ~spin_independent~ | ~int32_t~ | If 1, use same parameters for parallel and anti-parallel spins. Otherwise, 0. |
Computed data:
@ -371,7 +371,7 @@ typedef struct qmckl_jastrow_champ_struct{
uint64_t gl_date;
double rescale_factor_ee;
int32_t uninitialized;
bool spin_independent;
int32_t spin_independent;
bool provided;
} qmckl_jastrow_champ_struct;
@ -406,7 +406,7 @@ qmckl_exit_code qmckl_init_jastrow_champ(qmckl_context context) {
ctx->jastrow_champ.cord_num = -1;
ctx->jastrow_champ.dim_c_vector = -1;
ctx->jastrow_champ.type_nucl_num = -1;
ctx->jastrow_champ.spin_independent = false;
ctx->jastrow_champ.spin_independent = -1;
return QMCKL_SUCCESS;
}
@ -428,7 +428,7 @@ qmckl_exit_code qmckl_set_jastrow_champ_type_nucl_vector (qmckl_context contex
qmckl_exit_code qmckl_set_jastrow_champ_a_vector (qmckl_context context, const double * a_vector, const int64_t size_max);
qmckl_exit_code qmckl_set_jastrow_champ_b_vector (qmckl_context context, const double * b_vector, const int64_t size_max);
qmckl_exit_code qmckl_set_jastrow_champ_c_vector (qmckl_context context, const double * c_vector, const int64_t size_max);
qmckl_exit_code qmckl_set_jastrow_champ_spin_independent (qmckl_context context, const bool spin_independent);
qmckl_exit_code qmckl_set_jastrow_champ_spin_independent (qmckl_context context, const int32_t spin_independent);
#+end_src
#+NAME:pre2
@ -899,7 +899,7 @@ qmckl_set_jastrow_champ_rescale_factor_en(qmckl_context context,
}
qmckl_exit_code
qmckl_set_jastrow_champ_spin_independent(qmckl_context context, const bool spin_independent)
qmckl_set_jastrow_champ_spin_independent(qmckl_context context, const int32_t spin_independent)
{
int32_t mask = 1 << 10;
@ -1071,8 +1071,8 @@ qmckl_exit_code qmckl_finalize_jastrow_champ(qmckl_context context) {
use, intrinsic :: iso_c_binding
import
implicit none
integer (qmckl_context) , intent(in) , value :: context
logical(c_bool), intent(in), value :: spin_independent
integer(qmckl_context) , intent(in) , value :: context
integer(c_int32_t), intent(in), value :: spin_independent
end function qmckl_set_jastrow_champ_spin_independent
end interface
@ -1092,7 +1092,7 @@ qmckl_exit_code qmckl_get_jastrow_champ_c_vector (qmckl_context contex
qmckl_exit_code qmckl_get_jastrow_champ_rescale_factor_ee (const qmckl_context context, double* const rescale_factor_ee);
qmckl_exit_code qmckl_get_jastrow_champ_rescale_factor_en (const qmckl_context context, double* const rescale_factor_en, const int64_t size_max);
qmckl_exit_code qmckl_get_jastrow_champ_dim_c_vector (qmckl_context context, int64_t* const dim_c_vector);
qmckl_exit_code qmckl_get_jastrow_champ_spin_independent (qmckl_context context, bool* const spin_independent);
qmckl_exit_code qmckl_get_jastrow_champ_spin_independent (qmckl_context context, int32_t* const spin_independent);
#+end_src
@ -1500,7 +1500,7 @@ qmckl_exit_code qmckl_get_jastrow_champ_dim_c_vector(qmckl_context context, int6
return QMCKL_SUCCESS;
}
qmckl_exit_code qmckl_get_jastrow_champ_spin_independent(const qmckl_context context, bool* const spin_independent) {
qmckl_exit_code qmckl_get_jastrow_champ_spin_independent(const qmckl_context context, int32_t* const spin_independent) {
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
return qmckl_failwith( context,
@ -1635,8 +1635,8 @@ qmckl_exit_code qmckl_get_jastrow_champ_spin_independent(const qmckl_context con
use, intrinsic :: iso_c_binding
import
implicit none
integer (qmckl_context) , intent(in) , value :: context
logical(c_bool), intent(out) :: spin_independent
integer(qmckl_context) , intent(in) , value :: context
integer(c_int32_t), intent(out) :: spin_independent
end function qmckl_get_jastrow_champ_spin_independent
end interface
@ -1765,7 +1765,7 @@ assert(qmckl_nucleus_provided(context));
electron-electron jastrow factor \(J_{\text{ee}}\). Two values are
computed. The first one is for parallel spin pairs, and the
second one for antiparallel spin pairs.
If the ~spin_independent~ variable is set to ~true~, then
If the ~spin_independent~ variable is set to ~1~, then
$\delta^{\uparrow \downarrow}$ is always equal to one.
\[
@ -1896,7 +1896,7 @@ qmckl_exit_code qmckl_provide_jastrow_champ_asymp_jasb(qmckl_context context)
| ~bord_num~ | ~int64_t~ | in | Order of the polynomial |
| ~b_vector~ | ~double[bord_num+1]~ | in | Values of b |
| ~rescale_factor_ee~ | ~double~ | in | Electron coordinates |
| ~spin_independent~ | ~bool~ | in | If true, same parameters for parallel and anti-parallel pairs |
| ~spin_independent~ | ~int32_t~ | in | If 1, same parameters for parallel and anti-parallel pairs |
| ~asymp_jasb~ | ~double[2]~ | out | Asymptotic value |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
@ -1912,7 +1912,7 @@ function qmckl_compute_jastrow_champ_asymp_jasb_doc(context, &
integer (c_int64_t) , intent(in) , value :: bord_num
real (c_double ) , intent(in) :: b_vector(bord_num+1)
real (c_double ) , intent(in) , value :: rescale_factor_ee
logical (c_bool ) , intent(in) , value :: spin_independent
integer (c_int32_t) , intent(in) , value :: spin_independent
real (c_double ) , intent(out) :: asymp_jasb(2)
integer(qmckl_exit_code) :: info
@ -1933,7 +1933,7 @@ function qmckl_compute_jastrow_champ_asymp_jasb_doc(context, &
endif
asym_one = b_vector(1) * kappa_inv / (1.0d0 + b_vector(2) * kappa_inv)
if (spin_independent) then
if (spin_independent == 1) then
asymp_jasb(:) = (/asym_one, asym_one/)
else
asymp_jasb(:) = (/0.5d0*asym_one, asym_one/)
@ -1955,7 +1955,7 @@ qmckl_exit_code qmckl_compute_jastrow_champ_asymp_jasb_doc (const qmckl_context
const int64_t bord_num,
const double* b_vector,
const double rescale_factor_ee,
const bool spin_independent,
const int32_t spin_independent,
double* const asymp_jasb);
#+end_src
@ -1965,7 +1965,7 @@ qmckl_exit_code qmckl_compute_jastrow_champ_asymp_jasb_hpc (const qmckl_context
const int64_t bord_num,
const double* b_vector,
const double rescale_factor_ee,
const bool spin_independent,
const int32_t spin_independent,
double* const asymp_jasb );
#+end_src
@ -1975,7 +1975,7 @@ qmckl_compute_jastrow_champ_asymp_jasb_hpc (const qmckl_context context,
const int64_t bord_num,
const double* b_vector,
const double rescale_factor_ee,
const bool spin_independent,
const int32_t spin_independent,
double* const asymp_jasb )
{
@ -1997,7 +1997,7 @@ qmckl_compute_jastrow_champ_asymp_jasb_hpc (const qmckl_context context,
f += b_vector[k]*x;
}
asymp_jasb[0] = spin_independent ? asym_one + f : 0.5 * asym_one + f;
asymp_jasb[0] = spin_independent == 1 ? asym_one + f : 0.5 * asym_one + f;
asymp_jasb[1] = asym_one + f;
return QMCKL_SUCCESS;
@ -2010,7 +2010,7 @@ qmckl_exit_code qmckl_compute_jastrow_champ_asymp_jasb (const qmckl_context cont
const int64_t bord_num,
const double* b_vector,
const double rescale_factor_ee,
const bool spin_independent,
const int32_t spin_independent,
double* const asymp_jasb );
#+end_src
@ -2019,7 +2019,7 @@ qmckl_exit_code qmckl_compute_jastrow_champ_asymp_jasb (const qmckl_context cont
const int64_t bord_num,
const double* b_vector,
const double rescale_factor_ee,
const bool spin_independent,
const int32_t spin_independent,
double* const asymp_jasb )
{
#ifdef HAVE_HPC
@ -2073,7 +2073,7 @@ assert(!qmckl_jastrow_champ_provided(context));
/* Set the data */
rc = qmckl_check(context,
qmckl_set_jastrow_champ_spin_independent(context, false)
qmckl_set_jastrow_champ_spin_independent(context, 0)
);
rc = qmckl_check(context,
@ -2155,7 +2155,7 @@ assert(fabs(asymp_jasb[1]-1.043287918508297 ) < 1.e-12);
Calculate the electron-electron jastrow component ~factor_ee~ using the ~asymp_jasb~
component and the electron-electron rescaled distances ~ee_distance_rescaled~.
If the ~spin_independent~ variable is set to ~true~, then
If the ~spin_independent~ variable is set to ~1~, then
$\delta^{\uparrow \downarrow}$ is always equal to one.
\[
@ -2360,7 +2360,7 @@ function qmckl_compute_jastrow_champ_factor_ee_doc(context, &
real (c_double ) , intent(in) :: b_vector(bord_num+1)
real (c_double ) , intent(in) :: ee_distance_rescaled(elec_num,elec_num,walk_num)
real (c_double ) , intent(in) :: asymp_jasb(2)
logical (c_bool ) , intent(in), value :: spin_independent
integer (c_int32_t) , intent(in), value :: spin_independent
real (c_double ) , intent(out) :: factor_ee(walk_num)
integer(qmckl_exit_code) :: info
@ -2397,7 +2397,7 @@ function qmckl_compute_jastrow_champ_factor_ee_doc(context, &
do j=1,elec_num
do i=1,j-1
x = ee_distance_rescaled(i,j,nw)
if (spin_independent) then
if (spin_independent == 1) then
factor_ee(nw) = factor_ee(nw) + b_vector(1) * x / (1.d0 + b_vector(2) * x) - asymp_jasb(2)
else
if ( (j <= up_num).or.(i > up_num) ) then
@ -2430,7 +2430,7 @@ qmckl_compute_jastrow_champ_factor_ee_doc (const qmckl_context context,
const double* b_vector,
const double* ee_distance_rescaled,
const double* asymp_jasb,
const bool spin_independent,
const int32_t spin_independent,
double* const factor_ee );
#+end_src
@ -2444,7 +2444,7 @@ qmckl_compute_jastrow_champ_factor_ee_hpc (const qmckl_context context,
const double* b_vector,
const double* ee_distance_rescaled,
const double* asymp_jasb,
const bool spin_independent,
const int32_t spin_independent,
double* const factor_ee );
#+end_src
@ -2458,7 +2458,7 @@ qmckl_compute_jastrow_champ_factor_ee_hpc (const qmckl_context context,
const double* b_vector,
const double* ee_distance_rescaled,
const double* asymp_jasb,
const bool spin_independent,
const int32_t spin_independent,
double* const factor_ee )
{
@ -2486,7 +2486,7 @@ qmckl_compute_jastrow_champ_factor_ee_hpc (const qmckl_context context,
size_t ishift = nw * elec_num * elec_num;
if (spin_independent) {
if (spin_independent == 1) {
for (int j = 0; j < elec_num; ++j ) {
const double* xj = &(ee_distance_rescaled[j * elec_num + ishift]);
@ -2496,7 +2496,7 @@ qmckl_compute_jastrow_champ_factor_ee_hpc (const qmckl_context context,
}
} else {
for (int j = 0; j < up_num; ++j ) {
const double* xj = &(ee_distance_rescaled[j * elec_num + ishift]);
for (int i = 0; i < j ; ++i) {
@ -2548,7 +2548,7 @@ qmckl_compute_jastrow_champ_factor_ee (const qmckl_context context,
const double* b_vector,
const double* ee_distance_rescaled,
const double* asymp_jasb,
const bool spin_independent,
const int32_t spin_independent,
double* const factor_ee );
#+end_src
@ -2562,7 +2562,7 @@ qmckl_compute_jastrow_champ_factor_ee (const qmckl_context context,
const double* b_vector,
const double* ee_distance_rescaled,
const double* asymp_jasb,
const bool spin_independent,
const int32_t spin_independent,
double* const factor_ee )
{
@ -2772,14 +2772,15 @@ qmckl_exit_code qmckl_provide_jastrow_champ_factor_ee_gl(qmckl_context context)
}
rc = qmckl_compute_jastrow_champ_factor_ee_gl(context,
ctx->electron.walker.num,
ctx->electron.num,
ctx->electron.up_num,
ctx->jastrow_champ.bord_num,
ctx->jastrow_champ.b_vector,
ctx->jastrow_champ.ee_distance_rescaled,
ctx->jastrow_champ.ee_distance_rescaled_gl,
ctx->jastrow_champ.factor_ee_gl);
ctx->electron.walker.num,
ctx->electron.num,
ctx->electron.up_num,
ctx->jastrow_champ.bord_num,
ctx->jastrow_champ.b_vector,
ctx->jastrow_champ.ee_distance_rescaled,
ctx->jastrow_champ.ee_distance_rescaled_gl,
ctx->jastrow_champ.spin_independent,
ctx->jastrow_champ.factor_ee_gl);
if (rc != QMCKL_SUCCESS) {
return rc;
}
@ -2799,23 +2800,24 @@ qmckl_exit_code qmckl_provide_jastrow_champ_factor_ee_gl(qmckl_context context)
:END:
#+NAME: qmckl_factor_ee_gl_args
| Variable | Type | In/Out | Description |
|---------------------------+-------------------------------------------+--------+-----------------------------|
| ~context~ | ~qmckl_context~ | in | Global state |
| ~walk_num~ | ~int64_t~ | in | Number of walkers |
| ~elec_num~ | ~int64_t~ | in | Number of electrons |
| ~up_num~ | ~int64_t~ | in | Number of alpha electrons |
| ~bord_num~ | ~int64_t~ | in | Number of coefficients |
| ~b_vector~ | ~double[bord_num+1]~ | in | List of coefficients |
| ~ee_distance_rescaled~ | ~double[walk_num][elec_num][elec_num]~ | in | Electron-electron distances |
| ~ee_distance_rescaled_gl~ | ~double[walk_num][4][elec_num][elec_num]~ | in | Electron-electron distances |
| ~factor_ee_gl~ | ~double[walk_num][4][elec_num]~ | out | Electron-electron distances |
| Variable | Type | In/Out | Description |
|---------------------------+-------------------------------------------+--------+-----------------------------------------------------------|
| ~context~ | ~qmckl_context~ | in | Global state |
| ~walk_num~ | ~int64_t~ | in | Number of walkers |
| ~elec_num~ | ~int64_t~ | in | Number of electrons |
| ~up_num~ | ~int64_t~ | in | Number of alpha electrons |
| ~bord_num~ | ~int64_t~ | in | Number of coefficients |
| ~b_vector~ | ~double[bord_num+1]~ | in | List of coefficients |
| ~ee_distance_rescaled~ | ~double[walk_num][elec_num][elec_num]~ | in | Electron-electron distances |
| ~ee_distance_rescaled_gl~ | ~double[walk_num][4][elec_num][elec_num]~ | in | Electron-electron distances |
| ~spin_independent~ | ~int32_t~ | in | If 1, same parameters for parallel and antiparallel spins |
| ~factor_ee_gl~ | ~double[walk_num][4][elec_num]~ | out | Electron-electron distances |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
function qmckl_compute_jastrow_champ_factor_ee_gl_doc( &
context, walk_num, elec_num, up_num, bord_num, &
b_vector, ee_distance_rescaled, ee_distance_rescaled_gl, &
factor_ee_gl) &
spin_independent, factor_ee_gl) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
@ -2830,6 +2832,7 @@ function qmckl_compute_jastrow_champ_factor_ee_gl_doc( &
real (c_double ) , intent(in) :: b_vector(bord_num+1)
real (c_double ) , intent(in) :: ee_distance_rescaled(elec_num,elec_num,walk_num)
real (c_double ) , intent(in) :: ee_distance_rescaled_gl(4,elec_num,elec_num,walk_num)
integer (c_int32_t) , intent(in) , value :: spin_independent
real (c_double ) , intent(out) :: factor_ee_gl(elec_num,4,walk_num)
integer(qmckl_exit_code) :: info
@ -2861,6 +2864,11 @@ function qmckl_compute_jastrow_champ_factor_ee_gl_doc( &
return
endif
if ((spin_independent < 0).or.(spin_independent > 1)) then
info = QMCKL_INVALID_ARG_8
return
endif
do nw =1, walk_num
factor_ee_gl(:,:,nw) = 0.0d0
@ -2881,10 +2889,14 @@ function qmckl_compute_jastrow_champ_factor_ee_gl_doc( &
grad_c2 = dx(1)*dx(1) + dx(2)*dx(2) + dx(3)*dx(3)
if((i <= up_num .and. j <= up_num ) .or. (i > up_num .and. j > up_num)) then
f = 0.5d0 * b_vector(1) * invdenom2
else
if (spin_independent == 1) then
f = b_vector(1) * invdenom2
else
if((i <= up_num .and. j <= up_num ) .or. (i > up_num .and. j > up_num)) then
f = 0.5d0 * b_vector(1) * invdenom2
else
f = b_vector(1) * invdenom2
end if
end if
factor_ee_gl(i,1,nw) = factor_ee_gl(i,1,nw) + f * dx(1)
@ -2917,16 +2929,18 @@ end function qmckl_compute_jastrow_champ_factor_ee_gl_doc
#+end_src
#+begin_src c :comments org :tangle (eval c) :noweb yes
qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl_hpc(
const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
double* const factor_ee_gl ) {
qmckl_exit_code
qmckl_compute_jastrow_champ_factor_ee_gl_hpc(const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
const int32_t spin_independent,
double* const factor_ee_gl )
{
if (context == QMCKL_NULL_CONTEXT) return QMCKL_INVALID_CONTEXT;
if (walk_num <= 0) return QMCKL_INVALID_ARG_2;
@ -2936,14 +2950,15 @@ qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl_hpc(
if (b_vector == NULL) return QMCKL_INVALID_ARG_6;
if (ee_distance_rescaled == NULL) return QMCKL_INVALID_ARG_7;
if (ee_distance_rescaled_gl == NULL) return QMCKL_INVALID_ARG_8;
if (spin_independent & (int32_t) (-2)) return QMCKL_INVALID_ARG_8;
if (factor_ee_gl == NULL) return QMCKL_INVALID_ARG_9;
double kf[bord_num+1];
for (int k=0 ; k<=bord_num ; ++k) {
kf[k] = (double) k;
kf[k] = (double) k;
}
#pragma omp parallel for
#pragma omp parallel for
for (int nw = 0; nw < walk_num; ++nw) {
memset(&(factor_ee_gl[nw*4*elec_num]), 0, elec_num*4*sizeof(double));
@ -2969,9 +2984,12 @@ qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl_hpc(
const double grad_c2 = dx[0]*dx[0] + dx[1]*dx[1] + dx[2]*dx[2];
double f =
(i < up_num && j < up_num ) || (i >= up_num && j >= up_num) ?
0.5 * b_vector[0] * invdenom2 : b_vector[0] * invdenom2;
double f = b_vector[0] * invdenom2;
if ((spin_independent == 0) && (
((i < up_num) && (j < up_num)) ||
((i >= up_num) && (j >= up_num))) ) {
f *= 0.5;
}
factor_ee_gl_0[i] += f*dx[0];
factor_ee_gl_1[i] += f*dx[1];
@ -3007,66 +3025,70 @@ qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl_hpc(
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl (
const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
double* const factor_ee_gl );
qmckl_exit_code
qmckl_compute_jastrow_champ_factor_ee_gl (const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
const int32_t spin_independent,
double* const factor_ee_gl );
#+end_src
#+begin_src c :tangle (eval h_private_func) :comments org
qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl_hpc (
const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
double* const factor_ee_gl );
qmckl_exit_code
qmckl_compute_jastrow_champ_factor_ee_gl_hpc (const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
const int32_t spin_independent,
double* const factor_ee_gl );
#+end_src
#+begin_src c :tangle (eval h_private_func) :comments org
qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl_doc (
const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
double* const factor_ee_gl );
qmckl_exit_code
qmckl_compute_jastrow_champ_factor_ee_gl_doc (const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
const int32_t spin_independent,
double* const factor_ee_gl );
#+end_src
#+begin_src c :comments org :tangle (eval c) :noweb yes
qmckl_exit_code qmckl_compute_jastrow_champ_factor_ee_gl (
const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
double* const factor_ee_gl ) {
qmckl_exit_code
qmckl_compute_jastrow_champ_factor_ee_gl (const qmckl_context context,
const int64_t walk_num,
const int64_t elec_num,
const int64_t up_num,
const int64_t bord_num,
const double* b_vector,
const double* ee_distance_rescaled,
const double* ee_distance_rescaled_gl,
const int32_t spin_independent,
double* const factor_ee_gl )
{
#ifdef HAVE_HPC
return qmckl_compute_jastrow_champ_factor_ee_gl_hpc
return qmckl_compute_jastrow_champ_factor_ee_gl_hpc
#else
return qmckl_compute_jastrow_champ_factor_ee_gl_doc
return qmckl_compute_jastrow_champ_factor_ee_gl_doc
#endif
(context, walk_num, elec_num, up_num, bord_num, b_vector,
ee_distance_rescaled, ee_distance_rescaled_gl, factor_ee_gl );
ee_distance_rescaled, ee_distance_rescaled_gl, spin_independent, factor_ee_gl );
}
#+end_src
@ -3446,7 +3468,7 @@ qmckl_exit_code qmckl_compute_ee_distance_rescaled_hpc (
int64_t sze = elec_num*walk_num;
int64_t elec_num2= elec_num*elec_num;
qmckl_exit_code result = QMCKL_SUCCESS;
#pragma omp parallel
@ -3455,7 +3477,7 @@ qmckl_exit_code qmckl_compute_ee_distance_rescaled_hpc (
#pragma omp for
for (int64_t k=0 ; k<walk_num ; ++k)
{
rc |= qmckl_distance_rescaled(context, 'T', 'T', elec_num, elec_num,
rc |= qmckl_distance_rescaled(context, 'T', 'T', elec_num, elec_num,
&(coord[k*elec_num]), sze, &(coord[k*elec_num]), sze,
&(ee_distance_rescaled[k*elec_num2]), elec_num, rescale_factor_ee);
}
@ -3474,7 +3496,7 @@ qmckl_exit_code qmckl_compute_ee_distance_rescaled (
const int64_t walk_num,
const double* coord,
double* const ee_distance_rescaled )
{
{
#ifdef HAVE_HPC
return qmckl_compute_ee_distance_rescaled_hpc
#else
@ -3676,7 +3698,7 @@ function qmckl_compute_ee_distance_rescaled_gl_doc(context, &
integer (c_int64_t) , intent(in) , value :: walk_num
real (c_double ) , intent(in) :: coord(elec_num,walk_num,3)
real (c_double ) , intent(out) :: ee_distance_rescaled_gl(4,elec_num,elec_num,walk_num)
integer(qmckl_exit_code) :: info
integer(qmckl_exit_code) :: info
integer*8 :: k
@ -3762,7 +3784,7 @@ qmckl_exit_code qmckl_compute_ee_distance_rescaled_gl_hpc (
#pragma omp for
for (int64_t k=0 ; k<walk_num ; ++k)
{
rc |= qmckl_distance_rescaled_gl(context, 'T', 'T', elec_num, elec_num,
rc |= qmckl_distance_rescaled_gl(context, 'T', 'T', elec_num, elec_num,
&(coord[k*elec_num]), sze,
&(coord[k*elec_num]), sze,
&(ee_distance_rescaled_gl[k*elec_num2]), elec_num,
@ -4770,7 +4792,7 @@ function qmckl_compute_jastrow_champ_factor_en_gl_doc( &
f = a_vector(1, type_nucl_vector(a)+1) * invdenom2
grad_c2 = dx(1)*dx(1) + dx(2)*dx(2) + dx(3)*dx(3)
factor_en_gl(i,1,nw) = factor_en_gl(i,1,nw) + f * dx(1)
factor_en_gl(i,2,nw) = factor_en_gl(i,2,nw) + f * dx(2)
factor_en_gl(i,3,nw) = factor_en_gl(i,3,nw) + f * dx(3)
@ -4791,11 +4813,11 @@ function qmckl_compute_jastrow_champ_factor_en_gl_doc( &
x = x*x1
kf = kf + 1.d0
end do
end do
end do
end do
end function qmckl_compute_jastrow_champ_factor_en_gl_doc
#+end_src
@ -4949,7 +4971,7 @@ qmckl_compute_jastrow_champ_factor_en_gl (const qmckl_context context,
{
#ifdef HAVE_HPC
return qmckl_compute_jastrow_champ_factor_en_gl_hpc
#else
#else
return qmckl_compute_jastrow_champ_factor_en_gl_doc
#endif
(context, walk_num, elec_num, nucl_num, type_nucl_num, type_nucl_vector, aord_num,
@ -5292,7 +5314,7 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled_hpc (
if (en_distance_rescaled == NULL) return QMCKL_INVALID_ARG_10;
int64_t sze = elec_num*walk_num;
qmckl_exit_code result = QMCKL_SUCCESS;
#pragma omp parallel
{
@ -5302,7 +5324,7 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled_hpc (
{
for (int64_t a=0 ; a<nucl_num ; ++a) {
const double coord[3] = { nucl_coord[a], nucl_coord[a+nucl_num], nucl_coord[a+2*nucl_num] };
rc |= qmckl_distance_rescaled(context, 'T', 'N', elec_num, 1,
rc |= qmckl_distance_rescaled(context, 'T', 'N', elec_num, 1,
&(elec_coord[k*elec_num]), sze,
coord, 3,
&(en_distance_rescaled[elec_num*(a+nucl_num*k)]), elec_num,
@ -5680,7 +5702,7 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled_gl_hpc (
const double* elec_coord,
const double* nucl_coord,
double* const en_distance_rescaled_gl )
{
{
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) return QMCKL_INVALID_CONTEXT;
if (elec_num <= 0) return QMCKL_INVALID_ARG_2;
if (nucl_num <= 0) return QMCKL_INVALID_ARG_3;
@ -5703,7 +5725,7 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled_gl_hpc (
{
for (int64_t a=0 ; a<nucl_num ; ++a) {
const double coord[3] = { nucl_coord[a], nucl_coord[a+nucl_num], nucl_coord[a+2*nucl_num] };
rc |= qmckl_distance_rescaled_gl(context, 'T', 'T', elec_num, 1,
rc |= qmckl_distance_rescaled_gl(context, 'T', 'T', elec_num, 1,
&(elec_coord[k*elec_num]), sze,
coord, 1,
&(en_distance_rescaled_gl[4*elec_num*(a+nucl_num*k)]), elec_num,
@ -5729,7 +5751,7 @@ qmckl_exit_code qmckl_compute_en_distance_rescaled_gl (
const double* elec_coord,
const double* nucl_coord,
double* const en_distance_rescaled_gl )
{
{
#ifdef HAVE_HPC
return qmckl_compute_en_distance_rescaled_gl_hpc
#else
@ -8276,7 +8298,7 @@ qmckl_exit_code qmckl_compute_lkpm_combined_index_hpc (
}
}
}
return QMCKL_SUCCESS;
}
#+end_src
@ -8303,7 +8325,7 @@ qmckl_exit_code qmckl_compute_lkpm_combined_index_hpc (
end function qmckl_compute_lkpm_combined_index_doc
#+end_src
#+begin_src c :comments org :tangle (eval c) :noweb yes
qmckl_exit_code qmckl_compute_lkpm_combined_index (
const qmckl_context context,
@ -8322,7 +8344,7 @@ qmckl_exit_code qmckl_compute_lkpm_combined_index (
#+CALL: generate_c_header(table=lkpm_combined_index_args,rettyp=get_value("CRetType"),fname=get_value("Name"))
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
qmckl_exit_code qmckl_compute_lkpm_combined_index (
const qmckl_context context,
@ -8338,9 +8360,9 @@ qmckl_exit_code qmckl_compute_lkpm_combined_index (
const qmckl_context context,
const int64_t cord_num,
const int64_t dim_c_vector,
int64_t* const lkpm_combined_index );
int64_t* const lkpm_combined_index );
#+end_src
#+CALL: generate_c_header(table=lkpm_combined_index_args,rettyp=get_value("CRetType"),fname="qmckl_compute_lkpm_combined_index_hpc")
#+RESULTS:
@ -8349,7 +8371,7 @@ qmckl_exit_code qmckl_compute_lkpm_combined_index (
const qmckl_context context,
const int64_t cord_num,
const int64_t dim_c_vector,
int64_t* const lkpm_combined_index );
int64_t* const lkpm_combined_index );
#+end_src
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
@ -8360,7 +8382,7 @@ qmckl_exit_code qmckl_compute_lkpm_combined_index (
int64_t* const lkpm_combined_index );
#+end_src
#+CALL: generate_c_header(table=lkpm_combined_index_args,rettyp=get_value("CRetType"),fname=get_value("Name"))
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
qmckl_exit_code qmckl_compute_lkpm_combined_index (
const qmckl_context context,
@ -8662,7 +8684,7 @@ qmckl_compute_dtmp_c (const qmckl_context context,
double* const dtmp_c )
{
#ifdef HAVE_HPC
return qmckl_compute_dtmp_c_hpc
return qmckl_compute_dtmp_c_hpc
#else
return qmckl_compute_dtmp_c_doc
#endif
@ -9185,7 +9207,7 @@ integer function qmckl_compute_jastrow_champ_factor_een_naive_f( &
(een_rescaled_n(i,a,l,nw) + een_rescaled_n(j,a,l,nw)) * &
(een_rescaled_n(i,a,m,nw) * een_rescaled_n(j,a,m,nw))
end do
accu2 = accu2 + accu
accu2 = accu2 + accu
end do
factor_een(nw) = factor_een(nw) + accu2 * cn
end do
@ -9638,7 +9660,7 @@ qmckl_exit_code qmckl_provide_jastrow_champ_factor_een_gl(qmckl_context context)
}
/*
rc = qmckl_compute_jastrow_champ_factor_een_gl_naive (context,
rc = qmckl_compute_jastrow_champ_factor_een_gl_naive (context,
ctx->electron.walker.num,
ctx->electron.num,
ctx->nucleus.num,
@ -11074,9 +11096,9 @@ for (int64_t k=0 ; k< walk_num ; ++k) {
#+end_src
* End of files :noexport:
#+begin_src c :tangle (eval h_private_type)
#endif
#+end_src