mirror of
https://github.com/TREX-CoE/irpjast.git
synced 2024-12-22 20:36:08 +01:00
commit
71a78333d8
2
Makefile
2
Makefile
@ -1,5 +1,5 @@
|
|||||||
IRPF90 = irpf90 --codelet=factor_een:100000 #-a -d
|
IRPF90 = irpf90 --codelet=factor_een:100000 #-a -d
|
||||||
FC = ifort -xHost -g
|
FC = ifort -xHost -g -mkl=sequential
|
||||||
FCFLAGS= -O2 -ffree-line-length-none -I .
|
FCFLAGS= -O2 -ffree-line-length-none -I .
|
||||||
NINJA = ninja
|
NINJA = ninja
|
||||||
AR = ar
|
AR = ar
|
||||||
|
36
README.org
36
README.org
@ -1,3 +1,37 @@
|
|||||||
* IRPJAST
|
* IRPJAST
|
||||||
|
|
||||||
CHAMP's Jastrow factor computation using the IRPF90 method
|
CHAMP's Jastrow factor computation using the IRPF90 method
|
||||||
|
|
||||||
|
Original equation:
|
||||||
|
|
||||||
|
$$
|
||||||
|
\sum_{i=2}^{Ne} \sum_{j=1}^i \sum_{pkl} \sum_a^{Nn} c_{apkl}\, r_{ij}^k\, ( R_{ia}^l + R_{ja}^l) ( R_{ia} R_{ja})^m
|
||||||
|
$$
|
||||||
|
|
||||||
|
Expanding, one obtains:
|
||||||
|
|
||||||
|
$$
|
||||||
|
\sum_{i=2}^{Ne} \sum_{j=1}^i \sum_{pkl} \sum_a^{Nn} c_{apkl} R_{ia}^{p-k-l}\, r_{ij}^k\, R_{ja}^{p-k+l} + c_{apkl} R_{ia}^{p-k+l}\, r_{ij}^k\, R_{ja}^{p-k-l}
|
||||||
|
$$
|
||||||
|
|
||||||
|
The equation is symmetric if we exchange $i$ and $j$, so we can rewrite
|
||||||
|
|
||||||
|
$$
|
||||||
|
\sum_{i=1}^{Ne} \sum_{j=1}^{Ne} \sum_{pkl} \sum_a^{Nn} c_{apkl} R_{ia}^{p-k-l}\, r_{ij}^k\, R_{ja}^{p-k+l}
|
||||||
|
$$
|
||||||
|
|
||||||
|
If we reshape $R_{ja}^p$ as a matrix $R_{j,al}$ of size
|
||||||
|
$N_e \times N_n(N_c+1)$,
|
||||||
|
for every $k$, we can pre-compute the matrix product
|
||||||
|
|
||||||
|
$$
|
||||||
|
C_{i,al}^{k} = \sum_j r_{ij}^k\, R_{i,al}
|
||||||
|
$$
|
||||||
|
which can be computed efficiently with BLAS.
|
||||||
|
We can express the total Jastrow as:
|
||||||
|
|
||||||
|
$$
|
||||||
|
\sum_{i=1}^{Ne} \sum_{pkl} \sum_a^{Nn}
|
||||||
|
c_{apkl} R_{ia}^{p-k-l}\, C_{i,a(p-k+l)}^k
|
||||||
|
$$
|
||||||
|
|
||||||
|
@ -1,41 +1,50 @@
|
|||||||
BEGIN_PROVIDER [ double precision, factor_een ]
|
BEGIN_PROVIDER [ double precision, factor_een ]
|
||||||
implicit none
|
implicit none
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! ElectronE-electron-nuclei contribution to Jastrow factor
|
! ElectronE-electron-nuclei contribution to Jastrow factor
|
||||||
END_DOC
|
!
|
||||||
integer :: i, j, a, p, k, l, lmax, m
|
! 5436.20340250000
|
||||||
double precision :: rjam_cn
|
END_DOC
|
||||||
double precision :: cn
|
integer :: i, j, a, p, k, l, lmax, m, n
|
||||||
|
double precision :: cn, accu2, accu
|
||||||
|
double precision :: f(nnuc,0:ncord-2,0:ncord-2)
|
||||||
|
double precision :: tmp_c(nelec,nnuc,0:ncord,0:ncord-1)
|
||||||
|
|
||||||
|
! factor_een = factor_een_blas
|
||||||
|
! return
|
||||||
|
|
||||||
factor_een = 0.0d0
|
factor_een = 0.0d0
|
||||||
|
|
||||||
do p = 2, ncord
|
do p = 2, ncord
|
||||||
do k = 0, p - 1
|
do k = 0, p - 1
|
||||||
if (k /= 0) then
|
if (k /= 0) then
|
||||||
lmax = p - k
|
lmax = p - k
|
||||||
else
|
else
|
||||||
lmax = p - k - 2
|
lmax = p - k - 2
|
||||||
endif
|
endif
|
||||||
do l = 0, lmax
|
do l = 0, lmax
|
||||||
if ( iand(p - k - l, 1) == 1) cycle
|
if ( iand(p - k - l, 1) == 1) cycle
|
||||||
m = (p - k - l) / 2
|
|
||||||
do a = 1, nnuc
|
m = (p - k - l) / 2
|
||||||
cn = cord_vect_lkp(l, k, p, typenuc_arr(a))
|
|
||||||
rjam_cn = rescale_een_n(2, a, m) * cn
|
do a = 1, nnuc
|
||||||
factor_een = factor_een + rescale_een_e(1,2,k) * &
|
accu2 = 0.d0
|
||||||
(rescale_een_n(1,a,l) + rescale_een_n(2,a,l)) * &
|
cn = cord_vect_lkp(l, k, p, typenuc_arr(a))
|
||||||
rescale_een_n(1,a,m) * rjam_cn
|
do j = 1, nelec
|
||||||
do j = 3, nelec
|
accu = 0.d0
|
||||||
rjam_cn = rescale_een_n(j, a, m) * cn
|
do i = 1, nelec
|
||||||
do i = 1, j - 1
|
accu = accu + &
|
||||||
factor_een = factor_een + rescale_een_e(i,j,k) * &
|
rescale_een_e(i,j,k) * &
|
||||||
(rescale_een_n(i,a,l) + rescale_een_n(j,a,l)) * &
|
rescale_een_n(i,a,m) * &
|
||||||
rescale_een_n(i,a,m) * rjam_cn
|
rescale_een_n(j,a,m+l)
|
||||||
enddo
|
enddo
|
||||||
enddo
|
accu2 = accu2 + accu
|
||||||
enddo
|
enddo
|
||||||
|
factor_een = factor_een + accu2 * cn
|
||||||
enddo
|
enddo
|
||||||
enddo
|
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
END_PROVIDER
|
END_PROVIDER
|
||||||
@ -78,7 +87,7 @@ BEGIN_PROVIDER [ double precision, factor_een_deriv_e, (4, nelec) ]
|
|||||||
enddo
|
enddo
|
||||||
|
|
||||||
do i = 1, nelec
|
do i = 1, nelec
|
||||||
rial = rescale_een_n(i, a, l)
|
rial = rescale_een_n(i, a, l) + rjal
|
||||||
riam = rescale_een_n(i, a, m)
|
riam = rescale_een_n(i, a, m)
|
||||||
rijk = rescale_een_e(i, j, k)
|
rijk = rescale_een_e(i, j, k)
|
||||||
|
|
||||||
@ -86,24 +95,24 @@ BEGIN_PROVIDER [ double precision, factor_een_deriv_e, (4, nelec) ]
|
|||||||
drijk(ii) = rescale_een_e_deriv_e(ii, j, i, k)
|
drijk(ii) = rescale_een_e_deriv_e(ii, j, i, k)
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
v1 = rijk * (rial + rjal) ! v(x)
|
v1 = rijk * rial ! v(x)
|
||||||
v2 = rjam_cn * riam ! u(x)
|
v2 = rjam_cn * riam ! u(x)
|
||||||
|
|
||||||
lap1 = 0.0d0
|
lap1 = 0.0d0
|
||||||
lap2 = 0.0d0
|
lap2 = 0.0d0
|
||||||
do ii = 1, 3
|
do ii = 1, 3
|
||||||
d1 = drijk(ii) * (rial + rjal) + rijk * drjal(ii)
|
d1 = drijk(ii) * rial + rijk * drjal(ii)
|
||||||
d2 = drjam_cn(ii) * riam
|
d2 = drjam_cn(ii) * riam
|
||||||
lap1 = lap1 + d1 * d2
|
lap1 = lap1 + d1 * d2
|
||||||
lap2 = lap2 + drijk(ii) * drjal(ii)
|
lap2 = lap2 + drijk(ii) * drjal(ii)
|
||||||
factor_een_deriv_e(ii, j) += v1 * d2 + d1 * v2
|
factor_een_deriv_e(ii, j) = factor_een_deriv_e(ii, j) + v1 * d2 + d1 * v2
|
||||||
enddo
|
enddo
|
||||||
|
|
||||||
! v(x) u''(x) + 2 * u'(x) v'(x) + u(x) v''(x)
|
! v(x) u''(x) + 2 * u'(x) v'(x) + u(x) v''(x)
|
||||||
ii = 4
|
ii = 4
|
||||||
d1 = drijk(ii) * (rial + rjal) + rijk * drjal(ii) + 2.0d0 * lap2
|
d1 = drijk(ii) * rial + rijk * drjal(ii) + lap2 + lap2
|
||||||
d2 = drjam_cn(ii) * riam
|
d2 = drjam_cn(ii) * riam
|
||||||
factor_een_deriv_e(ii, j) += v1 * d2 + d1 * v2 + 2.0d0 * lap1
|
factor_een_deriv_e(ii, j) = factor_een_deriv_e(ii, j) + v1 * d2 + d1 * v2 + lap1 + lap1
|
||||||
|
|
||||||
enddo
|
enddo
|
||||||
enddo
|
enddo
|
||||||
|
126
el_nuc_el_blas.irp.f
Normal file
126
el_nuc_el_blas.irp.f
Normal file
@ -0,0 +1,126 @@
|
|||||||
|
BEGIN_PROVIDER [ double precision, factor_een_blas ]
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! ElectronE-electron-nuclei contribution to Jastrow factor
|
||||||
|
!
|
||||||
|
! 4124.84239750000
|
||||||
|
END_DOC
|
||||||
|
integer :: i, j, a, p, k, l, lmax, m, n
|
||||||
|
double precision :: cn(nnuc), accu
|
||||||
|
double precision :: f(nnuc,0:ncord-2,0:ncord-2)
|
||||||
|
double precision :: tmp_c(nelec,nnuc,0:ncord,0:ncord-1)
|
||||||
|
|
||||||
|
factor_een_blas = 0.0d0
|
||||||
|
|
||||||
|
! r_{ij}^k . R_{ja}^l -> tmp_c_{ia}^{kl}
|
||||||
|
do k=0,ncord-1
|
||||||
|
call dgemm('N','N', nelec, nnuc*(ncord+1), nelec, 1.d0, &
|
||||||
|
rescale_een_e(1,1,k), size(rescale_een_e,1), &
|
||||||
|
rescale_een_n(1,1,0), size(rescale_een_n,1), 0.d0, &
|
||||||
|
tmp_c(1,1,0,k), size(tmp_c,1))
|
||||||
|
enddo
|
||||||
|
|
||||||
|
do p = 2, ncord
|
||||||
|
do k = 0, p - 1
|
||||||
|
m = p-k
|
||||||
|
if (k > 0) then
|
||||||
|
lmax = m
|
||||||
|
else
|
||||||
|
lmax = m - 2
|
||||||
|
endif
|
||||||
|
|
||||||
|
n = shiftr(m,1)
|
||||||
|
do l = iand(m, 1), lmax, 2
|
||||||
|
|
||||||
|
do a = 1, nnuc
|
||||||
|
cn(a) = cord_vect_lkp(l, k, p, typenuc_arr(a))
|
||||||
|
enddo
|
||||||
|
|
||||||
|
do a = 1, nnuc
|
||||||
|
accu = 0.d0
|
||||||
|
do i=1,nelec
|
||||||
|
accu = accu + rescale_een_n(i,a,n) * tmp_c(i,a,n+l,k)
|
||||||
|
enddo
|
||||||
|
factor_een_blas = factor_een_blas + accu * cn(a)
|
||||||
|
enddo
|
||||||
|
n = n-1
|
||||||
|
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
|
||||||
|
END_PROVIDER
|
||||||
|
|
||||||
|
BEGIN_PROVIDER [ double precision, factor_een_deriv_e_blas, (4, nelec) ]
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! Dimensions 1-3 : dx, dy, dz
|
||||||
|
! Dimension 4 : d2x + d2y + d2z
|
||||||
|
END_DOC
|
||||||
|
integer :: i, ii, j, a, p, k, l, lmax, m
|
||||||
|
double precision :: riam, rjam_cn, rial, rjal, rijk
|
||||||
|
double precision, dimension(4) :: driam, drjam_cn, drial, drjal, drijk
|
||||||
|
double precision :: cn, v1, v2, d1, d2, lap1, lap2
|
||||||
|
|
||||||
|
factor_een_deriv_e_blas = 0.0d0
|
||||||
|
|
||||||
|
do p = 2, ncord
|
||||||
|
do k = 0 , p - 1
|
||||||
|
if (k /= 0) then
|
||||||
|
lmax = p - k
|
||||||
|
else
|
||||||
|
lmax = p - k - 2
|
||||||
|
endif
|
||||||
|
|
||||||
|
do l = 0, lmax
|
||||||
|
if ( iand(p - k - l, 1) == 1) cycle
|
||||||
|
m = (p - k - l) / 2
|
||||||
|
|
||||||
|
do a = 1, nnuc
|
||||||
|
cn = cord_vect_lkp(l, k, p, typenuc_arr(a))
|
||||||
|
|
||||||
|
do j = 1, nelec
|
||||||
|
rjal = rescale_een_n(j, a, l)
|
||||||
|
rjam_cn = rescale_een_n(j, a, m) * cn
|
||||||
|
|
||||||
|
do ii = 1, 4
|
||||||
|
drjal(ii) = rescale_een_n_deriv_e(ii, j, a, l)
|
||||||
|
drjam_cn(ii) = rescale_een_n_deriv_e(ii, j, a, m) * cn
|
||||||
|
enddo
|
||||||
|
|
||||||
|
do i = 1, nelec
|
||||||
|
rial = rescale_een_n(i, a, l) + rjal
|
||||||
|
riam = rescale_een_n(i, a, m)
|
||||||
|
rijk = rescale_een_e(i, j, k)
|
||||||
|
|
||||||
|
do ii = 1, 4
|
||||||
|
drijk(ii) = rescale_een_e_deriv_e(ii, j, i, k)
|
||||||
|
enddo
|
||||||
|
|
||||||
|
v1 = rijk * rial ! v(x)
|
||||||
|
v2 = rjam_cn * riam ! u(x)
|
||||||
|
|
||||||
|
lap1 = 0.0d0
|
||||||
|
lap2 = 0.0d0
|
||||||
|
do ii = 1, 3
|
||||||
|
d1 = drijk(ii) * rial + rijk * drjal(ii)
|
||||||
|
d2 = drjam_cn(ii) * riam
|
||||||
|
lap1 = lap1 + d1 * d2
|
||||||
|
lap2 = lap2 + drijk(ii) * drjal(ii)
|
||||||
|
factor_een_deriv_e_blas(ii, j) = factor_een_deriv_e_blas(ii, j) + v1 * d2 + d1 * v2
|
||||||
|
enddo
|
||||||
|
|
||||||
|
! v(x) u''(x) + 2 * u'(x) v'(x) + u(x) v''(x)
|
||||||
|
ii = 4
|
||||||
|
d1 = drijk(ii) * rial + rijk * drjal(ii) + lap2 + lap2
|
||||||
|
d2 = drjam_cn(ii) * riam
|
||||||
|
factor_een_deriv_e_blas(ii, j) = factor_een_deriv_e_blas(ii, j) + v1 * d2 + d1 * v2 + lap1 + lap1
|
||||||
|
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
|
||||||
|
END_PROVIDER
|
@ -107,24 +107,18 @@ BEGIN_PROVIDER [double precision, rescale_een_e, (nelec, nelec, 0:ncord)]
|
|||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
! R = exp(-kappa r) for electron-electron for $J_{een}$
|
! R = exp(-kappa r) for electron-electron for $J_{een}$
|
||||||
END_DOC
|
END_DOC
|
||||||
integer :: i, j, l
|
integer :: i, j, k, l
|
||||||
double precision :: x
|
double precision :: x
|
||||||
double precision, parameter :: f = dexp(1.d0)
|
double precision, parameter :: f = dexp(1.d0)
|
||||||
|
|
||||||
rescale_een_e(:, :, 0) = 1.d0
|
rescale_een_e(:, :, 0) = 1.d0
|
||||||
|
|
||||||
do j = 1, nelec
|
do l = 1, ncord
|
||||||
do i = 1, j-1
|
k=0
|
||||||
x = dexp(-kappa * elec_dist(i, j))
|
|
||||||
rescale_een_e(i, j, 1) = x
|
|
||||||
rescale_een_e(j, i, 1) = x
|
|
||||||
enddo
|
|
||||||
enddo
|
|
||||||
|
|
||||||
do l = 2, ncord
|
|
||||||
do j = 1, nelec
|
do j = 1, nelec
|
||||||
do i = 1, j-1
|
do i = 1, j-1
|
||||||
x = rescale_een_e(i, j, l-1) * rescale_een_e(i, j, 1)
|
k = k+1
|
||||||
|
x = rescale_een_e_ij(k,l)
|
||||||
rescale_een_e(i, j, l) = x
|
rescale_een_e(i, j, l) = x
|
||||||
rescale_een_e(j, i, l) = x
|
rescale_een_e(j, i, l) = x
|
||||||
enddo
|
enddo
|
||||||
@ -138,6 +132,33 @@ BEGIN_PROVIDER [double precision, rescale_een_e, (nelec, nelec, 0:ncord)]
|
|||||||
enddo
|
enddo
|
||||||
END_PROVIDER
|
END_PROVIDER
|
||||||
|
|
||||||
|
BEGIN_PROVIDER [double precision, rescale_een_e_ij, (nelec*(nelec-1)/2, 0:ncord)]
|
||||||
|
implicit none
|
||||||
|
BEGIN_DOC
|
||||||
|
! R = exp(-kappa r) for electron-electron for $J_{een}$
|
||||||
|
END_DOC
|
||||||
|
integer :: i, j, l,k
|
||||||
|
double precision :: x
|
||||||
|
double precision, parameter :: f = dexp(1.d0)
|
||||||
|
|
||||||
|
rescale_een_e_ij(:, 0) = 1.d0
|
||||||
|
|
||||||
|
k=0
|
||||||
|
do j = 1, nelec
|
||||||
|
do i = 1, j-1
|
||||||
|
k = k+1
|
||||||
|
rescale_een_e_ij(k, 1) = dexp(-kappa * elec_dist(i, j))
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
|
||||||
|
do l = 2, ncord
|
||||||
|
do k=1,(nelec*nelec-nelec)/2
|
||||||
|
rescale_een_e_ij(k, l) = rescale_een_e_ij(k, l-1) * rescale_een_e_ij(k, 1)
|
||||||
|
enddo
|
||||||
|
enddo
|
||||||
|
|
||||||
|
END_PROVIDER
|
||||||
|
|
||||||
BEGIN_PROVIDER [double precision, rescale_een_n, (nelec, nnuc, 0:ncord)]
|
BEGIN_PROVIDER [double precision, rescale_een_n, (nelec, nnuc, 0:ncord)]
|
||||||
implicit none
|
implicit none
|
||||||
BEGIN_DOC
|
BEGIN_DOC
|
||||||
|
Loading…
Reference in New Issue
Block a user