10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-06-28 08:02:33 +02:00

restructured mo_coef_complex provider; added mo_coef_real; maybe need to change ocaml?

This commit is contained in:
Kevin Gasperich 2020-01-29 16:23:00 -06:00
parent b0d27f8503
commit cc840cdbc1
4 changed files with 151 additions and 74 deletions

View File

@ -9,6 +9,12 @@ doc: Coefficient of the i-th |AO| on the j-th |MO|
interface: ezfio interface: ezfio
size: (ao_basis.ao_num,mo_basis.mo_num) size: (ao_basis.ao_num,mo_basis.mo_num)
[mo_coef_real]
type: double precision
doc: Imaginary part of the MO coefficient of the i-th |AO| on the j-th |MO|
interface: ezfio
size: (ao_basis.ao_num,mo_basis.mo_num)
[mo_coef_imag] [mo_coef_imag]
type: double precision type: double precision
doc: Imaginary part of the MO coefficient of the i-th |AO| on the j-th |MO| doc: Imaginary part of the MO coefficient of the i-th |AO| on the j-th |MO|

View File

@ -1,4 +1,4 @@
BEGIN_PROVIDER [ double precision, mo_coef_imag, (ao_num,mo_num) ] BEGIN_PROVIDER [ complex*16, mo_coef_complex, (ao_num,mo_num) ]
implicit none implicit none
BEGIN_DOC BEGIN_DOC
! Molecular orbital coefficients on |AO| basis set ! Molecular orbital coefficients on |AO| basis set
@ -8,14 +8,16 @@ BEGIN_PROVIDER [ double precision, mo_coef_imag, (ao_num,mo_num) ]
! mo_label : Label characterizing the |MOs| (local, canonical, natural, etc) ! mo_label : Label characterizing the |MOs| (local, canonical, natural, etc)
END_DOC END_DOC
integer :: i, j integer :: i, j
double precision, allocatable :: buffer(:,:) double precision, allocatable :: buffer_re(:,:),buffer_im(:,:)
logical :: exists logical :: exists_re,exists_im,exists
PROVIDE ezfio_filename PROVIDE ezfio_filename
if (mpi_master) then if (mpi_master) then
! Coefs ! Coefs
call ezfio_has_mo_basis_mo_coef_imag(exists) call ezfio_has_mo_basis_mo_coef_real(exists_re)
call ezfio_has_mo_basis_mo_coef_imag(exists_im)
exists = (exists_re.and.exists_im)
endif endif
IRP_IF MPI_DEBUG IRP_IF MPI_DEBUG
print *, irp_here, mpi_rank print *, irp_here, mpi_rank
@ -26,51 +28,106 @@ BEGIN_PROVIDER [ double precision, mo_coef_imag, (ao_num,mo_num) ]
integer :: ierr integer :: ierr
call MPI_BCAST(exists, 1, MPI_LOGICAL, 0, MPI_COMM_WORLD, ierr) call MPI_BCAST(exists, 1, MPI_LOGICAL, 0, MPI_COMM_WORLD, ierr)
if (ierr /= MPI_SUCCESS) then if (ierr /= MPI_SUCCESS) then
stop 'Unable to read mo_coef_imag with MPI' stop 'Unable to read mo_coef_real/imag with MPI'
endif endif
IRP_ENDIF IRP_ENDIF
if (exists) then if (exists) then
if (mpi_master) then if (mpi_master) then
call ezfio_get_mo_basis_mo_coef_imag(mo_coef_imag) allocate(buffer_re(ao_num,mo_num),buffer_im(ao_num,mo_num))
write(*,*) 'Read mo_coef_imag' call ezfio_get_mo_basis_mo_coef_real(buffer_re)
call ezfio_get_mo_basis_mo_coef_imag(buffer_im)
write(*,*) 'Read mo_coef_real/imag'
do i=1,mo_num
do j=1,ao_num
mo_coef_complex(j,i) = dcmplx(buffer_re(j,i),buffer_im(j,i))
enddo
enddo
deallocate(buffer_re,buffer_im)
endif endif
IRP_IF MPI IRP_IF MPI
call MPI_BCAST( mo_coef_imag, mo_num*ao_num, MPI_DOUBLE_PRECISION, 0, MPI_COMM_WORLD, ierr) call MPI_BCAST( mo_coef_complex, mo_num*ao_num, MPI_DOUBLE_COMPLEX, 0, MPI_COMM_WORLD, ierr)
if (ierr /= MPI_SUCCESS) then if (ierr /= MPI_SUCCESS) then
stop 'Unable to read mo_coef_imag with MPI' stop 'Unable to read mo_coef_real with MPI'
endif endif
IRP_ENDIF IRP_ENDIF
else else
! Orthonormalized AO basis ! Orthonormalized AO basis
do i=1,mo_num do i=1,mo_num
do j=1,ao_num do j=1,ao_num
mo_coef_imag(j,i) = 0.d0 mo_coef_complex(j,i) = ao_ortho_canonical_coef_complex(j,i)
enddo enddo
enddo enddo
endif endif
END_PROVIDER END_PROVIDER
BEGIN_PROVIDER [ complex*16, mo_coef_complex, (ao_num,mo_num) ] ! BEGIN_PROVIDER [ double precision, mo_coef_real, (ao_num,mo_num) ]
implicit none !&BEGIN_PROVIDER [ double precision, mo_coef_imag, (ao_num,mo_num) ]
BEGIN_DOC !&BEGIN_PROVIDER [ complex*16, mo_coef_complex, (ao_num,mo_num) ]
! Molecular orbital coefficients on |AO| basis set ! implicit none
! ! BEGIN_DOC
! mo_coef_complex(i,j) = coefficient of the i-th |AO| on the jth |MO| ! ! Molecular orbital coefficients on |AO| basis set
! ! !
! mo_label : Label characterizing the |MOs| (local, canonical, natural, etc) ! ! mo_coef_imag(i,j) = coefficient of the i-th |AO| on the jth |MO|
END_DOC ! !
integer :: i, j ! ! mo_label : Label characterizing the |MOs| (local, canonical, natural, etc)
PROVIDE ezfio_filename ! END_DOC
! integer :: i, j
! double precision, allocatable :: buffer_re(:,:),buffer_im(:,:)
! logical :: exists_re,exists_im,exists
! PROVIDE ezfio_filename
!
!
! if (mpi_master) then
! ! Coefs
! call ezfio_has_mo_basis_mo_coef_real(exists_re)
! call ezfio_has_mo_basis_mo_coef_imag(exists_im)
! exists = (exists_re.and.exists_im)
! endif
! IRP_IF MPI_DEBUG
! print *, irp_here, mpi_rank
! call MPI_BARRIER(MPI_COMM_WORLD, ierr)
! IRP_ENDIF
! IRP_IF MPI
! include 'mpif.h'
! integer :: ierr
! call MPI_BCAST(exists, 1, MPI_LOGICAL, 0, MPI_COMM_WORLD, ierr)
! if (ierr /= MPI_SUCCESS) then
! stop 'Unable to read mo_coef_real/imag with MPI'
! endif
! IRP_ENDIF
!
! if (exists) then
! if (mpi_master) then
! call ezfio_get_mo_basis_mo_coef_real(mo_coef_real)
! call ezfio_get_mo_basis_mo_coef_imag(mo_coef_imag)
! write(*,*) 'Read mo_coef_real/imag'
! endif
! IRP_IF MPI
! call MPI_BCAST( mo_coef_real, mo_num*ao_num, MPI_DOUBLE_PRECISION, 0, MPI_COMM_WORLD, ierr)
! if (ierr /= MPI_SUCCESS) then
! stop 'Unable to read mo_coef_real with MPI'
! endif
! call MPI_BCAST( mo_coef_imag, mo_num*ao_num, MPI_DOUBLE_PRECISION, 0, MPI_COMM_WORLD, ierr)
! if (ierr /= MPI_SUCCESS) then
! stop 'Unable to read mo_coef_imag with MPI'
! endif
! IRP_ENDIF
! do i=1,mo_num
! do j=1,ao_num
! mo_coef_complex(j,i) = dcmplx(mo_coef_real(j,i),mo_coef_imag(j,i))
! enddo
! enddo
! else
! ! Orthonormalized AO basis
! do i=1,mo_num
! do j=1,ao_num
! mo_coef_complex(j,i) = ao_ortho_canonical_coef_complex(j,i)
! enddo
! enddo
! endif
!END_PROVIDER
provide mo_coef mo_coef_imag
do i=1,mo_num
do j=1,ao_num
mo_coef_complex(j,i) = dcmplx(mo_coef(j,i),mo_coef_imag(j,i))
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [ complex*16, mo_coef_in_ao_ortho_basis_complex, (ao_num, mo_num) ] BEGIN_PROVIDER [ complex*16, mo_coef_in_ao_ortho_basis_complex, (ao_num, mo_num) ]
implicit none implicit none

