#+TITLE: Atomic Orbitals #+SETUPFILE: ../tools/theme.setup #+INCLUDE: ../tools/lib.org The atomic basis set is defined as a list of shells. Each shell $s$ is centered on a nucleus $A$, possesses a given angular momentum $l$ and a radial function $R_s$. The radial function is a linear combination of \emph{primitive} functions that can be of type Slater ($p=1$) or Gaussian ($p=2$): \[ R_s(\mathbf{r}) = \mathcal{N}_s |\mathbf{r}-\mathbf{R}_A|^{n_s} \sum_{k=1}^{N_{\text{prim}}} a_{ks}\, f_{ks} \exp \left( - \gamma_{ks} | \mathbf{r}-\mathbf{R}_A | ^p \right). \] In the case of Gaussian functions, $n_s$ is always zero. The normalization factor $\mathcal{N}_s$ ensures that all the functions of the shell are normalized to unity. Usually, basis sets are given a combination of normalized primitives, so the normalization coefficients of the primitives, $f_{ks}$, need also to be provided. Atomic orbitals (AOs) are defined as \[ \chi_i (\mathbf{r}) = \mathcal{M}_i\, P_{\eta(i)}(\mathbf{r})\, R_{\theta(i)} (\mathbf{r}) \] where $\theta(i)$ returns the shell on which the AO is expanded, and $\eta(i)$ denotes which angular function is chosen. Here, the parameter $\mathcal{M}_i$ is an extra parameter which allows the normalization of the different functions of the same shell to be different, as in GAMESS for example. In this section we describe first how the basis set is stored in the context, and then we present the kernels used to compute the values, gradients and Laplacian of the atomic basis functions. * Headers :noexport: #+begin_src elisp :noexport :results none (org-babel-lob-ingest "../tools/lib.org") #+end_src #+begin_src c :tangle (eval h_private_type) #ifndef QMCKL_AO_HPT #define QMCKL_AO_HPT #include #+end_src #+begin_src c :tangle (eval c_test) :noweb yes #include "qmckl.h" #include "assert.h" #ifdef HAVE_CONFIG_H #include "config.h" #endif #include #include #include "chbrclf.h" int main() { qmckl_context context; context = qmckl_context_create(); #+end_src #+begin_src c :tangle (eval c) #ifdef HAVE_CONFIG_H #include "config.h" #endif #ifdef HAVE_STDINT_H #include #elif HAVE_INTTYPES_H #include #endif #include #include #include #include #include "qmckl.h" #include "qmckl_context_private_type.h" #include "qmckl_memory_private_type.h" #include "qmckl_memory_private_func.h" #include "qmckl_ao_private_type.h" #include "qmckl_ao_private_func.h" #+end_src * Context The following arrays are stored in the context: |---------------------+---------------+----------------------------------------------------------------------| | ~type~ | | Gaussian (~'G'~) or Slater (~'S'~) | | ~shell_num~ | | Number of shells | | ~prim_num~ | | Total number of primitives | | ~nucleus_index~ | ~[nucl_num]~ | Index of the first shell of each nucleus | | ~nucleus_shell_num~ | ~[nucl_num]~ | Number of shells per nucleus | | ~shell_ang_mom~ | ~[shell_num]~ | Angular momentum of each shell | | ~shell_prim_num~ | ~[shell_num]~ | Number of primitives in each shell | | ~shell_prim_index~ | ~[shell_num]~ | Address of the first primitive of each shell in the ~EXPONENT~ array | | ~shell_factor~ | ~[shell_num]~ | Normalization factor for each shell | | ~exponent~ | ~[prim_num]~ | Array of exponents | | ~coefficient~ | ~[prim_num]~ | Array of coefficients | | ~prim_factor~ | ~[prim_num]~ | Normalization factors of the primtives | Computed data: |--------------------------+--------------------------------------+-----------------------------------------------------------------------------------------------| | ~coefficient_normalized~ | ~[prim_num]~ | Normalized primitive coefficients | | ~nucleus_prim_index~ | ~[nucl_num]~ | Index of the first primitive for each nucleus | | ~nucleus_max_ang_mom~ | ~[nucl_num]~ | Maximum angular momentum for each nucleus | | ~nucleus_range~ | ~[nucl_num]~ | Distance beyond which all the AOs are zero | |--------------------------+--------------------------------------+-----------------------------------------------------------------------------------------------| | ~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 | |--------------------------+--------------------------------------+-----------------------------------------------------------------------------------------------| | ~nucl_shell_index~ | ~[nucl_num]~ | Index of the first shell for each nucleus | | ~exponent_sorted~ | ~[prim_num]~ | Array of exponents for sorted primitives | | ~coeff_norm_sorted~ | ~[prim_num]~ | Array of normalized coefficients for sorted primitives | | ~prim_factor_sorted~ | ~[prim_num]~ | Normalization factors of the sorted primtives | For H_2 with the following basis set, #+BEGIN_EXAMPLE HYDROGEN S 5 1 3.387000E+01 6.068000E-03 2 5.095000E+00 4.530800E-02 3 1.159000E+00 2.028220E-01 4 3.258000E-01 5.039030E-01 5 1.027000E-01 3.834210E-01 S 1 1 3.258000E-01 1.000000E+00 S 1 1 1.027000E-01 1.000000E+00 P 1 1 1.407000E+00 1.000000E+00 P 1 1 3.880000E-01 1.000000E+00 D 1 1 1.057000E+00 1.0000000 #+END_EXAMPLE we have: #+BEGIN_EXAMPLE type = 'G' shell_num = 12 prim_num = 20 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.] shell_prim_num = [5, 1, 1, 1, 1, 1, 5, 1, 1, 1, 1, 1] shell_prim_index = [0 , 5 , 6 , 7 , 8 , 9 , 10, 15, 16, 17, 18, 19] exponent = [ 33.87, 5.095, 1.159, 0.3258, 0.1027, 0.3258, 0.1027, 1.407, 0.388, 1.057, 33.87, 5.095, 1.159, 0.3258, 0.1027, 0.3258, 0.1027, 1.407, 0.388, 1.057] coefficient = [ 0.006068, 0.045308, 0.202822, 0.503903, 0.383421, 1.0, 1.0, 1.0, 1.0, 1.0, 0.006068, 0.045308, 0.202822, 0.503903, 0.383421, 1.0, 1.0, 1.0, 1.0, 1.0] prim_factor = [ 1.0006253235944540e+01, 2.4169531573445120e+00, 7.9610924849766440e-01 3.0734305383061117e-01, 1.2929684417481876e-01, 3.0734305383061117e-01, 1.2929684417481876e-01, 2.1842769845268308e+00, 4.3649547399719840e-01, 1.8135965626177861e+00, 1.0006253235944540e+01, 2.4169531573445120e+00, 7.9610924849766440e-01, 3.0734305383061117e-01, 1.2929684417481876e-01, 3.0734305383061117e-01, 1.2929684417481876e-01, 2.1842769845268308e+00, 4.3649547399719840e-01, 1.8135965626177861e+00 ] #+END_EXAMPLE ** Data structure #+begin_src c :comments org :tangle (eval h_private_type) typedef struct qmckl_ao_basis_struct { int32_t uninitialized; int64_t shell_num; int64_t prim_num; int64_t * nucleus_index; int64_t * nucleus_shell_num; int32_t * shell_ang_mom; int64_t * shell_prim_num; int64_t * shell_prim_index; double * shell_factor; double * exponent ; double * coefficient ; double * prim_factor ; int64_t * nucleus_prim_index; double * coefficient_normalized; int32_t * nucleus_max_ang_mom; double * nucleus_range; double * primitive_vgl; int64_t primitive_vgl_date; double * shell_vgl; int64_t shell_vgl_date; bool provided; char type; } qmckl_ao_basis_struct; #+end_src The ~uninitialized~ integer contains one bit set to one for each initialization function which has not been called. It becomes equal to zero after all initialization functions have been called. The struct is then initialized and ~provided == true~. Some values are initialized by default, and are not concerned by this mechanism. #+begin_src c :comments org :tangle (eval h_private_func) qmckl_exit_code qmckl_init_ao_basis(qmckl_context context); #+end_src #+begin_src c :comments org :tangle (eval c) qmckl_exit_code qmckl_init_ao_basis(qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return false; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); ctx->ao_basis.uninitialized = (1 << 12) - 1; /* Default values */ /* ctx->ao_basis. ,*/ return QMCKL_SUCCESS; } #+end_src ** Access functions #+begin_src c :comments org :tangle (eval h_private_func) :exports none char qmckl_get_ao_basis_type (const qmckl_context context); 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_nucleus_index (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); double* qmckl_get_ao_basis_shell_factor (const qmckl_context context); double* qmckl_get_ao_basis_exponent (const qmckl_context context); double* qmckl_get_ao_basis_coefficient (const qmckl_context context); double* qmckl_get_ao_basis_prim_factor (const qmckl_context context); #+end_src When all the data for the AOs have been provided, the following function returns ~true~. #+begin_src c :comments org :tangle (eval h_func) bool qmckl_ao_basis_provided (const qmckl_context context); #+end_src #+NAME:post #+begin_src c :exports none if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } #+end_src #+begin_src c :comments org :tangle (eval c) :noweb yes :exports none char qmckl_get_ao_basis_type (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return (char) 0; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return (char) 0; } assert (ctx->ao_basis.type != (char) 0); return ctx->ao_basis.type; } int64_t qmckl_get_ao_basis_shell_num (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return (int64_t) 0; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 1; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return (int64_t) 0; } assert (ctx->ao_basis.shell_num > (int64_t) 0); return ctx->ao_basis.shell_num; } int64_t qmckl_get_ao_basis_prim_num (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return (int64_t) 0; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 2; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return (int64_t) 0; } assert (ctx->ao_basis.prim_num > (int64_t) 0); return ctx->ao_basis.prim_num; } int64_t* qmckl_get_ao_basis_nucleus_shell_num (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 3; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.nucleus_shell_num != NULL); return ctx->ao_basis.nucleus_shell_num ; } int64_t* qmckl_get_ao_basis_nucleus_index (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 4; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.nucleus_index != NULL); return ctx->ao_basis.nucleus_index ; } int32_t* qmckl_get_ao_basis_shell_ang_mom (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 5; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.shell_ang_mom != NULL); return ctx->ao_basis.shell_ang_mom; } int64_t* qmckl_get_ao_basis_shell_prim_num (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 6; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.shell_prim_num != NULL); return ctx->ao_basis.shell_prim_num; } int64_t* qmckl_get_ao_basis_shell_prim_index (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 7; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.shell_prim_index != NULL); return ctx->ao_basis.shell_prim_index; } double* qmckl_get_ao_basis_shell_factor (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 8; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.shell_factor != NULL); return ctx->ao_basis.shell_factor; } double* qmckl_get_ao_basis_exponent (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 9; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.exponent != NULL); return ctx->ao_basis.exponent; } double* qmckl_get_ao_basis_coefficient (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 10; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.coefficient != NULL); return ctx->ao_basis.coefficient; } double* qmckl_get_ao_basis_prim_factor (const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return NULL; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int32_t mask = 1 << 11; if ( (ctx->ao_basis.uninitialized & mask) != 0) { return NULL; } assert (ctx->ao_basis.prim_factor != NULL); return ctx->ao_basis.prim_factor; } bool qmckl_ao_basis_provided(const qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return false; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); return ctx->ao_basis.provided; } #+end_src ** Initialization functions To set the basis set, all the following functions need to be called. #+begin_src c :comments org :tangle (eval h_func) qmckl_exit_code qmckl_set_ao_basis_type (qmckl_context context, const char t); qmckl_exit_code qmckl_set_ao_basis_shell_num (qmckl_context context, const int64_t shell_num); qmckl_exit_code qmckl_set_ao_basis_prim_num (qmckl_context context, const int64_t prim_num); qmckl_exit_code qmckl_set_ao_basis_nucleus_index (qmckl_context context, const int64_t * nucleus_index); qmckl_exit_code qmckl_set_ao_basis_nucleus_shell_num(qmckl_context context, const int64_t * nucleus_shell_num); qmckl_exit_code qmckl_set_ao_basis_shell_ang_mom (qmckl_context context, const int32_t * shell_ang_mom); qmckl_exit_code qmckl_set_ao_basis_shell_prim_num (qmckl_context context, const int64_t * shell_prim_num); qmckl_exit_code qmckl_set_ao_basis_shell_prim_index (qmckl_context context, const int64_t * shell_prim_index); qmckl_exit_code qmckl_set_ao_basis_shell_factor (qmckl_context context, const double * shell_factor); qmckl_exit_code qmckl_set_ao_basis_exponent (qmckl_context context, const double * exponent); qmckl_exit_code qmckl_set_ao_basis_coefficient (qmckl_context context, const double * coefficient); qmckl_exit_code qmckl_set_ao_basis_prim_factor (qmckl_context context, const double * prim_factor); #+end_src #+NAME:pre2 #+begin_src c :exports none if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return QMCKL_NULL_CONTEXT; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; #+end_src #+NAME:post2 #+begin_src c :exports none ctx->ao_basis.uninitialized &= ~mask; ctx->ao_basis.provided = (ctx->ao_basis.uninitialized == 0); if (ctx->ao_basis.provided) { qmckl_exit_code rc_ = qmckl_finalize_basis(context); if (rc_ != QMCKL_SUCCESS) return rc_; } return QMCKL_SUCCESS; #+end_src #+begin_src c :comments org :tangle (eval c) :noweb yes :exports none qmckl_exit_code qmckl_set_ao_basis_type(qmckl_context context, const char t) { <> if (t != 'G' && t != 'S') { return qmckl_failwith( context, QMCKL_INVALID_ARG_2, "qmckl_set_ao_basis_type", NULL); } int32_t mask = 1; ctx->ao_basis.type = t; <> } qmckl_exit_code qmckl_set_ao_basis_shell_num(qmckl_context context, const int64_t shell_num) { <> if (shell_num <= 0) { return qmckl_failwith( context, QMCKL_INVALID_ARG_2, "qmckl_set_ao_basis_shell_num", "shell_num <= 0"); } int64_t prim_num = qmckl_get_ao_basis_prim_num(context); if (0L < prim_num && prim_num < shell_num) { return qmckl_failwith( context, QMCKL_INVALID_ARG_2, "qmckl_set_ao_basis_shell_num", "shell_num > prim_num"); } int32_t mask = 1 << 1; ctx->ao_basis.shell_num = shell_num; <> } qmckl_exit_code qmckl_set_ao_basis_prim_num(qmckl_context context, const int64_t prim_num) { <> if (prim_num <= 0) { return qmckl_failwith( context, QMCKL_INVALID_ARG_2, "qmckl_set_ao_basis_shell_num", "prim_num must be positive"); } int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (prim_num < shell_num) { return qmckl_failwith( context, QMCKL_INVALID_ARG_2, "qmckl_set_ao_basis_shell_num", "prim_num < shell_num"); } int32_t mask = 1 << 2; ctx->ao_basis.prim_num = prim_num; <> } qmckl_exit_code qmckl_set_ao_basis_nucleus_shell_num(qmckl_context context, const int64_t* nucleus_shell_num) { <> int32_t mask = 1 << 3; const int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (shell_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_nucleus_shell_num", "shell_num is not set"); } if (ctx->ao_basis.nucleus_shell_num != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.nucleus_shell_num); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_nucleus_shell_num", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = shell_num * sizeof(int64_t); int64_t* new_array = (int64_t*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_nucleus_shell_num", NULL); } memcpy(new_array, nucleus_shell_num, mem_info.size); ctx->ao_basis.nucleus_shell_num = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_nucleus_index(qmckl_context context, const int64_t* nucleus_index) { <> int32_t mask = 1 << 4; const int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (shell_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_nucleus_index", "shell_num is not set"); } if (ctx->ao_basis.nucleus_index != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.nucleus_index); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_nucleus_index", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = shell_num * sizeof(int64_t); int64_t* new_array = (int64_t*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_nucleus_index", NULL); } memcpy(new_array, nucleus_index, mem_info.size); ctx->ao_basis.nucleus_index = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_shell_ang_mom(qmckl_context context, const int32_t* shell_ang_mom) { <> int32_t mask = 1 << 5; const int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (shell_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_shell_ang_mom", "shell_num is not set"); } if (ctx->ao_basis.shell_ang_mom != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.shell_ang_mom); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_shell_ang_mom", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = shell_num * sizeof(char); int32_t * new_array = (int32_t*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_shell_ang_mom", NULL); } memcpy(new_array, shell_ang_mom, mem_info.size); ctx->ao_basis.shell_ang_mom = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_shell_prim_num(qmckl_context context, const int64_t* shell_prim_num) { <> int32_t mask = 1 << 6; const int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (shell_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_shell_prim_num", "shell_num is not set"); } if (ctx->ao_basis.shell_prim_num != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.shell_prim_num); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_shell_prim_num", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = shell_num * sizeof(int64_t); int64_t* new_array = (int64_t*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_shell_prim_num", NULL); } memcpy(new_array, shell_prim_num, mem_info.size); ctx->ao_basis.shell_prim_num = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_shell_prim_index(qmckl_context context, const int64_t* shell_prim_index) { <> int32_t mask = 1 << 7; const int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (shell_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_shell_prim_index", "shell_num is not set"); } if (ctx->ao_basis.shell_prim_index != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.shell_prim_index); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_shell_prim_index", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = shell_num * sizeof(int64_t); int64_t* new_array = (int64_t*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_shell_prim_index", NULL); } memcpy(new_array, shell_prim_index, mem_info.size); ctx->ao_basis.shell_prim_index = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_shell_factor(qmckl_context context, const double* shell_factor) { <> int32_t mask = 1 << 8; const int64_t shell_num = qmckl_get_ao_basis_shell_num(context); if (shell_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_shell_factor", "shell_num is not set"); } if (ctx->ao_basis.shell_factor != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.shell_factor); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_shell_factor", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = shell_num * sizeof(double); double* new_array = (double*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_shell_factor", NULL); } memcpy(new_array, shell_factor, mem_info.size); ctx->ao_basis.shell_factor = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_exponent(qmckl_context context, const double* exponent) { <> int32_t mask = 1 << 9; const int64_t prim_num = qmckl_get_ao_basis_prim_num(context); if (prim_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_exponent", "prim_num is not set"); } if (ctx->ao_basis.exponent != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.exponent); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_exponent", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = prim_num * sizeof(double); double* new_array = (double*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_exponent", NULL); } memcpy(new_array, exponent, mem_info.size); ctx->ao_basis.exponent = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_coefficient(qmckl_context context, const double* coefficient) { <> int32_t mask = 1 << 10; const int64_t prim_num = qmckl_get_ao_basis_prim_num(context); if (prim_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_coefficient", "prim_num is not set"); } if (ctx->ao_basis.coefficient != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.coefficient); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_coefficient", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = prim_num * sizeof(double); double* new_array = (double*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_coefficient", NULL); } memcpy(new_array, coefficient, mem_info.size); ctx->ao_basis.coefficient = new_array; <> } qmckl_exit_code qmckl_set_ao_basis_prim_factor(qmckl_context context, const double* prim_factor) { <> int32_t mask = 1 << 11; const int64_t prim_num = qmckl_get_ao_basis_prim_num(context); if (prim_num == 0L) { return qmckl_failwith( context, QMCKL_FAILURE, "qmckl_set_ao_basis_prim_factor", "prim_num is not set"); } if (ctx->ao_basis.prim_factor != NULL) { qmckl_exit_code rc = qmckl_free(context, ctx->ao_basis.prim_factor); if (rc != QMCKL_SUCCESS) { return qmckl_failwith( context, rc, "qmckl_set_ao_basis_prim_factor", NULL); } } qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = prim_num * sizeof(double); double* new_array = (double*) qmckl_malloc(context, mem_info); if (new_array == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_set_ao_basis_prim_factor", NULL); } memcpy(new_array, prim_factor, mem_info.size); ctx->ao_basis.prim_factor = new_array; <> } #+end_src When the basis set is completely entered, other data structures are computed to accelerate the calculations. The primitives within each contraction are sorted in ascending order of their exponents, such that as soon as a primitive is zero all the following functions vanish. Also, it is possible to compute a nuclear radius beyond which all the primitives are zero up to the numerical accuracy defined in the context. #+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none qmckl_exit_code qmckl_finalize_basis(qmckl_context context); #+end_src #+begin_src c :comments org :tangle (eval c) :noweb yes :exports none qmckl_exit_code qmckl_finalize_basis(qmckl_context context) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return QMCKL_INVALID_CONTEXT; } qmckl_context_struct* const ctx = (qmckl_context_struct* const) context; assert (ctx != NULL); int64_t nucl_num = 0; qmckl_exit_code rc = QMCKL_FAILURE; rc = qmckl_get_nucleus_num(context, &nucl_num); if (rc != QMCKL_SUCCESS) return rc; /* nucleus_prim_index */ { qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = (ctx->nucleus.num + (int64_t) 1) * sizeof(int64_t); ctx->ao_basis.nucleus_prim_index = (int64_t*) qmckl_malloc(context, mem_info); if (ctx->ao_basis.nucleus_prim_index == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "ao_basis.nucleus_prim_index", NULL); } for (int64_t i=0 ; iao_basis.nucleus_index[i]; ctx->ao_basis.nucleus_prim_index[i] = ctx->ao_basis.shell_prim_index[shell_idx]; } ctx->ao_basis.nucleus_prim_index[nucl_num] = ctx->ao_basis.prim_num; } /* Normalize coefficients */ { qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = ctx->ao_basis.prim_num * sizeof(double); ctx->ao_basis.coefficient_normalized = (double *) qmckl_malloc(context, mem_info); if (ctx->ao_basis.coefficient_normalized == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "ao_basis.coefficient_normalized", NULL); } for (int64_t ishell=0 ; ishell < ctx->ao_basis.shell_num ; ++ishell) { for (int64_t iprim=ctx->ao_basis.shell_prim_index[ishell] ; iprim < ctx->ao_basis.shell_prim_index[ishell]+ctx->ao_basis.shell_prim_num[ishell] ; ++iprim) { ctx->ao_basis.coefficient_normalized[iprim] = ctx->ao_basis.coefficient[iprim] * ctx->ao_basis.prim_factor[iprim] * ctx->ao_basis.shell_factor[ishell]; } } } /* Find max angular momentum on each nucleus */ { qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = ctx->nucleus.num * sizeof(int32_t); ctx->ao_basis.nucleus_max_ang_mom = (int32_t *) qmckl_malloc(context, mem_info); if (ctx->ao_basis.nucleus_max_ang_mom == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "ao_basis.nucleus_max_ang_mom", NULL); } for (int64_t inucl=0 ; inucl < ctx->nucleus.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] ; ++ishell) { ctx->ao_basis.nucleus_max_ang_mom[inucl] = ctx->ao_basis.nucleus_max_ang_mom[inucl] > ctx->ao_basis.shell_ang_mom[ishell] ? ctx->ao_basis.nucleus_max_ang_mom[inucl] : ctx->ao_basis.shell_ang_mom[ishell] ; } } } /* Find distance beyond which all AOs are zero. The distance is obtained by sqrt(log(cutoff)*range) */ { if (ctx->ao_basis.type == 'G') { qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = ctx->nucleus.num * sizeof(double); ctx->ao_basis.nucleus_range = (double *) qmckl_malloc(context, mem_info); if (ctx->ao_basis.nucleus_range == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "ao_basis.nucleus_range", NULL); } for (int64_t inucl=0 ; inucl < ctx->nucleus.num ; ++inucl) { ctx->ao_basis.nucleus_range[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] ; ++ishell) { for (int64_t iprim=ctx->ao_basis.shell_prim_index[ishell] ; iprim < ctx->ao_basis.shell_prim_index[ishell] + ctx->ao_basis.shell_prim_num[ishell] ; ++iprim) { double range = 1./ctx->ao_basis.exponent[iprim]; ctx->ao_basis.nucleus_range[inucl] = ctx->ao_basis.nucleus_range[inucl] > range ? ctx->ao_basis.nucleus_range[inucl] : range; } } } } } /* TODO : sort the basis set here */ return QMCKL_SUCCESS; } #+end_src ** Fortran interfaces #+begin_src f90 :tangle (eval fh_func) :comments org :exports none interface integer(c_int32_t) function qmckl_set_ao_basis_type (context, t) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context character(c_char) , intent(in) , value :: t end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_shell_num(context, num) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) , value :: num end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_prim_num(context, num) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) , value :: num end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_nucleus_index(context, idx) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) :: idx(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_nucleus_shell_num(context,shell_num) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) :: shell_num(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_shell_ang_mom(context,shell_ang_mom) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int32_t) , intent(in) :: shell_ang_mom(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_shell_prim_num(context,shell_prim_num) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) :: shell_prim_num(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_shell_prim_index(context,shell_prim_index) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) :: shell_prim_index(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_shell_factor(context,shell_factor) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double) , intent(in) :: shell_factor(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_exponent(context,exponent) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double) , intent(in) :: exponent(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_coefficient(context,coefficient) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double) , intent(in) :: coefficient(*) end function end interface interface integer(c_int32_t) function qmckl_set_ao_basis_prim_factor(context,prim_factor) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double) , intent(in) :: prim_factor(*) end function end interface #+end_src ** Test :noexport: #+begin_src c :tangle (eval c_test) :exports none :exports none const int64_t nucl_num = chbrclf_nucl_num; const double* nucl_charge = chbrclf_charge; const double* nucl_coord = &(chbrclf_nucl_coord[0][0]); qmckl_exit_code rc; rc = qmckl_set_nucleus_num (context, nucl_num); assert(rc == QMCKL_SUCCESS); rc = qmckl_set_nucleus_coord (context, 'T', &(nucl_coord[0])); assert(rc == QMCKL_SUCCESS); rc = qmckl_set_nucleus_charge(context, nucl_charge); assert(rc == QMCKL_SUCCESS); assert(qmckl_nucleus_provided(context)); const int64_t shell_num = chbrclf_shell_num; const int64_t prim_num = chbrclf_prim_num; const int64_t * nucleus_index = &(chbrclf_basis_nucleus_index[0]); const int64_t * nucleus_shell_num = &(chbrclf_basis_nucleus_shell_num[0]); const int32_t * shell_ang_mom = &(chbrclf_basis_shell_ang_mom[0]); const int64_t * shell_prim_num = &(chbrclf_basis_shell_prim_num[0]); const int64_t * shell_prim_index = &(chbrclf_basis_shell_prim_index[0]); const double * shell_factor = &(chbrclf_basis_shell_factor[0]); const double * exponent = &(chbrclf_basis_exponent[0]); const double * coefficient = &(chbrclf_basis_coefficient[0]); const double * prim_factor = &(chbrclf_basis_prim_factor[0]); char typ = 'G'; assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_type (context, typ); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_shell_num (context, shell_num); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_prim_num (context, prim_num); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_nucleus_index (context, nucleus_index); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_nucleus_shell_num (context, nucleus_shell_num); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_shell_ang_mom (context, shell_ang_mom); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_shell_factor (context, shell_factor); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_shell_prim_num (context, shell_prim_num); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_shell_prim_index (context, shell_prim_index); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_exponent (context, exponent); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_coefficient (context, coefficient); assert(rc == QMCKL_SUCCESS); assert(!qmckl_ao_basis_provided(context)); rc = qmckl_set_ao_basis_prim_factor (context, prim_factor); assert(rc == QMCKL_SUCCESS); assert(qmckl_ao_basis_provided(context)); #+end_src * Radial part ** TODO Helper functions to accelerate calculations ** General functions for Gaussian basis functions ~qmckl_ao_gaussian_vgl~ computes the values, gradients and Laplacians at a given point of ~n~ Gaussian functions centered at the same point: \[ v_i = \exp(-a_i |X-R|^2) \] \[ \nabla_x v_i = -2 a_i (X_x - R_x) v_i \] \[ \nabla_y v_i = -2 a_i (X_y - R_y) v_i \] \[ \nabla_z v_i = -2 a_i (X_z - R_z) v_i \] \[ \Delta v_i = a_i (4 |X-R|^2 a_i - 6) v_i \] |--------------+--------+------------------------------------------------------| | ~context~ | input | Global state | | ~X(3)~ | input | Array containing the coordinates of the points | | ~R(3)~ | input | Array containing the x,y,z coordinates of the center | | ~n~ | input | Number of computed Gaussians | | ~A(n)~ | input | Exponents of the Gaussians | | ~VGL(ldv,5)~ | output | Value, gradients and Laplacian of the Gaussians | | ~ldv~ | input | Leading dimension of array ~VGL~ | |--------------+--------+------------------------------------------------------| Requirements - ~context~ is not 0 - ~n~ > 0 - ~ldv~ >= 5 - ~A(i)~ > 0 for all ~i~ - ~X~ is allocated with at least $3 \times 8$ bytes - ~R~ is allocated with at least $3 \times 8$ bytes - ~A~ is allocated with at least $n \times 8$ bytes - ~VGL~ is allocated with at least $n \times 5 \times 8$ bytes #+begin_src c :tangle (eval h_func) qmckl_exit_code qmckl_ao_gaussian_vgl(const qmckl_context context, const double *X, const double *R, const int64_t *n, const int64_t *A, const double *VGL, const int64_t ldv); #+end_src #+begin_src f90 :tangle (eval f) integer function qmckl_ao_gaussian_vgl_f(context, X, R, n, A, VGL, ldv) result(info) use qmckl implicit none integer*8 , intent(in) :: context real*8 , intent(in) :: X(3), R(3) integer*8 , intent(in) :: n real*8 , intent(in) :: A(n) real*8 , intent(out) :: VGL(ldv,5) integer*8 , intent(in) :: ldv integer*8 :: i,j real*8 :: Y(3), r2, t, u, v info = QMCKL_SUCCESS if (context == QMCKL_NULL_CONTEXT) then info = QMCKL_INVALID_CONTEXT return endif if (n <= 0) then info = QMCKL_INVALID_ARG_4 return endif if (ldv < n) then info = QMCKL_INVALID_ARG_7 return endif do i=1,3 Y(i) = X(i) - R(i) end do r2 = Y(1)*Y(1) + Y(2)*Y(2) + Y(3)*Y(3) do i=1,n VGL(i,1) = dexp(-A(i) * r2) end do do i=1,n VGL(i,5) = A(i) * VGL(i,1) end do t = -2.d0 * ( X(1) - R(1) ) u = -2.d0 * ( X(2) - R(2) ) v = -2.d0 * ( X(3) - R(3) ) do i=1,n VGL(i,2) = t * VGL(i,5) VGL(i,3) = u * VGL(i,5) VGL(i,4) = v * VGL(i,5) end do t = 4.d0 * r2 do i=1,n VGL(i,5) = (t * A(i) - 6.d0) * VGL(i,5) end do end function qmckl_ao_gaussian_vgl_f #+end_src #+begin_src f90 :tangle (eval f) :exports none integer(c_int32_t) function qmckl_ao_gaussian_vgl(context, X, R, n, A, VGL, ldv) & bind(C) result(info) use, intrinsic :: iso_c_binding implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double) , intent(in) :: X(3), R(3) integer (c_int64_t) , intent(in) , value :: n real (c_double) , intent(in) :: A(n) real (c_double) , intent(out) :: VGL(ldv,5) integer (c_int64_t) , intent(in) , value :: ldv integer, external :: qmckl_ao_gaussian_vgl_f info = qmckl_ao_gaussian_vgl_f(context, X, R, n, A, VGL, ldv) end function qmckl_ao_gaussian_vgl #+end_src #+begin_src f90 :tangle (eval fh_func) :exports none interface integer(c_int32_t) function qmckl_ao_gaussian_vgl(context, X, R, n, A, VGL, ldv) & bind(C) use, intrinsic :: iso_c_binding integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) , value :: ldv integer (c_int64_t) , intent(in) , value :: n real (c_double) , intent(in) :: X(3), R(3), A(n) real (c_double) , intent(out) :: VGL(ldv,5) end function qmckl_ao_gaussian_vgl end interface #+end_src # Test #+begin_src f90 :tangle (eval f_test) integer(c_int32_t) function test_qmckl_ao_gaussian_vgl(context) bind(C) use qmckl implicit none integer(c_int64_t), intent(in), value :: context integer*8 :: n, ldv, j, i double precision :: X(3), R(3), Y(3), r2 double precision, allocatable :: VGL(:,:), A(:) double precision :: epsilon epsilon = qmckl_get_numprec_epsilon(context) X = (/ 1.1 , 2.2 , 3.3 /) R = (/ 0.1 , 1.2 , -2.3 /) Y(:) = X(:) - R(:) r2 = Y(1)**2 + Y(2)**2 + Y(3)**2 n = 10; ldv = 100; allocate (A(n), VGL(ldv,5)) do i=1,n A(i) = 0.0013 * dble(ishft(1,i)) end do test_qmckl_ao_gaussian_vgl = & qmckl_ao_gaussian_vgl(context, X, R, n, A, VGL, ldv) if (test_qmckl_ao_gaussian_vgl /= 0) return test_qmckl_ao_gaussian_vgl = -1 do i=1,n test_qmckl_ao_gaussian_vgl = -11 if (dabs(1.d0 - VGL(i,1) / (& dexp(-A(i) * r2) & )) > epsilon ) return test_qmckl_ao_gaussian_vgl = -12 if (dabs(1.d0 - VGL(i,2) / (& -2.d0 * A(i) * Y(1) * dexp(-A(i) * r2) & )) > epsilon ) return test_qmckl_ao_gaussian_vgl = -13 if (dabs(1.d0 - VGL(i,3) / (& -2.d0 * A(i) * Y(2) * dexp(-A(i) * r2) & )) > epsilon ) return test_qmckl_ao_gaussian_vgl = -14 if (dabs(1.d0 - VGL(i,4) / (& -2.d0 * A(i) * Y(3) * dexp(-A(i) * r2) & )) > epsilon ) return test_qmckl_ao_gaussian_vgl = -15 if (dabs(1.d0 - VGL(i,5) / (& A(i) * (4.d0*r2*A(i) - 6.d0) * dexp(-A(i) * r2) & )) > epsilon ) return end do test_qmckl_ao_gaussian_vgl = 0 deallocate(VGL) end function test_qmckl_ao_gaussian_vgl #+end_src #+begin_src c :tangle (eval c_test) :exports none int test_qmckl_ao_gaussian_vgl(qmckl_context context); assert(0 == test_qmckl_ao_gaussian_vgl(context)); #+end_src ** TODO General functions for Slater basis functions ** TODO General functions for Radial functions on a grid ** Computation of primitives *** Get #+begin_src c :comments org :tangle (eval h_func) :noweb yes qmckl_exit_code qmckl_get_ao_basis_primitive_vgl(qmckl_context context, double* const primitive_vgl); #+end_src #+begin_src c :comments org :tangle (eval c) :noweb yes :exports none qmckl_exit_code qmckl_get_ao_basis_primitive_vgl(qmckl_context context, double* const primitive_vgl) { if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return QMCKL_NULL_CONTEXT; } qmckl_exit_code rc; rc = qmckl_provide_ao_basis_primitive_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.prim_num * 5 * ctx->electron.num * ctx->electron.walk_num; memcpy(primitive_vgl, ctx->ao_basis.primitive_vgl, sze * sizeof(double)); return QMCKL_SUCCESS; } #+end_src *** Provide #+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none qmckl_exit_code qmckl_provide_ao_basis_primitive_vgl(qmckl_context context); #+end_src #+begin_src c :comments org :tangle (eval c) :noweb yes :exports none qmckl_exit_code qmckl_provide_ao_basis_primitive_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_basis_primitive_vgl", NULL); } /* Compute if necessary */ if (ctx->electron.coord_new_date > ctx->ao_basis.primitive_vgl_date) { /* Allocate array */ if (ctx->ao_basis.primitive_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 * ctx->electron.walk_num * sizeof(double); double* primitive_vgl = (double*) qmckl_malloc(context, mem_info); if (primitive_vgl == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_ao_basis_primitive_vgl", NULL); } ctx->ao_basis.primitive_vgl = primitive_vgl; } qmckl_exit_code rc; if (ctx->ao_basis.type == 'G') { rc = qmckl_compute_ao_basis_primitive_gaussian_vgl(context, ctx->ao_basis.prim_num, ctx->electron.num, ctx->nucleus.num, ctx->electron.walk_num, ctx->ao_basis.nucleus_prim_index, ctx->electron.coord_new, ctx->nucleus.coord, ctx->ao_basis.exponent, ctx->ao_basis.primitive_vgl); } else { return qmckl_failwith( context, QMCKL_FAILURE, "compute_ao_basis_primitive_vgl", "Not yet implemented"); } if (rc != QMCKL_SUCCESS) { return rc; } ctx->ao_basis.primitive_vgl_date = ctx->date; } return QMCKL_SUCCESS; } #+end_src *** Compute :PROPERTIES: :Name: qmckl_compute_ao_basis_primitive_gaussian_vgl :CRetType: qmckl_exit_code :FRetType: qmckl_exit_code :END: #+NAME: qmckl_ao_basis_primitive_gaussian_vgl_args | qmckl_context | context | in | Global state | | int64_t | prim_num | in | Number of primitives | | int64_t | elec_num | in | Number of electrons | | int64_t | nucl_num | in | Number of nuclei | | int64_t | walk_num | in | Number of walkers | | int64_t | nucleus_prim_index[nucl_num] | in | Index of the 1st primitive of each nucleus | | double | elec_coord[walk_num][3][elec_num] | in | Electron coordinates | | double | nucl_coord[3][elec_num] | in | Nuclear coordinates | | double | expo[prim_num] | in | Exponents of the primitives | | double | primitive_vgl[5][walk_num][elec_num][prim_num] | out | Value, gradients and Laplacian of the primitives | #+begin_src f90 :comments org :tangle (eval f) :noweb yes integer function qmckl_compute_ao_basis_primitive_gaussian_vgl_f(context, & prim_num, elec_num, nucl_num, walk_num, & nucleus_prim_index, elec_coord, nucl_coord, expo, primitive_vgl) & result(info) use qmckl implicit none integer(qmckl_context), intent(in) :: context integer*8 , intent(in) :: prim_num integer*8 , intent(in) :: nucl_num integer*8 , intent(in) :: elec_num integer*8 , intent(in) :: walk_num integer*8 , intent(in) :: nucleus_prim_index(nucl_num+1) double precision , intent(in) :: elec_coord(elec_num,3,walk_num) double precision , intent(in) :: nucl_coord(nucl_num,3) double precision , intent(in) :: expo(prim_num) double precision , intent(out) :: primitive_vgl(prim_num,elec_num,walk_num,5) integer*8 :: inucl, iprim, iwalk, ielec double precision :: x, y, z, two_a, ar2, r2, v, cutoff info = QMCKL_SUCCESS ! Don't compute exponentials when the result will be almost zero. cutoff = -dlog(1.d-15) do inucl=1,nucl_num ! C is zero-based, so shift bounds by one do iprim = nucleus_prim_index(inucl)+1, nucleus_prim_index(inucl+1) do iwalk = 1, walk_num do ielec = 1, elec_num x = elec_coord(ielec,1,iwalk) - nucl_coord(inucl,1) y = elec_coord(ielec,2,iwalk) - nucl_coord(inucl,2) z = elec_coord(ielec,3,iwalk) - nucl_coord(inucl,3) r2 = x*x + y*y + z*z ar2 = expo(iprim)*r2 if (ar2 > cutoff) cycle v = dexp(-ar2) two_a = -2.d0 * expo(iprim) * v primitive_vgl(iprim, ielec, iwalk, 1) = v primitive_vgl(iprim, ielec, iwalk, 2) = two_a * x primitive_vgl(iprim, ielec, iwalk, 3) = two_a * y primitive_vgl(iprim, ielec, iwalk, 4) = two_a * z primitive_vgl(iprim, ielec, iwalk, 5) = two_a * (3.d0 - 2.d0*ar2) end do end do end do end do end function qmckl_compute_ao_basis_primitive_gaussian_vgl_f #+end_src #+begin_src c :tangle (eval h_private_func) :comments org :exports none qmckl_exit_code qmckl_compute_ao_basis_primitive_gaussian_vgl( const qmckl_context context, const int64_t prim_num, const int64_t elec_num, const int64_t nucl_num, const int64_t walk_num, const int64_t* nucleus_prim_index, const double* elec_coord, const double* nucl_coord, const double* expo, double* const primitive_vgl); #+end_src #+CALL: generate_c_interface(table=qmckl_ao_basis_primitive_gaussian_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_basis_primitive_gaussian_vgl")) #+RESULTS: #+begin_src f90 :tangle (eval f) :comments org :exports none integer(c_int32_t) function qmckl_compute_ao_basis_primitive_gaussian_vgl & (context, prim_num, elec_num, nucl_num, walk_num, nucleus_prim_index, elec_coord, nucl_coord, expo, primitive_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 :: prim_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 integer (c_int64_t) , intent(in) :: nucleus_prim_index(nucl_num) real (c_double ) , intent(in) :: elec_coord(elec_num,3,walk_num) real (c_double ) , intent(in) :: nucl_coord(elec_num,3) real (c_double ) , intent(in) :: expo(prim_num) real (c_double ) , intent(out) :: primitive_vgl(elec_num,walk_num,5,prim_num) integer(c_int32_t), external :: qmckl_compute_ao_basis_primitive_gaussian_vgl_f info = qmckl_compute_ao_basis_primitive_gaussian_vgl_f & (context, prim_num, elec_num, nucl_num, walk_num, nucleus_prim_index, elec_coord, nucl_coord, expo, primitive_vgl) end function qmckl_compute_ao_basis_primitive_gaussian_vgl #+end_src #+begin_src python :results output :exports none import numpy as np def f(a,x,y): return np.exp( -a*(np.linalg.norm(x-y))**2 ) 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( [ 1.096243353458458e+00, 8.907054016973815e-01, 7.777092280258892e-01 ] ) nucl_2 = np.array( [ 1.168459237342663e+00, 1.125660720053393e+00, 2.833370314829343e+00 ] ) #double prim_vgl[prim_num][5][walk_num][elec_num]; a = 0.9059; x = elec_26_w1 ; y = nucl_1 print ( "[7][0][0][26] : %e"% f(a,x,y)) print ( "[7][1][0][26] : %e"% df(a,x,y,1)) print ( "[7][2][0][26] : %e"% df(a,x,y,2)) print ( "[7][3][0][26] : %e"% df(a,x,y,3)) print ( "[7][4][0][26] : %e"% lf(a,x,y)) a = 0.32578; x = elec_15_w2 ; y = nucl_2 print ( "[39][0][1][15] : %e"% f(a,x,y)) print ( "[39][1][1][15] : %e"% df(a,x,y,1)) print ( "[39][2][1][15] : %e"% df(a,x,y,2)) print ( "[39][3][1][15] : %e"% df(a,x,y,3)) print ( "[39][4][1][15] : %e"% lf(a,x,y)) #+end_src #+RESULTS: #+begin_example [7][0][0][26] : 1.050157e-03 [7][1][0][26] : -7.501497e-04 [7][2][0][26] : -3.825069e-03 [7][3][0][26] : 3.495056e-03 [7][4][0][26] : 2.040013e-02 [39][0][1][15] : 1.083038e-03 [39][1][1][15] : 2.378275e-03 [39][2][1][15] : 2.143086e-03 [39][3][1][15] : 4.332750e-04 [39][4][1][15] : 7.514605e-03 #+end_example *** Test #+begin_src c :tangle (eval c_test) :exports none { #define walk_num chbrclf_walk_num #define elec_num chbrclf_elec_num #define prim_num chbrclf_prim_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 prim_vgl[5][walk_num][elec_num][prim_num]; rc = qmckl_get_ao_basis_primitive_vgl(context, &(prim_vgl[0][0][0][0])); assert (rc == QMCKL_SUCCESS); assert( fabs(prim_vgl[0][0][26][7] - ( 1.0501570432064878E-003)) < 1.e-14 ); assert( fabs(prim_vgl[1][0][26][7] - (-7.5014974095310560E-004)) < 1.e-14 ); assert( fabs(prim_vgl[2][0][26][7] - (-3.8250692897610380E-003)) < 1.e-14 ); assert( fabs(prim_vgl[3][0][26][7] - ( 3.4950559194080275E-003)) < 1.e-14 ); assert( fabs(prim_vgl[4][0][26][7] - ( 2.0392163767356572E-002)) < 1.e-14 ); assert( fabs(prim_vgl[0][1][15][39] - ( 1.0825844173157661E-003)) < 1.e-14 ); assert( fabs(prim_vgl[1][1][15][39] - ( 2.3774237611651531E-003)) < 1.e-14 ); assert( fabs(prim_vgl[2][1][15][39] - ( 2.1423191526963063E-003)) < 1.e-14 ); assert( fabs(prim_vgl[3][1][15][39] - ( 4.3312003523048492E-004)) < 1.e-14 ); assert( fabs(prim_vgl[4][1][15][39] - ( 7.5174404780004771E-003)) < 1.e-14 ); } #+end_src *** Ideas for improvement #+begin_src c // m : walkers // j : electrons // l : primitives k=0; for (m=0 ; mao_basis.shell_num * 5 * ctx->electron.num * ctx->electron.walk_num; memcpy(shell_vgl, ctx->ao_basis.shell_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_basis_shell_vgl (context, shell_vgl) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context double precision, intent(out) :: shell_vgl(*) end function end interface #+end_src *** Provide #+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none qmckl_exit_code qmckl_provide_ao_basis_shell_vgl(qmckl_context context); #+end_src #+begin_src c :comments org :tangle (eval c) :noweb yes :exports none qmckl_exit_code qmckl_provide_ao_basis_shell_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_basis_shell_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.shell_vgl_date) { /* Allocate array */ 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 * ctx->electron.walk_num * sizeof(double); double* shell_vgl = (double*) qmckl_malloc(context, mem_info); if (shell_vgl == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_ao_basis_shell_vgl", NULL); } ctx->ao_basis.shell_vgl = shell_vgl; } qmckl_exit_code rc; if (ctx->ao_basis.type == 'G') { rc = qmckl_compute_ao_basis_shell_gaussian_vgl(context, ctx->ao_basis.prim_num, ctx->ao_basis.shell_num, ctx->electron.num, ctx->nucleus.num, ctx->electron.walk_num, ctx->ao_basis.nucleus_shell_num, ctx->ao_basis.nucleus_index, ctx->ao_basis.shell_prim_index, ctx->ao_basis.shell_prim_num, ctx->electron.coord_new, ctx->nucleus.coord, ctx->ao_basis.exponent, ctx->ao_basis.coefficient_normalized, ctx->ao_basis.shell_vgl); } else { return qmckl_failwith( context, QMCKL_FAILURE, "compute_ao_basis_shell_vgl", "Not yet implemented"); } if (rc != QMCKL_SUCCESS) { return rc; } ctx->ao_basis.shell_vgl_date = ctx->date; } return QMCKL_SUCCESS; } #+end_src *** Compute :PROPERTIES: :Name: qmckl_compute_ao_basis_shell_gaussian_vgl :CRetType: qmckl_exit_code :FRetType: qmckl_exit_code :END: #+NAME: qmckl_ao_basis_shell_gaussian_vgl_args | ~qmckl_context~ | ~context~ | in | Global state | | ~int64_t~ | ~prim_num~ | in | Number of primitives | | ~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 | | ~int64_t~ | ~nucleus_shell_num[nucl_num]~ | in | Number of shells for each nucleus | | ~int64_t~ | ~nucleus_index[nucl_num]~ | in | Index of the 1st shell of each nucleus | | ~int64_t~ | ~shell_prim_index[shell_num]~ | in | Index of the 1st primitive of each shell | | ~int64_t~ | ~shell_prim_num[shell_num]~ | in | Number of primitives per shell | | ~double~ | ~elec_coord[walk_num][3][elec_num]~ | in | Electron coordinates | | ~double~ | ~nucl_coord[3][elec_num]~ | in | Nuclear coordinates | | ~double~ | ~expo[prim_num]~ | in | Exponents of the primitives | | ~double~ | ~coef_normalized[prim_num]~ | in | Coefficients of the primitives | | ~double~ | ~shell_vgl[5][walk_num][elec_num][shell_num]~ | out | Value, gradients and Laplacian of the shells | #+begin_src f90 :comments org :tangle (eval f) :noweb yes integer function qmckl_compute_ao_basis_shell_gaussian_vgl_f(context, & prim_num, shell_num, elec_num, nucl_num, walk_num, & nucleus_shell_num, nucleus_index, shell_prim_index, shell_prim_num, & elec_coord, nucl_coord, expo, coef_normalized, shell_vgl) & result(info) use qmckl implicit none integer(qmckl_context), intent(in) :: context integer*8 , intent(in) :: prim_num integer*8 , intent(in) :: shell_num integer*8 , intent(in) :: nucl_num integer*8 , intent(in) :: elec_num integer*8 , intent(in) :: walk_num integer*8 , intent(in) :: nucleus_shell_num(nucl_num) integer*8 , intent(in) :: nucleus_index(nucl_num) integer*8 , intent(in) :: shell_prim_index(shell_num) integer*8 , intent(in) :: shell_prim_num(shell_num) double precision , intent(in) :: elec_coord(elec_num,3,walk_num) double precision , intent(in) :: nucl_coord(nucl_num,3) double precision , intent(in) :: expo(prim_num) double precision , intent(in) :: coef_normalized(prim_num) double precision , intent(out) :: shell_vgl(shell_num,elec_num,walk_num,5) integer*8 :: inucl, iprim, iwalk, ielec, ishell double precision :: x, y, z, two_a, ar2, r2, v, cutoff info = QMCKL_SUCCESS ! Don't compute exponentials when the result will be almost zero. ! TODO : Use numerical precision here cutoff = -dlog(1.d-15) do inucl=1,nucl_num do iwalk = 1, walk_num do ielec = 1, elec_num x = elec_coord(ielec,1,iwalk) - nucl_coord(inucl,1) y = elec_coord(ielec,2,iwalk) - nucl_coord(inucl,2) z = elec_coord(ielec,3,iwalk) - nucl_coord(inucl,3) r2 = x*x + y*y + z*z do ishell=nucleus_index(inucl)+1, nucleus_index(inucl)+nucleus_shell_num(inucl) ! C is zero-based, so shift bounds by one shell_vgl(ishell, ielec, iwalk, 1) = 0.d0 shell_vgl(ishell, ielec, iwalk, 2) = 0.d0 shell_vgl(ishell, ielec, iwalk, 3) = 0.d0 shell_vgl(ishell, ielec, iwalk, 4) = 0.d0 shell_vgl(ishell, ielec, iwalk, 5) = 0.d0 do iprim = shell_prim_index(ishell)+1, shell_prim_index(ishell)+shell_prim_num(ishell) ar2 = expo(iprim)*r2 if (ar2 > cutoff) then cycle end if v = coef_normalized(iprim) * dexp(-ar2) two_a = -2.d0 * expo(iprim) * v shell_vgl(ishell, ielec, iwalk, 1) = & shell_vgl(ishell, ielec, iwalk, 1) + v shell_vgl(ishell, ielec, iwalk, 2) = & shell_vgl(ishell, ielec, iwalk, 2) + two_a * x shell_vgl(ishell, ielec, iwalk, 3) = & shell_vgl(ishell, ielec, iwalk, 3) + two_a * y shell_vgl(ishell, ielec, iwalk, 4) = & shell_vgl(ishell, ielec, iwalk, 4) + two_a * z shell_vgl(ishell, ielec, iwalk, 5) = & shell_vgl(ishell, ielec, iwalk, 5) + two_a * (3.d0 - 2.d0*ar2) end do end do end do end do end do end function qmckl_compute_ao_basis_shell_gaussian_vgl_f #+end_src #+begin_src c :tangle (eval h_private_func) :comments org :exports none 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* nucleus_shell_num, const int64_t* shell_prim_index, const int64_t* nucleus_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); #+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")) #+RESULTS: #+begin_src f90 :tangle (eval f) :comments org :exports none integer(c_int32_t) function qmckl_compute_ao_basis_shell_gaussian_vgl & (context, & prim_num, & shell_num, & elec_num, & nucl_num, & walk_num, & nucleus_shell_num, & nucleus_index, & shell_prim_index, & shell_prim_num, & elec_coord, & nucl_coord, & expo, & coef_normalized, & shell_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 :: prim_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 integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num) integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num) integer (c_int64_t) , intent(in) :: shell_prim_index(shell_num) integer (c_int64_t) , intent(in) :: shell_prim_num(shell_num) real (c_double ) , intent(in) :: elec_coord(elec_num,3,walk_num) real (c_double ) , intent(in) :: nucl_coord(elec_num,3) real (c_double ) , intent(in) :: expo(prim_num) real (c_double ) , intent(in) :: coef_normalized(prim_num) real (c_double ) , intent(out) :: shell_vgl(shell_num,elec_num,walk_num,5) integer(c_int32_t), external :: qmckl_compute_ao_basis_shell_gaussian_vgl_f info = qmckl_compute_ao_basis_shell_gaussian_vgl_f & (context, & prim_num, & shell_num, & elec_num, & nucl_num, & walk_num, & nucleus_shell_num, & nucleus_index, & shell_prim_index, & shell_prim_num, & elec_coord, & nucl_coord, & expo, & coef_normalized, & shell_vgl) end function qmckl_compute_ao_basis_shell_gaussian_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( [ 1.096243353458458e+00, 8.907054016973815e-01, 7.777092280258892e-01 ] ) nucl_2 = np.array( [ 1.168459237342663e+00, 1.125660720053393e+00, 2.833370314829343e+00 ] ) #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 )] 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)) 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)] 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)) #+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 #+end_example *** 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 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 shell_vgl[5][walk_num][elec_num][shell_num]; rc = qmckl_get_ao_basis_shell_vgl(context, &(shell_vgl[0][0][0][0])); assert (rc == QMCKL_SUCCESS); printf(" shell_vgl[1][0][0][26] %25.15e\n", shell_vgl[0][0][26][1]); printf(" shell_vgl[1][1][0][26] %25.15e\n", shell_vgl[1][0][26][1]); printf(" shell_vgl[1][2][0][26] %25.15e\n", shell_vgl[2][0][26][1]); printf(" shell_vgl[1][3][0][26] %25.15e\n", shell_vgl[3][0][26][1]); printf(" shell_vgl[1][4][0][26] %25.15e\n", shell_vgl[4][0][26][1]); printf(" shell_vgl[14][0][1][15] %25.15e\n", shell_vgl[0][1][15][14]); printf(" shell_vgl[14][1][1][15] %25.15e\n", shell_vgl[1][1][15][14]); printf(" shell_vgl[14][2][1][15] %25.15e\n", shell_vgl[2][1][15][14]); printf(" shell_vgl[14][3][1][15] %25.15e\n", shell_vgl[3][1][15][14]); printf(" shell_vgl[14][4][1][15] %25.15e\n", shell_vgl[4][1][15][14]); assert( fabs(shell_vgl[0][0][26][1] - ( 3.564393437193868e-02)) < 1.e-14 ); assert( fabs(shell_vgl[1][0][26][1] - (-6.030177987072189e-03)) < 1.e-14 ); assert( fabs(shell_vgl[2][0][26][1] - (-3.074832579537582e-02)) < 1.e-14 ); assert( fabs(shell_vgl[3][0][26][1] - ( 2.809546963519935e-02)) < 1.e-14 ); assert( fabs(shell_vgl[4][0][26][1] - ( 1.896046117183968e-02)) < 1.e-14 ); assert( fabs(shell_vgl[0][1][15][14] - ( 5.928089771361000e-03)) < 1.e-14 ); assert( fabs(shell_vgl[1][1][15][14] - ( 4.355862296021654e-03)) < 1.e-14 ); assert( fabs(shell_vgl[2][1][15][14] - ( 3.925108924923650e-03)) < 1.e-14 ); assert( fabs(shell_vgl[3][1][15][14] - ( 7.935527784022099e-04)) < 1.e-14 ); assert( fabs(shell_vgl[4][1][15][14] - ( 2.708246573703548e-03)) < 1.e-14 ); } #+end_src * Polynomial part ** General functions for Powers of $x-X_i$ :PROPERTIES: :Name: qmckl_ao_power :CRetType: qmckl_exit_code :FRetType: qmckl_exit_code :END: The ~qmckl_ao_power~ function computes all the powers of the ~n~ input data up to the given maximum value given in input for each of the $n$ points: \[ P_{ik} = X_i^k \] #+NAME: qmckl_ao_power_args | qmckl_context | context | in | Global state | | int64_t | n | in | Number of values | | double | X[n] | in | Array containing the input values | | int32_t | LMAX[n] | in | Array containing the maximum power for each value | | double | P[n][ldp] | out | Array containing all the powers of ~X~ | | int64_t | ldp | in | Leading dimension of array ~P~ | *** Requirements - ~context~ is not ~QMCKL_NULL_CONTEXT~ - ~n~ > 0 - ~X~ is allocated with at least $n \times 8$ bytes - ~LMAX~ is allocated with at least $n \times 4$ bytes - ~P~ is allocated with at least $n \times \max_i \text{LMAX}_i \times 8$ bytes - ~LDP~ >= $\max_i$ ~LMAX[i]~ *** C Header #+CALL: generate_c_header(table=qmckl_ao_power_args,rettyp=get_value("CRetType"),fname="qmckl_ao_power") #+RESULTS: #+begin_src c :tangle (eval h_func) :comments org qmckl_exit_code qmckl_ao_power ( const qmckl_context context, const int64_t n, const double* X, const int32_t* LMAX, double* const P, const int64_t ldp ); #+end_src *** Source #+begin_src f90 :tangle (eval f) integer function qmckl_ao_power_f(context, n, X, LMAX, P, ldp) result(info) use qmckl implicit none integer*8 , intent(in) :: context integer*8 , intent(in) :: n real*8 , intent(in) :: X(n) integer , intent(in) :: LMAX(n) real*8 , intent(out) :: P(ldp,n) integer*8 , intent(in) :: ldp integer*8 :: i,k info = QMCKL_SUCCESS if (context == QMCKL_NULL_CONTEXT) then info = QMCKL_INVALID_CONTEXT return endif if (n <= ldp) then info = QMCKL_INVALID_ARG_2 return endif k = MAXVAL(LMAX) if (LDP < k) then info = QMCKL_INVALID_ARG_6 return endif if (k <= 0) then info = QMCKL_INVALID_ARG_4 return endif do i=1,n P(1,i) = X(i) do k=2,LMAX(i) P(k,i) = P(k-1,i) * X(i) end do end do end function qmckl_ao_power_f #+end_src *** C interface #+CALL: generate_c_interface(table=qmckl_ao_power_args,rettyp=get_value("CRetType"),fname="qmckl_ao_power") #+RESULTS: #+begin_src f90 :tangle (eval f) :comments org :exports none integer(c_int32_t) function qmckl_ao_power & (context, n, X, LMAX, P, ldp) & 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 :: n real (c_double ) , intent(in) :: X(n) integer (c_int32_t) , intent(in) :: LMAX(n) real (c_double ) , intent(out) :: P(ldp,n) integer (c_int64_t) , intent(in) , value :: ldp integer(c_int32_t), external :: qmckl_ao_power_f info = qmckl_ao_power_f & (context, n, X, LMAX, P, ldp) end function qmckl_ao_power #+end_src *** Fortran interface #+CALL: generate_f_interface(table=qmckl_ao_power_args,rettyp=get_value("CRetType"),fname="qmckl_ao_power") #+RESULTS: #+begin_src f90 :tangle (eval fh_func) :comments org :exports none interface integer(c_int32_t) function qmckl_ao_power & (context, n, X, LMAX, P, ldp) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context integer (c_int64_t) , intent(in) , value :: n real (c_double ) , intent(in) :: X(n) integer (c_int32_t) , intent(in) :: LMAX(n) real (c_double ) , intent(out) :: P(ldp,n) integer (c_int64_t) , intent(in) , value :: ldp end function qmckl_ao_power end interface #+end_src *** Test #+begin_src f90 :tangle (eval f_test) integer(c_int32_t) function test_qmckl_ao_power(context) bind(C) use qmckl implicit none integer(qmckl_context), intent(in), value :: context integer*8 :: n, LDP integer, allocatable :: LMAX(:) double precision, allocatable :: X(:), P(:,:) integer*8 :: i,j double precision :: epsilon epsilon = qmckl_get_numprec_epsilon(context) n = 100; LDP = 10; allocate(X(n), P(LDP,n), LMAX(n)) do j=1,n X(j) = -5.d0 + 0.1d0 * dble(j) LMAX(j) = 1 + int(mod(j, 5),4) end do test_qmckl_ao_power = qmckl_ao_power(context, n, X, LMAX, P, LDP) if (test_qmckl_ao_power /= QMCKL_SUCCESS) return test_qmckl_ao_power = QMCKL_FAILURE do j=1,n do i=1,LMAX(j) if ( X(j)**i == 0.d0 ) then if ( P(i,j) /= 0.d0) return else if ( dabs(1.d0 - P(i,j) / (X(j)**i)) > epsilon ) return end if end do end do test_qmckl_ao_power = QMCKL_SUCCESS deallocate(X,P,LMAX) end function test_qmckl_ao_power #+end_src #+begin_src c :tangle (eval c_test) :exports none int test_qmckl_ao_power(qmckl_context context); assert(0 == test_qmckl_ao_power(context)); #+end_src ** General functions for Value, Gradient and Laplacian of a polynomial :PROPERTIES: :Name: qmckl_ao_polynomial_vgl :CRetType: qmckl_exit_code :FRetType: qmckl_exit_code :END: A polynomial is centered on a nucleus $\mathbf{R}_i$ \[ P_l(\mathbf{r},\mathbf{R}_i) = (x-X_i)^a (y-Y_i)^b (z-Z_i)^c \] The gradients with respect to electron coordinates are \begin{eqnarray*} \frac{\partial }{\partial x} P_l\left(\mathbf{r},\mathbf{R}_i \right) & = & a (x-X_i)^{a-1} (y-Y_i)^b (z-Z_i)^c \\ \frac{\partial }{\partial y} P_l\left(\mathbf{r},\mathbf{R}_i \right) & = & b (x-X_i)^a (y-Y_i)^{b-1} (z-Z_i)^c \\ \frac{\partial }{\partial z} P_l\left(\mathbf{r},\mathbf{R}_i \right) & = & c (x-X_i)^a (y-Y_i)^b (z-Z_i)^{c-1} \\ \end{eqnarray*} and the Laplacian is \begin{eqnarray*} \left( \frac{\partial }{\partial x^2} + \frac{\partial }{\partial y^2} + \frac{\partial }{\partial z^2} \right) P_l \left(\mathbf{r},\mathbf{R}_i \right) & = & a(a-1) (x-X_i)^{a-2} (y-Y_i)^b (z-Z_i)^c + \\ && b(b-1) (x-X_i)^a (y-Y_i)^{b-1} (z-Z_i)^c + \\ && c(c-1) (x-X_i)^a (y-Y_i)^b (z-Z_i)^{c-1}. \end{eqnarray*} ~qmckl_ao_polynomial_vgl~ computes the values, gradients and Laplacians at a given point in space, of all polynomials with an angular momentum up to ~lmax~. #+NAME: qmckl_ao_polynomial_vgl_args | qmckl_context | context | in | Global state | | double | X[3] | in | Array containing the coordinates of the points | | double | R[3] | in | Array containing the x,y,z coordinates of the center | | int32_t | lmax | in | Maximum angular momentum | | int64_t | n | inout | Number of computed polynomials | | int32_t | L[n][ldl] | out | Contains a,b,c for all ~n~ results | | int64_t | ldl | in | Leading dimension of ~L~ | | double | VGL[n][ldv] | out | Value, gradients and Laplacian of the polynomials | | int64_t | ldv | in | Leading dimension of array ~VGL~ | *** Requirements - ~context~ is not ~QMCKL_NULL_CONTEXT~ - ~n~ > 0 - ~lmax~ >= 0 - ~ldl~ >= 3 - ~ldv~ >= 5 - ~X~ is allocated with at least $3 \times 8$ bytes - ~R~ is allocated with at least $3 \times 8$ bytes - ~n~ >= ~(lmax+1)(lmax+2)(lmax+3)/6~ - ~L~ is allocated with at least $3 \times n \times 4$ bytes - ~VGL~ is allocated with at least $5 \times n \times 8$ bytes - On output, ~n~ should be equal to ~(lmax+1)(lmax+2)(lmax+3)/6~ - On output, the powers are given in the following order (l=a+b+c): - Increasing values of ~l~ - Within a given value of ~l~, alphabetical order of the string made by a*"x" + b*"y" + c*"z" (in Python notation). For example, with a=0, b=2 and c=1 the string is "yyz" *** C Header #+CALL: generate_c_header(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("CRetType"),fname=get_value("Name")) #+RESULTS: #+begin_src c :tangle (eval h_func) :comments org qmckl_exit_code qmckl_ao_polynomial_vgl ( const qmckl_context context, const double* X, const double* R, const int32_t lmax, int64_t* n, int32_t* const L, const int64_t ldl, double* const VGL, const int64_t ldv ); #+end_src *** Source #+begin_src f90 :tangle (eval f) integer function qmckl_ao_polynomial_vgl_f(context, X, R, lmax, n, L, ldl, VGL, ldv) result(info) use qmckl implicit none integer*8 , intent(in) :: context real*8 , intent(in) :: X(3), R(3) integer , intent(in) :: lmax integer*8 , intent(out) :: n integer , intent(out) :: L(ldl,(lmax+1)*(lmax+2)*(lmax+3)/6) integer*8 , intent(in) :: ldl real*8 , intent(out) :: VGL(ldv,(lmax+1)*(lmax+2)*(lmax+3)/6) integer*8 , intent(in) :: ldv integer*8 :: i,j integer :: a,b,c,d 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 info = 0 if (context == QMCKL_NULL_CONTEXT) then info = QMCKL_INVALID_CONTEXT return endif if (lmax < 0) then info = QMCKL_INVALID_ARG_4 return endif if (ldl < 3) then info = QMCKL_INVALID_ARG_7 return endif if (ldv < 5) then info = QMCKL_INVALID_ARG_9 return endif do i=1,3 Y(i) = X(i) - R(i) end do lmax_array(1:3) = lmax if (lmax == 0) then VGL(1,1) = 1.d0 vgL(2:5,1) = 0.d0 l(1:3,1) = 0 n=1 else if (lmax > 0) then pows(-2:0,1:3) = 1.d0 do i=1,lmax pows(i,1) = pows(i-1,1) * Y(1) pows(i,2) = pows(i-1,2) * Y(2) pows(i,3) = pows(i-1,3) * Y(3) end do VGL(1:5,1:4) = 0.d0 l (1:3,1:4) = 0 VGL(1 ,1 ) = 1.d0 vgl(1:5,2:4) = 0.d0 l (1,2) = 1 vgl(1,2) = pows(1,1) vgL(2,2) = 1.d0 l (2,3) = 1 vgl(1,3) = pows(1,2) vgL(3,3) = 1.d0 l (3,4) = 1 vgl(1,4) = pows(1,3) vgL(4,4) = 1.d0 n=4 endif ! l>=2 dd = 2.d0 do d=2,lmax da = dd do a=d,0,-1 db = dd-da do b=d-a,0,-1 c = d - a - b dc = dd - da - db n = n+1 l(1,n) = a l(2,n) = b l(3,n) = c xy = pows(a,1) * pows(b,2) yz = pows(b,2) * pows(c,3) xz = pows(a,1) * pows(c,3) vgl(1,n) = xy * pows(c,3) xy = dc * xy xz = db * xz yz = da * yz vgl(2,n) = pows(a-1,1) * yz vgl(3,n) = pows(b-1,2) * xz vgl(4,n) = pows(c-1,3) * xy vgl(5,n) = & (da-1.d0) * pows(a-2,1) * yz + & (db-1.d0) * pows(b-2,2) * xz + & (dc-1.d0) * pows(c-2,3) * xy db = db - 1.d0 end do da = da - 1.d0 end do dd = dd + 1.d0 end do info = QMCKL_SUCCESS end function qmckl_ao_polynomial_vgl_f #+end_src *** C interface #+CALL: generate_c_interface(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("CRetType"),fname=get_value("Name")) #+RESULTS: #+begin_src f90 :tangle (eval f) :comments org :exports none integer(c_int32_t) function qmckl_ao_polynomial_vgl & (context, X, R, lmax, n, L, ldl, VGL, ldv) & bind(C) result(info) use, intrinsic :: iso_c_binding implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double ) , intent(in) :: X(3) real (c_double ) , intent(in) :: R(3) integer (c_int32_t) , intent(in) , value :: lmax integer (c_int64_t) , intent(inout) :: n integer (c_int32_t) , intent(out) :: L(ldl,n) integer (c_int64_t) , intent(in) , value :: ldl real (c_double ) , intent(out) :: VGL(ldv,n) integer (c_int64_t) , intent(in) , value :: ldv integer(c_int32_t), external :: qmckl_ao_polynomial_vgl_f info = qmckl_ao_polynomial_vgl_f & (context, X, R, lmax, n, L, ldl, VGL, ldv) end function qmckl_ao_polynomial_vgl #+end_src *** Fortran interface #+CALL: generate_f_interface(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("FRetType"),fname=get_value("Name")) #+RESULTS: #+begin_src f90 :tangle (eval fh_func) :comments org :exports none interface integer(c_int32_t) function qmckl_ao_polynomial_vgl & (context, X, R, lmax, n, L, ldl, VGL, ldv) & bind(C) use, intrinsic :: iso_c_binding import implicit none integer (c_int64_t) , intent(in) , value :: context real (c_double ) , intent(in) :: X(3) real (c_double ) , intent(in) :: R(3) integer (c_int32_t) , intent(in) , value :: lmax integer (c_int64_t) , intent(inout) :: n integer (c_int32_t) , intent(out) :: L(ldl,n) integer (c_int64_t) , intent(in) , value :: ldl real (c_double ) , intent(out) :: VGL(ldv,n) integer (c_int64_t) , intent(in) , value :: ldv end function qmckl_ao_polynomial_vgl end interface #+end_src *** Test #+begin_src f90 :tangle (eval f_test) integer(c_int32_t) function test_qmckl_ao_polynomial_vgl(context) bind(C) use qmckl implicit none integer(c_int64_t), intent(in), value :: context integer :: lmax, d, i integer, allocatable :: L(:,:) integer*8 :: n, ldl, ldv, j double precision :: X(3), R(3), Y(3) double precision, allocatable :: VGL(:,:) double precision :: w double precision :: epsilon epsilon = qmckl_get_numprec_epsilon(context) X = (/ 1.1 , 2.2 , 3.3 /) R = (/ 0.1 , 1.2 , -2.3 /) Y(:) = X(:) - R(:) lmax = 4; ldl = 3; ldv = 100; d = (lmax+1)*(lmax+2)*(lmax+3)/6 allocate (L(ldl,d), VGL(ldv,d)) test_qmckl_ao_polynomial_vgl = & qmckl_ao_polynomial_vgl(context, X, R, lmax, n, L, ldl, VGL, ldv) if (test_qmckl_ao_polynomial_vgl /= QMCKL_SUCCESS) return if (n /= d) return do j=1,n test_qmckl_ao_polynomial_vgl = QMCKL_FAILURE do i=1,3 if (L(i,j) < 0) return end do test_qmckl_ao_polynomial_vgl = QMCKL_FAILURE if (dabs(1.d0 - VGL(1,j) / (& Y(1)**L(1,j) * Y(2)**L(2,j) * Y(3)**L(3,j) & )) > epsilon ) return test_qmckl_ao_polynomial_vgl = QMCKL_FAILURE if (L(1,j) < 1) then if (VGL(2,j) /= 0.d0) return else if (dabs(1.d0 - VGL(2,j) / (& L(1,j) * Y(1)**(L(1,j)-1) * Y(2)**L(2,j) * Y(3)**L(3,j) & )) > epsilon ) return end if test_qmckl_ao_polynomial_vgl = QMCKL_FAILURE if (L(2,j) < 1) then if (VGL(3,j) /= 0.d0) return else if (dabs(1.d0 - VGL(3,j) / (& L(2,j) * Y(1)**L(1,j) * Y(2)**(L(2,j)-1) * Y(3)**L(3,j) & )) > epsilon ) return end if test_qmckl_ao_polynomial_vgl = QMCKL_FAILURE if (L(3,j) < 1) then if (VGL(4,j) /= 0.d0) return else if (dabs(1.d0 - VGL(4,j) / (& L(3,j) * Y(1)**L(1,j) * Y(2)**L(2,j) * Y(3)**(L(3,j)-1) & )) > epsilon ) return end if test_qmckl_ao_polynomial_vgl = QMCKL_FAILURE w = 0.d0 if (L(1,j) > 1) then w = w + L(1,j) * (L(1,j)-1) * Y(1)**(L(1,j)-2) * Y(2)**L(2,j) * Y(3)**L(3,j) end if if (L(2,j) > 1) then w = w + L(2,j) * (L(2,j)-1) * Y(1)**L(1,j) * Y(2)**(L(2,j)-2) * Y(3)**L(3,j) end if if (L(3,j) > 1) then w = w + L(3,j) * (L(3,j)-1) * Y(1)**L(1,j) * Y(2)**L(2,j) * Y(3)**(L(3,j)-2) end if if (dabs(1.d0 - VGL(5,j) / w) > epsilon ) return end do test_qmckl_ao_polynomial_vgl = QMCKL_SUCCESS deallocate(L,VGL) end function test_qmckl_ao_polynomial_vgl #+end_src #+begin_src c :tangle (eval c_test) int test_qmckl_ao_polynomial_vgl(qmckl_context context); assert(0 == test_qmckl_ao_polynomial_vgl(context)); #+end_src * Combining radial and polynomial parts * End of files :noexport: #+begin_src c :tangle (eval h_private_type) #endif #+end_src *** Test #+begin_src c :tangle (eval c_test) rc = qmckl_context_destroy(context); assert (rc == QMCKL_SUCCESS); return 0; } #+end_src **✸ Compute file names #+begin_src emacs-lisp ; The following is required to compute the file names (setq pwd (file-name-directory buffer-file-name)) (setq name (file-name-nondirectory (substring buffer-file-name 0 -4))) (setq f (concat pwd name "_f.f90")) (setq fh (concat pwd name "_fh.f90")) (setq c (concat pwd name ".c")) (setq h (concat name ".h")) (setq h_private (concat name "_private.h")) (setq c_test (concat pwd "test_" name ".c")) (setq f_test (concat pwd "test_" name "_f.f90")) ; Minted (require 'ox-latex) (setq org-latex-listings 'minted) (add-to-list 'org-latex-packages-alist '("" "listings")) (add-to-list 'org-latex-packages-alist '("" "color")) #+end_src #+RESULTS: | | color | | | listings | # -*- mode: org -*- # vim: syntax=c