10
0
mirror of https://github.com/LCPQ/quantum_package synced 2025-01-03 10:05:57 +01:00

Fixing travis

This commit is contained in:
Anthony Scemama 2017-06-07 21:56:46 +02:00
parent c49720fdae
commit df92bd52e4
4 changed files with 88 additions and 97 deletions

View File

@ -471,6 +471,7 @@ end subroutine
end if end if
norm_left -= pt2_weight(i) norm_left -= pt2_weight(i)
end do end do
first_det_of_comb = max(1,first_det_of_comb)
call write_int(6, first_det_of_comb-1, 'Size of deterministic set') call write_int(6, first_det_of_comb-1, 'Size of deterministic set')
comb_step = (1d0 - pt2_cweight(first_det_of_comb-1)) * comb_step comb_step = (1d0 - pt2_cweight(first_det_of_comb-1)) * comb_step

View File

@ -2,7 +2,7 @@
type: Threshold type: Threshold
doc: Threshold on the magnitude of the smallest eigenvalues of the overlap matrix in the AO basis doc: Threshold on the magnitude of the smallest eigenvalues of the overlap matrix in the AO basis
interface: ezfio,provider,ocaml interface: ezfio,provider,ocaml
default: 1.e-6 default: 1.e-5
[max_dim_diis] [max_dim_diis]
type: integer type: integer

View File

