10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-10-09 17:43:18 +02:00
quantum_package/src/MO_Basis/utils.irp.f

286 lines
7.7 KiB
Fortran
Raw Normal View History

2014-04-03 01:50:22 +02:00
subroutine save_mos
implicit none
double precision, allocatable :: buffer(:,:)
integer :: i,j
2015-06-08 14:49:10 +02:00
call system('$QP_ROOT/scripts/save_current_mos.sh '//trim(ezfio_filename))
2014-04-03 01:50:22 +02:00
2014-09-18 17:01:43 +02:00
call ezfio_set_mo_basis_mo_tot_num(mo_tot_num)
2014-04-03 01:50:22 +02:00
call ezfio_set_mo_basis_mo_label(mo_label)
2014-12-25 23:53:29 +01:00
call ezfio_set_mo_basis_ao_md5(ao_md5)
2014-04-03 01:50:22 +02:00
allocate ( buffer(ao_num,mo_tot_num) )
buffer = 0.d0
do j = 1, mo_tot_num
do i = 1, ao_num
buffer(i,j) = mo_coef(i,j)
enddo
enddo
call ezfio_set_mo_basis_mo_coef(buffer)
call ezfio_set_mo_basis_mo_occ(mo_occ)
2014-04-03 01:50:22 +02:00
deallocate (buffer)
end
subroutine save_mos_truncated(n)
implicit none
double precision, allocatable :: buffer(:,:)
integer :: i,j,n
call system('$QP_ROOT/scripts/save_current_mos.sh '//trim(ezfio_filename))
call ezfio_set_mo_basis_mo_tot_num(n)
call ezfio_set_mo_basis_mo_label(mo_label)
call ezfio_set_mo_basis_ao_md5(ao_md5)
allocate ( buffer(ao_num,n) )
buffer = 0.d0
do j = 1, n
do i = 1, ao_num
buffer(i,j) = mo_coef(i,j)
enddo
enddo
call ezfio_set_mo_basis_mo_coef(buffer)
call ezfio_set_mo_basis_mo_occ(mo_occ)
deallocate (buffer)
end
subroutine mo_as_eigvectors_of_mo_matrix(matrix,n,m,label,sign)
2014-04-03 01:50:22 +02:00
implicit none
integer,intent(in) :: n,m, sign
2014-04-03 01:50:22 +02:00
character*(64), intent(in) :: label
double precision, intent(in) :: matrix(n,m)
integer :: i,j
double precision, allocatable :: mo_coef_new(:,:), R(:,:),eigvalues(:), A(:,:)
2014-04-03 01:50:22 +02:00
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: mo_coef_new, R
2015-06-24 09:37:38 +02:00
call write_time(output_mo_basis)
2014-04-03 01:50:22 +02:00
if (m /= mo_tot_num) then
print *, irp_here, ': Error : m/= mo_tot_num'
2014-06-20 18:35:26 +02:00
stop 1
2014-04-03 01:50:22 +02:00
endif
allocate(A(n,m),R(n,m),mo_coef_new(ao_num_align,m),eigvalues(m))
if (sign == -1) then
do j=1,m
do i=1,n
A(i,j) = -matrix(i,j)
enddo
enddo
else
do j=1,m
do i=1,n
A(i,j) = matrix(i,j)
enddo
enddo
endif
2014-04-03 01:50:22 +02:00
mo_coef_new = mo_coef
call lapack_diag(eigvalues,R,A,n,m)
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') 'MOs are now **'//trim(label)//'**'
write (output_mo_basis,'(A)') ''
write (output_mo_basis,'(A)') 'Eigenvalues'
write (output_mo_basis,'(A)') '-----------'
write (output_mo_basis,'(A)') ''
write (output_mo_basis,'(A)') '======== ================'
if (sign == -1) then
do i=1,m
eigvalues(i) = -eigvalues(i)
enddo
endif
do i=1,m
2017-04-12 19:29:21 +02:00
write (output_mo_basis,'(I8,1X,F16.10)') i,eigvalues(i)
2014-04-03 01:50:22 +02:00
enddo
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') '======== ================'
write (output_mo_basis,'(A)') ''
2014-04-03 01:50:22 +02:00
call dgemm('N','N',ao_num,m,m,1.d0,mo_coef_new,size(mo_coef_new,1),R,size(R,1),0.d0,mo_coef,size(mo_coef,1))
deallocate(A,mo_coef_new,R,eigvalues)
2015-06-24 09:37:38 +02:00
call write_time(output_mo_basis)
2014-04-03 01:50:22 +02:00
mo_label = label
end
2015-12-08 13:24:43 +01:00
subroutine mo_as_svd_vectors_of_mo_matrix(matrix,lda,m,n,label)
implicit none
integer,intent(in) :: lda,m,n
character*(64), intent(in) :: label
double precision, intent(in) :: matrix(lda,n)
integer :: i,j
double precision, allocatable :: mo_coef_new(:,:), U(:,:),D(:), A(:,:), Vt(:,:), work(:)
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: mo_coef_new, U, Vt, A
call write_time(output_mo_basis)
if (m /= mo_tot_num) then
print *, irp_here, ': Error : m/= mo_tot_num'
stop 1
endif
allocate(A(lda,n),U(lda,n),mo_coef_new(ao_num_align,m),D(m),Vt(lda,n))
do j=1,n
do i=1,m
A(i,j) = matrix(i,j)
enddo
enddo
mo_coef_new = mo_coef
call svd(A,lda,U,lda,D,Vt,lda,m,n)
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') 'MOs are now **'//trim(label)//'**'
write (output_mo_basis,'(A)') ''
write (output_mo_basis,'(A)') 'Eigenvalues'
write (output_mo_basis,'(A)') '-----------'
write (output_mo_basis,'(A)') ''
write (output_mo_basis,'(A)') '======== ================'
2015-12-08 13:24:43 +01:00
do i=1,m
2017-04-12 19:29:21 +02:00
write (output_mo_basis,'(I8,1X,F16.10)') i,D(i)
2015-12-08 13:24:43 +01:00
enddo
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') '======== ================'
write (output_mo_basis,'(A)') ''
2015-12-08 13:24:43 +01:00
call dgemm('N','N',ao_num,m,m,1.d0,mo_coef_new,size(mo_coef_new,1),U,size(U,1),0.d0,mo_coef,size(mo_coef,1))
deallocate(A,mo_coef_new,U,Vt,D)
call write_time(output_mo_basis)
mo_label = label
end
2015-01-07 17:59:31 +01:00
subroutine mo_as_eigvectors_of_mo_matrix_sort_by_observable(matrix,observable,n,m,label)
implicit none
integer,intent(in) :: n,m
character*(64), intent(in) :: label
double precision, intent(in) :: matrix(n,m),observable(n,n)
double precision, allocatable :: mo_coef_new(:,:), R(:,:),eigvalues(:),value(:)
integer,allocatable :: iorder(:)
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: mo_coef_new, R
2015-06-24 09:37:38 +02:00
call write_time(output_mo_basis)
2015-01-07 17:59:31 +01:00
if (m /= mo_tot_num) then
print *, irp_here, ': Error : m/= mo_tot_num'
stop 1
endif
allocate(R(n,m),mo_coef_new(ao_num_align,m),eigvalues(m),value(m),iorder(m))
mo_coef_new = mo_coef
call lapack_diag(eigvalues,R,matrix,size(matrix,1),size(matrix,2))
do i = 1, mo_tot_num
value(i) = 0.d0
iorder(i) = i
do j = 1, mo_tot_num
do k = 1, mo_tot_num
value(i) += R(k,i) * R(j,i) * observable(k,j)
enddo
enddo
! print*,'value(i) = ',i,value(i)
enddo
integer :: i,j,k,index
double precision :: R_tmp(m,m),obs(m)
print*,'sort ....'
call dsort(value,iorder,n)
do i = 1, mo_tot_num
index = iorder(i)
! print*,'iorder(i) = ',iorder(i)
obs(i) = eigvalues(iorder(i))
do j = 1, mo_tot_num
R_tmp(j,i) = R(j,index)
enddo
enddo
do i = 1, mo_tot_num
do j = 1, mo_tot_num
R(j,i) = R_tmp(j,i)
enddo
enddo
do i = 1, mo_tot_num
value(i) = 0.d0
do j = 1, mo_tot_num
do k = 1, mo_tot_num
value(i) += R(k,i) * R(j,i) * observable(k,j)
enddo
enddo
print*,'i = ',i
print*,'value(i) = ',value(i)
print*,'obs(i) = ',obs(i)
print*,''
enddo
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') 'MOs are now **'//trim(label)//'**'
write (output_mo_basis,'(A)') ''
write (output_mo_basis,'(A)') 'Eigenvalues'
write (output_mo_basis,'(A)') '-----------'
write (output_mo_basis,'(A)') ''
write (output_mo_basis,'(A)') '======== ================'
2015-01-07 17:59:31 +01:00
do i = 1, m
2017-04-12 19:29:21 +02:00
write (output_mo_basis,'(I8,1X,F16.10)') i,eigvalues(i)
2015-01-07 17:59:31 +01:00
enddo
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') '======== ================'
write (output_mo_basis,'(A)') ''
2015-01-07 17:59:31 +01:00
call dgemm('N','N',ao_num,m,m,1.d0,mo_coef_new,size(mo_coef_new,1),R,size(R,1),0.d0,mo_coef,size(mo_coef,1))
deallocate(mo_coef_new,R,eigvalues)
2015-06-24 09:37:38 +02:00
call write_time(output_mo_basis)
2015-01-07 17:59:31 +01:00
mo_label = label
SOFT_TOUCH mo_coef mo_label
end
2015-03-19 21:14:52 +01:00
subroutine mo_sort_by_observable(observable,label)
implicit none
character*(64), intent(in) :: label
double precision, intent(in) :: observable(mo_tot_num)
double precision, allocatable :: mo_coef_new(:,:),value(:)
integer,allocatable :: iorder(:)
allocate(mo_coef_new(ao_num_align,mo_tot_num),value(mo_tot_num),iorder(mo_tot_num))
print*,'allocate !'
mo_coef_new = mo_coef
do i = 1, mo_tot_num
iorder(i) = i
value(i) = observable(i)
enddo
integer :: i,j,k,index
print*,'sort ....'
call dsort(value,iorder,mo_tot_num)
do i = 1, mo_tot_num
index = iorder(i)
do j = 1, mo_tot_num
mo_coef(j,i) = mo_coef_new(j,index)
enddo
enddo
2016-02-19 00:20:28 +01:00
write (output_mo_basis,'(A)') 'MOs are now **'//trim(label)//'**'
write (output_mo_basis,'(A)') ''
2015-03-19 21:14:52 +01:00
deallocate(mo_coef_new,value)
2015-06-24 09:37:38 +02:00
! call write_time(output_mo_basis)
2015-03-19 21:14:52 +01:00
mo_label = label
SOFT_TOUCH mo_coef mo_label
end
2016-04-17 22:25:25 +02:00
subroutine give_all_mos_at_r(r,mos_array)
implicit none
double precision, intent(in) :: r(3)
double precision, intent(out) :: mos_array(mo_tot_num)
double precision :: aos_array(ao_num),accu
integer :: i,j
call give_all_aos_at_r(r,aos_array)
do i = 1, mo_tot_num
accu = 0.d0
do j = 1, ao_num
Bugs to fix (#50) * Add config for knl * Add mising readme * Add .gitignore * Add pseudo to qp_convert * Working pseudo * Dressed matrix for pt2 works for one state * now eigenfunction of S^2 * minor modifs in printing * Fixed the perturbation with psi_ref instead of psi_det * Trying do really fo sin free multiple excitations * Beginning to merge MRCC and MRPT * final version of MRPT, at least I hope * Fix 404: Update Zlib Url. * Delete ifort_knl.cfg * Update module_handler.py * Update pot_ao_pseudo_ints.irp.f * Update map_module.f90 * Restaure map_module.f90 * Update configure * Update configure * Update sort.irp.f * Update sort.irp.f * Update selection.irp.f * Update selection.irp.f * Update dressing.irp.f * TApplencourt IRPF90 -> LCPQ * Remove `irpf90.make` in dependency * Update configure * Missing PROVIDE * Missing PROVIDE * Missing PROVIDE * Missing PROVIDE * Update configure * pouet * density based mrpt2 * debugging FOBOCI * Added SCF_density * New version of FOBOCI * added density.irp.f * minor changes in plugins/FOBOCI/SC2_1h1p.irp.f * added track_orb.irp.f * minor changes * minor modifs in FOBOCI * med * Minor changes * minor changes * strange things in MRPT * minor modifs mend * Fix #185 (Graphviz API / Python 2.6) * beginning to debug dft * fixed the factor 2 in lebedev * DFT integration works for non overlapping densities * DFT begins to work with lda * KS LDA is okay * added core integrals * mend * Beginning logn range integrals * Trying to handle two sets of integrals * beginning to clean erf integrals * Handling of two different mo and ao integrals map
2017-04-20 08:36:11 +02:00
accu += mo_coef(j,i) * aos_array(j)
2016-04-17 22:25:25 +02:00
enddo
mos_array(i) = accu
enddo
end