2021-09-15 18:30:41 +02:00
|
|
|
#+TITLE: Molecular Orbitals
|
|
|
|
#+SETUPFILE: ../tools/theme.setup
|
|
|
|
#+INCLUDE: ../tools/lib.org
|
|
|
|
|
2022-01-17 16:09:41 +01:00
|
|
|
The molecular orbitals (MOs) are defined in the basis of AOs along with a AO to MO
|
2021-09-15 18:30:41 +02:00
|
|
|
coefficient matrix \[C\]. Using these coefficients (e.g. from Hartree Fock SCF method)
|
|
|
|
the MOs are defined as follows:
|
|
|
|
|
|
|
|
\[
|
|
|
|
\phi_i(\mathbf{r}) = C_i * \chi_i (\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
|
|
|
|
|
|
|
In this section we demonstrate how to use the QMCkl specific DGEMM
|
|
|
|
function to calculate the MOs.
|
|
|
|
|
|
|
|
|
|
|
|
* Headers :noexport:
|
|
|
|
#+begin_src elisp :noexport :results none
|
|
|
|
(org-babel-lob-ingest "../tools/lib.org")
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
|
2021-10-14 21:40:14 +02:00
|
|
|
#+begin_src c :tangle (eval h_private_func)
|
|
|
|
#ifndef QMCKL_MO_HPF
|
|
|
|
#define QMCKL_MO_HPF
|
|
|
|
#+end_src
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
#+begin_src c :tangle (eval h_private_type)
|
|
|
|
#ifndef QMCKL_MO_HPT
|
|
|
|
#define QMCKL_MO_HPT
|
|
|
|
|
|
|
|
#include <stdbool.h>
|
|
|
|
#+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 <stdio.h>
|
|
|
|
#include <math.h>
|
|
|
|
#include "chbrclf.h"
|
|
|
|
#include "qmckl_ao_private_func.h"
|
|
|
|
#include "qmckl_mo_private_func.h"
|
|
|
|
|
|
|
|
int main() {
|
|
|
|
qmckl_context context;
|
|
|
|
context = qmckl_context_create();
|
|
|
|
|
|
|
|
qmckl_exit_code rc;
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :tangle (eval c)
|
|
|
|
#ifdef HAVE_CONFIG_H
|
|
|
|
#include "config.h"
|
|
|
|
#endif
|
|
|
|
|
|
|
|
#ifdef HAVE_STDINT_H
|
|
|
|
#include <stdint.h>
|
|
|
|
#elif HAVE_INTTYPES_H
|
|
|
|
#include <inttypes.h>
|
|
|
|
#endif
|
|
|
|
|
|
|
|
#include <stdlib.h>
|
|
|
|
#include <string.h>
|
|
|
|
#include <stdbool.h>
|
|
|
|
#include <assert.h>
|
|
|
|
|
|
|
|
#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"
|
|
|
|
#include "qmckl_mo_private_type.h"
|
|
|
|
#include "qmckl_mo_private_func.h"
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
* Context
|
|
|
|
|
|
|
|
The following arrays are stored in the context:
|
|
|
|
|
2021-10-14 17:20:23 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
|-----------------+--------------------+----------------------------------------|
|
|
|
|
| ~mo_num~ | | Number of MOs |
|
|
|
|
| ~coefficient~ | ~[mo_num][ao_num]~ | Orbital coefficients |
|
|
|
|
| ~coefficient_t~ | ~[ao_num][mo_num]~ | Transposed of the Orbital coefficients |
|
|
|
|
|-----------------+--------------------+----------------------------------------|
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
Computed data:
|
2022-01-17 16:09:41 +01:00
|
|
|
|
2022-02-11 16:07:25 +01:00
|
|
|
|---------------+--------------------------+-------------------------------------------------------------------------------------|
|
|
|
|
| ~mo_vgl~ | ~[point_num][5][mo_num]~ | Value, gradients, Laplacian of the MOs at point positions |
|
|
|
|
| ~mo_vgl_date~ | ~uint64_t~ | Late modification date of Value, gradients, Laplacian of the MOs at point positions |
|
|
|
|
|---------------+--------------------------+-------------------------------------------------------------------------------------|
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
** Data structure
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval h_private_type)
|
|
|
|
typedef struct qmckl_mo_basis_struct {
|
2022-01-17 16:09:41 +01:00
|
|
|
int64_t mo_num;
|
2022-03-28 16:53:36 +02:00
|
|
|
double * restrict coefficient;
|
|
|
|
double * restrict coefficient_t;
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
double * restrict mo_vgl;
|
2022-02-16 15:14:41 +01:00
|
|
|
uint64_t mo_vgl_date;
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
int32_t uninitialized;
|
|
|
|
bool provided;
|
|
|
|
} qmckl_mo_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.
|
|
|
|
|
2022-01-17 16:09:41 +01:00
|
|
|
#+begin_src c :comments org :tangle (eval h_private_func)
|
2021-10-25 18:45:06 +02:00
|
|
|
qmckl_exit_code qmckl_init_mo_basis(qmckl_context context);
|
|
|
|
#+end_src
|
2022-01-17 16:09:41 +01:00
|
|
|
|
2021-10-25 18:45:06 +02:00
|
|
|
#+begin_src c :comments org :tangle (eval c)
|
|
|
|
qmckl_exit_code qmckl_init_mo_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->mo_basis.uninitialized = (1 << 2) - 1;
|
|
|
|
|
|
|
|
return QMCKL_SUCCESS;
|
|
|
|
}
|
|
|
|
#+end_src
|
2022-01-17 16:09:41 +01:00
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
** Access functions
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
#+begin_src c :comments org :tangle (eval h_func) :exports none
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_mo_num (const qmckl_context context,
|
|
|
|
int64_t* mo_num);
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
|
2021-10-12 00:10:01 +02:00
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_mo_num (const qmckl_context context,
|
|
|
|
int64_t* mo_num)
|
|
|
|
{
|
|
|
|
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_CONTEXT,
|
|
|
|
"qmckl_get_mo_basis_mo_num",
|
|
|
|
NULL);
|
|
|
|
return (int64_t) 0;
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
|
|
|
|
assert (ctx != NULL);
|
|
|
|
|
|
|
|
int32_t mask = 1;
|
|
|
|
|
|
|
|
if ( (ctx->mo_basis.uninitialized & mask) != 0) {
|
2021-10-12 00:10:01 +02:00
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_NOT_PROVIDED,
|
|
|
|
"qmckl_get_mo_basis_mo_num",
|
|
|
|
NULL);
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
assert (ctx->mo_basis.mo_num > (int64_t) 0);
|
|
|
|
,*mo_num = ctx->mo_basis.mo_num;
|
|
|
|
return QMCKL_SUCCESS;
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval h_func) :exports none
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_coefficient (const qmckl_context context,
|
|
|
|
double* const coefficient,
|
|
|
|
const int64_t size_max);
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :exports none
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_coefficient (const qmckl_context context,
|
|
|
|
double* const coefficient,
|
|
|
|
const int64_t size_max)
|
|
|
|
{
|
2021-09-15 18:30:41 +02:00
|
|
|
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
|
2021-10-12 00:10:01 +02:00
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_CONTEXT,
|
|
|
|
"qmckl_get_mo_basis_coefficient",
|
|
|
|
NULL);
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
|
|
|
|
assert (ctx != NULL);
|
|
|
|
|
|
|
|
int32_t mask = 1 << 1;
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
if ( (ctx->ao_basis.uninitialized & mask) != 0) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_NOT_PROVIDED,
|
|
|
|
"qmckl_get_mo_basis_coefficient",
|
|
|
|
NULL);
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
if (coefficient == NULL) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_ARG_2,
|
|
|
|
"qmckl_get_mo_basis_coefficient",
|
|
|
|
"NULL pointer");
|
|
|
|
}
|
|
|
|
|
|
|
|
if (size_max < ctx->ao_basis.ao_num * ctx->mo_basis.mo_num) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_ARG_3,
|
|
|
|
"qmckl_get_mo_basis_coefficient",
|
|
|
|
"Array too small. Expected mo_num * ao_num");
|
|
|
|
}
|
|
|
|
|
|
|
|
assert (ctx->mo_basis.coefficient != NULL);
|
|
|
|
memcpy(coefficient, ctx->mo_basis.coefficient,
|
2022-01-06 18:54:31 +01:00
|
|
|
ctx->ao_basis.ao_num * ctx->mo_basis.mo_num * sizeof(double));
|
2021-10-12 00:10:01 +02:00
|
|
|
|
|
|
|
return QMCKL_SUCCESS;
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
2021-10-12 00:10:01 +02:00
|
|
|
#+end_src
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
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_mo_basis_provided (const qmckl_context context);
|
|
|
|
#+end_src
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
#+begin_src c :comments org :tangle (eval c) :exports none
|
2021-09-27 18:30:59 +02:00
|
|
|
bool qmckl_mo_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->mo_basis.provided;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-04-20 15:55:59 +02:00
|
|
|
*** Fortran interfaces
|
|
|
|
|
2022-04-20 16:11:06 +02:00
|
|
|
#+begin_src f90 :tangle (eval fh_func) :comments org
|
2022-04-20 15:55:59 +02:00
|
|
|
interface
|
|
|
|
integer(c_int32_t) function qmckl_get_mo_basis_mo_num (context, &
|
|
|
|
mo_num) bind(C)
|
|
|
|
use, intrinsic :: iso_c_binding
|
|
|
|
import
|
|
|
|
implicit none
|
|
|
|
integer (c_int64_t) , intent(in) , value :: context
|
|
|
|
integer (c_int64_t) , intent(out) :: mo_num
|
|
|
|
end function qmckl_get_mo_basis_mo_num
|
|
|
|
end interface
|
|
|
|
|
|
|
|
interface
|
|
|
|
integer(c_int32_t) function qmckl_get_mo_basis_coefficient(context, &
|
|
|
|
coefficient, size_max) bind(C)
|
|
|
|
use, intrinsic :: iso_c_binding
|
|
|
|
import
|
|
|
|
implicit none
|
|
|
|
integer (c_int64_t) , intent(in) , value :: context
|
|
|
|
double precision, intent(out) :: coefficient(*)
|
2022-04-20 16:11:06 +02:00
|
|
|
integer (c_int64_t) , intent(in) , value :: size_max
|
2022-04-20 15:55:59 +02:00
|
|
|
end function qmckl_get_mo_basis_coefficient
|
|
|
|
end interface
|
|
|
|
|
|
|
|
#+end_src
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
** 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_mo_basis_mo_num (qmckl_context context, const int64_t mo_num);
|
|
|
|
qmckl_exit_code qmckl_set_mo_basis_coefficient (qmckl_context context, const double * coefficient);
|
|
|
|
#+end_src
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
#+NAME:pre
|
2021-09-15 18:30:41 +02:00
|
|
|
#+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
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
#+NAME:post
|
2021-09-15 18:30:41 +02:00
|
|
|
#+begin_src c :exports none
|
|
|
|
ctx->mo_basis.uninitialized &= ~mask;
|
|
|
|
ctx->mo_basis.provided = (ctx->mo_basis.uninitialized == 0);
|
|
|
|
if (ctx->mo_basis.provided) {
|
2021-10-25 18:45:06 +02:00
|
|
|
qmckl_exit_code rc_ = qmckl_finalize_mo_basis(context);
|
2021-09-15 18:30:41 +02:00
|
|
|
if (rc_ != QMCKL_SUCCESS) return rc_;
|
2021-10-25 18:45:06 +02:00
|
|
|
}
|
2021-09-15 18:30:41 +02:00
|
|
|
return QMCKL_SUCCESS;
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
|
|
|
|
qmckl_exit_code qmckl_set_mo_basis_mo_num(qmckl_context context, const int64_t mo_num) {
|
2021-10-12 00:10:01 +02:00
|
|
|
<<pre>>
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
if (mo_num <= 0) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_ARG_2,
|
|
|
|
"qmckl_set_mo_basis_mo_num",
|
|
|
|
"mo_num <= 0");
|
|
|
|
}
|
|
|
|
|
2021-10-14 16:39:05 +02:00
|
|
|
int32_t mask = 1 ;
|
2021-09-15 18:30:41 +02:00
|
|
|
ctx->mo_basis.mo_num = mo_num;
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
<<post>>
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
qmckl_exit_code qmckl_set_mo_basis_coefficient(qmckl_context context, const double* coefficient) {
|
2021-10-12 00:10:01 +02:00
|
|
|
<<pre>>
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2021-10-14 16:39:05 +02:00
|
|
|
int32_t mask = 1 << 1;
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
if (ctx->mo_basis.coefficient != NULL) {
|
|
|
|
qmckl_exit_code rc = qmckl_free(context, ctx->mo_basis.coefficient);
|
|
|
|
if (rc != QMCKL_SUCCESS) {
|
|
|
|
return qmckl_failwith( context, rc,
|
|
|
|
"qmckl_set_mo_basis_coefficient",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
|
|
|
|
mem_info.size = ctx->ao_basis.ao_num * ctx->mo_basis.mo_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_mo_basis_coefficient",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
|
|
|
memcpy(new_array, coefficient, mem_info.size);
|
|
|
|
|
|
|
|
ctx->mo_basis.coefficient = new_array;
|
|
|
|
|
2021-10-12 00:10:01 +02:00
|
|
|
<<post>>
|
2021-09-15 18:30:41 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
When the basis set is completely entered, other data structures are
|
|
|
|
computed to accelerate the calculations.
|
|
|
|
|
2021-10-25 18:45:06 +02:00
|
|
|
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
|
|
|
|
qmckl_exit_code qmckl_finalize_mo_basis(qmckl_context context);
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
|
|
|
|
qmckl_exit_code qmckl_finalize_mo_basis(qmckl_context context) {
|
2022-01-17 15:54:21 +01:00
|
|
|
|
2021-10-25 18:45:06 +02:00
|
|
|
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_CONTEXT,
|
|
|
|
"qmckl_finalize_mo_basis",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
|
|
|
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
|
|
|
|
assert (ctx != NULL);
|
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
|
|
|
|
mem_info.size = ctx->ao_basis.ao_num * ctx->mo_basis.mo_num * sizeof(double);
|
|
|
|
double* new_array = (double*) qmckl_malloc(context, mem_info);
|
|
|
|
if (new_array == NULL) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_ALLOCATION_FAILED,
|
|
|
|
"qmckl_finalize_mo_basis",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
|
|
|
assert (ctx->mo_basis.coefficient != NULL);
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
if (ctx->mo_basis.coefficient_t != NULL) {
|
|
|
|
qmckl_exit_code rc = qmckl_free(context, ctx->mo_basis.coefficient);
|
|
|
|
if (rc != QMCKL_SUCCESS) {
|
|
|
|
return qmckl_failwith( context, rc,
|
|
|
|
"qmckl_finalize_mo_basis",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
for (int64_t i=0 ; i<ctx->ao_basis.ao_num ; ++i) {
|
|
|
|
for (int64_t j=0 ; j<ctx->mo_basis.mo_num ; ++j) {
|
|
|
|
new_array[i*ctx->mo_basis.mo_num + j] = ctx->mo_basis.coefficient[j*ctx->ao_basis.ao_num + i];
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
ctx->mo_basis.coefficient_t = new_array;
|
2022-01-17 15:54:21 +01:00
|
|
|
qmckl_exit_code rc = QMCKL_SUCCESS;
|
2021-10-25 18:45:06 +02:00
|
|
|
return rc;
|
|
|
|
}
|
|
|
|
#+end_src
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
* Computation
|
|
|
|
|
|
|
|
** Computation of MOs
|
|
|
|
|
|
|
|
*** Get
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval h_func) :noweb yes
|
2022-03-28 16:53:36 +02:00
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_mo_vgl(qmckl_context context,
|
|
|
|
double* const mo_vgl,
|
|
|
|
const int64_t size_max);
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
|
2022-03-28 16:53:36 +02:00
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_mo_vgl(qmckl_context context,
|
|
|
|
double* const mo_vgl,
|
|
|
|
const int64_t size_max)
|
|
|
|
{
|
2022-01-17 16:09:41 +01:00
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
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;
|
|
|
|
|
2021-09-27 18:17:21 +02:00
|
|
|
rc = qmckl_provide_mo_vgl(context);
|
|
|
|
if (rc != QMCKL_SUCCESS) return rc;
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
|
|
|
|
assert (ctx != NULL);
|
|
|
|
|
2022-03-28 18:26:20 +02:00
|
|
|
const int64_t sze = ctx->point.num * 5 * ctx->mo_basis.mo_num;
|
2022-03-28 16:53:36 +02:00
|
|
|
if (size_max < sze) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_ARG_3,
|
|
|
|
"qmckl_get_mo_basis_mo_vgl",
|
|
|
|
"input array too small");
|
|
|
|
}
|
2021-09-15 18:30:41 +02:00
|
|
|
memcpy(mo_vgl, ctx->mo_basis.mo_vgl, sze * sizeof(double));
|
|
|
|
|
|
|
|
return QMCKL_SUCCESS;
|
|
|
|
}
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
|
2022-03-28 16:53:36 +02:00
|
|
|
interface
|
|
|
|
integer(c_int32_t) function qmckl_get_mo_basis_mo_vgl (context, &
|
|
|
|
mo_vgl, size_max) bind(C)
|
|
|
|
use, intrinsic :: iso_c_binding
|
|
|
|
import
|
|
|
|
implicit none
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
integer (c_int64_t) , intent(in) , value :: context
|
|
|
|
double precision, intent(out) :: mo_vgl(*)
|
|
|
|
integer (c_int64_t) , intent(in) , value :: size_max
|
|
|
|
end function qmckl_get_mo_basis_mo_vgl
|
|
|
|
end interface
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
Uses the given array to compute the VGL.
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval h_func) :noweb yes
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_mo_vgl_inplace (qmckl_context context,
|
|
|
|
double* const mo_vgl,
|
|
|
|
const int64_t size_max);
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_get_mo_basis_mo_vgl_inplace (qmckl_context context,
|
|
|
|
double* const mo_vgl,
|
|
|
|
const int64_t size_max)
|
|
|
|
{
|
|
|
|
|
|
|
|
if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_CONTEXT,
|
|
|
|
"qmckl_get_mo_basis_mo_vgl",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
|
|
|
qmckl_exit_code rc;
|
|
|
|
|
|
|
|
qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
|
|
|
|
assert (ctx != NULL);
|
|
|
|
|
2022-03-28 18:26:20 +02:00
|
|
|
const int64_t sze = ctx->mo_basis.mo_num * 5 * ctx->point.num;
|
2022-03-28 16:53:36 +02:00
|
|
|
if (size_max < sze) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_INVALID_ARG_3,
|
|
|
|
"qmckl_get_mo_basis_mo_vgl",
|
|
|
|
"input array too small");
|
|
|
|
}
|
|
|
|
|
|
|
|
rc = qmckl_context_touch(context);
|
|
|
|
if (rc != QMCKL_SUCCESS) return rc;
|
|
|
|
|
|
|
|
double* old_array = ctx->mo_basis.mo_vgl;
|
|
|
|
|
|
|
|
ctx->mo_basis.mo_vgl = mo_vgl;
|
|
|
|
|
|
|
|
rc = qmckl_provide_mo_vgl(context);
|
|
|
|
if (rc != QMCKL_SUCCESS) return rc;
|
|
|
|
|
|
|
|
ctx->mo_basis.mo_vgl = old_array;
|
|
|
|
|
|
|
|
return QMCKL_SUCCESS;
|
|
|
|
}
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
|
|
|
|
interface
|
|
|
|
integer(c_int32_t) function qmckl_get_mo_basis_mo_vgl_inplace (context, &
|
|
|
|
mo_vgl, size_max) bind(C)
|
|
|
|
use, intrinsic :: iso_c_binding
|
|
|
|
import
|
|
|
|
implicit none
|
|
|
|
integer (c_int64_t) , intent(in) , value :: context
|
|
|
|
double precision, intent(out) :: mo_vgl(*)
|
|
|
|
integer (c_int64_t) , intent(in) , value :: size_max
|
|
|
|
end function qmckl_get_mo_basis_mo_vgl_inplace
|
|
|
|
end interface
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
*** Provide
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval h_private_func) :noweb yes :exports none
|
2021-09-27 18:17:21 +02:00
|
|
|
qmckl_exit_code qmckl_provide_mo_vgl(qmckl_context context);
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+begin_src c :comments org :tangle (eval c) :noweb yes :exports none
|
2021-09-27 18:17:21 +02:00
|
|
|
qmckl_exit_code qmckl_provide_mo_vgl(qmckl_context context)
|
2021-09-15 18:30:41 +02:00
|
|
|
{
|
|
|
|
|
2022-01-17 16:09:41 +01:00
|
|
|
qmckl_exit_code rc;
|
2021-09-15 18:30:41 +02:00
|
|
|
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",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
2021-10-25 18:45:06 +02:00
|
|
|
rc = qmckl_provide_ao_vgl(context);
|
|
|
|
if (rc != QMCKL_SUCCESS) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_NOT_PROVIDED,
|
|
|
|
"qmckl_ao_basis",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
2021-09-27 18:39:45 +02:00
|
|
|
if (!ctx->mo_basis.provided) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_NOT_PROVIDED,
|
|
|
|
"qmckl_mo_basis",
|
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
/* Compute if necessary */
|
2022-02-11 16:07:25 +01:00
|
|
|
if (ctx->point.date > ctx->mo_basis.mo_vgl_date) {
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
/* Allocate array */
|
|
|
|
if (ctx->mo_basis.mo_vgl == NULL) {
|
|
|
|
|
|
|
|
qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero;
|
2022-02-11 16:07:25 +01:00
|
|
|
mem_info.size = 5 * ctx->point.num * ctx->mo_basis.mo_num * sizeof(double);
|
2021-09-15 18:30:41 +02:00
|
|
|
double* mo_vgl = (double*) qmckl_malloc(context, mem_info);
|
|
|
|
|
|
|
|
if (mo_vgl == NULL) {
|
|
|
|
return qmckl_failwith( context,
|
|
|
|
QMCKL_ALLOCATION_FAILED,
|
2022-03-28 16:53:36 +02:00
|
|
|
"qmckl_mo_basis_mo_vgl",
|
2021-09-15 18:30:41 +02:00
|
|
|
NULL);
|
|
|
|
}
|
|
|
|
ctx->mo_basis.mo_vgl = mo_vgl;
|
|
|
|
}
|
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
rc = qmckl_compute_mo_basis_mo_vgl(context,
|
|
|
|
ctx->ao_basis.ao_num,
|
|
|
|
ctx->mo_basis.mo_num,
|
|
|
|
ctx->point.num,
|
|
|
|
ctx->mo_basis.coefficient_t,
|
|
|
|
ctx->ao_basis.ao_vgl,
|
|
|
|
ctx->mo_basis.mo_vgl);
|
2021-09-27 22:06:49 +02:00
|
|
|
if (rc != QMCKL_SUCCESS) {
|
|
|
|
return rc;
|
|
|
|
}
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
ctx->mo_basis.mo_vgl_date = ctx->date;
|
|
|
|
}
|
|
|
|
|
|
|
|
return QMCKL_SUCCESS;
|
|
|
|
}
|
|
|
|
#+end_src
|
2022-01-17 16:09:41 +01:00
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
*** Compute
|
|
|
|
:PROPERTIES:
|
2022-03-28 16:53:36 +02:00
|
|
|
:Name: qmckl_compute_mo_basis_mo_vgl
|
2021-09-15 18:30:41 +02:00
|
|
|
:CRetType: qmckl_exit_code
|
|
|
|
:FRetType: qmckl_exit_code
|
|
|
|
:END:
|
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
#+NAME: qmckl_mo_basis_mo_vgl_args
|
|
|
|
| Variable | Type | In/Out | Description |
|
|
|
|
|---------------------+--------------------------------+--------+-------------------------------------------------|
|
|
|
|
| ~context~ | ~qmckl_context~ | in | Global state |
|
|
|
|
| ~ao_num~ | ~int64_t~ | in | Number of AOs |
|
|
|
|
| ~mo_num~ | ~int64_t~ | in | Number of MOs |
|
|
|
|
| ~point_num~ | ~int64_t~ | in | Number of points |
|
|
|
|
| ~coef_normalized_t~ | ~double[mo_num][ao_num]~ | in | Transpose of the AO to MO transformation matrix |
|
|
|
|
| ~ao_vgl~ | ~double[point_num][5][ao_num]~ | in | Value, gradients and Laplacian of the AOs |
|
|
|
|
| ~mo_vgl~ | ~double[point_num][5][mo_num]~ | out | Value, gradients and Laplacian of the MOs |
|
|
|
|
|
|
|
|
|
|
|
|
The matrix of AO values is very sparse, so we use a sparse-dense
|
|
|
|
matrix multiplication instead of a dgemm, as exposed in
|
|
|
|
https://dx.doi.org/10.1007/978-3-642-38718-0_14.
|
2021-10-11 23:02:09 +02:00
|
|
|
|
2022-04-20 16:11:06 +02:00
|
|
|
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
|
2022-03-28 16:53:36 +02:00
|
|
|
integer function qmckl_compute_mo_basis_mo_vgl_doc_f(context, &
|
2022-02-11 16:07:25 +01:00
|
|
|
ao_num, mo_num, point_num, &
|
2022-03-28 16:53:36 +02:00
|
|
|
coef_normalized_t, ao_vgl, mo_vgl) &
|
2021-09-15 18:30:41 +02:00
|
|
|
result(info)
|
|
|
|
use qmckl
|
|
|
|
implicit none
|
|
|
|
integer(qmckl_context), intent(in) :: context
|
|
|
|
integer*8 , intent(in) :: ao_num, mo_num
|
2022-02-11 16:07:25 +01:00
|
|
|
integer*8 , intent(in) :: point_num
|
|
|
|
double precision , intent(in) :: ao_vgl(ao_num,5,point_num)
|
2022-03-28 18:29:29 +02:00
|
|
|
double precision , intent(in) :: coef_normalized_t(mo_num,ao_num)
|
2022-02-11 16:07:25 +01:00
|
|
|
double precision , intent(out) :: mo_vgl(mo_num,5,point_num)
|
2022-03-28 18:29:29 +02:00
|
|
|
integer*8 :: i,j,k
|
2022-03-28 16:53:36 +02:00
|
|
|
double precision :: c1, c2, c3, c4, c5
|
|
|
|
|
|
|
|
do j=1,point_num
|
|
|
|
mo_vgl(:,:,j) = 0.d0
|
|
|
|
do k=1,ao_num
|
|
|
|
if (ao_vgl(k,1,j) /= 0.d0) then
|
|
|
|
c1 = ao_vgl(k,1,j)
|
|
|
|
c2 = ao_vgl(k,2,j)
|
|
|
|
c3 = ao_vgl(k,3,j)
|
|
|
|
c4 = ao_vgl(k,4,j)
|
|
|
|
c5 = ao_vgl(k,5,j)
|
|
|
|
do i=1,mo_num
|
|
|
|
mo_vgl(i,1,j) = mo_vgl(i,1,j) + coef_normalized_t(i,k) * c1
|
|
|
|
mo_vgl(i,2,j) = mo_vgl(i,2,j) + coef_normalized_t(i,k) * c2
|
|
|
|
mo_vgl(i,3,j) = mo_vgl(i,3,j) + coef_normalized_t(i,k) * c3
|
|
|
|
mo_vgl(i,4,j) = mo_vgl(i,4,j) + coef_normalized_t(i,k) * c4
|
|
|
|
mo_vgl(i,5,j) = mo_vgl(i,5,j) + coef_normalized_t(i,k) * c5
|
|
|
|
end do
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
end function qmckl_compute_mo_basis_mo_vgl_doc_f
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
#+CALL: generate_c_header(table=qmckl_mo_basis_mo_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_mo_basis_mo_vgl"))
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
#+RESULTS:
|
2021-09-27 22:06:49 +02:00
|
|
|
#+begin_src c :tangle (eval h_func) :comments org
|
2022-03-28 16:53:36 +02:00
|
|
|
qmckl_exit_code qmckl_compute_mo_basis_mo_vgl (
|
|
|
|
const qmckl_context context,
|
|
|
|
const int64_t ao_num,
|
|
|
|
const int64_t mo_num,
|
|
|
|
const int64_t point_num,
|
|
|
|
const double* coef_normalized_t,
|
|
|
|
const double* ao_vgl,
|
2022-04-20 16:11:06 +02:00
|
|
|
double* const mo_vgl );
|
2021-09-27 22:06:49 +02:00
|
|
|
#+end_src
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
#+CALL: generate_c_header(table=qmckl_mo_basis_mo_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_mo_basis_mo_vgl_doc"))
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
#+RESULTS:
|
|
|
|
#+begin_src c :tangle (eval h_func) :comments org
|
|
|
|
qmckl_exit_code qmckl_compute_mo_basis_mo_vgl_doc (
|
|
|
|
const qmckl_context context,
|
|
|
|
const int64_t ao_num,
|
|
|
|
const int64_t mo_num,
|
|
|
|
const int64_t point_num,
|
|
|
|
const double* coef_normalized_t,
|
|
|
|
const double* ao_vgl,
|
2022-04-20 16:11:06 +02:00
|
|
|
double* const mo_vgl );
|
2022-03-28 16:53:36 +02:00
|
|
|
#+end_src
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
#+CALL: generate_c_interface(table=qmckl_mo_basis_mo_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_compute_mo_basis_mo_vgl_doc"))
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
#+RESULTS:
|
2021-09-27 22:06:49 +02:00
|
|
|
#+begin_src f90 :tangle (eval f) :comments org :exports none
|
2022-03-28 16:53:36 +02:00
|
|
|
integer(c_int32_t) function qmckl_compute_mo_basis_mo_vgl_doc &
|
|
|
|
(context, ao_num, mo_num, point_num, coef_normalized_t, ao_vgl, mo_vgl) &
|
|
|
|
bind(C) result(info)
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
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 :: mo_num
|
2022-02-11 16:07:25 +01:00
|
|
|
integer (c_int64_t) , intent(in) , value :: point_num
|
2022-03-28 16:53:36 +02:00
|
|
|
real (c_double ) , intent(in) :: coef_normalized_t(ao_num,mo_num)
|
2022-02-11 16:07:25 +01:00
|
|
|
real (c_double ) , intent(in) :: ao_vgl(ao_num,5,point_num)
|
|
|
|
real (c_double ) , intent(out) :: mo_vgl(mo_num,5,point_num)
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
integer(c_int32_t), external :: qmckl_compute_mo_basis_mo_vgl_doc_f
|
|
|
|
info = qmckl_compute_mo_basis_mo_vgl_doc_f &
|
|
|
|
(context, ao_num, mo_num, point_num, coef_normalized_t, ao_vgl, mo_vgl)
|
|
|
|
|
|
|
|
end function qmckl_compute_mo_basis_mo_vgl_doc
|
|
|
|
#+end_src
|
|
|
|
|
2022-04-20 16:11:06 +02:00
|
|
|
#+begin_src c :tangle (eval c) :comments org
|
2022-03-28 16:53:36 +02:00
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_compute_mo_basis_mo_vgl (const qmckl_context context,
|
|
|
|
const int64_t ao_num,
|
|
|
|
const int64_t mo_num,
|
|
|
|
const int64_t point_num,
|
|
|
|
const double* coef_normalized_t,
|
|
|
|
const double* ao_vgl,
|
|
|
|
double* const mo_vgl )
|
|
|
|
{
|
|
|
|
#ifdef HAVE_HPC
|
|
|
|
return qmckl_compute_mo_basis_mo_vgl_hpc (context, ao_num, mo_num, point_num, coef_normalized_t, ao_vgl, mo_vgl);
|
|
|
|
#else
|
|
|
|
return qmckl_compute_mo_basis_mo_vgl_doc (context, ao_num, mo_num, point_num, coef_normalized_t, ao_vgl, mo_vgl);
|
|
|
|
#endif
|
|
|
|
}
|
|
|
|
#+end_src
|
2022-04-20 16:11:06 +02:00
|
|
|
|
|
|
|
|
2022-03-28 16:53:36 +02:00
|
|
|
*** HPC version
|
|
|
|
|
2022-04-20 16:11:06 +02:00
|
|
|
|
|
|
|
#+begin_src c :tangle (eval h_func) :comments org
|
2022-03-28 16:53:36 +02:00
|
|
|
#ifdef HAVE_HPC
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_compute_mo_basis_mo_vgl_hpc (const qmckl_context context,
|
|
|
|
const int64_t ao_num,
|
|
|
|
const int64_t mo_num,
|
|
|
|
const int64_t point_num,
|
|
|
|
const double* coef_normalized_t,
|
|
|
|
const double* ao_vgl,
|
|
|
|
double* const mo_vgl );
|
|
|
|
#endif
|
|
|
|
#+end_src
|
2021-09-15 18:30:41 +02:00
|
|
|
|
2022-04-20 16:11:06 +02:00
|
|
|
#+begin_src c :tangle (eval c) :comments org
|
2022-03-28 16:53:36 +02:00
|
|
|
#ifdef HAVE_HPC
|
|
|
|
qmckl_exit_code
|
|
|
|
qmckl_compute_mo_basis_mo_vgl_hpc (const qmckl_context context,
|
|
|
|
const int64_t ao_num,
|
|
|
|
const int64_t mo_num,
|
|
|
|
const int64_t point_num,
|
|
|
|
const double* restrict coef_normalized_t,
|
|
|
|
const double* restrict ao_vgl,
|
|
|
|
double* restrict const mo_vgl )
|
|
|
|
{
|
2022-03-28 17:37:50 +02:00
|
|
|
#ifdef HAVE_OPENMP
|
|
|
|
#pragma omp parallel for
|
|
|
|
#endif
|
|
|
|
for (int64_t ipoint=0 ; ipoint < point_num ; ++ipoint) {
|
|
|
|
double* restrict const vgl1 = &(mo_vgl[ipoint*5*mo_num]);
|
|
|
|
const double* restrict avgl1 = &(ao_vgl[ipoint*5*ao_num]);
|
|
|
|
double* restrict const vgl2 = vgl1 + mo_num;
|
|
|
|
double* restrict const vgl3 = vgl1 + (mo_num << 1);
|
|
|
|
double* restrict const vgl4 = vgl1 + (mo_num << 1) + mo_num;
|
|
|
|
double* restrict const vgl5 = vgl1 + (mo_num << 2);
|
|
|
|
const double* restrict avgl2 = avgl1 + ao_num;
|
|
|
|
const double* restrict avgl3 = avgl1 + (ao_num << 1);
|
|
|
|
const double* restrict avgl4 = avgl1 + (ao_num << 1) + ao_num;
|
|
|
|
const double* restrict avgl5 = avgl1 + (ao_num << 2);
|
|
|
|
|
|
|
|
for (int64_t i=0 ; i<mo_num ; ++i) {
|
|
|
|
vgl1[i] = 0.;
|
|
|
|
vgl2[i] = 0.;
|
|
|
|
vgl3[i] = 0.;
|
|
|
|
vgl4[i] = 0.;
|
|
|
|
vgl5[i] = 0.;
|
|
|
|
}
|
|
|
|
|
|
|
|
int64_t nidx=0;
|
|
|
|
int64_t idx[ao_num];
|
|
|
|
double av1[ao_num];
|
|
|
|
double av2[ao_num];
|
|
|
|
double av3[ao_num];
|
|
|
|
double av4[ao_num];
|
|
|
|
double av5[ao_num];
|
|
|
|
for (int64_t k=0 ; k<ao_num ; ++k) {
|
|
|
|
const double* restrict ck1 = coef_normalized_t + k*mo_num;
|
|
|
|
if (avgl1[k] != 0.) {
|
|
|
|
idx[nidx] = k;
|
|
|
|
av1[nidx] = avgl1[k];
|
|
|
|
av2[nidx] = avgl2[k];
|
|
|
|
av3[nidx] = avgl3[k];
|
|
|
|
av4[nidx] = avgl4[k];
|
|
|
|
av5[nidx] = avgl5[k];
|
|
|
|
++nidx;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-03-28 17:58:03 +02:00
|
|
|
int64_t n;
|
|
|
|
for (n=0 ; n < nidx-4 ; n+=4) {
|
2022-03-28 17:37:50 +02:00
|
|
|
int64_t k = idx[n];
|
2022-03-28 17:58:03 +02:00
|
|
|
const double* restrict ck1 = coef_normalized_t + idx[n ]*mo_num;
|
|
|
|
const double* restrict ck2 = coef_normalized_t + idx[n+1]*mo_num;
|
|
|
|
const double* restrict ck3 = coef_normalized_t + idx[n+2]*mo_num;
|
|
|
|
const double* restrict ck4 = coef_normalized_t + idx[n+3]*mo_num;
|
|
|
|
|
|
|
|
const double a11 = av1[n ];
|
|
|
|
const double a21 = av1[n+1];
|
|
|
|
const double a31 = av1[n+2];
|
|
|
|
const double a41 = av1[n+3];
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-03-28 17:58:03 +02:00
|
|
|
const double a12 = av2[n ];
|
|
|
|
const double a22 = av2[n+1];
|
|
|
|
const double a32 = av2[n+2];
|
|
|
|
const double a42 = av2[n+3];
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-03-28 17:58:03 +02:00
|
|
|
const double a13 = av3[n ];
|
|
|
|
const double a23 = av3[n+1];
|
|
|
|
const double a33 = av3[n+2];
|
|
|
|
const double a43 = av3[n+3];
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-03-28 17:58:03 +02:00
|
|
|
const double a14 = av4[n ];
|
|
|
|
const double a24 = av4[n+1];
|
|
|
|
const double a34 = av4[n+2];
|
|
|
|
const double a44 = av4[n+3];
|
2022-04-20 16:11:06 +02:00
|
|
|
|
2022-03-28 17:58:03 +02:00
|
|
|
const double a15 = av5[n ];
|
|
|
|
const double a25 = av5[n+1];
|
|
|
|
const double a35 = av5[n+2];
|
|
|
|
const double a45 = av5[n+3];
|
|
|
|
|
|
|
|
#ifdef HAVE_OPENMP
|
|
|
|
#pragma omp simd
|
|
|
|
#endif
|
|
|
|
for (int64_t i=0 ; i<mo_num ; ++i) {
|
|
|
|
vgl1[i] = vgl1[i] + ck1[i] * a11 + ck2[i] * a21 + ck3[i] * a31 + ck4[i] * a41;
|
|
|
|
vgl2[i] = vgl2[i] + ck1[i] * a12 + ck2[i] * a22 + ck3[i] * a32 + ck4[i] * a42;
|
|
|
|
vgl3[i] = vgl3[i] + ck1[i] * a13 + ck2[i] * a23 + ck3[i] * a33 + ck4[i] * a43;
|
|
|
|
vgl4[i] = vgl4[i] + ck1[i] * a14 + ck2[i] * a24 + ck3[i] * a34 + ck4[i] * a44;
|
|
|
|
vgl5[i] = vgl5[i] + ck1[i] * a15 + ck2[i] * a25 + ck3[i] * a35 + ck4[i] * a45;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
int64_t n0 = nidx-4;
|
|
|
|
n0 = n0 < 0 ? 0 : n0;
|
|
|
|
for (int64_t n=n0 ; n < nidx ; n+=1) {
|
|
|
|
const double* restrict ck = coef_normalized_t + idx[n]*mo_num;
|
|
|
|
const double a1 = av1[n];
|
|
|
|
const double a2 = av2[n];
|
|
|
|
const double a3 = av3[n];
|
|
|
|
const double a4 = av4[n];
|
|
|
|
const double a5 = av5[n];
|
|
|
|
|
2022-03-28 17:37:50 +02:00
|
|
|
#ifdef HAVE_OPENMP
|
|
|
|
#pragma omp simd
|
|
|
|
#endif
|
2022-03-28 17:58:03 +02:00
|
|
|
for (int64_t i=0 ; i<mo_num ; ++i) {
|
|
|
|
vgl1[i] += ck[i] * a1;
|
|
|
|
vgl2[i] += ck[i] * a2;
|
|
|
|
vgl3[i] += ck[i] * a3;
|
|
|
|
vgl4[i] += ck[i] * a4;
|
|
|
|
vgl5[i] += ck[i] * a5;
|
2022-03-28 17:37:50 +02:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
return QMCKL_SUCCESS;
|
2022-03-28 16:53:36 +02:00
|
|
|
}
|
|
|
|
#endif
|
2021-09-27 22:06:49 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
*** Test
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
#+begin_src python :results output :exports none
|
|
|
|
import numpy as np
|
|
|
|
|
|
|
|
def f(a,x,y):
|
2022-01-17 16:09:41 +01:00
|
|
|
return np.sum( [c * np.exp( -b*(np.linalg.norm(x-y))**2) for b,c in a] )
|
2021-09-15 18:30:41 +02:00
|
|
|
|
|
|
|
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 * 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] : %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 *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 ( "[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
|
|
|
|
|
|
|
|
#+begin_src c :tangle (eval c_test) :exports none
|
2021-09-27 17:18:43 +02:00
|
|
|
{
|
|
|
|
#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]);
|
|
|
|
const int64_t nucl_num = chbrclf_nucl_num;
|
|
|
|
const double* nucl_charge = chbrclf_charge;
|
|
|
|
const double* nucl_coord = &(chbrclf_nucl_coord[0][0]);
|
|
|
|
|
2022-01-17 16:09:41 +01:00
|
|
|
rc = qmckl_set_electron_num (context, elec_up_num, elec_dn_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert (rc == QMCKL_SUCCESS);
|
|
|
|
|
|
|
|
rc = qmckl_set_electron_walk_num (context, walk_num);
|
|
|
|
assert (rc == QMCKL_SUCCESS);
|
|
|
|
|
|
|
|
assert(qmckl_electron_provided(context));
|
|
|
|
|
2022-01-17 16:09:41 +01:00
|
|
|
rc = qmckl_set_electron_coord (context, 'N', elec_coord, walk_num*elec_num*3);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
|
|
|
rc = qmckl_set_nucleus_num (context, nucl_num);
|
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
2022-01-23 16:18:46 +01:00
|
|
|
rc = qmckl_set_nucleus_coord (context, 'T', &(nucl_coord[0]), nucl_num*3);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
2022-01-23 16:18:46 +01:00
|
|
|
rc = qmckl_set_nucleus_charge(context, nucl_charge, nucl_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
|
|
|
assert(qmckl_nucleus_provided(context));
|
|
|
|
|
|
|
|
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]);
|
|
|
|
const double * ao_factor = &(chbrclf_basis_ao_factor[0]);
|
|
|
|
|
|
|
|
const 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, chbrclf_shell_num);
|
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
|
|
|
rc = qmckl_set_ao_basis_prim_num (context, chbrclf_prim_num);
|
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_nucleus_index (context, nucleus_index, nucl_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_nucleus_shell_num (context, nucleus_shell_num, nucl_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_shell_ang_mom (context, shell_ang_mom, chbrclf_shell_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_shell_factor (context, shell_factor, chbrclf_shell_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_shell_prim_num (context, shell_prim_num, chbrclf_shell_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_shell_prim_index (context, shell_prim_index, chbrclf_shell_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_exponent (context, exponent, chbrclf_prim_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_coefficient (context, coefficient, chbrclf_prim_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
assert(!qmckl_ao_basis_provided(context));
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_prim_factor (context, prim_factor, chbrclf_prim_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
|
|
|
rc = qmckl_set_ao_basis_ao_num(context, chbrclf_ao_num);
|
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
2022-01-06 02:28:13 +01:00
|
|
|
rc = qmckl_set_ao_basis_ao_factor (context, ao_factor, chbrclf_ao_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert(rc == QMCKL_SUCCESS);
|
|
|
|
|
|
|
|
assert(qmckl_ao_basis_provided(context));
|
|
|
|
|
|
|
|
|
2022-02-11 16:07:25 +01:00
|
|
|
double ao_vgl[walk_num*elec_num][5][chbrclf_ao_num];
|
2021-09-27 17:18:43 +02:00
|
|
|
|
2022-02-11 16:07:25 +01:00
|
|
|
rc = qmckl_get_ao_basis_ao_vgl(context, &(ao_vgl[0][0][0]),
|
2022-01-06 02:28:13 +01:00
|
|
|
(int64_t) 5*walk_num*elec_num*chbrclf_ao_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert (rc == QMCKL_SUCCESS);
|
|
|
|
|
2021-09-27 18:17:21 +02:00
|
|
|
/* Set up MO data */
|
2021-09-29 18:12:41 +02:00
|
|
|
const int64_t mo_num = chbrclf_mo_num;
|
2021-09-27 18:17:21 +02:00
|
|
|
rc = qmckl_set_mo_basis_mo_num(context, mo_num);
|
2021-09-27 17:18:43 +02:00
|
|
|
assert (rc == QMCKL_SUCCESS);
|
|
|
|
|
2021-09-29 18:12:41 +02:00
|
|
|
const double * mo_coefficient = &(chbrclf_mo_coef[0]);
|
2021-09-27 18:17:21 +02:00
|
|
|
|
2021-09-29 18:12:41 +02:00
|
|
|
rc = qmckl_set_mo_basis_coefficient(context, mo_coefficient);
|
2021-09-27 18:17:21 +02:00
|
|
|
assert (rc == QMCKL_SUCCESS);
|
|
|
|
|
2021-09-27 18:30:59 +02:00
|
|
|
assert(qmckl_mo_basis_provided(context));
|
2021-09-27 18:17:21 +02:00
|
|
|
|
2022-02-11 16:07:25 +01:00
|
|
|
double mo_vgl[walk_num*elec_num][5][chbrclf_mo_num];
|
2022-03-28 16:53:36 +02:00
|
|
|
rc = qmckl_get_mo_basis_mo_vgl(context, &(mo_vgl[0][0][0]), walk_num*elec_num*5*chbrclf_mo_num);
|
2021-09-27 23:36:11 +02:00
|
|
|
assert (rc == QMCKL_SUCCESS);
|
2021-09-27 18:17:21 +02:00
|
|
|
|
2021-10-04 09:08:38 +02:00
|
|
|
// Test overlap of MO
|
2021-10-07 14:01:40 +02:00
|
|
|
//double point_x[10];
|
|
|
|
//double point_y[10];
|
|
|
|
//double point_z[10];
|
|
|
|
//int32_t npoints=10;
|
2021-10-04 09:08:38 +02:00
|
|
|
//// obtain points
|
|
|
|
//double dr = 20./(npoints-1);
|
|
|
|
//double dr3 = dr*dr*dr;
|
2021-09-29 18:55:02 +02:00
|
|
|
//
|
2021-10-04 09:08:38 +02:00
|
|
|
//for (int i=0;i<npoints;++i) {
|
|
|
|
// point_x[i] = -10. + dr*i;
|
|
|
|
// point_y[i] = -10. + dr*i;
|
|
|
|
// point_z[i] = -10. + dr*i;
|
2021-09-29 18:12:41 +02:00
|
|
|
//}
|
2021-10-04 09:08:38 +02:00
|
|
|
//
|
|
|
|
//double ovlmo1 = 0.0;
|
|
|
|
//// Calculate overlap
|
|
|
|
//for (int i=0;i<npoints;++i) {
|
|
|
|
// fflush(stdout);
|
|
|
|
// for (int j=0;j<npoints;++j) {
|
2021-10-07 14:01:40 +02:00
|
|
|
// printf(" .. ");
|
2021-10-04 09:08:38 +02:00
|
|
|
// for (int k=0;k<npoints;++k) {
|
2021-10-07 14:01:40 +02:00
|
|
|
// printf(" . ");
|
2021-10-04 09:08:38 +02:00
|
|
|
// // Set point
|
|
|
|
// elec_coord[0] = point_x[i];
|
|
|
|
// elec_coord[1] = point_y[j];
|
|
|
|
// elec_coord[2] = point_z[k];
|
2022-01-17 16:09:41 +01:00
|
|
|
// rc = qmckl_set_electron_coord (context, 'N', elec_coord);
|
2021-10-04 09:08:38 +02:00
|
|
|
// assert(rc == QMCKL_SUCCESS);
|
|
|
|
//
|
|
|
|
// // Calculate value of MO (1st electron)
|
|
|
|
// double mo_vgl[5][walk_num][elec_num][chbrclf_mo_num];
|
2022-03-28 16:53:36 +02:00
|
|
|
// rc = qmckl_get_mo_basis_mo_vgl(context, &(mo_vgl[0][0][0][0]));
|
2021-10-04 09:08:38 +02:00
|
|
|
// assert (rc == QMCKL_SUCCESS);
|
2021-10-04 09:50:16 +02:00
|
|
|
// ovlmo1 += mo_vgl[0][0][0][0]*mo_vgl[0][0][0][0]*dr3;
|
2021-10-04 09:08:38 +02:00
|
|
|
// }
|
|
|
|
// }
|
2021-09-29 18:55:02 +02:00
|
|
|
//}
|
2021-10-04 09:08:38 +02:00
|
|
|
//printf("OVL MO1 = %10.15f\n",ovlmo1);
|
2021-09-29 18:55:02 +02:00
|
|
|
|
2021-09-29 18:12:41 +02:00
|
|
|
|
2021-09-27 22:01:13 +02:00
|
|
|
printf("\n");
|
2022-02-11 16:07:25 +01:00
|
|
|
printf(" mo_vgl mo_vgl[0][26][219] %25.15e\n", mo_vgl[2][0][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[1][26][219] %25.15e\n", mo_vgl[2][1][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[0][26][220] %25.15e\n", mo_vgl[2][0][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[1][26][220] %25.15e\n", mo_vgl[2][1][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[0][26][221] %25.15e\n", mo_vgl[2][0][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[1][26][221] %25.15e\n", mo_vgl[2][1][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[0][26][222] %25.15e\n", mo_vgl[2][0][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[1][26][222] %25.15e\n", mo_vgl[2][1][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[0][26][223] %25.15e\n", mo_vgl[2][0][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[1][26][223] %25.15e\n", mo_vgl[2][1][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[0][26][224] %25.15e\n", mo_vgl[2][0][3]);
|
|
|
|
printf(" mo_vgl mo_vgl[1][26][224] %25.15e\n", mo_vgl[2][1][3]);
|
2021-09-27 22:01:13 +02:00
|
|
|
printf("\n");
|
2022-01-17 16:09:41 +01:00
|
|
|
}
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
#+end_src
|
|
|
|
|
|
|
|
* End of files :noexport:
|
|
|
|
|
|
|
|
#+begin_src c :tangle (eval h_private_type)
|
|
|
|
#endif
|
|
|
|
#+end_src
|
|
|
|
|
2021-10-14 21:40:14 +02:00
|
|
|
#+begin_src c :tangle (eval h_private_func)
|
|
|
|
#endif
|
|
|
|
#+end_src
|
|
|
|
|
2021-09-15 18:30:41 +02:00
|
|
|
*** 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
|
|
|
|
|
|
|
|
|
|
|
|
# -*- mode: org -*-
|
|
|
|
# vim: syntax=c
|
|
|
|
|
|
|
|
|