guess_huckel.F90 Source File


This file depends on

sourcefile~~guess_huckel.f90~~EfferentGraph sourcefile~guess_huckel.f90 guess_huckel.F90 sourcefile~basis_tools.f90 basis_tools.F90 sourcefile~guess_huckel.f90->sourcefile~basis_tools.f90 sourcefile~c_interop.f90 c_interop.F90 sourcefile~guess_huckel.f90->sourcefile~c_interop.f90 sourcefile~constants_io.f90 constants_io.F90 sourcefile~guess_huckel.f90->sourcefile~constants_io.f90 sourcefile~guess.f90 guess.F90 sourcefile~guess_huckel.f90->sourcefile~guess.f90 sourcefile~huckel.f90 huckel.F90 sourcefile~guess_huckel.f90->sourcefile~huckel.f90 sourcefile~mathlib.f90 mathlib.F90 sourcefile~guess_huckel.f90->sourcefile~mathlib.f90 sourcefile~messages.f90 messages.F90 sourcefile~guess_huckel.f90->sourcefile~messages.f90 sourcefile~parallel.f90 parallel.F90 sourcefile~guess_huckel.f90->sourcefile~parallel.f90 sourcefile~precision.f90 precision.F90 sourcefile~guess_huckel.f90->sourcefile~precision.f90 sourcefile~printing.f90 printing.F90 sourcefile~guess_huckel.f90->sourcefile~printing.f90 sourcefile~strings.f90 strings.F90 sourcefile~guess_huckel.f90->sourcefile~strings.f90 sourcefile~tagarray_driver.f90 tagarray_driver.F90 sourcefile~guess_huckel.f90->sourcefile~tagarray_driver.f90 sourcefile~types.f90 types.F90 sourcefile~guess_huckel.f90->sourcefile~types.f90 sourcefile~util.f90 util.F90 sourcefile~guess_huckel.f90->sourcefile~util.f90 sourcefile~basis_tools.f90->sourcefile~constants_io.f90 sourcefile~basis_tools.f90->sourcefile~messages.f90 sourcefile~basis_tools.f90->sourcefile~parallel.f90 sourcefile~basis_tools.f90->sourcefile~precision.f90 sourcefile~atomic_structure.f90 atomic_structure.F90 sourcefile~basis_tools.f90->sourcefile~atomic_structure.f90 sourcefile~basis_library.f90 basis_library.F90 sourcefile~basis_tools.f90->sourcefile~basis_library.f90 sourcefile~constants.f90 constants.F90 sourcefile~basis_tools.f90->sourcefile~constants.f90 sourcefile~elements.f90 elements.F90 sourcefile~basis_tools.f90->sourcefile~elements.f90 sourcefile~c_interop.f90->sourcefile~basis_tools.f90 sourcefile~c_interop.f90->sourcefile~messages.f90 sourcefile~c_interop.f90->sourcefile~strings.f90 sourcefile~c_interop.f90->sourcefile~tagarray_driver.f90 sourcefile~c_interop.f90->sourcefile~types.f90 sourcefile~guess.f90->sourcefile~basis_tools.f90 sourcefile~guess.f90->sourcefile~mathlib.f90 sourcefile~guess.f90->sourcefile~messages.f90 sourcefile~guess.f90->sourcefile~precision.f90 sourcefile~guess.f90->sourcefile~types.f90 sourcefile~eigen.f90 eigen.F90 sourcefile~guess.f90->sourcefile~eigen.f90 sourcefile~oqp_linalg.f90 oqp_linalg.F90 sourcefile~guess.f90->sourcefile~oqp_linalg.f90 sourcefile~huckel.f90->sourcefile~basis_tools.f90 sourcefile~huckel.f90->sourcefile~guess.f90 sourcefile~huckel.f90->sourcefile~mathlib.f90 sourcefile~huckel.f90->sourcefile~messages.f90 sourcefile~huckel.f90->sourcefile~precision.f90 sourcefile~huckel.f90->sourcefile~types.f90 sourcefile~huckel.f90->sourcefile~atomic_structure.f90 sourcefile~huckel.f90->sourcefile~constants.f90 sourcefile~huckel.f90->sourcefile~eigen.f90 sourcefile~huckel_lut.f90 huckel_lut.F90 sourcefile~huckel.f90->sourcefile~huckel_lut.f90 sourcefile~int1.f90 int1.F90 sourcefile~huckel.f90->sourcefile~int1.f90 sourcefile~huckel.f90->sourcefile~oqp_linalg.f90 sourcefile~mathlib.f90->sourcefile~messages.f90 sourcefile~mathlib.f90->sourcefile~precision.f90 sourcefile~mathlib.f90->sourcefile~eigen.f90 sourcefile~mathlib.f90->sourcefile~oqp_linalg.f90 sourcefile~messages.f90->sourcefile~constants_io.f90 sourcefile~messages.f90->sourcefile~precision.f90 sourcefile~parallel.f90->sourcefile~precision.f90 sourcefile~printing.f90->sourcefile~basis_tools.f90 sourcefile~printing.f90->sourcefile~constants_io.f90 sourcefile~printing.f90->sourcefile~messages.f90 sourcefile~printing.f90->sourcefile~precision.f90 sourcefile~printing.f90->sourcefile~tagarray_driver.f90 sourcefile~printing.f90->sourcefile~types.f90 sourcefile~tagarray_driver.f90->sourcefile~messages.f90 sourcefile~types.f90->sourcefile~basis_tools.f90 sourcefile~types.f90->sourcefile~parallel.f90 sourcefile~types.f90->sourcefile~precision.f90 sourcefile~types.f90->sourcefile~atomic_structure.f90 sourcefile~functionals.f90 functionals.F90 sourcefile~types.f90->sourcefile~functionals.f90 sourcefile~util.f90->sourcefile~precision.f90 sourcefile~atomic_structure.f90->sourcefile~strings.f90 sourcefile~basis_library.f90->sourcefile~constants_io.f90 sourcefile~basis_library.f90->sourcefile~strings.f90 sourcefile~basis_library.f90->sourcefile~constants.f90 sourcefile~basis_library.f90->sourcefile~elements.f90 sourcefile~constants.f90->sourcefile~precision.f90 sourcefile~eigen.f90->sourcefile~messages.f90 sourcefile~eigen.f90->sourcefile~precision.f90 sourcefile~eigen.f90->sourcefile~oqp_linalg.f90 sourcefile~mathlib_types.f90 mathlib_types.F90 sourcefile~eigen.f90->sourcefile~mathlib_types.f90 sourcefile~elements.f90->sourcefile~strings.f90 sourcefile~physical_constants.f90 physical_constants.F90 sourcefile~elements.f90->sourcefile~physical_constants.f90 sourcefile~functionals.f90->sourcefile~messages.f90 sourcefile~functionals.f90->sourcefile~precision.f90 sourcefile~int1.f90->sourcefile~basis_tools.f90 sourcefile~int1.f90->sourcefile~constants_io.f90 sourcefile~int1.f90->sourcefile~messages.f90 sourcefile~int1.f90->sourcefile~parallel.f90 sourcefile~int1.f90->sourcefile~precision.f90 sourcefile~int1.f90->sourcefile~printing.f90 sourcefile~ecp.f90 ecp.F90 sourcefile~int1.f90->sourcefile~ecp.f90 sourcefile~mod_1e_primitives.f90 mod_1e_primitives.F90 sourcefile~int1.f90->sourcefile~mod_1e_primitives.f90 sourcefile~mod_shell_tools.f90 mod_shell_tools.F90 sourcefile~int1.f90->sourcefile~mod_shell_tools.f90 sourcefile~blas_wrap.f90 blas_wrap.F90 sourcefile~oqp_linalg.f90->sourcefile~blas_wrap.f90 sourcefile~lapack_wrap.f90 lapack_wrap.F90 sourcefile~oqp_linalg.f90->sourcefile~lapack_wrap.f90 sourcefile~blas_wrap.f90->sourcefile~messages.f90 sourcefile~blas_wrap.f90->sourcefile~precision.f90 sourcefile~blas_wrap.f90->sourcefile~mathlib_types.f90 sourcefile~ecp.f90->sourcefile~basis_tools.f90 sourcefile~ecp.f90->sourcefile~precision.f90 sourcefile~ecp.f90->sourcefile~constants.f90 sourcefile~ecpint.f90 ecpint.F90 sourcefile~ecp.f90->sourcefile~ecpint.f90 sourcefile~lapack_wrap.f90->sourcefile~messages.f90 sourcefile~lapack_wrap.f90->sourcefile~precision.f90 sourcefile~lapack_wrap.f90->sourcefile~mathlib_types.f90 sourcefile~mod_1e_primitives.f90->sourcefile~constants.f90 sourcefile~mod_1e_primitives.f90->sourcefile~mod_shell_tools.f90 sourcefile~mod_gauss_hermite.f90 mod_gauss_hermite.F90 sourcefile~mod_1e_primitives.f90->sourcefile~mod_gauss_hermite.f90 sourcefile~rys.f90 rys.F90 sourcefile~mod_1e_primitives.f90->sourcefile~rys.f90 sourcefile~xyzorder.f90 xyzorder.F90 sourcefile~mod_1e_primitives.f90->sourcefile~xyzorder.f90 sourcefile~mod_shell_tools.f90->sourcefile~basis_tools.f90 sourcefile~mod_shell_tools.f90->sourcefile~precision.f90 sourcefile~mod_gauss_hermite.f90->sourcefile~precision.f90 sourcefile~rys.f90->sourcefile~precision.f90 sourcefile~rys.f90->sourcefile~constants.f90 sourcefile~rys_lut.f90 rys_lut.F90 sourcefile~rys.f90->sourcefile~rys_lut.f90

