10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-07-04 10:25:58 +02:00
quantum_package/plugins/shiftedbk/shifted_bk_routines.irp.f

320 lines
10 KiB
Fortran
Raw Normal View History

use selection_types
2018-03-21 12:06:26 +01:00
2018-04-10 14:25:28 +02:00
BEGIN_PROVIDER [ double precision, global_sum_alpha2, (N_states) ]
&BEGIN_PROVIDER [ double precision, slave_sum_alpha2, (N_states, Nproc) ]
global_sum_alpha2 = 0d0
slave_sum_alpha2 = 0d0
END_PROVIDER
2018-03-26 12:06:18 +02:00
BEGIN_PROVIDER [ double precision, fock_diag_tmp_, (2,mo_tot_num+1,Nproc) ]
2018-04-04 11:32:27 +02:00
&BEGIN_PROVIDER [ integer, n_det_add ]
2018-03-26 12:06:18 +02:00
&BEGIN_PROVIDER [ double precision, a_h_i, (N_det, Nproc) ]
&BEGIN_PROVIDER [ double precision, a_s2_i, (N_det, Nproc) ]
&BEGIN_PROVIDER [ type(selection_buffer), sb, (Nproc) ]
2018-04-04 11:32:27 +02:00
&BEGIN_PROVIDER [ type(selection_buffer), global_sb ]
&BEGIN_PROVIDER [ type(selection_buffer), mini_sb ]
&BEGIN_PROVIDER [ double precision, N_det_increase_factor ]
implicit none
fock_diag_tmp_(:,:,:) = 0.d0
integer :: i
2018-05-23 18:16:33 +02:00
N_det_increase_factor = dble(N_states)
n_det_add = max(1, int(float(N_det) * N_det_increase_factor))
call create_selection_buffer(n_det_add, n_det_add*2, global_sb)
call create_selection_buffer(n_det_add, n_det_add*2, mini_sb)
do i=1,Nproc
call create_selection_buffer(n_det_add, n_det_add*2, sb(i))
end do
a_h_i = 0d0
a_s2_i = 0d0
END_PROVIDER
2018-03-21 12:06:26 +01:00
2018-04-04 11:32:27 +02:00
BEGIN_PROVIDER [ integer, N_dress_int_buffer ]
&BEGIN_PROVIDER [ integer, N_dress_double_buffer ]
&BEGIN_PROVIDER [ integer, N_dress_det_buffer ]
implicit none
N_dress_int_buffer = 1
2018-04-10 14:25:28 +02:00
N_dress_double_buffer = n_det_add+N_states
2018-04-04 11:32:27 +02:00
N_dress_det_buffer = n_det_add
END_PROVIDER
subroutine generator_done(i_gen, int_buf, double_buf, det_buf, N_buf, iproc)
2018-03-26 12:06:18 +02:00
implicit none
2018-04-04 11:32:27 +02:00
integer, intent(in) :: i_gen, iproc
integer, intent(out) :: int_buf(N_dress_int_buffer), N_buf(3)
double precision, intent(out) :: double_buf(N_dress_double_buffer)
integer(bit_kind), intent(out) :: det_buf(N_int, 2, N_dress_det_buffer)
integer :: i
2018-05-22 19:11:10 +02:00
int_buf(:) = 0
2018-03-30 18:16:00 +02:00
2018-04-04 11:32:27 +02:00
call sort_selection_buffer(sb(iproc))
2018-04-04 11:32:27 +02:00
if(sb(iproc)%cur > 0) then
!$OMP CRITICAL
call merge_selection_buffers(sb(iproc), mini_sb)
2018-04-27 12:41:39 +02:00
!call sort_selection_buffer(mini_sb)
2018-04-04 11:32:27 +02:00
do i=1,Nproc
mini_sb%mini = min(sb(i)%mini, mini_sb%mini)
end do
do i=1,Nproc
sb(i)%mini = mini_sb%mini
2018-04-04 11:32:27 +02:00
end do
!$OMP END CRITICAL
end if
call truncate_to_mini(sb(iproc))
det_buf(:,:,:sb(iproc)%cur) = sb(iproc)%det(:,:,:sb(iproc)%cur)
double_buf(:sb(iproc)%cur) = sb(iproc)%val(:sb(iproc)%cur)
double_buf(sb(iproc)%cur+1:sb(iproc)%cur+N_states) = slave_sum_alpha2(:,iproc)
N_buf(1) = 1
N_buf(2) = sb(iproc)%cur+N_states
N_buf(3) = sb(iproc)%cur
2018-04-10 14:25:28 +02:00
2018-04-04 11:32:27 +02:00
sb(iproc)%cur = 0
2018-04-10 14:25:28 +02:00
slave_sum_alpha2(:,iproc) = 0d0
2018-03-30 18:16:00 +02:00
end subroutine
2018-03-26 12:06:18 +02:00
2018-04-04 11:32:27 +02:00
subroutine generator_start(i_gen, iproc)
implicit none
integer, intent(in) :: i_gen, iproc
integer :: i
call build_fock_tmp(fock_diag_tmp_(1,1,iproc),psi_det_generators(1,1,i_gen),N_int)
end subroutine
2018-03-26 12:06:18 +02:00
2018-04-04 11:32:27 +02:00
subroutine dress_pulled(ind, int_buf, double_buf, det_buf, N_buf)
2018-03-30 18:16:00 +02:00
use bitmasks
implicit none
2018-04-04 11:32:27 +02:00
integer, intent(in) :: ind, N_buf(3)
2018-03-30 18:16:00 +02:00
integer, intent(in) :: int_buf(*)
double precision, intent(in) :: double_buf(*)
integer(bit_kind), intent(in) :: det_buf(N_int,2,*)
2018-04-04 11:32:27 +02:00
integer :: i
2018-04-10 14:25:28 +02:00
do i=1,N_buf(3)
2018-04-04 11:32:27 +02:00
call add_to_selection_buffer(global_sb, det_buf(1,1,i), double_buf(i))
end do
2018-04-10 14:25:28 +02:00
if(N_buf(3) + N_states /= N_buf(2)) stop "buf size"
!$OMP CRITICAL
global_sum_alpha2(:) += double_buf(N_buf(3)+1:N_buf(2))
!$OMP END CRITICAL
2018-03-30 18:16:00 +02:00
end subroutine
2018-03-21 12:06:26 +01:00
subroutine delta_ij_done()
2018-03-23 11:06:52 +01:00
use bitmasks
implicit none
2018-04-04 11:32:27 +02:00
integer :: i, old_det_gen
2018-03-26 13:08:26 +02:00
integer(bit_kind), allocatable :: old_generators(:,:,:)
2018-03-26 13:08:26 +02:00
allocate(old_generators(N_int, 2, N_det_generators))
old_generators(:,:,:) = psi_det_generators(:,:,:N_det_generators)
old_det_gen = N_det_generators
2018-05-23 18:16:33 +02:00
! Add buffer only when the last state is computed
call unique_selection_buffer(global_sb)
call sort_selection_buffer(global_sb)
call fill_H_apply_buffer_no_selection(global_sb%cur,global_sb%det,N_int,0)
call copy_H_apply_buffer_to_wf()
if (s2_eig.or.(N_states > 1) ) then
call make_s2_eigenfunction
endif
2018-05-23 18:16:33 +02:00
call undress_with_alpha(old_generators, old_det_gen, psi_det(1,1,N_det_delta_ij+1), N_det-N_det_delta_ij)
call save_wavefunction
2018-05-22 19:11:10 +02:00
end subroutine
2018-03-26 13:08:26 +02:00
subroutine undress_with_alpha(old_generators, old_det_gen, alpha, n_alpha)
2018-03-23 11:06:52 +01:00
use bitmasks
implicit none
integer(bit_kind), intent(in) :: alpha(N_int,2,n_alpha)
integer, intent(in) :: n_alpha
integer, allocatable :: minilist(:)
integer(bit_kind), allocatable :: det_minilist(:,:,:)
double precision, allocatable :: delta_ij_loc(:,:,:,:)
2018-03-26 13:08:26 +02:00
integer :: exc(0:2,2,2), h1, h2, p1, p2, s1, s2
integer :: i, j, k, ex, n_minilist, iproc, degree
2018-04-10 14:25:28 +02:00
double precision :: haa, contrib, phase, c_alpha(N_states,Nproc), s_c_alpha(N_states)
2018-03-26 13:08:26 +02:00
logical :: ok
2018-03-23 11:06:52 +01:00
integer, external :: omp_get_thread_num
2018-03-26 13:08:26 +02:00
integer,intent(in) :: old_det_gen
integer(bit_kind), intent(in) :: old_generators(N_int, 2, old_det_gen)
allocate(minilist(N_det_delta_ij), det_minilist(N_int, 2, N_det_delta_ij), delta_ij_loc(N_states, N_det_delta_ij, 2, Nproc))
2018-04-10 14:25:28 +02:00
c_alpha = 0d0
2018-03-23 11:06:52 +01:00
delta_ij_loc = 0d0
2018-03-26 13:08:26 +02:00
!$OMP PARALLEL DO DEFAULT(SHARED) SCHEDULE(STATIC) PRIVATE(i, j, iproc, n_minilist, ex) &
2018-04-10 14:25:28 +02:00
!$OMP PRIVATE(det_minilist, minilist, haa, contrib, s_c_alpha) &
2018-03-26 13:08:26 +02:00
!$OMP PRIVATE(exc, h1, h2, p1, p2, s1, s2, phase, degree, ok)
do i=n_alpha,1,-1
2018-03-23 11:06:52 +01:00
iproc = omp_get_thread_num()+1
if(mod(i,10000) == 0) print *, "UNDRESSING", i, "/", n_alpha, iproc
n_minilist = 0
2018-03-26 13:08:26 +02:00
ok = .false.
do j=1, old_det_gen
call get_excitation_degree(alpha(1,1,i), old_generators(1,1,j), ex, N_int)
if(ex <= 2) then
call get_excitation(old_generators(1,1,j), alpha(1,1,i), exc,degree,phase,N_int)
call decode_exc(exc,degree,h1,p1,h2,p2,s1,s2)
ok = (mo_class(h1)(1:1) == 'A' .or. mo_class(h1)(1:1) == 'I') .and. &
(mo_class(p1)(1:1) == 'A' .or. mo_class(p1)(1:1) == 'V')
if(ok .and. degree == 2) then
ok = (mo_class(h2)(1:1) == 'A' .or. mo_class(h2)(1:1) == 'I') .and. &
(mo_class(p2)(1:1) == 'A' .or. mo_class(p2)(1:1) == 'V')
end if
if(ok) exit
end if
end do
if(.not. ok) cycle
do j=1, N_det_delta_ij
2018-03-23 11:06:52 +01:00
call get_excitation_degree(alpha(1,1,i), psi_det(1,1,j), ex, N_int)
if(ex <= 2) then
n_minilist += 1
det_minilist(:,:,n_minilist) = psi_det(:,:,j)
minilist(n_minilist) = j
end if
end do
2018-03-26 13:08:26 +02:00
call i_h_j(alpha(1,1,i), alpha(1,1,i), N_int, haa)
call dress_with_alpha_(N_states, N_det_delta_ij, N_int, delta_ij_loc(1,1,1,iproc), &
2018-04-10 14:25:28 +02:00
minilist, det_minilist, n_minilist, alpha(1,1,i), haa, contrib, s_c_alpha, iproc)
c_alpha(:,iproc) += s_c_alpha(:)**2
2018-03-23 11:06:52 +01:00
end do
!$OMP END PARALLEL DO
do i=2,Nproc
delta_ij_loc(:,:,:,1) += delta_ij_loc(:,:,:,i)
2018-04-10 14:25:28 +02:00
c_alpha(:,1) += c_alpha(:,i)
2018-03-23 11:06:52 +01:00
end do
2018-05-17 16:02:51 +02:00
delta_ij_tmp(:,:,1) -= delta_ij_loc(:,:,1,1)
delta_ij_tmp(:,:,2) -= delta_ij_loc(:,:,2,1)
2018-04-10 14:25:28 +02:00
2018-05-14 13:00:04 +02:00
!print *, "SUM ALPHA2 PRE", global_sum_alpha2
2018-04-10 14:25:28 +02:00
!global_sum_alpha2(:) -= c_alpha(:,1)
2018-05-17 16:02:51 +02:00
print *, "SUM C_ALPHA^2 =", global_sum_alpha2(:)
!print *, "*** DRESSINS DIVIDED BY 1+SUM C_ALPHA^2 ***"
!do i=1,N_states
! delta_ij_tmp(i,:,:) = delta_ij_tmp(i,:,:) / (1d0 + global_sum_alpha2(i))
!end do
2018-04-10 14:25:28 +02:00
global_sum_alpha2 = 0d0
2018-03-23 11:06:52 +01:00
end subroutine
2018-04-10 14:25:28 +02:00
subroutine dress_with_alpha_(Nstates,Ndet,Nint,delta_ij_loc,minilist, det_minilist, n_minilist, alpha, haa, contrib, c_alpha, iproc)
2018-03-23 11:06:52 +01:00
use bitmasks
implicit none
2018-03-26 12:06:18 +02:00
BEGIN_DOC
!delta_ij_loc(:,:,1) : dressing column for H
!delta_ij_loc(:,:,2) : dressing column for S2
2018-05-25 17:40:27 +02:00
!minilist : indices of determinants connected to alpha ( in psi_det )
2018-03-26 12:06:18 +02:00
!n_minilist : size of minilist
!alpha : alpha determinant
END_DOC
2018-03-23 11:06:52 +01:00
integer, intent(in) :: Nint, Ndet, Nstates, n_minilist, iproc
2018-03-26 12:06:18 +02:00
integer(bit_kind), intent(in) :: alpha(Nint,2), det_minilist(Nint, 2, n_minilist)
integer,intent(in) :: minilist(n_minilist)
double precision, intent(inout) :: delta_ij_loc(Nstates,N_det,2)
2018-05-09 15:03:31 +02:00
double precision, intent(out) :: contrib, c_alpha(N_states)
double precision,intent(in) :: haa
2018-03-23 11:06:52 +01:00
double precision :: hij, sij
2018-03-26 12:06:18 +02:00
double precision, external :: diag_H_mat_elem_fock
integer :: i,j,k,l,m, l_sd
double precision :: hdress, sdress
2018-05-09 15:03:31 +02:00
double precision :: de, a_h_psi(Nstates)
2018-03-26 12:06:18 +02:00
a_h_psi = 0d0
do l_sd=1,n_minilist
call i_h_j_s2(alpha,det_minilist(1,1,l_sd),N_int,hij, sij)
a_h_i(l_sd, iproc) = hij
a_s2_i(l_sd, iproc) = sij
do i=1,Nstates
a_h_psi(i) += hij * psi_coef(minilist(l_sd), i)
end do
end do
contrib = 0d0
2018-03-26 12:06:18 +02:00
do i=1,Nstates
2018-05-07 15:21:28 +02:00
de = dress_E0_denominator(i) - haa
2018-03-26 12:06:18 +02:00
if(DABS(de) < 1D-5) cycle
2018-04-10 14:25:28 +02:00
c_alpha(i) = a_h_psi(i) / de
contrib = min(contrib, c_alpha(i) * a_h_psi(i))
2018-03-26 12:06:18 +02:00
do l_sd=1,n_minilist
2018-04-10 14:25:28 +02:00
hdress = c_alpha(i) * a_h_i(l_sd, iproc)
sdress = c_alpha(i) * a_s2_i(l_sd, iproc)
2018-05-14 13:00:04 +02:00
!if(c_alpha(i) * a_s2_i(l_sd, iproc) > 1d-1) then
! call debug_det(det_minilist(1,1,l_sd), N_int)
! call debug_det(alpha,N_int)
!end if
2018-03-26 12:06:18 +02:00
delta_ij_loc(i, minilist(l_sd), 1) += hdress
delta_ij_loc(i, minilist(l_sd), 2) += sdress
end do
end do
end subroutine
2018-03-23 11:06:52 +01:00
subroutine dress_with_alpha_buffer(Nstates,Ndet,Nint,delta_ij_loc, i_gen, minilist, det_minilist, n_minilist, alpha, iproc)
use bitmasks
implicit none
BEGIN_DOC
!delta_ij_loc(:,:,1) : dressing column for H
!delta_ij_loc(:,:,2) : dressing column for S2
!i_gen : generator index in psi_det_generators
2018-05-25 17:40:27 +02:00
!minilist : indices of determinants connected to alpha ( in psi_det )
2018-03-23 11:06:52 +01:00
!n_minilist : size of minilist
!alpha : alpha determinant
END_DOC
integer, intent(in) :: Nint, Ndet, Nstates, n_minilist, iproc, i_gen
integer(bit_kind), intent(in) :: alpha(Nint,2), det_minilist(Nint, 2, n_minilist)
integer,intent(in) :: minilist(n_minilist)
double precision, intent(inout) :: delta_ij_loc(Nstates,N_det,2)
double precision, external :: diag_H_mat_elem_fock
2018-04-10 14:25:28 +02:00
double precision :: haa, contrib, c_alpha(N_states)
2018-03-23 11:06:52 +01:00
2018-05-02 14:32:41 +02:00
haa = diag_H_mat_elem_fock(psi_det_generators(1,1,i_gen),alpha,fock_diag_tmp_(1,1,iproc),N_int)
Squashed commit of the following: commit 96715abd7bc0645b994fc4fff1c764e7bec3042e Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Thu Sep 6 12:08:34 2018 +0200 Tasks commit e43b1e2faff5ad1ec4ecd2eb99e8d86e0000a9a2 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Thu Sep 6 11:47:02 2018 +0200 Fixed print commit c498c8944b5695b953909711a2e9a542e1bc6015 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Thu Sep 6 11:11:51 2018 +0200 PT2 and shiftedBk fixed commit 965cf0361d54df096c9bfca93f9d50ecd946c198 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 18:48:59 2018 +0200 Shifted Bk multistate broken commit 87ef641b65a122fa1256605cbe098c1a0de04bc0 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 17:23:38 2018 +0200 PT2 fixed commit a2adb533bcaf96191a24ff5fcef86cd14ac00697 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 16:55:05 2018 +0200 Working on PT2 (broken) commit 33f52991b65ac42ced5521ef0e713df765268860 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 12:13:23 2018 +0200 Fixed missing argument commit 712bf75f76421880299dc65b30acfda3d531f709 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 11:42:31 2018 +0200 Fixed floating invalid in PT2 commit cf2412ebd99f9acf573a5180603611f1c3e35155 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 11:34:37 2018 +0200 n_states_diag >= n_states commit bb415435e4d9be72a285ec789c3d63046a9173e4 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Wed Sep 5 11:23:41 2018 +0200 Fixed final print commit f2d339cf7b3535359bfc9d3057225b5957cf9e16 Merge: cfa8e1dc 52ca18c1 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 20:51:14 2018 +0200 Merge branch 'thesis' of ssh://github.com/garniron/quantum_package into thesis commit cfa8e1dc34c78f666aad677d2362a0b23e738e55 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 20:50:09 2018 +0200 restored relative_error commit bccad69c7785d79eb079bd15cd8d7143503ba9b5 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 20:47:11 2018 +0200 uninitialized variable commit 52ca18c1525cbe51f290813a4706928f227af142 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 20:47:11 2018 +0200 uninitialized variable commit ba0094f5f88153da295f479568bea9c6f494098c Merge: 093e3fd0 68458296 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 20:07:15 2018 +0200 merge with garniroy commit 68458296dcc86c57bab25582e138ae057a66f019 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Tue Sep 4 18:43:39 2018 +0200 Almost working but still broken commit 9ebb88cbf32eed90ed2c1ff985128f271bcc1c15 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Tue Sep 4 18:05:00 2018 +0200 Cleaning commit 873035e01635a1a0576cc0e8bf834beb3c60fc80 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Tue Sep 4 17:31:45 2018 +0200 Squashed commit of the following: commit 4b9c435dce0f3b3078d573e66fd32b40fca26497 Merge: 74e559c8 093e3fd0 Author: Anthony Scemama <scemama@irsamc.ups-tlse.fr> Date: Tue Sep 4 16:58:51 2018 +0200 Merge branch 'thesis' of git://github.com/garniron/quantum_package into garniron-thesis commit 093e3fd021ca5ed73fe3be530fa18e79ceb0085e Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 16:13:00 2018 +0200 removed ungodly hack commit 8529a0f3f66d2f1ba9d717b0de60b268a4a821a7 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:57:19 2018 +0200 reduced prints in pt2_stoch commit 03b8f353bd94b837c02e570a8d3fb68689003135 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:41:46 2018 +0200 teeth building check for pt2_stoch commit 0d91b9310a7d580e8bc97a6a0af7a62a8506f8c0 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:35:04 2018 +0200 timestamp of first pull commit 34d9fa01657ed5d159a1c62f8f24cdb266d48153 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:27:10 2018 +0200 potential numerical precision bug commit 9a0f900d8c57bcbca36a09c11616e234528e1643 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:09:51 2018 +0200 tests if teeth can be built commit dda0dc34df82a59fac3c3abb2b14943fdb23c0f8 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 17:48:04 2018 +0200 corrected pt2_find_sample commit a521f0cb82bd333ce9879b4dfd87ae1abcf2baaa Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 16:08:02 2018 +0200 tasks get by batches of Nproc commit 997a5a1265951394753d3154a6136ca40177ba8d Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 14:18:04 2018 +0200 buffered task_id send commit 99ea7948e0d58023f7ce02386014c8c3d493deb4 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 12:29:12 2018 +0200 unbalanced fragmentation commit abb3b7e08bebdcaf84208ebbe404441d3571b1f9 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sun Sep 2 17:18:44 2018 +0200 overflow of pt2_J commit 8df49f394b588ca304ce2ee7d8f8095455573784 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sun Sep 2 15:58:48 2018 +0200 removed useless computation of intermediate checkpoints commit 4ba5b79eb391819a4966c7cd19b2e2fdf2ed5129 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sun Sep 2 15:50:14 2018 +0200 dressing only sent for chosen checkpoint commit a4a6a69459321b1b210c2f51ffbef4f792943919 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sat Sep 1 17:01:56 2018 +0200 cumulative dot_F commit 6a7f04cb79ec921d5829a011aa8d8b2ea224760d Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sat Sep 1 16:58:07 2018 +0200 simpler purge commit 168ca2f2e29b902ab37afc58672b678998f0212a Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 21:07:01 2018 +0200 task list optimized commit de4a0d0caf0e1d441e0244a9ba7c65ae535d58b1 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 18:57:03 2018 +0200 removed print commit fee31d4e3e70c2edf1aeda202d3f067eddc17532 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 18:56:23 2018 +0200 dress fragmentation commit 02893a419de07c922708fe8ee2b828d9ebfd2e9b Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 15:52:16 2018 +0200 bug in blocked search - replaced with thesis version commit bb6e073cf10039f8a333d120b88dda109169528b Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 21:24:45 2018 +0200 ungodly hack to prevent double providing commit 0609e8c627ea7bf88d617e45dffeb64c9aea37b6 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 20:52:05 2018 +0200 debugging commit a254fdd7cffd8348846c910f240965104c969a09 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 15:24:07 2018 +0200 parallel bug commit 2a6c1941d45be78083ec8c1dccb453f51a12c338 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 11:43:11 2018 +0200 corrected when relative_error=0d0 commit bac039bdf1bd00673c0795b1cb23adb593e9d426 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 10:58:17 2018 +0200 relative error 1d-5 commit aae9d203ecbf8d53accc6bc35c2ab56aeeae78a6 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 10:07:02 2018 +0200 potential fragmentation bug commit ad69f39f99d0b0dd73f556fb13d1d55337c5b066 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Wed Aug 29 20:54:58 2018 +0200 dress_zmq re-implemented commit d78f64732a5493d7f10c7c80b564005e63a133fc Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Wed Aug 29 11:30:19 2018 +0200 pt2_stoch re-implemented commit 4b9b54e19ac7459589681e5ff7aa358dde9f5fd5 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Aug 28 10:24:38 2018 +0200 removed test for phase_mask_bit commit 3abccca5e35948e54a659cacccea42fbfcf4c296 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 3 23:44:05 2018 +0200 phasemask_bit commit 093e3fd021ca5ed73fe3be530fa18e79ceb0085e Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 16:13:00 2018 +0200 removed ungodly hack commit 8529a0f3f66d2f1ba9d717b0de60b268a4a821a7 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:57:19 2018 +0200 reduced prints in pt2_stoch commit 03b8f353bd94b837c02e570a8d3fb68689003135 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:41:46 2018 +0200 teeth building check for pt2_stoch commit 0d91b9310a7d580e8bc97a6a0af7a62a8506f8c0 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:35:04 2018 +0200 timestamp of first pull commit 34d9fa01657ed5d159a1c62f8f24cdb266d48153 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:27:10 2018 +0200 potential numerical precision bug commit 9a0f900d8c57bcbca36a09c11616e234528e1643 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Sep 4 14:09:51 2018 +0200 tests if teeth can be built commit dda0dc34df82a59fac3c3abb2b14943fdb23c0f8 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 17:48:04 2018 +0200 corrected pt2_find_sample commit a521f0cb82bd333ce9879b4dfd87ae1abcf2baaa Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 16:08:02 2018 +0200 tasks get by batches of Nproc commit 997a5a1265951394753d3154a6136ca40177ba8d Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 14:18:04 2018 +0200 buffered task_id send commit 99ea7948e0d58023f7ce02386014c8c3d493deb4 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Mon Sep 3 12:29:12 2018 +0200 unbalanced fragmentation commit abb3b7e08bebdcaf84208ebbe404441d3571b1f9 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sun Sep 2 17:18:44 2018 +0200 overflow of pt2_J commit 8df49f394b588ca304ce2ee7d8f8095455573784 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sun Sep 2 15:58:48 2018 +0200 removed useless computation of intermediate checkpoints commit 4ba5b79eb391819a4966c7cd19b2e2fdf2ed5129 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sun Sep 2 15:50:14 2018 +0200 dressing only sent for chosen checkpoint commit a4a6a69459321b1b210c2f51ffbef4f792943919 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sat Sep 1 17:01:56 2018 +0200 cumulative dot_F commit 6a7f04cb79ec921d5829a011aa8d8b2ea224760d Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Sat Sep 1 16:58:07 2018 +0200 simpler purge commit 168ca2f2e29b902ab37afc58672b678998f0212a Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 21:07:01 2018 +0200 task list optimized commit de4a0d0caf0e1d441e0244a9ba7c65ae535d58b1 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 18:57:03 2018 +0200 removed print commit fee31d4e3e70c2edf1aeda202d3f067eddc17532 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 18:56:23 2018 +0200 dress fragmentation commit 02893a419de07c922708fe8ee2b828d9ebfd2e9b Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 31 15:52:16 2018 +0200 bug in blocked search - replaced with thesis version commit bb6e073cf10039f8a333d120b88dda109169528b Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 21:24:45 2018 +0200 ungodly hack to prevent double providing commit 0609e8c627ea7bf88d617e45dffeb64c9aea37b6 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 20:52:05 2018 +0200 debugging commit a254fdd7cffd8348846c910f240965104c969a09 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 15:24:07 2018 +0200 parallel bug commit 2a6c1941d45be78083ec8c1dccb453f51a12c338 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 11:43:11 2018 +0200 corrected when relative_error=0d0 commit bac039bdf1bd00673c0795b1cb23adb593e9d426 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 10:58:17 2018 +0200 relative error 1d-5 commit aae9d203ecbf8d53accc6bc35c2ab56aeeae78a6 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Thu Aug 30 10:07:02 2018 +0200 potential fragmentation bug commit ad69f39f99d0b0dd73f556fb13d1d55337c5b066 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Wed Aug 29 20:54:58 2018 +0200 dress_zmq re-implemented commit d78f64732a5493d7f10c7c80b564005e63a133fc Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Wed Aug 29 11:30:19 2018 +0200 pt2_stoch re-implemented commit 4b9b54e19ac7459589681e5ff7aa358dde9f5fd5 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Tue Aug 28 10:24:38 2018 +0200 removed test for phase_mask_bit commit 3abccca5e35948e54a659cacccea42fbfcf4c296 Author: Yann Garniron <yann.garniron@yahoo.fr> Date: Fri Aug 3 23:44:05 2018 +0200 phasemask_bit
2018-09-06 15:19:51 +02:00
2018-04-10 14:25:28 +02:00
call dress_with_alpha_(Nstates, Ndet, Nint, delta_ij_loc, minilist, det_minilist, n_minilist, alpha, haa, contrib, c_alpha, iproc)
slave_sum_alpha2(:,iproc) += c_alpha(:)**2
2018-04-04 11:32:27 +02:00
if(contrib < sb(iproc)%mini) then
call add_to_selection_buffer(sb(iproc), alpha, contrib)
end if
2018-03-21 12:06:26 +01:00
end subroutine
2018-03-26 12:06:18 +02:00
BEGIN_PROVIDER [ logical, initialize_E0_denominator ]
implicit none
BEGIN_DOC
! If true, initialize pt2_E0_denominator
END_DOC
initialize_E0_denominator = .True.
END_PROVIDER