View File

@ -1,6 +1,6 @@
subroutine save_mos subroutine save_mos
implicit none implicit none
double precision, allocatable :: buffer(:,:) double precision, allocatable :: buffer(:,:),buffer_im(:,:)
integer :: i,j integer :: i,j
!TODO: change this for periodic? !TODO: change this for periodic?
! save real/imag parts of mo_coef_complex ! save real/imag parts of mo_coef_complex
@ -10,62 +10,75 @@ subroutine save_mos
call ezfio_set_mo_basis_mo_num(mo_num) call ezfio_set_mo_basis_mo_num(mo_num)
call ezfio_set_mo_basis_mo_label(mo_label) call ezfio_set_mo_basis_mo_label(mo_label)
call ezfio_set_mo_basis_ao_md5(ao_md5) call ezfio_set_mo_basis_ao_md5(ao_md5)
allocate ( buffer(ao_num,mo_num) )
buffer = 0.d0
do j = 1, mo_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)
call ezfio_set_mo_basis_mo_class(mo_class)
if (is_periodic) then if (is_periodic) then
allocate ( buffer(ao_num,mo_num),buffer_im(ao_num,mo_num))
buffer = 0.d0
buffer_im = 0.d0
do j = 1, mo_num
do i = 1, ao_num
buffer(i,j) = dble(mo_coef_complex(i,j))
buffer_im(i,j) = dimag(mo_coef_complex(i,j))
enddo
enddo
call ezfio_set_mo_basis_mo_coef_real(buffer)
call ezfio_set_mo_basis_mo_coef_imag(buffer_im)
deallocate (buffer,buffer_im)
else
allocate ( buffer(ao_num,mo_num) )
buffer = 0.d0 buffer = 0.d0
do j = 1, mo_num do j = 1, mo_num
do i = 1, ao_num do i = 1, ao_num
buffer(i,j) = mo_coef_imag(i,j) buffer(i,j) = mo_coef(i,j)
enddo enddo
enddo enddo
call ezfio_set_mo_basis_mo_coef_imag(buffer) call ezfio_set_mo_basis_mo_coef(buffer)
deallocate (buffer)
endif endif
deallocate (buffer) call ezfio_set_mo_basis_mo_occ(mo_occ)
call ezfio_set_mo_basis_mo_class(mo_class)
end end
subroutine save_mos_no_occ subroutine save_mos_no_occ
implicit none implicit none
double precision, allocatable :: buffer(:,:) double precision, allocatable :: buffer(:,:),buffer_im(:,:)
integer :: i,j integer :: i,j
call system('$QP_ROOT/scripts/save_current_mos.sh '//trim(ezfio_filename)) call system('$QP_ROOT/scripts/save_current_mos.sh '//trim(ezfio_filename))
!call ezfio_set_mo_basis_mo_num(mo_num) !call ezfio_set_mo_basis_mo_num(mo_num)
!call ezfio_set_mo_basis_mo_label(mo_label) !call ezfio_set_mo_basis_mo_label(mo_label)
!call ezfio_set_mo_basis_ao_md5(ao_md5) !call ezfio_set_mo_basis_ao_md5(ao_md5)
allocate ( buffer(ao_num,mo_num) )
buffer = 0.d0
do j = 1, mo_num
do i = 1, ao_num
buffer(i,j) = mo_coef(i,j)
enddo
enddo
call ezfio_set_mo_basis_mo_coef(buffer)
if (is_periodic) then if (is_periodic) then
allocate ( buffer(ao_num,mo_num),buffer_im(ao_num,mo_num))
buffer = 0.d0
buffer_im = 0.d0
do j = 1, mo_num
do i = 1, ao_num
buffer(i,j) = dble(mo_coef_complex(i,j))
buffer_im(i,j) = dimag(mo_coef_complex(i,j))
enddo
enddo
call ezfio_set_mo_basis_mo_coef_real(buffer)
call ezfio_set_mo_basis_mo_coef_imag(buffer_im)
deallocate (buffer,buffer_im)
else
allocate ( buffer(ao_num,mo_num) )
buffer = 0.d0 buffer = 0.d0
do j = 1, mo_num do j = 1, mo_num
do i = 1, ao_num do i = 1, ao_num
buffer(i,j) = mo_coef_imag(i,j) buffer(i,j) = mo_coef(i,j)
enddo enddo
enddo enddo
call ezfio_set_mo_basis_mo_coef_imag(buffer) call ezfio_set_mo_basis_mo_coef(buffer)
deallocate (buffer)
endif endif
deallocate (buffer)
end end
subroutine save_mos_truncated(n) subroutine save_mos_truncated(n)
implicit none implicit none
double precision, allocatable :: buffer(:,:),buffer_im(:,:)
double precision, allocatable :: buffer(:,:) double precision, allocatable :: buffer(:,:)
integer :: i,j,n integer :: i,j,n
@ -74,26 +87,32 @@ subroutine save_mos_truncated(n)
call ezfio_set_mo_basis_mo_num(n) call ezfio_set_mo_basis_mo_num(n)
call ezfio_set_mo_basis_mo_label(mo_label) call ezfio_set_mo_basis_mo_label(mo_label)
call ezfio_set_mo_basis_ao_md5(ao_md5) 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)
call ezfio_set_mo_basis_mo_class(mo_class)
if (is_periodic) then if (is_periodic) then
allocate ( buffer(ao_num,n),buffer_im(ao_num,n))
buffer = 0.d0
buffer_im = 0.d0
do j = 1, n
do i = 1, ao_num
buffer(i,j) = dble(mo_coef_complex(i,j))
buffer_im(i,j) = dimag(mo_coef_complex(i,j))
enddo
enddo
call ezfio_set_mo_basis_mo_coef_real(buffer)
call ezfio_set_mo_basis_mo_coef_imag(buffer_im)
deallocate (buffer,buffer_im)
else
allocate ( buffer(ao_num,n) )
buffer = 0.d0 buffer = 0.d0
do j = 1, n do j = 1, n
do i = 1, ao_num do i = 1, ao_num
buffer(i,j) = mo_coef_imag(i,j) buffer(i,j) = mo_coef(i,j)
enddo enddo
enddo enddo
call ezfio_set_mo_basis_mo_coef_imag(buffer) call ezfio_set_mo_basis_mo_coef(buffer)
deallocate (buffer)
endif endif
deallocate (buffer) call ezfio_set_mo_basis_mo_occ(mo_occ)
call ezfio_set_mo_basis_mo_class(mo_class)
end end

View File

@ -5,7 +5,6 @@ compare master-features_periodic
TODO: TODO:
ao_ints ao_ints
reverse index reverse index
s_half_inv_complex
ao_overlap_abs for complex ao_overlap_abs for complex
ao_integrals_n_e_per_atom_complex? ao_integrals_n_e_per_atom_complex?
not implemented for periodic: not implemented for periodic:
@ -16,14 +15,10 @@ ao_ints
compute_ao_integrals_jl compute_ao_integrals_jl
mo_basis mo_basis
decide how to handle real/imag/complex parts of mo_coef (maybe just need to chage save_mos?)
reorder_core_orb: implement for periodic reorder_core_orb: implement for periodic
save_mos_no_occ: implement for periodic
scf scf
finish complex DIIS
finish ao_two_e_integral_{alpha,beta}_complex (need reverse index?) finish ao_two_e_integral_{alpha,beta}_complex (need reverse index?)
finish eigenvectors_Fock_matrix_AO_complex
mo_two_e_ints mo_two_e_ints
not started not started