Source Code

module guess_huckel_mod

  implicit none

  character(len=*), parameter :: module_name = "guess_huckel_mod"

contains

  subroutine guess_huckel_C(c_handle) bind(C, name="guess_huckel")
    use c_interop, only: oqp_handle_t, oqp_handle_get_info
    use types, only: information
    type(oqp_handle_t) :: c_handle
    type(information), pointer :: inf
    inf => oqp_handle_get_info(c_handle)
    call guess_huckel(inf)
  end subroutine guess_huckel_C

  subroutine guess_huckel(infos)
    use precision, only: dp
    use types, only: information
    use io_constants, only: IW
    use oqp_tagarray_driver
    use basis_tools, only: basis_set
    use guess, only: get_ab_initio_density
    use mathlib, only: matrix_invsqrt
    use huckel, only: huckel_guess
    use util, only: measure_time
    use messages, only: show_message, WITH_ABORT
    use strings, only: Cstring, fstring
    use printing, only: print_module_info
    use oqp_tagarray_driver
    use iso_c_binding, only: c_char
    use parallel, only: par_env_t

    implicit none

    character(len=*), parameter :: subroutine_name = "guess_huckel"

    type(information), target, intent(inout) :: infos
    integer :: i, nbf, nbf2

    type(basis_set), pointer :: basis
    type(basis_set) :: huckel_basis
    character(len=:), allocatable :: basis_file
    logical :: err
    integer , parameter :: root = 0
    type(par_env_t) :: pe
  ! tagarray
    real(kind=dp), contiguous, pointer :: &
      Smat(:), &
      dmat_a(:), mo_a(:,:), mo_energy_a(:), &
      dmat_b(:), mo_b(:,:), mo_energy_b(:)
    character(len=*), parameter :: tags_alpha(3) = (/ character(len=80) :: &
      OQP_DM_A, OQP_E_MO_A, OQP_VEC_MO_A /)
    character(len=*), parameter :: tags_beta(3) = (/ character(len=80) :: &
      OQP_DM_B, OQP_E_MO_B, OQP_VEC_MO_B /)
    character(len=*), parameter :: tags_general(2) = (/ character(len=80) :: &
      OQP_SM, OQP_hbasis_filename /)
    character(len=1,kind=c_char), contiguous, pointer :: basis_filename(:)
  !
  ! Section of Tagarray for the log filename
  ! We are getting lot file name from Python via tagarray
  !

    call data_has_tags(infos%dat, tags_general, module_name, subroutine_name, with_abort)
    call tagarray_get_data(infos%dat, OQP_hbasis_filename, basis_filename)
    allocate(character(ubound(basis_filename,1)) :: basis_file)
    do i = 1, ubound(basis_filename,1)
       basis_file(i:i) = basis_filename(i)
    end do
   !
  ! Files open
  ! 1. XYZ: Read : Geometric data, ATOMS
  ! 3. LOG: Read Write: Main output file
  !
    open (unit=IW, file=infos%log_filename, position="append")
  !
    call print_module_info('Guess_Huckel','Initial guess using Huckel theory')
  !
  ! Readings
  ! load basis set
    basis => infos%basis
    call pe%init(infos%mpiinfo%comm, infos%mpiinfo%usempi)

    if (pe%rank == root) then
      call huckel_basis%from_file(basis_file, infos%atoms, err)
    end if

  ! Checking error of basis set reading..
    infos%control%basis_set_issue = err
    call pe%bcast(infos%control%basis_set_issue, 1)

    basis%atoms => infos%atoms
  !  Allocate H, S ,T and D matrices
    nbf = basis%nbf
    nbf2 =nbf*(nbf+1)/2

    ! clean data
    call infos%dat%remove_records(tags_alpha)
    call infos%dat%remove_records(tags_beta)

    ! load general data
    call data_has_tags(infos%dat, tags_general, module_name, subroutine_name, WITH_ABORT)
    call tagarray_get_data(infos%dat, OQP_SM, smat)

    ! allocate alpha
    call infos%dat%reserve_data(OQP_DM_A, TA_TYPE_REAL64, nbf2, comment=OQP_DM_A_comment)
    call infos%dat%reserve_data(OQP_E_MO_A, TA_TYPE_REAL64, nbf, comment=OQP_E_MO_A_comment)
    call infos%dat%reserve_data(OQP_VEC_MO_A, TA_TYPE_REAL64, nbf*nbf, (/ nbf, nbf /), comment=OQP_VEC_MO_A_comment)

    ! load alpha data
    call data_has_tags(infos%dat, tags_alpha, module_name, subroutine_name, WITH_ABORT)
    call tagarray_get_data(infos%dat, OQP_DM_A, dmat_a)
    call tagarray_get_data(infos%dat, OQP_E_MO_A, mo_energy_a)
    call tagarray_get_data(infos%dat, OQP_VEC_MO_A, mo_a)

  ! UHF/ROHF
    if (infos%control%scftype >= 2) then
      ! allocate beta
      call infos%dat%reserve_data(OQP_DM_B, TA_TYPE_REAL64, nbf2, comment=OQP_DM_B_comment)
      call infos%dat%reserve_data(OQP_E_MO_B, TA_TYPE_REAL64, nbf, comment=OQP_E_MO_B_comment)
      call infos%dat%reserve_data(OQP_VEC_MO_B, TA_TYPE_REAL64, nbf*nbf, (/ nbf, nbf /), comment=OQP_VEC_MO_B_comment)

      ! load beta
      call data_has_tags(infos%dat, tags_beta, module_name, subroutine_name, WITH_ABORT)
      call tagarray_get_data(infos%dat, OQP_DM_B, dmat_b)
      call tagarray_get_data(infos%dat, OQP_E_MO_B, mo_energy_b)
      call tagarray_get_data(infos%dat, OQP_VEC_MO_B, mo_b)
    end if

  ! Calculate Huckel MO
    if (pe%rank == root) then
      call huckel_guess(Smat, MO_A, infos, basis, huckel_basis)
    endif
  !  For ROHF/UHF
    if (INFOS%control%scftype >= 2) MO_B = MO_A
  !
  ! Calculate Density Matrix
    if (pe%rank == root) then
  ! RHF
      if (infos%control%scftype == 1) then
        call get_ab_initio_density(Dmat_A, MO_A, infos=infos, basis=basis)
  ! ROHF/UHF
      else
        call get_ab_initio_density(Dmat_A, MO_A, Dmat_B, MO_B, infos, basis)
      endif
    endif
    ! Broadcast MO and density matrices to all processes
    call pe%bcast(MO_A, nbf*nbf)
    if (infos%control%scftype >= 2) then
      call pe%bcast(MO_B, nbf*nbf)
    endif
    ! Broadcast the density matrices to all processes
    if (infos%control%scftype == 1) then
      call pe%bcast(Dmat_A, nbf2)
    else
      call pe%bcast(Dmat_A, nbf2)
      call pe%bcast(Dmat_B, nbf2)
    endif
    call pe%barrier()
    write (iw, '(/x,a,/)') '...... End of initial orbital guess ......'
    call measure_time(print_total=1, log_unit=iw)
    close(iw)
  end subroutine guess_huckel

end module guess_huckel_mod