1
0
mirror of https://github.com/TREX-CoE/qmckl.git synced 2024-07-11 22:04:04 +02:00

Added AO polynomials

This commit is contained in:
Anthony Scemama 2021-07-11 23:20:40 +02:00
parent fba34e4982
commit 75a93d12c6
3 changed files with 650 additions and 58 deletions

View File

@ -59,6 +59,7 @@ gradients and Laplacian of the atomic basis functions.
#include <stdio.h>
#include <math.h>
#include "chbrclf.h"
#include "qmckl_ao_private_func.h"
int main() {
qmckl_context context;
@ -107,8 +108,9 @@ int main() {
| ~coefficient~ | ~[prim_num]~ | Array of coefficients |
| ~prim_factor~ | ~[prim_num]~ | Normalization factors of the primtives |
| ~ao_num~ | | Number of AOs |
| ~ao_factor~ | ~[ao_num]~ | Normalization factor of the AO |
| ~ao_cartesian~ | | If true, use polynomials. Otherwise, use spherical harmonics |
| ~ao_factor~ | ~[ao_num]~ | Normalization factor of the AO |
| ~ao_shell~ | ~[ao_num]~ | For each AO, specify to which shell it belongs |
Computed data:
@ -121,7 +123,9 @@ int main() {
| ~primitive_vgl~ | ~[5][walk_num][elec_num][prim_num]~ | Value, gradients, Laplacian of the primitives at electron positions |
| ~primitive_vgl_date~ | ~uint64_t~ | Late modification date of Value, gradients, Laplacian of the primitives at electron positions |
| ~shell_vgl~ | ~[5][walk_num][elec_num][shell_num]~ | Value, gradients, Laplacian of the primitives at electron positions |
| ~shell_vgl_date~ | ~uint64_t~ | Late modification date of Value, gradients, Laplacian of the shells at electron positions |
| ~shell_vgl_date~ | ~uint64_t~ | Late modification date of Value, gradients, Laplacian of the AOs at electron positions |
| ~ao_vgl~ | ~[5][walk_num][elec_num][ao_num]~ | Value, gradients, Laplacian of the primitives at electron positions |
| ~ao_vgl_date~ | ~uint64_t~ | Late modification date of Value, gradients, Laplacian of the AOs at electron positions |
|--------------------------+--------------------------------------+-----------------------------------------------------------------------------------------------|
| ~nucl_shell_index~ | ~[nucl_num]~ | Index of the first shell for each nucleus |
| ~exponent_sorted~ | ~[prim_num]~ | Array of exponents for sorted primitives |
@ -156,6 +160,7 @@ D 1
type = 'G'
shell_num = 12
prim_num = 20
ao_num = 38
nucleus_index = [0 , 6]
shell_ang_mom = [0, 0, 0, 1, 1, 2, 0, 0, 0, 1, 1, 2]
shell_factor = [ 1., 1., 1., 1., 1., 1., 1., 1., 1., 1., 1., 1.]
@ -202,6 +207,8 @@ typedef struct qmckl_ao_basis_struct {
int64_t primitive_vgl_date;
double * shell_vgl;
int64_t shell_vgl_date;
double * ao_vgl;
int64_t ao_vgl_date;
int32_t uninitialized;
bool provided;
@ -248,6 +255,7 @@ int64_t qmckl_get_ao_basis_shell_num (const qmckl_context context);
int64_t qmckl_get_ao_basis_prim_num (const qmckl_context context);
int64_t qmckl_get_ao_basis_ao_num (const qmckl_context context);
int64_t* qmckl_get_ao_basis_nucleus_index (const qmckl_context context);
int64_t* qmckl_get_ao_basis_nucleus_shell_num(const qmckl_context context);
int32_t* qmckl_get_ao_basis_shell_ang_mom (const qmckl_context context);
int64_t* qmckl_get_ao_basis_shell_prim_num (const qmckl_context context);
int64_t* qmckl_get_ao_basis_shell_prim_index (const qmckl_context context);
@ -775,7 +783,7 @@ qmckl_exit_code qmckl_set_ao_basis_shell_ang_mom(qmckl_context context, const i
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
mem_info.size = shell_num * sizeof(char);
mem_info.size = shell_num * sizeof(int32_t);
int32_t * new_array = (int32_t*) qmckl_malloc(context, mem_info);
if (new_array == NULL) {
@ -1208,7 +1216,7 @@ qmckl_exit_code qmckl_finalize_basis(qmckl_context context) {
NULL);
}
for (int64_t inucl=0 ; inucl < ctx->nucleus.num ; ++inucl) {
for (int64_t inucl=0 ; inucl < nucl_num ; ++inucl) {
ctx->ao_basis.nucleus_max_ang_mom[inucl] = 0;
for (int64_t ishell=ctx->ao_basis.nucleus_index[inucl] ;
ishell < ctx->ao_basis.nucleus_index[inucl] + ctx->ao_basis.nucleus_shell_num[inucl] ;
@ -1417,7 +1425,7 @@ const double * coefficient = &(chbrclf_basis_coefficient[0]);
const double * prim_factor = &(chbrclf_basis_prim_factor[0]);
const double * ao_factor = &(chbrclf_basis_ao_factor[0]);
char typ = 'G';
const char typ = 'G';
assert(!qmckl_ao_basis_provided(context));
@ -1476,6 +1484,84 @@ assert(rc == QMCKL_SUCCESS);
assert(qmckl_ao_basis_provided(context));
int64_t shell_num_test ;
int64_t prim_num_test ;
int64_t ao_num_test ;
int64_t * nucleus_index_test ;
int64_t * nucleus_shell_num_test;
int32_t * shell_ang_mom_test ;
int64_t * shell_prim_num_test ;
int64_t * shell_prim_index_test ;
double * shell_factor_test ;
double * exponent_test ;
double * coefficient_test ;
double * prim_factor_test ;
double * ao_factor_test ;
char typ_test ;
typ_test = qmckl_get_ao_basis_type (context);
assert(typ == typ_test);
shell_num_test = qmckl_get_ao_basis_shell_num (context);
assert(shell_num == shell_num_test);
prim_num_test = qmckl_get_ao_basis_prim_num (context);
assert(prim_num == prim_num_test);
nucleus_index_test = qmckl_get_ao_basis_nucleus_index (context);
for (int64_t i=0 ; i < nucl_num ; ++i) {
assert(nucleus_index_test[i] == nucleus_index[i]);
}
nucleus_shell_num_test = qmckl_get_ao_basis_nucleus_shell_num (context);
for (int64_t i=0 ; i < nucl_num ; ++i) {
assert(nucleus_shell_num_test[i] == nucleus_shell_num[i]);
}
shell_ang_mom_test = qmckl_get_ao_basis_shell_ang_mom (context);
for (int64_t i=0 ; i < shell_num ; ++i) {
assert(shell_ang_mom_test[i] == shell_ang_mom[i]);
}
shell_factor_test = qmckl_get_ao_basis_shell_factor (context);
for (int64_t i=0 ; i < shell_num ; ++i) {
assert(shell_factor_test[i] == shell_factor[i]);
}
shell_prim_num_test = qmckl_get_ao_basis_shell_prim_num (context);
for (int64_t i=0 ; i < shell_num ; ++i) {
assert(shell_prim_num_test[i] == shell_prim_num[i]);
}
shell_prim_index_test = qmckl_get_ao_basis_shell_prim_index (context);
for (int64_t i=0 ; i < shell_num ; ++i) {
assert(shell_prim_index_test[i] == shell_prim_index[i]);
}
exponent_test = qmckl_get_ao_basis_exponent(context);
for (int64_t i=0 ; i < prim_num ; ++i) {
assert(exponent_test[i] == exponent[i]);
}
coefficient_test = qmckl_get_ao_basis_coefficient(context);
for (int64_t i=0 ; i < prim_num ; ++i) {
assert(coefficient_test[i] == coefficient[i]);
}
prim_factor_test = qmckl_get_ao_basis_prim_factor (context);
for (int64_t i=0 ; i < prim_num ; ++i) {
assert(prim_factor_test[i] == prim_factor[i]);
}
ao_num_test = qmckl_get_ao_basis_ao_num(context);
assert(ao_num == ao_num_test);
ao_factor_test = qmckl_get_ao_basis_ao_factor (context);
for (int64_t i=0 ; i < ao_num ; ++i) {
assert(ao_factor_test[i] == ao_factor[i]);
}
#+end_src
* Radial part
@ -2136,7 +2222,7 @@ qmckl_exit_code qmckl_provide_ao_basis_shell_vgl(qmckl_context context)
if (ctx->ao_basis.shell_vgl == NULL) {
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
mem_info.size = ctx->ao_basis.prim_num * 5 * ctx->electron.num *
mem_info.size = ctx->ao_basis.shell_num * 5 * ctx->electron.num *
ctx->electron.walk_num * sizeof(double);
double* shell_vgl = (double*) qmckl_malloc(context, mem_info);
@ -2296,24 +2382,26 @@ integer function qmckl_compute_ao_basis_shell_gaussian_vgl_f(context, &
end function qmckl_compute_ao_basis_shell_gaussian_vgl_f
#+end_src
# #+CALL: generate_c_header(table=qmckl_ao_basis_shell_gaussian_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_basis_shell_gaussian_vgl"))
#+begin_src c :tangle (eval h_private_func) :comments org :exports none
qmckl_exit_code qmckl_compute_ao_basis_shell_gaussian_vgl(
#+RESULTS:
#+begin_src c :tangle (eval h_private_func) :comments org
qmckl_exit_code qmckl_compute_ao_basis_shell_gaussian_vgl (
const qmckl_context context,
const int64_t prim_num,
const int64_t shell_num,
const int64_t elec_num,
const int64_t nucl_num,
const int64_t walk_num,
const int64_t prim_num,
const int64_t shell_num,
const int64_t elec_num,
const int64_t nucl_num,
const int64_t walk_num,
const int64_t* nucleus_shell_num,
const int64_t* shell_prim_index,
const int64_t* nucleus_index,
const int64_t* shell_prim_index,
const int64_t* shell_prim_num,
const double* elec_coord,
const double* nucl_coord,
const double* expo,
const double* coef_normalized,
double* const shell_vgl);
const double* elec_coord,
const double* nucl_coord,
const double* expo,
const double* coef_normalized,
double* const shell_vgl );
#+end_src
#+CALL: generate_c_interface(table=qmckl_ao_basis_shell_gaussian_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_basis_shell_gaussian_vgl"))
@ -2408,50 +2496,50 @@ nucl_2 = np.array( [ 1.168459237342663e+00, 1.125660720053393e+00, 2.83337031
#double prim_vgl[prim_num][5][walk_num][elec_num];
x = elec_26_w1 ; y = nucl_1
a = [( 8.236000E+03, -1.130000E-04 ),
( 1.235000E+03, -8.780000E-04 ),
( 2.808000E+02, -4.540000E-03 ),
( 7.927000E+01, -1.813300E-02 ),
( 2.559000E+01, -5.576000E-02 ),
( 8.997000E+00, -1.268950E-01 ),
( 3.319000E+00, -1.703520E-01 ),
( 9.059000E-01, 1.403820E-01 ),
( 3.643000E-01, 5.986840E-01 ),
( 1.285000E-01, 3.953890E-01 )]
a = [( 8.236000E+03, -1.130000E-04 * 6.1616545431994848e+02 ),
( 1.235000E+03, -8.780000E-04 * 1.4847738511079908e+02 ),
( 2.808000E+02, -4.540000E-03 * 4.8888635917437597e+01 ),
( 7.927000E+01, -1.813300E-02 * 1.8933972232608955e+01 ),
( 2.559000E+01, -5.576000E-02 * 8.1089160941724145e+00 ),
( 8.997000E+00, -1.268950E-01 * 3.7024003863155635e+00 ),
( 3.319000E+00, -1.703520E-01 * 1.7525302846177560e+00 ),
( 9.059000E-01, 1.403820E-01 * 6.6179013183966806e-01 ),
( 3.643000E-01, 5.986840E-01 * 3.3419848027174592e-01 ),
( 1.285000E-01, 3.953890E-01 * 1.5296336817449557e-01 )]
print ( "[1][0][0][26] : %e"% f(a,x,y))
print ( "[1][1][0][26] : %e"% df(a,x,y,1))
print ( "[1][2][0][26] : %e"% df(a,x,y,2))
print ( "[1][3][0][26] : %e"% df(a,x,y,3))
print ( "[1][4][0][26] : %e"% lf(a,x,y))
print ( "[1][0][0][26] : %25.15e"% f(a,x,y))
print ( "[1][1][0][26] : %25.15e"% df(a,x,y,1))
print ( "[1][2][0][26] : %25.15e"% df(a,x,y,2))
print ( "[1][3][0][26] : %25.15e"% df(a,x,y,3))
print ( "[1][4][0][26] : %25.15e"% lf(a,x,y))
x = elec_15_w2 ; y = nucl_2
a = [(3.387000E+01, 6.068000E-03),
(5.095000E+00, 4.530800E-02),
(1.159000E+00, 2.028220E-01),
(3.258000E-01, 5.039030E-01),
(1.027000E-01, 3.834210E-01)]
a = [(3.387000E+01, 6.068000E-03 *1.0006253235944540e+01),
(5.095000E+00, 4.530800E-02 *2.4169531573445120e+00),
(1.159000E+00, 2.028220E-01 *7.9610924849766440e-01),
(3.258000E-01, 5.039030E-01 *3.0734305383061117e-01),
(1.027000E-01, 3.834210E-01 *1.2929684417481876e-01)]
print ( "[14][0][1][15] : %e"% f(a,x,y))
print ( "[14][1][1][15] : %e"% df(a,x,y,1))
print ( "[14][2][1][15] : %e"% df(a,x,y,2))
print ( "[14][3][1][15] : %e"% df(a,x,y,3))
print ( "[14][4][1][15] : %e"% lf(a,x,y))
print ( "[0][1][15][14] : %25.15e"% f(a,x,y))
print ( "[1][1][15][14] : %25.15e"% df(a,x,y,1))
print ( "[2][1][15][14] : %25.15e"% df(a,x,y,2))
print ( "[3][1][15][14] : %25.15e"% df(a,x,y,3))
print ( "[4][1][15][14] : %25.15e"% lf(a,x,y))
#+end_src
#+RESULTS:
#+begin_example
[1][0][0][26] : 1.875569e-01
[1][1][0][26] : -2.615250e-02
[1][2][0][26] : -1.333535e-01
[1][3][0][26] : 1.218483e-01
[1][4][0][26] : 3.227973e-02
[14][0][1][15] : 4.509748e-02
[14][1][1][15] : 3.203918e-02
[14][2][1][15] : 2.887081e-02
[14][3][1][15] : 5.836910e-03
[14][4][1][15] : 1.564721e-02
[1][0][0][26] : 3.564393437193867e-02
[1][1][0][26] : -6.030177988891605e-03
[1][2][0][26] : -3.074832579871845e-02
[1][3][0][26] : 2.809546963133958e-02
[1][4][0][26] : 1.903338597841753e-02
[0][1][15][14] : 5.928089771361000e-03
[1][1][15][14] : 4.355862298893037e-03
[2][1][15][14] : 3.925108924950765e-03
[3][1][15][14] : 7.935527764416084e-04
[4][1][15][14] : 2.697495005143935e-03
#+end_example
*** Test
@ -2511,6 +2599,22 @@ assert( fabs(shell_vgl[4][1][15][14] - ( 2.708246573703548e-03)) < 1.e-14 );
#+end_src
* Polynomial part
Going from the atomic basis set to AOs implies a systematic
construction of all the angular functions of each shell. We
consider two cases for the angular functions: the real-valued
spherical harmonics, and the polynomials in Cartesian coordinates.
In the case of spherical harmonics, the AOs are ordered in
increasing magnetic quantum number ($-l \le m \le l$), and in the
case of polynomials we choose the canonical ordering, i.e
\begin{eqnarray}
p & : & p_x, p_y, p_z \nonumber \\
d & : & d_{xx}, d_{xy}, d_{xz}, d_{yy}, d_{yz}, d_{zz} \nonumber \\
f & : & f_{xxx}, f_{xxy}, f_{xxz}, f_{xyy}, f_{xyz}, f_{xzz}, f_{yyy}, f_{yyz}, f_{yzz}, f_{zzz} \nonumber \\
{\rm etc.} \nonumber
\end{eqnarray}
** General functions for Powers of $x-X_i$
:PROPERTIES:
:Name: qmckl_ao_power
@ -2814,7 +2918,6 @@ integer function qmckl_ao_polynomial_vgl_f(context, X, R, lmax, n, L, ldl, VGL,
real*8 :: Y(3)
integer :: lmax_array(3)
real*8 :: pows(-2:lmax,3)
integer, external :: qmckl_ao_power_f
double precision :: xy, yz, xz
double precision :: da, db, dc, dd
@ -3084,6 +3187,497 @@ end function test_qmckl_ao_polynomial_vgl
#+end_src
* Combining radial and polynomial parts
*** TODO Get
#+begin_src c :comments org :tangle (eval h_func) :noweb yes
qmckl_exit_code qmckl_get_ao_vgl(qmckl_context context, double* const ao_vgl);
#+end_src
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
qmckl_exit_code qmckl_get_ao_vgl(qmckl_context context, double* const ao_vgl) {
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
return QMCKL_NULL_CONTEXT;
}
qmckl_exit_code rc;
rc = qmckl_provide_ao_vgl(context);
if (rc != QMCKL_SUCCESS) return rc;
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
assert (ctx != NULL);
size_t sze = ctx->ao_basis.ao_num * 5 * ctx->electron.num * ctx->electron.walk_num;
memcpy(ao_vgl, ctx->ao_basis.ao_vgl, sze * sizeof(double));
return QMCKL_SUCCESS;
}
#+end_src
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_get_ao_vgl (context, ao_vgl) &
bind(C)
use, intrinsic :: iso_c_binding
import
implicit none
integer (c_int64_t) , intent(in) , value :: context
double precision, intent(out) :: ao_vgl(*)
end function
end interface
#+end_src
*** TODO Provide
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
qmckl_exit_code qmckl_provide_ao_vgl(qmckl_context context);
#+end_src
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
qmckl_exit_code qmckl_provide_ao_vgl(qmckl_context context)
{
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
return QMCKL_NULL_CONTEXT;
}
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
assert (ctx != NULL);
if (!ctx->ao_basis.provided) {
return qmckl_failwith( context,
QMCKL_NOT_PROVIDED,
"qmckl_ao_vgl",
NULL);
}
if(!(ctx->electron.provided)) {
return qmckl_failwith( context,
QMCKL_NOT_PROVIDED,
"qmckl_electron",
NULL);
}
/* Compute if necessary */
if (ctx->electron.coord_new_date > ctx->ao_basis.ao_vgl_date) {
qmckl_exit_code rc;
/* Provide required data */
rc = qmckl_provide_ao_basis_shell_vgl(context);
if (rc != QMCKL_SUCCESS) {
return qmckl_failwith( context, rc, "qmckl_provide_ao_basis_shell_vgl", NULL);
}
/* Allocate array */
if (ctx->ao_basis.ao_vgl == NULL) {
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
mem_info.size = ctx->ao_basis.ao_num * 5 * ctx->electron.num *
ctx->electron.walk_num * sizeof(double);
double* ao_vgl = (double*) qmckl_malloc(context, mem_info);
if (ao_vgl == NULL) {
return qmckl_failwith( context,
QMCKL_ALLOCATION_FAILED,
"qmckl_ao_basis_ao_vgl",
NULL);
}
ctx->ao_basis.ao_vgl = ao_vgl;
}
rc = qmckl_compute_ao_vgl(context,
ctx->ao_basis.ao_num,
ctx->ao_basis.shell_num,
ctx->electron.num,
ctx->nucleus.num,
ctx->electron.walk_num,
ctx->electron.coord_new,
ctx->nucleus.coord,
ctx->ao_basis.nucleus_index,
ctx->ao_basis.nucleus_shell_num,
ctx->ao_basis.nucleus_range,
ctx->ao_basis.nucleus_max_ang_mom,
ctx->ao_basis.shell_ang_mom,
ctx->ao_basis.shell_vgl,
ctx->ao_basis.ao_vgl);
if (rc != QMCKL_SUCCESS) {
return rc;
}
ctx->ao_basis.ao_vgl_date = ctx->date;
}
return QMCKL_SUCCESS;
}
#+end_src
*** TODO Compute
:PROPERTIES:
:Name: qmckl_compute_ao_vgl
:CRetType: qmckl_exit_code
:FRetType: qmckl_exit_code
:END:
#+NAME: qmckl_ao_vgl_args
| ~qmckl_context~ | ~context~ | in | Global state |
| ~int64_t~ | ~ao_num~ | in | Number of AOs |
| ~int64_t~ | ~shell_num~ | in | Number of shells |
| ~int64_t~ | ~elec_num~ | in | Number of electrons |
| ~int64_t~ | ~nucl_num~ | in | Number of nuclei |
| ~int64_t~ | ~walk_num~ | in | Number of walkers |
| ~double~ | ~elec_coord[walk_num][3][elec_num]~ | in | Electron coordinates |
| ~double~ | ~nucl_coord[3][nucl_num]~ | in | Nuclear coordinates |
| ~int64_t~ | ~nucleus_index[nucl_num]~ | in | Index of the 1st shell of each nucleus |
| ~int64_t~ | ~nucleus_shell_num[nucl_num]~ | in | Number of shells per nucleus |
| ~double~ | ~nucleus_range[nucl_num]~ | in | Range beyond which all is zero |
| ~int32_t~ | ~nucleus_max_ang_mom[nucl_num]~ | in | Maximum angular momentum per nucleus |
| ~int32_t~ | ~shell_ang_mom[shell_num]~ | in | Angular momentum of each shell |
| ~double~ | ~shell_vgl[5][walk_num][elec_num][shell_num]~ | in | Value, gradients and Laplacian of the shells |
| ~double~ | ~ao_vgl[5][walk_num][elec_num][ao_num]~ | out | Value, gradients and Laplacian of the AOs |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_ao_vgl_f(context, &
ao_num, shell_num, elec_num, nucl_num, walk_num, &
elec_coord, nucl_coord, nucleus_index, nucleus_shell_num, &
nucleus_range, nucleus_max_ang_mom, shell_ang_mom, &
shell_vgl, ao_vgl) &
result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in) :: context
integer*8 , intent(in) :: ao_num
integer*8 , intent(in) :: shell_num
integer*8 , intent(in) :: elec_num
integer*8 , intent(in) :: nucl_num
integer*8 , intent(in) :: walk_num
double precision , intent(in) :: elec_coord(elec_num,3,walk_num)
double precision , intent(in) :: nucl_coord(nucl_num,3)
integer*8 , intent(in) :: nucleus_index(nucl_num)
integer*8 , intent(in) :: nucleus_shell_num(nucl_num)
double precision , intent(in) :: nucleus_range(nucl_num)
integer , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer , intent(in) :: shell_ang_mom(shell_num)
double precision , intent(in) :: shell_vgl(shell_num,elec_num,walk_num,5)
double precision , intent(out) :: ao_vgl(ao_num,elec_num,walk_num,5)
double precision :: e_coord(3), n_coord(3)
integer*8 :: n_poly
integer :: l, il, k
integer*8 :: ielec, inucl, ishell, iwalk
integer :: lstart(0:20)
double precision :: x, y, z, r2
double precision :: cutoff
integer, external :: qmckl_ao_polynomial_vgl_f
double precision, allocatable :: poly_vgl(:,:)
integer , allocatable :: powers(:,:)
allocate(poly_vgl(5,ao_num), powers(3,ao_num))
! Pre-computed data
do l=0,20
lstart(l) = l*(l+1)*(l+2)/6 +1
end do
info = QMCKL_SUCCESS
! Don't compute polynomials when the radial part is zero.
! TODO : Use numerical precision here
cutoff = -dlog(1.d-15)
do iwalk = 1,walk_num
do ielec = 1, elec_num
e_coord(1) = elec_coord(ielec,1,iwalk)
e_coord(2) = elec_coord(ielec,2,iwalk)
e_coord(3) = elec_coord(ielec,3,iwalk)
k=1
do inucl=1,nucl_num
n_coord(1) = nucl_coord(inucl,1)
n_coord(2) = nucl_coord(inucl,2)
n_coord(3) = nucl_coord(inucl,3)
! Test if the electron is in the range of the nucleus
x = e_coord(1) - n_coord(1)
y = e_coord(2) - n_coord(2)
z = e_coord(3) - n_coord(3)
r2 = x*x + z*z + z*z
if (r2 > cutoff*nucleus_range(inucl)) then
cycle
end if
! Compute polynomials
info = qmckl_ao_polynomial_vgl_f(context, e_coord, n_coord, &
nucleus_max_ang_mom(inucl), n_poly, powers, 3_8, &
poly_vgl, 5_8)
! Loop over shells
do ishell = nucleus_index(inucl)+1, nucleus_index(inucl)+nucleus_shell_num(inucl)
l = shell_ang_mom(ishell)
do il = lstart(l), lstart(l+1)-1
! Value
ao_vgl(k,ielec,iwalk,1) = &
poly_vgl(1,il) * shell_vgl(ishell,ielec,iwalk,1)
! Grad_x
ao_vgl(k,ielec,iwalk,2) = &
poly_vgl(2,il) * shell_vgl(ishell,ielec,iwalk,1) + &
poly_vgl(1,il) * shell_vgl(ishell,ielec,iwalk,2)
! Grad_y
ao_vgl(k,ielec,iwalk,3) = &
poly_vgl(3,il) * shell_vgl(ishell,ielec,iwalk,1) + &
poly_vgl(1,il) * shell_vgl(ishell,ielec,iwalk,3)
! Grad_z
ao_vgl(k,ielec,iwalk,4) = &
poly_vgl(4,il) * shell_vgl(ishell,ielec,iwalk,1) + &
poly_vgl(1,il) * shell_vgl(ishell,ielec,iwalk,4)
! Lapl_z
ao_vgl(k,ielec,iwalk,5) = &
poly_vgl(5,il) * shell_vgl(ishell,ielec,iwalk,1) + &
poly_vgl(1,il) * shell_vgl(ishell,ielec,iwalk,5) + &
2.d0 * ( &
poly_vgl(2,il) * shell_vgl(ishell,ielec,iwalk,2) + &
poly_vgl(3,il) * shell_vgl(ishell,ielec,iwalk,3) + &
poly_vgl(4,il) * shell_vgl(ishell,ielec,iwalk,4) )
k = k+1
end do
end do
end do
end do
end do
deallocate(poly_vgl, powers)
end function qmckl_compute_ao_vgl_f
#+end_src
# #+CALL: generate_c_header(table=qmckl_ao_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_vgl"))
#+RESULTS:
#+begin_src c :tangle (eval h_private_func) :comments org
qmckl_exit_code qmckl_compute_ao_vgl (
const qmckl_context context,
const int64_t ao_num,
const int64_t shell_num,
const int64_t elec_num,
const int64_t nucl_num,
const int64_t walk_num,
const double* elec_coord,
const double* nucl_coord,
const int64_t* nucleus_index,
const int64_t* nucleus_shell_num,
const double* nucleus_range,
const int32_t* nucleus_max_ang_mom,
const int32_t* shell_ang_mom,
const double* shell_vgl,
double* const ao_vgl );
#+end_src
#+CALL: generate_c_interface(table=qmckl_ao_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_vgl"))
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_compute_ao_vgl &
(context, &
ao_num, &
shell_num, &
elec_num, &
nucl_num, &
walk_num, &
elec_coord, &
nucl_coord, &
nucleus_index, &
nucleus_shell_num, &
nucleus_range, &
nucleus_max_ang_mom, &
shell_ang_mom, &
shell_vgl, &
ao_vgl) &
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 :: ao_num
integer (c_int64_t) , intent(in) , value :: shell_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 :: walk_num
real (c_double ) , intent(in) :: elec_coord(elec_num,3,walk_num)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(nucl_num)
integer (c_int32_t) , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer (c_int32_t) , intent(in) :: shell_ang_mom(shell_num)
real (c_double ) , intent(in) :: shell_vgl(shell_num,elec_num,walk_num,5)
real (c_double ) , intent(out) :: ao_vgl(ao_num,elec_num,walk_num,5)
integer(c_int32_t), external :: qmckl_compute_ao_vgl_f
info = qmckl_compute_ao_vgl_f &
(context, &
ao_num, &
shell_num, &
elec_num, &
nucl_num, &
walk_num, &
elec_coord, &
nucl_coord, &
nucleus_index, &
nucleus_shell_num, &
nucleus_range, &
nucleus_max_ang_mom, &
shell_ang_mom, &
shell_vgl, &
ao_vgl)
end function qmckl_compute_ao_vgl
#+end_src
#+begin_src python :results output :exports none
import numpy as np
def f(a,x,y):
return np.sum( [c * np.exp( -b*(np.linalg.norm(x-y))**2) for b,c in a] )
def df(a,x,y,n):
h0 = 1.e-6
if n == 1: h = np.array([h0,0.,0.])
elif n == 2: h = np.array([0.,h0,0.])
elif n == 3: h = np.array([0.,0.,h0])
return ( f(a,x+h,y) - f(a,x-h,y) ) / (2.*h0)
def d2f(a,x,y,n):
h0 = 1.e-6
if n == 1: h = np.array([h0,0.,0.])
elif n == 2: h = np.array([0.,h0,0.])
elif n == 3: h = np.array([0.,0.,h0])
return ( f(a,x+h,y) - 2.*f(a,x,y) + f(a,x-h,y) ) / h0**2
def lf(a,x,y):
return d2f(a,x,y,1) + d2f(a,x,y,2) + d2f(a,x,y,3)
elec_26_w1 = np.array( [ 1.49050402641, 2.90106987953, -1.05920815468 ] )
elec_15_w2 = np.array( [ -2.20180344582,-1.9113150239, 2.2193744778600002 ] )
nucl_1 = np.array( [ -2.302574592081335e+00, -3.542027060505035e-01, -5.334129934317614e-02] )
#double prim_vgl[prim_num][5][walk_num][elec_num];
x = elec_26_w1 ; y = nucl_1
a = [( 403.830000, 0.001473 * 5.9876577632594533e+04),
( 121.170000, 0.012672 * 7.2836806319891484e+03),
( 46.345000, 0.058045 * 1.3549226646722386e+03),
( 19.721000, 0.170510 * 3.0376315094739988e+02),
( 8.862400, 0.318596 * 7.4924579607137730e+01),
( 3.996200, 0.384502 * 1.8590543353806009e+01),
( 1.763600, 0.273774 * 4.4423176930919421e+00),
( 0.706190, 0.074397 * 8.9541051939952665e-01)]
print ( "[0][0][26][219] : %25.15e"%(f(a,x,y) * (x[0] - y[0])**2) )
print ( "[1][0][26][219] : %25.15e"%(df(a,x,y,1)* (x[0] - y[0]) * (x[1] - y[1]) + 2.*f(a,x,y) * (x[0] - y[0])) )
print ( "[0][0][26][220] : %25.15e"%(f(a,x,y) * (x[0] - y[0]) * (x[1] - y[1])) )
print ( "[1][0][26][220] : %25.15e"%(df(a,x,y,1)* (x[0] - y[0]) * (x[1] - y[1]) + f(a,x,y) * (x[1] - y[1])) )
print ( "[0][0][26][221] : %25.15e"%(f(a,x,y) * (x[0] - y[0]) * (x[2] - y[2])) )
print ( "[1][0][26][221] : %25.15e"%(df(a,x,y,1)* (x[0] - y[0]) * (x[2] - y[2]) + f(a,x,y) * (x[2] - y[2])) )
print ( "[0][0][26][222] : %25.15e"%(f(a,x,y) * (x[1] - y[1]) * (x[1] - y[1])) )
print ( "[1][0][26][222] : %25.15e"%(df(a,x,y,1)* (x[1] - y[1]) * (x[1] - y[1])) )
print ( "[0][0][26][223] : %25.15e"%(f(a,x,y) * (x[1] - y[1]) * (x[2] - y[2])) )
print ( "[1][0][26][223] : %25.15e"%(df(a,x,y,1)* (x[1] - y[1]) * (x[2] - y[2])) )
print ( "[0][0][26][224] : %25.15e"%(f(a,x,y) * (x[2] - y[2]) * (x[2] - y[2])) )
print ( "[1][0][26][224] : %25.15e"%(df(a,x,y,1)* (x[2] - y[2]) * (x[2] - y[2])) )
#+end_src
#+RESULTS:
#+begin_example
[0][0][26][219] : 1.020302912653649e-08
[1][0][26][219] : -4.153046808203204e-08
[0][0][26][220] : 8.756380857379661e-09
[1][0][26][220] : -4.460176677299534e-08
[0][0][26][221] : -2.705688401075445e-09
[1][0][26][221] : 1.378177639720419e-08
[0][0][26][222] : 7.514847283937212e-09
[1][0][26][222] : -4.025905373647693e-08
[0][0][26][223] : -2.322059246071533e-09
[1][0][26][223] : 1.243989457599443e-08
[0][0][26][224] : 7.175074806631758e-10
[1][0][26][224] : -3.843880138733679e-09
#+end_example
*** TODO Test
#+begin_src c :tangle (eval c_test) :exports none
{
#define walk_num chbrclf_walk_num
#define elec_num chbrclf_elec_num
#define shell_num chbrclf_shell_num
#define ao_num chbrclf_ao_num
int64_t elec_up_num = chbrclf_elec_up_num;
int64_t elec_dn_num = chbrclf_elec_dn_num;
double* elec_coord = &(chbrclf_elec_coord[0][0][0]);
rc = qmckl_set_electron_num (context, elec_up_num, elec_dn_num);
assert (rc == QMCKL_SUCCESS);
rc = qmckl_set_electron_walk_num (context, walk_num);
assert (rc == QMCKL_SUCCESS);
assert(qmckl_electron_provided(context));
rc = qmckl_set_electron_coord (context, 'N', elec_coord);
assert(rc == QMCKL_SUCCESS);
double ao_vgl[5][walk_num][elec_num][ao_num];
rc = qmckl_get_ao_vgl(context, &(ao_vgl[0][0][0][0]));
assert (rc == QMCKL_SUCCESS);
printf("\n");
printf(" ao_vgl ao_vgl[0][0][26][219] %25.15e\n", ao_vgl[0][0][26][219]);
printf(" ao_vgl ao_vgl[1][0][26][219] %25.15e\n", ao_vgl[1][0][26][219]);
printf(" ao_vgl ao_vgl[0][0][26][220] %25.15e\n", ao_vgl[0][0][26][220]);
printf(" ao_vgl ao_vgl[1][0][26][220] %25.15e\n", ao_vgl[1][0][26][220]);
printf(" ao_vgl ao_vgl[0][0][26][221] %25.15e\n", ao_vgl[0][0][26][221]);
printf(" ao_vgl ao_vgl[1][0][26][221] %25.15e\n", ao_vgl[1][0][26][221]);
printf(" ao_vgl ao_vgl[0][0][26][222] %25.15e\n", ao_vgl[0][0][26][222]);
printf(" ao_vgl ao_vgl[1][0][26][222] %25.15e\n", ao_vgl[1][0][26][222]);
printf(" ao_vgl ao_vgl[0][0][26][223] %25.15e\n", ao_vgl[0][0][26][223]);
printf(" ao_vgl ao_vgl[1][0][26][223] %25.15e\n", ao_vgl[1][0][26][223]);
printf(" ao_vgl ao_vgl[0][0][26][224] %25.15e\n", ao_vgl[0][0][26][224]);
printf(" ao_vgl ao_vgl[1][0][26][224] %25.15e\n", ao_vgl[1][0][26][224]);
printf("\n");
assert( fabs(ao_vgl[0][0][26][219] - ( 1.020298798341620e-08)) < 1.e-14 );
assert( fabs(ao_vgl[1][0][26][219] - ( -4.928035238010602e-08)) < 1.e-14 );
assert( fabs(ao_vgl[0][0][26][220] - ( 8.756345547784206e-09)) < 1.e-14 );
assert( fabs(ao_vgl[1][0][26][220] - ( -4.460158690983819e-08)) < 1.e-14 );
assert( fabs(ao_vgl[0][0][26][221] - ( -2.705677490544664e-09)) < 1.e-14 );
assert( fabs(ao_vgl[1][0][26][221] - ( 1.378172082017231e-08)) < 1.e-14 );
assert( fabs(ao_vgl[0][0][26][222] - ( 7.514816980753531e-09)) < 1.e-14 );
assert( fabs(ao_vgl[1][0][26][222] - ( -4.025889138635182e-08)) < 1.e-14 );
assert( fabs(ao_vgl[0][0][26][223] - ( -2.322049882502961e-09)) < 1.e-14 );
assert( fabs(ao_vgl[1][0][26][223] - ( 1.243984441042288e-08)) < 1.e-14 );
assert( fabs(ao_vgl[0][0][26][224] - ( 7.175045873560788e-10)) < 1.e-14 );
assert( fabs(ao_vgl[1][0][26][224] - ( -3.843864637762753e-09)) < 1.e-14 );
}
#+end_src
* End of files :noexport:
#+begin_src c :tangle (eval h_private_type)

View File

@ -528,7 +528,7 @@ F 1
#define chbrclf_prim_num 297
#define chbrclf_ao_num 263
int64_t chbrclf_basis_nucleus_index[chbrclf_nucl_num] = {0, 14, 23, 27, 53};
int64_t chbrclf_basis_nucleus_index[chbrclf_nucl_num] = {0, 14, 23, 37, 53};
int64_t chbrclf_basis_nucleus_shell_num[chbrclf_nucl_num] = {14, 9, 14, 16, 19};

View File

@ -96,7 +96,6 @@ def parse_table(table):
return result
#+END_SRC
*** Generates a C header
#+NAME: generate_c_header
@ -129,7 +128,6 @@ return template
#+END_SRC
*** Generates a C interface to the Fortran function
#+NAME: generate_c_interface