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

Cleaning Fortran

This commit is contained in:
Anthony Scemama 2023-09-22 16:41:43 +02:00
parent 50fa3aa754
commit 0d9af3c497
3 changed files with 281 additions and 677 deletions

View File

@ -3346,15 +3346,18 @@ qmckl_ao_gaussian_vgl(const qmckl_context context,
#+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
function qmckl_ao_gaussian_vgl(context, X, R, n, A, VGL, ldv) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
implicit none
integer*8 , intent(in) :: context
double precision , intent(in) :: X(3), R(3)
integer*8 , intent(in) :: n
double precision , intent(in) :: A(n)
double precision , intent(out) :: VGL(ldv,5)
integer*8 , intent(in) :: ldv
integer (qmckl_context) , 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 (qmckl_exit_code) :: info
integer*8 :: i,j
double precision :: Y(3), r2, t, u, v
@ -3405,36 +3408,21 @@ integer function qmckl_ao_gaussian_vgl_f(context, X, R, n, A, VGL, ldv) result(i
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)
function qmckl_ao_gaussian_vgl(context, &
X, R, n, A, VGL, ldv) bind(C) result(info)
use, intrinsic :: iso_c_binding
integer (c_int64_t) , intent(in) , value :: context
use qmckl_constants
integer (qmckl_context) , 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)
integer(qmckl_exit_code) :: info
end function qmckl_ao_gaussian_vgl
end interface
#+end_src
@ -3442,11 +3430,12 @@ end function qmckl_ao_gaussian_vgl
*** Test :noexport:
#+begin_src f90 :tangle (eval f_test)
integer(c_int32_t) function test_qmckl_ao_gaussian_vgl(context) bind(C)
function test_qmckl_ao_gaussian_vgl(context) bind(C)
use qmckl
implicit none
integer(c_int64_t), intent(in), value :: context
integer(qmckl_context), intent(in), value :: context
integer(qmckl_exit_code) :: test_qmckl_ao_gaussian_vgl
integer*8 :: n, ldv, j, i
double precision :: X(3), R(3), Y(3), r2, z
@ -3558,7 +3547,7 @@ assert(0 == test_qmckl_ao_gaussian_vgl(context));
| ~prim_num~ | ~int64_t~ | in | Number of primitives |
| ~point_num~ | ~int64_t~ | in | Number of points considered |
| ~nucl_num~ | ~int64_t~ | in | Number of nuclei |
| ~nucleus_prim_index~ | ~int64_t[nucl_num]~ | in | Index of the 1st primitive of each nucleus |
| ~nucleus_prim_index~ | ~int64_t[nucl_num+1]~ | in | Index of the 1st primitive of each nucleus |
| ~coord~ | ~double[3][point_num]~ | in | Coordinates |
| ~nucl_coord~ | ~double[3][nucl_num]~ | in | Nuclear coordinates |
| ~expo~ | ~double[prim_num]~ | in | Exponents of the primitives |
@ -3582,23 +3571,24 @@ assert(0 == test_qmckl_ao_gaussian_vgl(context));
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_ao_basis_primitive_gaussian_vgl_f( &
context, prim_num, point_num, nucl_num, &
nucleus_prim_index, coord, nucl_coord, &
expo, primitive_vgl) &
result(info)
function qmckl_compute_ao_basis_primitive_gaussian_vgl &
(context, prim_num, point_num, nucl_num, nucleus_prim_index, coord, nucl_coord, expo, primitive_vgl) &
bind(C) result(info)
use qmckl
use, intrinsic :: iso_c_binding
use qmckl_constants
implicit none
integer(qmckl_context), intent(in) :: context
integer*8 , intent(in) :: prim_num
integer*8 , intent(in) :: nucl_num
integer*8 , intent(in) :: point_num
integer*8 , intent(in) :: nucleus_prim_index(nucl_num+1)
double precision , intent(in) :: coord(point_num,3)
double precision , intent(in) :: nucl_coord(nucl_num,3)
double precision , intent(in) :: expo(prim_num)
double precision , intent(out) :: primitive_vgl(prim_num,5,point_num)
integer (qmckl_context), intent(in) , value :: context
integer (c_int64_t) , intent(in) , value :: prim_num
integer (c_int64_t) , intent(in) , value :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
integer (c_int64_t) , intent(in) :: nucleus_prim_index(nucl_num+1)
real (c_double ) , intent(in) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
real (c_double ) , intent(in) :: expo(prim_num)
real (c_double ) , intent(out) :: primitive_vgl(prim_num,5,point_num)
integer(qmckl_exit_code) :: info
integer*8 :: inucl, iprim, ipoint
double precision :: x, y, z, two_a, ar2, r2, v, cutoff
@ -3633,53 +3623,10 @@ integer function qmckl_compute_ao_basis_primitive_gaussian_vgl_f( &
end do
end do
end function qmckl_compute_ao_basis_primitive_gaussian_vgl_f
#+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, &
point_num, &
nucl_num, &
nucleus_prim_index, &
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 :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
integer (c_int64_t) , intent(in) :: nucleus_prim_index(nucl_num)
real (c_double ) , intent(in) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
real (c_double ) , intent(in) :: expo(prim_num)
real (c_double ) , intent(out) :: primitive_vgl(prim_num,5,point_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, &
point_num, &
nucl_num, &
nucleus_prim_index, &
coord, &
nucl_coord, &
expo, &
primitive_vgl)
end function qmckl_compute_ao_basis_primitive_gaussian_vgl
#+end_src
*** Provide :noexport:
#+CALL: write_provider_header( group="ao_basis", data="primitive_vgl" )
@ -3924,29 +3871,33 @@ print ( "[7][4][26] : %e"% lf(a,x,y))
#+end_src
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_ao_basis_shell_gaussian_vgl_f( &
function qmckl_compute_ao_basis_shell_gaussian_vgl( &
context, prim_num, shell_num, point_num, nucl_num, &
nucleus_shell_num, nucleus_index, nucleus_range, &
shell_prim_index, shell_prim_num, coord, nucl_coord, &
expo, coef_normalized, shell_vgl) &
result(info)
use qmckl
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
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) :: point_num
integer*8 , intent(in) :: nucleus_shell_num(nucl_num)
integer*8 , intent(in) :: nucleus_index(nucl_num)
double precision , intent(in) :: nucleus_range(nucl_num)
integer*8 , intent(in) :: shell_prim_index(shell_num)
integer*8 , intent(in) :: shell_prim_num(shell_num)
double precision , intent(in) :: coord(point_num,3)
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,5,point_num)
integer (qmckl_context), 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 :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(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) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_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,5,point_num)
integer(qmckl_exit_code) :: info
integer*8 :: inucl, iprim, ipoint, ishell
integer*8 :: ishell_start, ishell_end
@ -4020,71 +3971,10 @@ integer function qmckl_compute_ao_basis_shell_gaussian_vgl_f( &
end do
end function qmckl_compute_ao_basis_shell_gaussian_vgl_f
#+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, &
point_num, &
nucl_num, &
nucleus_shell_num, &
nucleus_index, &
nucleus_range, &
shell_prim_index, &
shell_prim_num, &
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 :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(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) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_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,5,point_num)
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, &
point_num, &
nucl_num, &
nucleus_shell_num, &
nucleus_index, &
nucleus_range, &
shell_prim_index, &
shell_prim_num, &
coord, &
nucl_coord, &
expo, &
coef_normalized, &
shell_vgl)
end function qmckl_compute_ao_basis_shell_gaussian_vgl
#+end_src
*** Provide :noexport:
#+CALL: write_provider_header( group="ao_basis", data="shell_vgl" )
@ -4351,17 +4241,22 @@ print ( "[1][4][26] : %25.15e"% lf(a,x,y))
#+end_src
#+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
function qmckl_ao_power(context, n, X, LMAX, P, ldp) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
integer*8 :: i,k
use qmckl_constants
implicit none
integer (qmckl_context), 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(qmckl_exit_code) :: info
integer(c_int64_t) :: i,k
info = QMCKL_SUCCESS
@ -4393,7 +4288,7 @@ integer function qmckl_ao_power_f(context, n, X, LMAX, P, ldp) result(info)
end do
end do
end function qmckl_ao_power_f
end function qmckl_ao_power
#+end_src
#+CALL: generate_f_interface(table=qmckl_ao_power_args,rettyp=get_value("CRetType"),fname="qmckl_ao_power")
@ -4401,14 +4296,14 @@ end function qmckl_ao_power_f
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_ao_power &
integer(qmckl_exit_code) 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 (qmckl_context), 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)
@ -4419,39 +4314,15 @@ end function qmckl_ao_power_f
end interface
#+end_src
#+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
*** Test :noexport:
#+begin_src f90 :tangle (eval f_test)
integer(c_int32_t) function test_qmckl_ao_power(context) bind(C)
function test_qmckl_ao_power(context) bind(C)
use qmckl
implicit none
integer(qmckl_context), intent(in), value :: context
integer(qmckl_exit_code) :: test_qmckl_ao_power
integer*8 :: n, LDP
integer, allocatable :: LMAX(:)
@ -4614,26 +4485,33 @@ qmckl_ao_polynomial_vgl (const qmckl_context context,
{
#ifdef HAVE_HPC
//return qmckl_ao_polynomial_vgl_hpc (context, X, R, lmax, n, L, ldl, VGL, ldv);
return qmckl_ao_polynomial_vgl_doc (context, X, R, lmax, n, L, ldl, VGL, ldv);
return qmckl_ao_polynomial_vgl_doc
#else
return qmckl_ao_polynomial_vgl_doc (context, X, R, lmax, n, L, ldl, VGL, ldv);
return qmckl_ao_polynomial_vgl_doc
#endif
(context, X, R, lmax, n, L, ldl, VGL, ldv);
}
#+end_src
#+begin_src f90 :tangle (eval f)
integer function qmckl_ao_polynomial_vgl_doc_f (context, &
X, R, lmax, n, L, ldl, VGL, ldv) result(info)
use qmckl
function qmckl_ao_polynomial_vgl_doc (context, &
X, R, lmax, n, L, ldl, VGL, ldv) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
implicit none
integer*8, intent(in) :: context
double precision, 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
double precision, intent(out) :: VGL(ldv,(lmax+1)*(lmax+2)*(lmax+3)/6)
integer*8, intent(in) :: ldv
integer (qmckl_context), 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,(lmax+1)*(lmax+2)*(lmax+3)/6)
integer (c_int64_t) , intent(in) , value :: ldl
real (c_double ) , intent(out) :: VGL(ldv,(lmax+1)*(lmax+2)*(lmax+3)/6)
integer (c_int64_t) , intent(in) , value :: ldv
integer(qmckl_exit_code) :: info
integer*8 :: i,j
integer :: a,b,c,d
@ -4747,34 +4625,6 @@ integer function qmckl_ao_polynomial_vgl_doc_f (context, &
info = QMCKL_SUCCESS
end function qmckl_ao_polynomial_vgl_doc_f
#+end_src
#+CALL: generate_c_interface(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_ao_polynomial_vgl_doc" )
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_ao_polynomial_vgl_doc &
(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_doc_f
info = qmckl_ao_polynomial_vgl_doc_f &
(context, X, R, lmax, n, L, ldl, VGL, ldv)
end function qmckl_ao_polynomial_vgl_doc
#+end_src
@ -4783,14 +4633,14 @@ end function qmckl_ao_polynomial_vgl_doc_f
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_ao_polynomial_vgl_doc &
integer(qmckl_exit_code) function qmckl_ao_polynomial_vgl_doc &
(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
integer (qmckl_context), 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
@ -4809,14 +4659,14 @@ end function qmckl_ao_polynomial_vgl_doc_f
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_ao_polynomial_vgl &
integer(qmckl_exit_code) 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
integer (qmckl_context), 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
@ -4891,26 +4741,34 @@ qmckl_ao_polynomial_transp_vgl (const qmckl_context context,
const int64_t ldv )
{
#ifdef HAVE_HPC
return qmckl_ao_polynomial_transp_vgl_hpc (context, X, R, lmax, n, L, ldl, VGL, ldv);
return qmckl_ao_polynomial_transp_vgl_hpc
#else
return qmckl_ao_polynomial_transp_vgl_doc (context, X, R, lmax, n, L, ldl, VGL, ldv);
return qmckl_ao_polynomial_transp_vgl_doc
#endif
(context, X, R, lmax, n, L, ldl, VGL, ldv);
}
#+end_src
#+begin_src f90 :tangle (eval f)
integer function qmckl_ao_polynomial_transp_vgl_doc_f (context, &
X, R, lmax, n, L, ldl, VGL, ldv) result(info)
use qmckl
function qmckl_ao_polynomial_transp_vgl_doc (context, &
X, R, lmax, n, L, ldl, VGL, ldv) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
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,5)
integer*8 , intent(in) :: ldv
integer (qmckl_context), 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,(lmax+1)*(lmax+2)*(lmax+3)/6)
integer (c_int64_t) , intent(in) , value :: ldl
real (c_double ) , intent(out) :: VGL(ldv,5)
integer (c_int64_t) , intent(in) , value :: ldv
integer(qmckl_exit_code) :: info
integer*8 :: i,j
integer :: a,b,c,d
@ -5023,7 +4881,7 @@ integer function qmckl_ao_polynomial_transp_vgl_doc_f (context, &
info = QMCKL_SUCCESS
end function qmckl_ao_polynomial_transp_vgl_doc_f
end function qmckl_ao_polynomial_transp_vgl_doc
#+end_src
#+begin_src c :tangle (eval c) :comments org
@ -5319,80 +5177,24 @@ qmckl_ao_polynomial_transp_vgl_hpc (const qmckl_context context,
}
#+end_src
#+CALL: generate_c_interface(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("CRetType"),fname="qmckl_ao_polynomial_transp_vgl_doc")
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_ao_polynomial_transp_vgl_doc &
(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_transp_vgl_doc_f
info = qmckl_ao_polynomial_transp_vgl_doc_f &
(context, X, R, lmax, n, L, ldl, VGL, ldv)
end function qmckl_ao_polynomial_transp_vgl_doc
#+end_src
#+CALL: generate_f_interface(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("FRetType"),fname="qmckl_ao_polynomial_transp_vgl_doc")
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_ao_polynomial_transp_vgl_doc &
integer(qmckl_exit_code) function qmckl_ao_polynomial_transp_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
integer (qmckl_context), 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_transp_vgl_doc
end interface
#+end_src
#+CALL: generate_f_interface(table=qmckl_ao_polynomial_vgl_args,rettyp=get_value("FRetType"),fname="qmckl_ao_polynomial_transp_vgl")
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_ao_polynomial_transp_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)
real (c_double ) , intent(out) :: VGL(ldv,5)
integer (c_int64_t) , intent(in) , value :: ldv
end function qmckl_ao_polynomial_transp_vgl
@ -5402,11 +5204,11 @@ qmckl_ao_polynomial_transp_vgl_hpc (const qmckl_context context,
*** Test :noexport:
#+begin_src f90 :tangle (eval f_test)
integer(c_int32_t) function test_qmckl_ao_polynomial_vgl(context) bind(C)
function test_qmckl_ao_polynomial_vgl(context) bind(C)
use qmckl
implicit none
integer(c_int64_t), intent(in), value :: context
integer(qmckl_context), intent(in), value :: context
integer :: lmax, d, i
integer, allocatable :: L(:,:)
@ -5415,6 +5217,7 @@ integer(c_int32_t) function test_qmckl_ao_polynomial_vgl(context) bind(C)
double precision, allocatable :: VGL(:,:)
double precision :: w
double precision :: epsilon
integer(qmckl_exit_code) :: test_qmckl_ao_polynomial_vgl
epsilon = qmckl_get_numprec_epsilon(context)
@ -5565,29 +5368,33 @@ for (int32_t ldl=3 ; ldl<=5 ; ++ldl) {
| ~ao_value~ | ~double[point_num][ao_num]~ | out | Values of the AOs |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_ao_value_doc_f(context, &
function qmckl_compute_ao_value_doc(context, &
ao_num, shell_num, point_num, nucl_num, &
coord, nucl_coord, nucleus_index, nucleus_shell_num, &
nucleus_range, nucleus_max_ang_mom, shell_ang_mom, &
ao_factor, shell_vgl, ao_value) &
result(info)
use qmckl
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
use qmckl, only: qmckl_ao_polynomial_vgl
implicit none
integer(qmckl_context), intent(in) :: context
integer*8 , intent(in) :: ao_num
integer*8 , intent(in) :: shell_num
integer*8 , intent(in) :: point_num
integer*8 , intent(in) :: nucl_num
double precision , intent(in) :: coord(point_num,3)
double precision , intent(in) :: nucl_coord(nucl_num,3)
integer*8 , intent(in) :: nucleus_index(nucl_num)
integer*8 , intent(in) :: nucleus_shell_num(nucl_num)
double precision , intent(in) :: nucleus_range(nucl_num)
integer , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer , intent(in) :: shell_ang_mom(shell_num)
double precision , intent(in) :: ao_factor(ao_num)
double precision , intent(in) :: shell_vgl(shell_num,5,point_num)
double precision , intent(out) :: ao_value(ao_num,point_num)
integer (qmckl_context), intent(in) , value :: context
integer (c_int64_t) , intent(in) , value :: ao_num
integer (c_int64_t) , intent(in) , value :: shell_num
integer (c_int64_t) , intent(in) , value :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
real (c_double ) , intent(in) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(nucl_num)
integer (c_int32_t) , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer (c_int32_t) , intent(in) :: shell_ang_mom(shell_num)
real (c_double ) , intent(in) :: ao_factor(ao_num)
real (c_double ) , intent(in) :: shell_vgl(shell_num,5,point_num)
real (c_double ) , intent(out) :: ao_value(ao_num,point_num)
integer(qmckl_exit_code) :: info
double precision :: e_coord(3), n_coord(3)
integer*8 :: n_poly
@ -5597,7 +5404,6 @@ integer function qmckl_compute_ao_value_doc_f(context, &
integer :: lstart(0:20)
double precision :: x, y, z, r2
double precision :: cutoff
integer, external :: qmckl_ao_polynomial_vgl_doc_f
double precision, allocatable :: poly_vgl(:,:)
integer , allocatable :: powers(:,:), ao_index(:)
@ -5646,7 +5452,7 @@ integer function qmckl_compute_ao_value_doc_f(context, &
end if
! Compute polynomials
info = qmckl_ao_polynomial_vgl_doc_f(context, e_coord, n_coord, &
info = qmckl_ao_polynomial_vgl(context, e_coord, n_coord, &
nucleus_max_ang_mom(inucl), n_poly, powers, 3_8, &
poly_vgl, 5_8)
@ -5667,7 +5473,7 @@ integer function qmckl_compute_ao_value_doc_f(context, &
end do
deallocate(poly_vgl, powers)
end function qmckl_compute_ao_value_doc_f
end function qmckl_compute_ao_value_doc
#+end_src
*** HPC version
@ -5994,68 +5800,6 @@ qmckl_compute_ao_value_hpc_gaussian (const qmckl_context context,
#endif
#+end_src
#+CALL: generate_c_interface(table=qmckl_ao_value_args_doc,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_value_doc"))
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_compute_ao_value_doc &
(context, &
ao_num, &
shell_num, &
point_num, &
nucl_num, &
coord, &
nucl_coord, &
nucleus_index, &
nucleus_shell_num, &
nucleus_range, &
nucleus_max_ang_mom, &
shell_ang_mom, &
ao_factor, &
shell_vgl, &
ao_value) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
implicit none
integer (c_int64_t) , intent(in) , value :: context
integer (c_int64_t) , intent(in) , value :: ao_num
integer (c_int64_t) , intent(in) , value :: shell_num
integer (c_int64_t) , intent(in) , value :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
real (c_double ) , intent(in) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(nucl_num)
integer (c_int32_t) , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer (c_int32_t) , intent(in) :: shell_ang_mom(shell_num)
real (c_double ) , intent(in) :: ao_factor(ao_num)
real (c_double ) , intent(in) :: shell_vgl(shell_num,5,point_num)
real (c_double ) , intent(out) :: ao_value(ao_num,point_num)
integer(c_int32_t), external :: qmckl_compute_ao_value_doc_f
info = qmckl_compute_ao_value_doc_f &
(context, &
ao_num, &
shell_num, &
point_num, &
nucl_num, &
coord, &
nucl_coord, &
nucleus_index, &
nucleus_shell_num, &
nucleus_range, &
nucleus_max_ang_mom, &
shell_ang_mom, &
ao_factor, &
shell_vgl, &
ao_value)
end function qmckl_compute_ao_value_doc
#+end_src
**** Provide :noexport:
#+CALL: write_provider_header( group="ao_basis", data="ao_value" )
@ -6369,29 +6113,32 @@ assert( fabs(ao_value[26][224] - ( 7.175045873560788e-10)) < 1.e-14 );
| ~ao_vgl~ | ~double[point_num][5][ao_num]~ | out | Value, gradients and Laplacian of the AOs |
#+begin_src f90 :comments org :tangle (eval f) :noweb yes
integer function qmckl_compute_ao_vgl_doc_f(context, &
function qmckl_compute_ao_vgl_doc(context, &
ao_num, shell_num, point_num, nucl_num, &
coord, nucl_coord, nucleus_index, nucleus_shell_num, &
nucleus_range, nucleus_max_ang_mom, shell_ang_mom, &
ao_factor, shell_vgl, ao_vgl) &
result(info)
use qmckl
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
use qmckl, only : qmckl_ao_polynomial_vgl
implicit none
integer(qmckl_context), intent(in) :: context
integer*8 , intent(in) :: ao_num
integer*8 , intent(in) :: shell_num
integer*8 , intent(in) :: point_num
integer*8 , intent(in) :: nucl_num
double precision , intent(in) :: coord(point_num,3)
double precision , intent(in) :: nucl_coord(nucl_num,3)
integer*8 , intent(in) :: nucleus_index(nucl_num)
integer*8 , intent(in) :: nucleus_shell_num(nucl_num)
double precision , intent(in) :: nucleus_range(nucl_num)
integer , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer , intent(in) :: shell_ang_mom(shell_num)
double precision , intent(in) :: ao_factor(ao_num)
double precision , intent(in) :: shell_vgl(shell_num,5,point_num)
double precision , intent(out) :: ao_vgl(ao_num,5,point_num)
integer (qmckl_context), intent(in) , value :: context
integer (c_int64_t) , intent(in) , value :: ao_num
integer (c_int64_t) , intent(in) , value :: shell_num
integer (c_int64_t) , intent(in) , value :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
real (c_double ) , intent(in) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(nucl_num)
integer (c_int32_t) , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer (c_int32_t) , intent(in) :: shell_ang_mom(shell_num)
real (c_double ) , intent(in) :: ao_factor(ao_num)
real (c_double ) , intent(in) :: shell_vgl(shell_num,5,point_num)
real (c_double ) , intent(out) :: ao_vgl(ao_num,5,point_num)
integer(qmckl_exit_code) :: info
double precision :: e_coord(3), n_coord(3)
integer*8 :: n_poly
@ -6401,7 +6148,6 @@ integer function qmckl_compute_ao_vgl_doc_f(context, &
integer :: lstart(0:20)
double precision :: x, y, z, r2
double precision :: cutoff
integer, external :: qmckl_ao_polynomial_vgl_doc_f
double precision, allocatable :: poly_vgl(:,:)
integer , allocatable :: powers(:,:), ao_index(:)
@ -6450,7 +6196,7 @@ integer function qmckl_compute_ao_vgl_doc_f(context, &
end if
! Compute polynomials
info = qmckl_ao_polynomial_vgl_doc_f(context, e_coord, n_coord, &
info = qmckl_ao_polynomial_vgl(context, e_coord, n_coord, &
nucleus_max_ang_mom(inucl), n_poly, powers, 3_8, &
poly_vgl, 5_8)
@ -6500,7 +6246,7 @@ integer function qmckl_compute_ao_vgl_doc_f(context, &
end do
deallocate(poly_vgl, powers)
end function qmckl_compute_ao_vgl_doc_f
end function qmckl_compute_ao_vgl_doc
#+end_src
*** HPC version
@ -6905,68 +6651,6 @@ qmckl_compute_ao_vgl_hpc_gaussian (
#endif
#+end_src
#+CALL: generate_c_interface(table=qmckl_ao_vgl_args_doc,rettyp=get_value("CRetType"),fname="qmckl_compute_ao_vgl_doc"))
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_compute_ao_vgl_doc &
(context, &
ao_num, &
shell_num, &
point_num, &
nucl_num, &
coord, &
nucl_coord, &
nucleus_index, &
nucleus_shell_num, &
nucleus_range, &
nucleus_max_ang_mom, &
shell_ang_mom, &
ao_factor, &
shell_vgl, &
ao_vgl) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
implicit none
integer (c_int64_t) , intent(in) , value :: context
integer (c_int64_t) , intent(in) , value :: ao_num
integer (c_int64_t) , intent(in) , value :: shell_num
integer (c_int64_t) , intent(in) , value :: point_num
integer (c_int64_t) , intent(in) , value :: nucl_num
real (c_double ) , intent(in) :: coord(point_num,3)
real (c_double ) , intent(in) :: nucl_coord(nucl_num,3)
integer (c_int64_t) , intent(in) :: nucleus_index(nucl_num)
integer (c_int64_t) , intent(in) :: nucleus_shell_num(nucl_num)
real (c_double ) , intent(in) :: nucleus_range(nucl_num)
integer (c_int32_t) , intent(in) :: nucleus_max_ang_mom(nucl_num)
integer (c_int32_t) , intent(in) :: shell_ang_mom(shell_num)
real (c_double ) , intent(in) :: ao_factor(ao_num)
real (c_double ) , intent(in) :: shell_vgl(shell_num,5,point_num)
real (c_double ) , intent(out) :: ao_vgl(ao_num,5,point_num)
integer(c_int32_t), external :: qmckl_compute_ao_vgl_doc_f
info = qmckl_compute_ao_vgl_doc_f &
(context, &
ao_num, &
shell_num, &
point_num, &
nucl_num, &
coord, &
nucl_coord, &
nucleus_index, &
nucleus_shell_num, &
nucleus_range, &
nucleus_max_ang_mom, &
shell_ang_mom, &
ao_factor, &
shell_vgl, &
ao_vgl)
end function qmckl_compute_ao_vgl_doc
#+end_src
**** Provide :noexport:
#+CALL: write_provider_header( group="ao_basis", data="ao_vgl" )

View File

@ -231,19 +231,19 @@ end function qmckl_distance_sq
This function is more efficient when ~A~ and ~B~ are
transposed.
#+CALL: generate_f_interface(table=qmckl_distance_sq_args,rettyp=get_value("FRetType"),fname=get_value("Name"))
#+CALL: generate_f_interface(table=qmckl_distance_sq_args,fname=get_value("Name"))
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_distance_sq &
integer(qmckl_exit_code) function qmckl_distance_sq &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc) &
bind(C)
use, intrinsic :: iso_c_binding
import
implicit none
integer (c_int64_t) , intent(in) , value :: context
integer (qmckl_context), intent(in) , value :: context
character(c_char ) , intent(in) , value :: transa
character(c_char ) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
@ -467,20 +467,24 @@ end function test_qmckl_distance_sq
*** Source
#+begin_src f90 :tangle (eval f)
integer function qmckl_distance_f(context, transa, transb, m, n, &
function qmckl_distance(context, transa, transb, m, n, &
A, LDA, B, LDB, C, LDC) &
result(info)
use qmckl
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
implicit none
integer(qmckl_context) , intent(in) :: context
character , intent(in) :: transa, transb
integer*8 , intent(in) :: m, n
integer*8 , intent(in) :: lda
real*8 , intent(in) :: A(lda,*)
integer*8 , intent(in) :: ldb
real*8 , intent(in) :: B(ldb,*)
integer*8 , intent(in) :: ldc
real*8 , intent(out) :: C(ldc,*)
integer(qmckl_context), intent(in), value :: context
character(c_char) , intent(in) , value :: transa
character(c_char) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
integer (c_int64_t) , intent(in) , value :: n
real (c_double ) , intent(in) :: A(lda,*)
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(ldc,n)
integer (c_int64_t) , intent(in) , value :: ldc
integer (qmckl_exit_code) :: info
integer*8 :: i,j
real*8 :: x, y, z
@ -605,58 +609,22 @@ integer function qmckl_distance_f(context, transa, transb, m, n, &
end select
end function qmckl_distance_f
#+end_src
*** Performance
This function is more efficient when ~A~ and ~B~ are transposed.
** C interface :noexport:
#+CALL: generate_c_interface(table=qmckl_distance_args,rettyp=get_value("FRetType"),fname=get_value("Name"))
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_distance &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
implicit none
integer (c_int64_t) , intent(in) , value :: context
character(c_char) , intent(in) , value :: transa
character(c_char) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
integer (c_int64_t) , intent(in) , value :: n
real (c_double ) , intent(in) :: A(lda,*)
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(ldc,n)
integer (c_int64_t) , intent(in) , value :: ldc
integer(c_int32_t), external :: qmckl_distance_f
info = qmckl_distance_f &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc)
end function qmckl_distance
#+end_src
#+CALL: generate_f_interface(table=qmckl_distance_args,rettyp=get_value("FRetType"),fname=get_value("Name"))
#+CALL: generate_f_interface(table=qmckl_distance_args,fname="qmckl_distance")
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_distance &
integer(qmckl_exit_code) function qmckl_distance &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc) &
bind(C)
use, intrinsic :: iso_c_binding
import
implicit none
integer (c_int64_t) , intent(in) , value :: context
integer (qmckl_context), intent(in) , value :: context
character(c_char ) , intent(in) , value :: transa
character(c_char ) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
@ -672,6 +640,10 @@ end function qmckl_distance_f
end interface
#+end_src
*** Performance
This function is more efficient when ~A~ and ~B~ are transposed.
*** Test :noexport:
#+begin_src f90 :tangle (eval f_test)
@ -891,21 +863,25 @@ end function test_qmckl_dist
*** Source
#+begin_src f90 :tangle (eval f)
integer function qmckl_distance_rescaled_f(context, transa, transb, m, n, &
function qmckl_distance_rescaled(context, transa, transb, m, n, &
A, LDA, B, LDB, C, LDC, rescale_factor_kappa) &
result(info)
use qmckl
bind(C) result(info)
use, intrinsic :: iso_c_binding
use qmckl_constants
implicit none
integer(qmckl_context) , intent(in) :: context
character , intent(in) :: transa, transb
integer*8 , intent(in) :: m, n
integer*8 , intent(in) :: lda
real*8 , intent(in) :: A(lda,*)
integer*8 , intent(in) :: ldb
real*8 , intent(in) :: B(ldb,*)
integer*8 , intent(in) :: ldc
real*8 , intent(out) :: C(ldc,*)
real*8 , intent(in) :: rescale_factor_kappa
integer (qmckl_context), intent(in) , value :: context
character(c_char ) , intent(in) , value :: transa
character(c_char ) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
integer (c_int64_t) , intent(in) , value :: n
real (c_double ) , intent(in) :: A(lda,*)
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(ldc,n)
integer (c_int64_t) , intent(in) , value :: ldc
real (c_double ) , intent(in) , value :: rescale_factor_kappa
integer(qmckl_exit_code) :: info
integer*8 :: i,j
real*8 :: x, y, z, dist, rescale_factor_kappa_inv
@ -1032,7 +1008,7 @@ integer function qmckl_distance_rescaled_f(context, transa, transb, m, n, &
end select
end function qmckl_distance_rescaled_f
end function qmckl_distance_rescaled
#+end_src
*** Performance
@ -1041,50 +1017,19 @@ end function qmckl_distance_rescaled_f
** C interface :noexport:
#+CALL: generate_c_interface(table=qmckl_distance_rescaled_args,rettyp=get_value("FRetType"),fname=get_value("Name"))
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_distance_rescaled &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc, rescale_factor_kappa) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
implicit none
integer (c_int64_t) , intent(in) , value :: context
character(c_char) , intent(in) , value :: transa
character(c_char) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
integer (c_int64_t) , intent(in) , value :: n
real (c_double ) , intent(in) :: A(lda,*)
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(ldc,n)
integer (c_int64_t) , intent(in) , value :: ldc
real (c_double ) , intent(in) , value :: rescale_factor_kappa
integer(c_int32_t), external :: qmckl_distance_rescaled_f
info = qmckl_distance_rescaled_f &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc, rescale_factor_kappa)
end function qmckl_distance_rescaled
#+end_src
#+CALL: generate_f_interface(table=qmckl_distance_rescaled_args,rettyp=get_value("FRetType"),fname=get_value("Name"))
#+CALL: generate_f_interface(table=qmckl_distance_rescaled_args,fname="qmckl_distance_rescaled")
#+RESULTS:
#+begin_src f90 :tangle (eval fh_func) :comments org :exports none
interface
integer(c_int32_t) function qmckl_distance_rescaled &
integer(qmckl_exit_code) function qmckl_distance_rescaled &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc, rescale_factor_kappa) &
bind(C)
use, intrinsic :: iso_c_binding
import
implicit none
integer (c_int64_t) , intent(in) , value :: context
integer (qmckl_context), intent(in) , value :: context
character(c_char ) , intent(in) , value :: transa
character(c_char ) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
@ -1345,7 +1290,7 @@ end function test_qmckl_dist_rescaled
| ~lda~ | ~int64_t~ | in | Leading dimension of array ~A~ |
| ~B~ | ~double[][ldb]~ | in | Array containing the $n \times 3$ matrix $B$ |
| ~ldb~ | ~int64_t~ | in | Leading dimension of array ~B~ |
| ~C~ | ~double[4][n][ldc]~ | out | Array containing the $4 \times m \times n$ matrix $C$ |
| ~C~ | ~double[n][ldc][4]~ | out | Array containing the $4 \times m \times n$ matrix $C$ |
| ~ldc~ | ~int64_t~ | in | Leading dimension of array ~C~ |
| ~rescale_factor_kappa~ | ~double~ | in | Factor for calculating rescaled distances derivatives |
@ -1383,21 +1328,26 @@ end function test_qmckl_dist_rescaled
#+end_src
#+begin_src f90 :tangle (eval f)
integer function qmckl_distance_rescaled_gl_f(context, transa, transb, m, n, &
function qmckl_distance_rescaled_gl(context, transa, transb, m, n, &
A, LDA, B, LDB, C, LDC, rescale_factor_kappa) &
result(info)
use qmckl
bind(C) result(info)
use qmckl_constants
use, intrinsic :: iso_c_binding
implicit none
integer(qmckl_context) , intent(in) :: context
character , intent(in) :: transa, transb
integer*8 , intent(in) :: m, n
integer*8 , intent(in) :: lda
real*8 , intent(in) :: A(lda,*)
integer*8 , intent(in) :: ldb
real*8 , intent(in) :: B(ldb,*)
integer*8 , intent(in) :: ldc
real*8 , intent(out) :: C(4,ldc,*)
real*8 , intent(in) :: rescale_factor_kappa
integer(qmckl_exit_code) :: info
integer (qmckl_context), intent(in) , value :: context
character(c_char ) , intent(in) , value :: transa
character(c_char ) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
integer (c_int64_t) , intent(in) , value :: n
real (c_double ) , intent(in) :: A(lda,*)
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(4,ldc,n)
integer (c_int64_t) , intent(in) , value :: ldc
real (c_double ) , intent(in) , value :: rescale_factor_kappa
integer*8 :: i,j
real*8 :: x, y, z, dist, dist_inv
@ -1559,55 +1509,25 @@ integer function qmckl_distance_rescaled_gl_f(context, transa, transb, m, n, &
end select
end function qmckl_distance_rescaled_gl_f
end function qmckl_distance_rescaled_gl
#+end_src
This function is more efficient when ~A~ and ~B~ are transposed.
#+CALL: generate_c_interface(table=qmckl_distance_rescaled_gl_args,fname=get_value("Name"))
#+RESULTS:
#+begin_src f90 :tangle (eval f) :comments org :exports none
integer(c_int32_t) function qmckl_distance_rescaled_gl &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc, rescale_factor_kappa) &
bind(C) result(info)
use, intrinsic :: iso_c_binding
implicit none
integer (c_int64_t) , intent(in) , value :: context
character(c_char) , intent(in) , value :: transa
character(c_char) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
integer (c_int64_t) , intent(in) , value :: n
real (c_double ) , intent(in) :: A(lda,*)
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(ldc,n,4)
integer (c_int64_t) , intent(in) , value :: ldc
real (c_double ) , intent(in) , value :: rescale_factor_kappa
integer(c_int32_t), external :: qmckl_distance_rescaled_gl_f
info = qmckl_distance_rescaled_gl_f &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc, rescale_factor_kappa)
end function qmckl_distance_rescaled_gl
#+end_src
#+CALL: generate_f_interface(table=qmckl_distance_rescaled_gl_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_distance_rescaled_gl &
integer(qmckl_exit_code) function qmckl_distance_rescaled_gl &
(context, transa, transb, m, n, A, lda, B, ldb, C, ldc, rescale_factor_kappa) &
bind(C)
use, intrinsic :: iso_c_binding
import
implicit none
integer (c_int64_t) , intent(in) , value :: context
integer (qmckl_context), intent(in) , value :: context
character(c_char ) , intent(in) , value :: transa
character(c_char ) , intent(in) , value :: transb
integer (c_int64_t) , intent(in) , value :: m
@ -1616,7 +1536,7 @@ end function qmckl_distance_rescaled_gl_f
integer (c_int64_t) , intent(in) , value :: lda
real (c_double ) , intent(in) :: B(ldb,*)
integer (c_int64_t) , intent(in) , value :: ldb
real (c_double ) , intent(out) :: C(ldc,n,4)
real (c_double ) , intent(out) :: C(4,ldc,n)
integer (c_int64_t) , intent(in) , value :: ldc
real (c_double ) , intent(in) , value :: rescale_factor_kappa

View File

@ -35,10 +35,10 @@
*** Fortran-C type conversions
#+NAME:f_of_c
#+BEGIN_SRC python :var table=test :var rettyp="integer" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval f) :comments org :exports none"
#+BEGIN_SRC python :var table=test :var rettyp="qmckl_exit_code" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval f) :comments org :exports none"
f_of_c_d = { '' : ''
, 'qmckl_context' : 'integer (c_int64_t)'
, 'qmckl_exit_code' : 'integer (c_int32_t)'
, 'qmckl_context' : 'integer (qmckl_context)'
, 'qmckl_exit_code' : 'integer (qmckl_exit_code)'
, 'bool' : 'logical*8'
, 'int32_t' : 'integer (c_int32_t)'
, 'int64_t' : 'integer (c_int64_t)'
@ -53,8 +53,8 @@ f_of_c_d = { '' : ''
#+NAME:c_of_f
#+BEGIN_SRC python :var table=test :var rettyp="integer" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval f) :comments org :exports none"
ctypeid_d = { '' : ''
, 'qmckl_context' : 'integer(c_int64_t)'
, 'qmckl_exit_code' : 'integer(c_int32_t)'
, 'qmckl_context' : 'integer(qmckl_context)'
, 'qmckl_exit_code' : 'integer(qmckl_exit_code)'
, 'integer' : 'integer(c_int32_t)'
, 'integer*8' : 'integer(c_int64_t)'
, 'real' : 'real(c_float)'
@ -164,7 +164,7 @@ return template
*** Generates a C interface to the Fortran function
#+NAME: generate_c_interface
#+BEGIN_SRC python :var table=[] :var rettyp="integer" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval f) :comments org :exports none"
#+BEGIN_SRC python :var table=[] :var rettyp="qmckl_exit_code" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval f) :comments org :exports none"
<<c_of_f>>
<<f_of_c>>
<<parse_table>>
@ -220,7 +220,7 @@ return results
*** Generates a Fortran interface to the C function
#+NAME: generate_f_interface
#+BEGIN_SRC python :var table=test :var rettyp="integer" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval fh_func) :comments org :exports none"
#+BEGIN_SRC python :var table=test :var rettyp="qmckl_exit_code" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval fh_func) :comments org :exports none"
<<c_of_f>>
<<f_of_c>>
<<parse_table>>
@ -269,7 +269,7 @@ return results
#+END_SRC
#+NAME: generate_private_f_interface
#+BEGIN_SRC python :var table=test :var rettyp="integer" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval fh_private_func) :comments org :exports none"
#+BEGIN_SRC python :var table=test :var rettyp="qmckl_exit_code" :var fname=[] :results value :noweb yes :wrap "src f90 :tangle (eval fh_private_func) :comments org :exports none"
<<c_of_f>>
<<f_of_c>>
<<parse_table>>