@ -10,100 +10,94 @@
integer, allocatable :: iwork(:) integer, allocatable :: iwork(:)
double precision, allocatable :: work(:), F(:,:), S(:,:) double precision, allocatable :: work(:), F(:,:), S(:,:)
allocate( F(mo_tot_num_align,mo_tot_num) ) allocate( F(mo_tot_num,mo_tot_num) )
do j=1,mo_tot_num do j=1,mo_tot_num
do i=1,mo_tot_num do i=1,mo_tot_num
F(i,j) = Fock_matrix_mo(i,j) F(i,j) = Fock_matrix_mo(i,j)
enddo enddo
enddo
if(no_oa_or_av_opt)then
integer :: iorb,jorb
do i = 1, n_act_orb
iorb = list_act(i)
do j = 1, n_inact_orb
jorb = list_inact(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo enddo
if(no_oa_or_av_opt)then do j = 1, n_virt_orb
integer :: iorb,jorb jorb = list_virt(j)
do i = 1, n_act_orb F(iorb,jorb) = 0.d0
iorb = list_act(i) F(jorb,iorb) = 0.d0
do j = 1, n_inact_orb
jorb = list_inact(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo
do j = 1, n_virt_orb
jorb = list_virt(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo
do j = 1, n_core_orb
jorb = list_core(j)
F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo
enddo
endif
! Insert level shift here
do i = elec_beta_num+1, elec_alpha_num
F(i,i) += 0.5d0*level_shift
enddo enddo
do j = 1, n_core_orb
do i = elec_alpha_num+1, mo_tot_num jorb = list_core(j)
F(i,i) += level_shift F(iorb,jorb) = 0.d0
F(jorb,iorb) = 0.d0
enddo enddo
enddo
n = mo_tot_num endif
lwork = 1+6*n + 2*n*n
liwork = 3 + 5*n
allocate(work(lwork))
allocate(iwork(liwork) ) ! Insert level shift here
do i = elec_beta_num+1, elec_alpha_num
lwork = -1 F(i,i) += 0.5d0*level_shift
liwork = -1 enddo
call dsyevd( 'V', 'U', mo_tot_num, F, & do i = elec_alpha_num+1, mo_tot_num
size(F,1), diagonal_Fock_matrix_mo, & F(i,i) += level_shift
work, lwork, iwork, liwork, info) enddo
if (info /= 0) then n = mo_tot_num
print *, irp_here//' DSYEVD failed : ', info lwork = 1+6*n + 2*n*n
stop 1 liwork = 3 + 5*n
endif
lwork = int(work(1)) allocate(work(lwork))
liwork = iwork(1) allocate(iwork(liwork) )
deallocate(iwork)
deallocate(work) lwork = -1
liwork = -1
allocate(work(lwork))
allocate(iwork(liwork) ) call dsyevd( 'V', 'U', mo_tot_num, F, &
call dsyevd( 'V', 'U', mo_tot_num, F, & size(F,1), diagonal_Fock_matrix_mo, &
size(F,1), diagonal_Fock_matrix_mo, & work, lwork, iwork, liwork, info)
work, lwork, iwork, liwork, info)
deallocate(iwork) if (info /= 0) then
print *, irp_here//' DSYEVD failed : ', info
stop 1
if (info /= 0) then endif
info = 0 lwork = int(work(1))
do j=1,mo_tot_num liwork = iwork(1)
do i=1,mo_tot_num deallocate(iwork)
F(i,j) = Fock_matrix_mo(i,j) deallocate(work)
enddo
enddo allocate(work(lwork))
call dsyev( 'V', 'L', mo_tot_num, F, & allocate(iwork(liwork) )
size(F,1), diagonal_Fock_matrix_mo, & call dsyevd( 'V', 'U', mo_tot_num, F, &
work, lwork, info) size(F,1), diagonal_Fock_matrix_mo, &
work, lwork, iwork, liwork, info)
if (info /= 0) then deallocate(iwork)
print *, irp_here//' DSYEV failed : ', info
stop 1
endif if (info /= 0) then
endif call dsyev( 'V', 'L', mo_tot_num, F, &
size(F,1), diagonal_Fock_matrix_mo, &
call dgemm('N','N',ao_num,mo_tot_num,mo_tot_num, 1.d0, & work, lwork, info)
mo_coef, size(mo_coef,1), F, size(F,1), &
0.d0, eigenvectors_Fock_matrix_mo, size(eigenvectors_Fock_matrix_mo,1)) if (info /= 0) then
deallocate(work, F) print *, irp_here//' DSYEV failed : ', info
stop 1
endif
endif
call dgemm('N','N',ao_num,mo_tot_num,mo_tot_num, 1.d0, &
mo_coef, size(mo_coef,1), F, size(F,1), &
0.d0, eigenvectors_Fock_matrix_mo, size(eigenvectors_Fock_matrix_mo,1))
deallocate(work, F)
END_PROVIDER END_PROVIDER

View File

@ -110,10 +110,6 @@ subroutine mrcc_part_dress(delta_ij_, delta_ii_,delta_ij_s2_, delta_ii_s2_,i_gen
allocate(idx_alpha(0:psi_det_size), degree_alpha(psi_det_size)) allocate(idx_alpha(0:psi_det_size), degree_alpha(psi_det_size))
call create_minilist_find_previous(key_mask, psi_det_generators, miniList, i_generator-1, N_miniList, fullMatch, Nint) call create_minilist_find_previous(key_mask, psi_det_generators, miniList, i_generator-1, N_miniList, fullMatch, Nint)
! if(fullMatch) then
! return
! end if
allocate(ptr_microlist(0:mo_tot_num*2+1), & allocate(ptr_microlist(0:mo_tot_num*2+1), &
N_microlist(0:mo_tot_num*2) ) N_microlist(0:mo_tot_num*2) )
allocate( microlist(Nint,2,N_minilist*4), & allocate( microlist(Nint,2,N_minilist*4), &
@ -141,7 +137,7 @@ subroutine mrcc_part_dress(delta_ij_, delta_ii_,delta_ij_s2_, delta_ii_s2_,i_gen
if(N_minilist == 0) return if(N_minilist == 0) return
if(key_mask(1,1) /= 0) then !!!!!!!!!!! PAS GENERAL !!!!!!!!! if(sum(abs(key_mask(1:N_int,1))) /= 0)
allocate(microlist_zero(Nint,2,N_minilist), idx_microlist_zero(N_minilist)) allocate(microlist_zero(Nint,2,N_minilist), idx_microlist_zero(N_minilist))
allocate( microlist(Nint,2,N_minilist*4), & allocate( microlist(Nint,2,N_minilist*4), &