10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-26 15:12:14 +02:00

Merge pull request #207 from scemama/master

Bug fix in mmap
This commit is contained in:
Thomas Applencourt 2017-07-12 13:39:45 -05:00 committed by GitHub
commit bc89110eaf
19 changed files with 141 additions and 122 deletions

2
configure vendored
View File

@ -77,7 +77,7 @@ path_github = {"head": "http://github.com", "tail": "archive/master.tar.gz"}
ocaml = Info(
url='http://raw.github.com/ocaml/opam/master/shell/opam_installer.sh',
description=' Ocaml, Opam and the Core library (it will take some time roughly 20min)',
description=' OCaml, Opam and the Core library (it will take some time roughly 20min)',
default_path=join(QP_ROOT_BIN, "opam"))
m4 = Info(

View File

@ -150,7 +150,7 @@ let of_xyz_file
in
let result =
try
int_of_string x > 0
(int_of_string @@ String.strip x) > 0
with
| Failure "int_of_string" -> false
in

View File

@ -1,5 +1,5 @@
===============
Ocaml scripts
OCaml scripts
===============
This directory contains all the scripts that control the input/output

View File

@ -19,7 +19,7 @@ BEGIN_PROVIDER [double precision, FPS_SPF_Matrix_AO, (AO_num, AO_num)]
END_DOC
double precision, allocatable :: scratch(:,:)
allocate( &
scratch(AO_num_align, AO_num) &
scratch(AO_num, AO_num) &
)
! Compute FP
@ -71,7 +71,7 @@ END_PROVIDER
BEGIN_PROVIDER [ double precision, eigenvalues_Fock_matrix_AO, (AO_num) ]
&BEGIN_PROVIDER [ double precision, eigenvectors_Fock_matrix_AO, (AO_num_align,AO_num) ]
&BEGIN_PROVIDER [ double precision, eigenvectors_Fock_matrix_AO, (AO_num,AO_num) ]
BEGIN_DOC
! Eigenvalues and eigenvectors of the Fock matrix over the AO basis
@ -85,7 +85,7 @@ END_PROVIDER
lwork = 3*AO_num - 1
allocate( &
scratch(AO_num_align,AO_num), &
scratch(AO_num,AO_num), &
work(lwork), &
Xt(AO_num,AO_num) &
)
@ -137,7 +137,7 @@ END_PROVIDER
END_PROVIDER
BEGIN_PROVIDER [ double precision, X_matrix_AO, (AO_num_align,AO_num) ]
BEGIN_PROVIDER [ double precision, X_matrix_AO, (AO_num,AO_num) ]
BEGIN_DOC
! Matrix X = S^{-1/2} obtained by SVD

View File

@ -1,4 +1,4 @@
BEGIN_PROVIDER [ double precision, Fock_matrix_mo, (mo_tot_num_align,mo_tot_num) ]
BEGIN_PROVIDER [ double precision, Fock_matrix_mo, (mo_tot_num,mo_tot_num) ]
&BEGIN_PROVIDER [ double precision, Fock_matrix_diag_mo, (mo_tot_num)]
implicit none
BEGIN_DOC
@ -81,8 +81,8 @@ END_PROVIDER
BEGIN_PROVIDER [ double precision, Fock_matrix_ao_alpha, (ao_num_align, ao_num) ]
&BEGIN_PROVIDER [ double precision, Fock_matrix_ao_beta, (ao_num_align, ao_num) ]
BEGIN_PROVIDER [ double precision, Fock_matrix_ao_alpha, (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, Fock_matrix_ao_beta, (ao_num, ao_num) ]
implicit none
BEGIN_DOC
! Alpha Fock matrix in AO basis set
@ -90,7 +90,6 @@ END_PROVIDER
integer :: i,j
do j=1,ao_num
!DIR$ VECTOR ALIGNED
do i=1,ao_num
Fock_matrix_ao_alpha(i,j) = ao_mono_elec_integral(i,j) + ao_bi_elec_integral_alpha(i,j)
Fock_matrix_ao_beta (i,j) = ao_mono_elec_integral(i,j) + ao_bi_elec_integral_beta (i,j)
@ -100,8 +99,8 @@ END_PROVIDER
END_PROVIDER
BEGIN_PROVIDER [ double precision, ao_bi_elec_integral_alpha, (ao_num_align, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_bi_elec_integral_beta , (ao_num_align, ao_num) ]
BEGIN_PROVIDER [ double precision, ao_bi_elec_integral_alpha, (ao_num, ao_num) ]
&BEGIN_PROVIDER [ double precision, ao_bi_elec_integral_beta , (ao_num, ao_num) ]
use map_module
implicit none
BEGIN_DOC
@ -115,8 +114,6 @@ END_PROVIDER
double precision :: ao_bielec_integral, local_threshold
double precision, allocatable :: ao_bi_elec_integral_alpha_tmp(:,:)
double precision, allocatable :: ao_bi_elec_integral_beta_tmp(:,:)
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: ao_bi_elec_integral_beta_tmp
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: ao_bi_elec_integral_alpha_tmp
ao_bi_elec_integral_alpha = 0.d0
ao_bi_elec_integral_beta = 0.d0
@ -126,13 +123,13 @@ END_PROVIDER
!$OMP PRIVATE(i,j,l,k1,k,integral,ii,jj,kk,ll,i8,keys,values,p,q,r,s,i0,j0,k0,l0, &
!$OMP ao_bi_elec_integral_alpha_tmp,ao_bi_elec_integral_beta_tmp, c0, c1, c2, &
!$OMP local_threshold)&
!$OMP SHARED(ao_num,ao_num_align,HF_density_matrix_ao_alpha,HF_density_matrix_ao_beta,&
!$OMP SHARED(ao_num,HF_density_matrix_ao_alpha,HF_density_matrix_ao_beta,&
!$OMP ao_integrals_map,ao_integrals_threshold, ao_bielec_integral_schwartz, &
!$OMP ao_overlap_abs, ao_bi_elec_integral_alpha, ao_bi_elec_integral_beta)
allocate(keys(1), values(1))
allocate(ao_bi_elec_integral_alpha_tmp(ao_num_align,ao_num), &
ao_bi_elec_integral_beta_tmp(ao_num_align,ao_num))
allocate(ao_bi_elec_integral_alpha_tmp(ao_num,ao_num), &
ao_bi_elec_integral_beta_tmp(ao_num,ao_num))
ao_bi_elec_integral_alpha_tmp = 0.d0
ao_bi_elec_integral_beta_tmp = 0.d0
@ -212,13 +209,13 @@ END_PROVIDER
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE(i,j,l,k1,k,integral,ii,jj,kk,ll,i8,keys,values,n_elements_max, &
!$OMP n_elements,ao_bi_elec_integral_alpha_tmp,ao_bi_elec_integral_beta_tmp)&
!$OMP SHARED(ao_num,ao_num_align,HF_density_matrix_ao_alpha,HF_density_matrix_ao_beta,&
!$OMP SHARED(ao_num,HF_density_matrix_ao_alpha,HF_density_matrix_ao_beta,&
!$OMP ao_integrals_map, ao_bi_elec_integral_alpha, ao_bi_elec_integral_beta)
call get_cache_map_n_elements_max(ao_integrals_map,n_elements_max)
allocate(keys(n_elements_max), values(n_elements_max))
allocate(ao_bi_elec_integral_alpha_tmp(ao_num_align,ao_num), &
ao_bi_elec_integral_beta_tmp(ao_num_align,ao_num))
allocate(ao_bi_elec_integral_alpha_tmp(ao_num,ao_num), &
ao_bi_elec_integral_beta_tmp(ao_num,ao_num))
ao_bi_elec_integral_alpha_tmp = 0.d0
ao_bi_elec_integral_beta_tmp = 0.d0
@ -261,42 +258,40 @@ END_PROVIDER
END_PROVIDER
BEGIN_PROVIDER [ double precision, Fock_matrix_mo_alpha, (mo_tot_num_align,mo_tot_num) ]
BEGIN_PROVIDER [ double precision, Fock_matrix_mo_alpha, (mo_tot_num,mo_tot_num) ]
implicit none
BEGIN_DOC
! Fock matrix on the MO basis
END_DOC
double precision, allocatable :: T(:,:)
allocate ( T(ao_num_align,mo_tot_num) )
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: T
allocate ( T(ao_num,mo_tot_num) )
call dgemm('N','N', ao_num, mo_tot_num, ao_num, &
1.d0, Fock_matrix_ao_alpha,size(Fock_matrix_ao_alpha,1), &
mo_coef, size(mo_coef,1), &
0.d0, T, ao_num_align)
0.d0, T, size(T,1))
call dgemm('T','N', mo_tot_num, mo_tot_num, ao_num, &
1.d0, mo_coef,size(mo_coef,1), &
T, size(T,1), &
0.d0, Fock_matrix_mo_alpha, mo_tot_num_align)
0.d0, Fock_matrix_mo_alpha, size(Fock_matrix_mo_alpha,1))
deallocate(T)
END_PROVIDER
BEGIN_PROVIDER [ double precision, Fock_matrix_mo_beta, (mo_tot_num_align,mo_tot_num) ]
BEGIN_PROVIDER [ double precision, Fock_matrix_mo_beta, (mo_tot_num,mo_tot_num) ]
implicit none
BEGIN_DOC
! Fock matrix on the MO basis
END_DOC
double precision, allocatable :: T(:,:)
allocate ( T(ao_num_align,mo_tot_num) )
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: T
allocate ( T(ao_num,mo_tot_num) )
call dgemm('N','N', ao_num, mo_tot_num, ao_num, &
1.d0, Fock_matrix_ao_beta,size(Fock_matrix_ao_beta,1), &
mo_coef, size(mo_coef,1), &
0.d0, T, ao_num_align)
0.d0, T, size(T,1))
call dgemm('T','N', mo_tot_num, mo_tot_num, ao_num, &
1.d0, mo_coef,size(mo_coef,1), &
T, size(T,1), &
0.d0, Fock_matrix_mo_beta, mo_tot_num_align)
0.d0, Fock_matrix_mo_beta, size(Fock_matrix_mo_beta,1))
deallocate(T)
END_PROVIDER
@ -319,7 +314,7 @@ BEGIN_PROVIDER [ double precision, HF_energy ]
END_PROVIDER
BEGIN_PROVIDER [ double precision, Fock_matrix_ao, (ao_num_align, ao_num) ]
BEGIN_PROVIDER [ double precision, Fock_matrix_ao, (ao_num, ao_num) ]
implicit none
BEGIN_DOC
! Fock matrix in AO basis set
@ -330,8 +325,7 @@ BEGIN_PROVIDER [ double precision, Fock_matrix_ao, (ao_num_align, ao_num) ]
then
integer :: i,j
do j=1,ao_num
!DIR$ VECTOR ALIGNED
do i=1,ao_num_align
do i=1,ao_num
Fock_matrix_ao(i,j) = Fock_matrix_ao_alpha(i,j)
enddo
enddo
@ -339,7 +333,7 @@ BEGIN_PROVIDER [ double precision, Fock_matrix_ao, (ao_num_align, ao_num) ]
double precision, allocatable :: T(:,:), M(:,:)
integer :: ierr
! F_ao = S C F_mo C^t S
allocate (T(ao_num_align,ao_num),M(ao_num_align,ao_num),stat=ierr)
allocate (T(ao_num,ao_num),M(ao_num,ao_num),stat=ierr)
if (ierr /=0 ) then
print *, irp_here, ' : allocation failed'
endif
@ -391,7 +385,7 @@ subroutine Fock_mo_to_ao(FMO,LDFMO,FAO,LDFAO)
double precision, allocatable :: T(:,:), M(:,:)
integer :: ierr
! F_ao = S C F_mo C^t S
allocate (T(ao_num_align,ao_num),M(ao_num_align,ao_num),stat=ierr)
allocate (T(ao_num,ao_num),M(ao_num,ao_num),stat=ierr)
if (ierr /=0 ) then
print *, irp_here, ' : allocation failed'
endif

View File

@ -1,4 +1,4 @@
BEGIN_PROVIDER [double precision, HF_density_matrix_ao_alpha, (ao_num_align,ao_num) ]
BEGIN_PROVIDER [double precision, HF_density_matrix_ao_alpha, (ao_num,ao_num) ]
implicit none
BEGIN_DOC
! S^-1 x Alpha density matrix in the AO basis x S^-1
@ -11,7 +11,7 @@ BEGIN_PROVIDER [double precision, HF_density_matrix_ao_alpha, (ao_num_align,ao_n
END_PROVIDER
BEGIN_PROVIDER [ double precision, HF_density_matrix_ao_beta, (ao_num_align,ao_num) ]
BEGIN_PROVIDER [ double precision, HF_density_matrix_ao_beta, (ao_num,ao_num) ]
implicit none
BEGIN_DOC
! S^-1 Beta density matrix in the AO basis x S^-1
@ -24,7 +24,7 @@ BEGIN_PROVIDER [ double precision, HF_density_matrix_ao_beta, (ao_num_align,ao_
END_PROVIDER
BEGIN_PROVIDER [ double precision, HF_density_matrix_ao, (ao_num_align,ao_num) ]
BEGIN_PROVIDER [ double precision, HF_density_matrix_ao, (ao_num,ao_num) ]
implicit none
BEGIN_DOC
! S^-1 Density matrix in the AO basis S^-1

View File

@ -12,13 +12,13 @@ subroutine damping_SCF
character :: save_char
allocate( &
D_alpha( ao_num_align, ao_num ), &
D_beta( ao_num_align, ao_num ), &
F_new( ao_num_align, ao_num ), &
D_new_alpha( ao_num_align, ao_num ), &
D_new_beta( ao_num_align, ao_num ), &
delta_alpha( ao_num_align, ao_num ), &
delta_beta( ao_num_align, ao_num ))
D_alpha( ao_num, ao_num ), &
D_beta( ao_num, ao_num ), &
F_new( ao_num, ao_num ), &
D_new_alpha( ao_num, ao_num ), &
D_new_beta( ao_num, ao_num ), &
delta_alpha( ao_num, ao_num ), &
delta_beta( ao_num, ao_num ))
do j=1,ao_num
do i=1,ao_num

View File

@ -1,5 +1,5 @@
BEGIN_PROVIDER [ double precision, diagonal_Fock_matrix_mo, (ao_num) ]
&BEGIN_PROVIDER [ double precision, eigenvectors_Fock_matrix_mo, (ao_num_align,mo_tot_num) ]
&BEGIN_PROVIDER [ double precision, eigenvectors_Fock_matrix_mo, (ao_num,mo_tot_num) ]
implicit none
BEGIN_DOC
! Diagonal Fock matrix in the MO basis

View File

@ -17,7 +17,6 @@ subroutine huckel_guess
c = 0.5d0 * 1.75d0
do j=1,ao_num
!DIR$ VECTOR ALIGNED
do i=1,ao_num
Fock_matrix_ao(i,j) = c*ao_overlap(i,j)*(ao_mono_elec_integral_diag(i) + &
ao_mono_elec_integral_diag(j))

View File

@ -8,6 +8,7 @@ BEGIN_PROVIDER [ double precision , ao_value_at_nucl, (ao_num,nucl_num) ]
do k=1,nucl_num
do i=1,ao_num
ao_value_at_nucl(i,k) = 0.d0
x = nucl_coord(ao_nucl(i),1) - nucl_coord(k,1)
y = nucl_coord(ao_nucl(i),2) - nucl_coord(k,2)
z = nucl_coord(ao_nucl(i),3) - nucl_coord(k,3)
@ -15,7 +16,6 @@ BEGIN_PROVIDER [ double precision , ao_value_at_nucl, (ao_num,nucl_num) ]
if (poly == 0.d0) cycle
r2 = (x*x) + (y*y) + (z*z)
ao_value_at_nucl(i,k) = 0.d0
do j=1,ao_prim_num(i)
expo = ao_expo_ordered_transp(j,i)*r2
if (expo > 40.d0) cycle

View File

@ -452,19 +452,19 @@ BEGIN_PROVIDER [ double precision, GauSla$X_matrix, (ao_num, nucl_num) ]
integer :: aGau(3)
!TODO
logical :: read
integer :: iunit
integer :: getunitandopen
inquire(FILE=trim(ezfio_filename)//'/work/GauSla$X.dat',EXIST=read)
if (read) then
print *, 'READ $X'
iunit = getunitandopen(trim(ezfio_filename)//'/work/GauSla$X.dat','r')
else
print *, 'WRITE $X'
iunit = getunitandopen(trim(ezfio_filename)//'/work/GauSla$X.inp','w')
write(iunit,*) '{'
endif
! logical :: read
! integer :: iunit
! integer :: getunitandopen
!
! inquire(FILE=trim(ezfio_filename)//'/work/GauSla$X.dat',EXIST=read)
! if (read) then
! print *, 'READ $X'
! iunit = getunitandopen(trim(ezfio_filename)//'/work/GauSla$X.dat','r')
! else
! print *, 'WRITE $X'
! iunit = getunitandopen(trim(ezfio_filename)//'/work/GauSla$X.inp','w')
! write(iunit,*) '{'
! endif
!TODO
do k=1,nucl_num
@ -478,22 +478,22 @@ BEGIN_PROVIDER [ double precision, GauSla$X_matrix, (ao_num, nucl_num) ]
do j=1,ao_prim_num(i)
expGau = ao_expo_ordered_transp(j,i)
! call GauSla$X(expGau,cGau,aGau,expSla,cSla,res)
if (read) then
call GauSla$X_read(expGau,cGau,aGau,expSla,cSla,res,iunit)
else
call GauSla$X_write(expGau,cGau,aGau,expSla,cSla,res,iunit)
endif
call GauSla$X(expGau,cGau,aGau,expSla,cSla,res)
! if (read) then
! call GauSla$X_read(expGau,cGau,aGau,expSla,cSla,res,iunit)
! else
! call GauSla$X_write(expGau,cGau,aGau,expSla,cSla,res,iunit)
! endif
GauSla$X_matrix(i,k) += ao_coef_normalized_ordered_transp(j,i) * res
enddo
enddo
enddo
if (.not.read) then
write(iunit,*) '0.}'
endif
close(iunit)
! if (.not.read) then
! write(iunit,*) '0.}'
! endif
! close(iunit)
END_PROVIDER

View File

@ -0,0 +1,8 @@
program create_wf
read_wf = .true.
SOFT_TOUCH read_wf
PROVIDE H_apply_buffer_allocated
call generate_all_alpha_beta_det_products
call diagonalize_ci
call save_wavefunction
end

View File

@ -5,11 +5,18 @@ program e_curve
double precision :: norm, E, hij, num, ci, cj
integer, allocatable :: iorder(:)
double precision , allocatable :: norm_sort(:)
double precision :: e_0(N_states)
PROVIDE mo_bielec_integrals_in_map
nab = n_det_alpha_unique+n_det_beta_unique
allocate ( norm_sort(0:nab), iorder(0:nab) )
double precision, allocatable :: u_t(:,:), v_t(:,:), s_t(:,:)
double precision, allocatable :: u_0(:,:), v_0(:,:)
allocate(u_t(N_states,N_det),v_t(N_states,N_det),s_t(N_states,N_det))
allocate(u_0(N_states,N_det),v_0(N_states,N_det))
norm_sort(0) = 0.d0
iorder(0) = 0
@ -59,47 +66,45 @@ program e_curve
if (thresh > norm_sort(j)) then
cycle
endif
num = 0.d0
norm = 0.d0
m = 0
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(k,kk,l,det_i,det_j,ci,cj,hij) REDUCTION(+:norm,m,num)
allocate( det_i(N_int,2), det_j(N_int,2))
!$OMP DO SCHEDULE(guided)
do k=1,n_det
if (psi_bilinear_matrix_values(k,1) == 0.d0) then
cycle
endif
ci = psi_bilinear_matrix_values(k,1)
do kk=1,N_int
det_i(kk,1) = psi_det_alpha_unique(kk,psi_bilinear_matrix_rows(k))
det_i(kk,2) = psi_det_beta_unique(kk,psi_bilinear_matrix_columns(k))
enddo
do l=1,n_det
if (psi_bilinear_matrix_values(l,1) == 0.d0) then
cycle
endif
cj = psi_bilinear_matrix_values(l,1)
do kk=1,N_int
det_j(kk,1) = psi_det_alpha_unique(kk,psi_bilinear_matrix_rows(l))
det_j(kk,2) = psi_det_beta_unique(kk,psi_bilinear_matrix_columns(l))
enddo
call i_h_j(det_i, det_j, N_int, hij)
num = num + ci*cj*hij
enddo
norm = norm + ci*ci
m = m+1
u_0 = psi_bilinear_matrix_values(1:N_det,1:N_states)
v_t = 0.d0
s_t = 0.d0
call dtranspose( &
u_0, &
size(u_0, 1), &
u_t, &
size(u_t, 1), &
N_det, N_states)
call H_S2_u_0_nstates_openmp_work(v_t,s_t,u_t,N_states,N_det,1,N_det,0,1)
call dtranspose( &
v_t, &
size(v_t, 1), &
v_0, &
size(v_0, 1), &
N_states, N_det)
double precision, external :: u_dot_u, u_dot_v
do i=1,N_states
e_0(i) = u_dot_v(v_t(1,i),u_0(1,i),N_det)/u_dot_u(u_0(1,i),N_det)
enddo
!$OMP END DO
deallocate (det_i,det_j)
!$OMP END PARALLEL
m = 0
do k=1,n_det
if (psi_bilinear_matrix_values(k,1) /= 0.d0) then
m = m+1
endif
enddo
if (m == 0) then
exit
endif
E = num / norm + nuclear_repulsion
E = E_0(1) + nuclear_repulsion
norm = u_dot_u(u_0(1,1),N_det)
print '(E9.1,2X,I8,2X,F10.2,2X,F10.8,2X,F12.6)', thresh, m, &
dble( elec_alpha_num**3 + elec_alpha_num**2 * (nab-1) ) / &
dble( elec_alpha_num**3 + elec_alpha_num**2 * (j-1)), norm, E
thresh = thresh * 2.d0
thresh = thresh * dsqrt(10.d0)
enddo
print *, '=========================================================='

View File

@ -106,7 +106,7 @@ def is_bool(str_):
def get_type_dict():
"""
This function makes the correspondance between the type of value read in
ezfio.cfg into the f90 and Ocaml Type
ezfio.cfg into the f90 and OCaml Type
return fancy_type[fancy_type] = namedtuple('Type', 'ocaml fortran')
For example fancy_type['Ndet'].fortran = interger
.ocaml = int
@ -617,7 +617,7 @@ def save_ocaml_input(module_lower, str_ocaml_input):
def get_l_module_with_auto_generate_ocaml_lower():
"""
Get all modules which have EZFIO.cfg with Ocaml data
Get all modules which have EZFIO.cfg with OCaml data
(NB `search` in all the lines and `match` only in one)
"""

View File

@ -304,3 +304,4 @@ subroutine make_s2_eigenfunction
end

View File

@ -3,5 +3,7 @@ program save_natorb
touch read_wf
call save_natural_mos
call save_ref_determinant
call ezfio_set_integrals_bielec_disk_access_mo_integrals('None')
call ezfio_set_integrals_monoelec_disk_access_mo_one_integrals('None')
end

View File

@ -625,9 +625,16 @@ subroutine create_wf_of_psi_bilinear_matrix(truncate)
integer :: idx
integer, external :: get_index_in_psi_det_sorted_bit
double precision :: norm(N_states)
PROVIDE psi_bilinear_matrix
call generate_all_alpha_beta_det_products
norm = 0.d0
!$OMP PARALLEL DO DEFAULT(NONE) &
!$OMP PRIVATE(i,j,k,idx,tmp_det) &
!$OMP SHARED(N_det_alpha_unique, N_det_beta_unique, N_det, &
!$OMP N_int, N_states, norm, psi_det_beta_unique, &
!$OMP psi_det_alpha_unique, psi_bilinear_matrix, &
!$OMP psi_coef_sorted_bit)
do j=1,N_det_beta_unique
do k=1,N_int
tmp_det(k,2) = psi_det_beta_unique(k,j)
@ -640,11 +647,14 @@ subroutine create_wf_of_psi_bilinear_matrix(truncate)
if (idx > 0) then
do k=1,N_states
psi_coef_sorted_bit(idx,k) = psi_bilinear_matrix(i,j,k)
!$OMP ATOMIC
norm(k) += psi_bilinear_matrix(i,j,k)
enddo
endif
enddo
enddo
!$OMP END PARALLEL DO
do k=1,N_states
norm(k) = 1.d0/dsqrt(norm(k))
do i=1,N_det
@ -688,7 +698,7 @@ subroutine generate_all_alpha_beta_det_products
!$OMP PRIVATE(i,j,k,l,tmp_det,iproc)
!$ iproc = omp_get_thread_num()
allocate (tmp_det(N_int,2,N_det_alpha_unique))
!$OMP DO
!$OMP DO SCHEDULE(static,1)
do j=1,N_det_beta_unique
l = 1
do i=1,N_det_alpha_unique

View File

@ -53,17 +53,17 @@ subroutine map_save_to_disk(filename,map)
map % consolidated = .True.
! call munmap( (/ map % map_size + 2_8 /), 8, fd(1), c_pointer(1))
! call mmap(trim(filename)//'_consolidated_idx', (/ map % map_size + 2_8 /), 8, fd(1), .True., c_pointer(1))
! call c_f_pointer(c_pointer(1),map % consolidated_idx, (/ map % map_size +2_8/))
!
! call munmap( (/ map % n_elements /), cache_key_kind, fd(2), c_pointer(2))
! call mmap(trim(filename)//'_consolidated_key', (/ map % n_elements /), cache_key_kind, fd(2), .True., c_pointer(2))
! call c_f_pointer(c_pointer(2),map % consolidated_key, (/ map % n_elements /))
!
! call munmap( (/ map % n_elements /), integral_kind, fd(3), c_pointer(3))
! call mmap(trim(filename)//'_consolidated_value', (/ map % n_elements /), integral_kind, fd(3), .True., c_pointer(3))
! call c_f_pointer(c_pointer(3),map % consolidated_value, (/ map % n_elements /))
call munmap( (/ map % map_size + 2_8 /), 8, fd(1), c_pointer(1))
call mmap(trim(filename)//'_consolidated_idx', (/ map % map_size + 2_8 /), 8, fd(1), .True., c_pointer(1))
call c_f_pointer(c_pointer(1),map % consolidated_idx, (/ map % map_size +2_8/))
call munmap( (/ map % n_elements /), cache_key_kind, fd(2), c_pointer(2))
call mmap(trim(filename)//'_consolidated_key', (/ map % n_elements /), cache_key_kind, fd(2), .True., c_pointer(2))
call c_f_pointer(c_pointer(2),map % consolidated_key, (/ map % n_elements /))
call munmap( (/ map % n_elements /), integral_kind, fd(3), c_pointer(3))
call mmap(trim(filename)//'_consolidated_value', (/ map % n_elements /), integral_kind, fd(3), .True., c_pointer(3))
call c_f_pointer(c_pointer(3),map % consolidated_value, (/ map % n_elements /))
end

View File

@ -36,7 +36,7 @@ module mmap_module
integer, intent(out) :: fd ! File descriptor
type(c_ptr), intent(out) :: map ! C Pointer
integer(c_long) :: length
integer(c_size_t) :: length
integer(c_int) :: fd_
length = PRODUCT( shape(:) ) * bytes
@ -56,7 +56,7 @@ module mmap_module
integer, intent(in) :: fd ! File descriptor
type(c_ptr), intent(in) :: map ! C pointer
integer(c_long) :: length
integer(c_size_t) :: length
integer(c_int) :: fd_
length = PRODUCT( shape(:) ) * bytes