hf_energy.f90 Source File


This file depends on

sourcefile~~hf_energy.f90~~EfferentGraph sourcefile~hf_energy.f90 hf_energy.f90 sourcefile~basis_tools.f90 basis_tools.F90 sourcefile~hf_energy.f90->sourcefile~basis_tools.f90 sourcefile~c_interop.f90 c_interop.F90 sourcefile~hf_energy.f90->sourcefile~c_interop.f90 sourcefile~constants_io.f90 constants_io.F90 sourcefile~hf_energy.f90->sourcefile~constants_io.f90 sourcefile~dft.f90 dft.F90 sourcefile~hf_energy.f90->sourcefile~dft.f90 sourcefile~dft_molgrid.f90 dft_molgrid.F90 sourcefile~hf_energy.f90->sourcefile~dft_molgrid.f90 sourcefile~messages.f90 messages.F90 sourcefile~hf_energy.f90->sourcefile~messages.f90 sourcefile~printing.f90 printing.F90 sourcefile~hf_energy.f90->sourcefile~printing.f90 sourcefile~scf.f90 scf.F90 sourcefile~hf_energy.f90->sourcefile~scf.f90 sourcefile~strings.f90 strings.F90 sourcefile~hf_energy.f90->sourcefile~strings.f90 sourcefile~tagarray_driver.f90 tagarray_driver.F90 sourcefile~hf_energy.f90->sourcefile~tagarray_driver.f90 sourcefile~types.f90 types.F90 sourcefile~hf_energy.f90->sourcefile~types.f90 sourcefile~basis_tools.f90->sourcefile~constants_io.f90 sourcefile~basis_tools.f90->sourcefile~messages.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~parallel.f90 parallel.F90 sourcefile~basis_tools.f90->sourcefile~parallel.f90 sourcefile~precision.f90 precision.F90 sourcefile~basis_tools.f90->sourcefile~precision.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~dft.f90->sourcefile~basis_tools.f90 sourcefile~dft.f90->sourcefile~constants_io.f90 sourcefile~dft.f90->sourcefile~dft_molgrid.f90 sourcefile~dft.f90->sourcefile~messages.f90 sourcefile~dft.f90->sourcefile~strings.f90 sourcefile~dft.f90->sourcefile~tagarray_driver.f90 sourcefile~dft.f90->sourcefile~types.f90 sourcefile~bragg_slater.f90 bragg_slater.F90 sourcefile~dft.f90->sourcefile~bragg_slater.f90 sourcefile~dft.f90->sourcefile~constants.f90 sourcefile~dft_fuzzycell.f90 dft_fuzzycell.F90 sourcefile~dft.f90->sourcefile~dft_fuzzycell.f90 sourcefile~dft_gridint_energy.f90 dft_gridint_energy.F90 sourcefile~dft.f90->sourcefile~dft_gridint_energy.f90 sourcefile~dft_gridint_grad.f90 dft_gridint_grad.F90 sourcefile~dft.f90->sourcefile~dft_gridint_grad.f90 sourcefile~grid_storage.f90 grid_storage.F90 sourcefile~dft.f90->sourcefile~grid_storage.f90 sourcefile~libxc.f90 libxc.F90 sourcefile~dft.f90->sourcefile~libxc.f90 sourcefile~mathlib.f90 mathlib.F90 sourcefile~dft.f90->sourcefile~mathlib.f90 sourcefile~physical_constants.f90 physical_constants.F90 sourcefile~dft.f90->sourcefile~physical_constants.f90 sourcefile~dft.f90->sourcefile~precision.f90 sourcefile~radial_grid_types.f90 radial_grid_types.F90 sourcefile~dft.f90->sourcefile~radial_grid_types.f90 sourcefile~dft_molgrid.f90->sourcefile~bragg_slater.f90 sourcefile~dft_molgrid.f90->sourcefile~constants.f90 sourcefile~dft_partfunc.f90 dft_partfunc.F90 sourcefile~dft_molgrid.f90->sourcefile~dft_partfunc.f90 sourcefile~dft_molgrid.f90->sourcefile~grid_storage.f90 sourcefile~lebedev.f90 lebedev.F90 sourcefile~dft_molgrid.f90->sourcefile~lebedev.f90 sourcefile~dft_molgrid.f90->sourcefile~precision.f90 sourcefile~messages.f90->sourcefile~constants_io.f90 sourcefile~messages.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~tagarray_driver.f90 sourcefile~printing.f90->sourcefile~types.f90 sourcefile~printing.f90->sourcefile~precision.f90 sourcefile~scf.f90->sourcefile~basis_tools.f90 sourcefile~scf.f90->sourcefile~constants_io.f90 sourcefile~scf.f90->sourcefile~dft.f90 sourcefile~scf.f90->sourcefile~dft_molgrid.f90 sourcefile~scf.f90->sourcefile~messages.f90 sourcefile~scf.f90->sourcefile~printing.f90 sourcefile~scf.f90->sourcefile~tagarray_driver.f90 sourcefile~scf.f90->sourcefile~types.f90 sourcefile~scf.f90->sourcefile~constants.f90 sourcefile~guess.f90 guess.F90 sourcefile~scf.f90->sourcefile~guess.f90 sourcefile~int2.f90 int2.F90 sourcefile~scf.f90->sourcefile~int2.f90 sourcefile~scf.f90->sourcefile~mathlib.f90 sourcefile~oqp_linalg.f90 oqp_linalg.F90 sourcefile~scf.f90->sourcefile~oqp_linalg.f90 sourcefile~scf.f90->sourcefile~precision.f90 sourcefile~scf_converger.f90 scf_converger.F90 sourcefile~scf.f90->sourcefile~scf_converger.f90 sourcefile~util.f90 util.F90 sourcefile~scf.f90->sourcefile~util.f90 sourcefile~tagarray_driver.f90->sourcefile~messages.f90 sourcefile~types.f90->sourcefile~basis_tools.f90 sourcefile~types.f90->sourcefile~atomic_structure.f90 sourcefile~functionals.f90 functionals.F90 sourcefile~types.f90->sourcefile~functionals.f90 sourcefile~types.f90->sourcefile~parallel.f90 sourcefile~types.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~bragg_slater.f90->sourcefile~precision.f90 sourcefile~constants.f90->sourcefile~precision.f90 sourcefile~dft_fuzzycell.f90->sourcefile~basis_tools.f90 sourcefile~dft_fuzzycell.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_fuzzycell.f90->sourcefile~atomic_structure.f90 sourcefile~dft_fuzzycell.f90->sourcefile~constants.f90 sourcefile~dft_fuzzycell.f90->sourcefile~dft_partfunc.f90 sourcefile~dft_fuzzycell.f90->sourcefile~grid_storage.f90 sourcefile~dft_fuzzycell.f90->sourcefile~precision.f90 sourcefile~dft_gridint_energy.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_energy.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_energy.f90->sourcefile~types.f90 sourcefile~dft_gridint_energy.f90->sourcefile~oqp_linalg.f90 sourcefile~dft_gridint_energy.f90->sourcefile~precision.f90 sourcefile~dft_gridint.f90 dft_gridint.F90 sourcefile~dft_gridint_energy.f90->sourcefile~dft_gridint.f90 sourcefile~dft_gridint_grad.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint_grad.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint_grad.f90->sourcefile~types.f90 sourcefile~dft_gridint_grad.f90->sourcefile~precision.f90 sourcefile~dft_gridint_grad.f90->sourcefile~dft_gridint.f90 sourcefile~dft_partfunc.f90->sourcefile~precision.f90 sourcefile~elements.f90->sourcefile~strings.f90 sourcefile~elements.f90->sourcefile~physical_constants.f90 sourcefile~functionals.f90->sourcefile~messages.f90 sourcefile~functionals.f90->sourcefile~precision.f90 sourcefile~grid_storage.f90->sourcefile~precision.f90 sourcefile~guess.f90->sourcefile~basis_tools.f90 sourcefile~guess.f90->sourcefile~messages.f90 sourcefile~guess.f90->sourcefile~types.f90 sourcefile~guess.f90->sourcefile~mathlib.f90 sourcefile~guess.f90->sourcefile~oqp_linalg.f90 sourcefile~guess.f90->sourcefile~precision.f90 sourcefile~eigen.f90 eigen.F90 sourcefile~guess.f90->sourcefile~eigen.f90 sourcefile~int2.f90->sourcefile~basis_tools.f90 sourcefile~int2.f90->sourcefile~constants_io.f90 sourcefile~int2.f90->sourcefile~messages.f90 sourcefile~int2.f90->sourcefile~tagarray_driver.f90 sourcefile~int2.f90->sourcefile~types.f90 sourcefile~int2.f90->sourcefile~atomic_structure.f90 sourcefile~int2.f90->sourcefile~constants.f90 sourcefile~int2.f90->sourcefile~parallel.f90 sourcefile~int2.f90->sourcefile~precision.f90 sourcefile~int2_pairs.f90 int2_pairs.F90 sourcefile~int2.f90->sourcefile~int2_pairs.f90 sourcefile~int_libint.f90 int_libint.F90 sourcefile~int2.f90->sourcefile~int_libint.f90 sourcefile~int_rotaxis.f90 int_rotaxis.F90 sourcefile~int2.f90->sourcefile~int_rotaxis.f90 sourcefile~int_rys.f90 int_rys.F90 sourcefile~int2.f90->sourcefile~int_rys.f90 sourcefile~lebedev.f90->sourcefile~messages.f90 sourcefile~lebedev.f90->sourcefile~precision.f90 sourcefile~libxc.f90->sourcefile~messages.f90 sourcefile~libxc.f90->sourcefile~types.f90 sourcefile~libxc.f90->sourcefile~functionals.f90 sourcefile~libxc.f90->sourcefile~precision.f90 sourcefile~mathlib.f90->sourcefile~messages.f90 sourcefile~mathlib.f90->sourcefile~oqp_linalg.f90 sourcefile~mathlib.f90->sourcefile~precision.f90 sourcefile~mathlib.f90->sourcefile~eigen.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~parallel.f90->sourcefile~precision.f90 sourcefile~radial_grid_types.f90->sourcefile~precision.f90 sourcefile~scf_converger.f90->sourcefile~constants_io.f90 sourcefile~scf_converger.f90->sourcefile~printing.f90 sourcefile~scf_converger.f90->sourcefile~mathlib.f90 sourcefile~scf_converger.f90->sourcefile~oqp_linalg.f90 sourcefile~scf_converger.f90->sourcefile~precision.f90 sourcefile~util.f90->sourcefile~precision.f90 sourcefile~blas_wrap.f90->sourcefile~messages.f90 sourcefile~blas_wrap.f90->sourcefile~precision.f90 sourcefile~mathlib_types.f90 mathlib_types.F90 sourcefile~blas_wrap.f90->sourcefile~mathlib_types.f90 sourcefile~dft_gridint.f90->sourcefile~basis_tools.f90 sourcefile~dft_gridint.f90->sourcefile~constants_io.f90 sourcefile~dft_gridint.f90->sourcefile~dft_molgrid.f90 sourcefile~dft_gridint.f90->sourcefile~functionals.f90 sourcefile~dft_gridint.f90->sourcefile~oqp_linalg.f90 sourcefile~dft_gridint.f90->sourcefile~parallel.f90 sourcefile~dft_gridint.f90->sourcefile~precision.f90 sourcefile~dft_xc_libxc.f90 dft_xc_libxc.F90 sourcefile~dft_gridint.f90->sourcefile~dft_xc_libxc.f90 sourcefile~eigen.f90->sourcefile~messages.f90 sourcefile~eigen.f90->sourcefile~oqp_linalg.f90 sourcefile~eigen.f90->sourcefile~precision.f90 sourcefile~eigen.f90->sourcefile~mathlib_types.f90 sourcefile~int2_pairs.f90->sourcefile~basis_tools.f90 sourcefile~int2_pairs.f90->sourcefile~precision.f90 sourcefile~int_libint.f90->sourcefile~basis_tools.f90 sourcefile~int_libint.f90->sourcefile~constants.f90 sourcefile~int_libint.f90->sourcefile~precision.f90 sourcefile~int_libint.f90->sourcefile~int2_pairs.f90 sourcefile~libint_f.f90 libint_f.F90 sourcefile~int_libint.f90->sourcefile~libint_f.f90 sourcefile~int_rotaxis.f90->sourcefile~basis_tools.f90 sourcefile~int_rotaxis.f90->sourcefile~constants.f90 sourcefile~int_rotaxis.f90->sourcefile~precision.f90 sourcefile~int_rotaxis.f90->sourcefile~int2_pairs.f90 sourcefile~int_rys.f90->sourcefile~basis_tools.f90 sourcefile~int_rys.f90->sourcefile~constants.f90 sourcefile~int_rys.f90->sourcefile~precision.f90 sourcefile~int_rys.f90->sourcefile~int2_pairs.f90 sourcefile~rys.f90 rys.F90 sourcefile~int_rys.f90->sourcefile~rys.f90 sourcefile~lapack_wrap.f90->sourcefile~messages.f90 sourcefile~lapack_wrap.f90->sourcefile~precision.f90 sourcefile~lapack_wrap.f90->sourcefile~mathlib_types.f90 sourcefile~dft_xc_libxc.f90->sourcefile~functionals.f90 sourcefile~dft_xc_libxc.f90->sourcefile~precision.f90 sourcefile~dft_xclib.f90 dft_xclib.F90 sourcefile~dft_xc_libxc.f90->sourcefile~dft_xclib.f90 sourcefile~rys.f90->sourcefile~constants.f90 sourcefile~rys.f90->sourcefile~precision.f90 sourcefile~rys_lut.f90 rys_lut.F90 sourcefile~rys.f90->sourcefile~rys_lut.f90 sourcefile~dft_xclib.f90->sourcefile~functionals.f90 sourcefile~dft_xclib.f90->sourcefile~precision.f90

Source Code

module hf_energy_mod

  implicit none

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

  private

  public hf_energy

contains

  subroutine hf_energy_C(c_handle) bind(C, name="hf_energy")
    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 hf_energy(inf)
  end subroutine hf_energy_C

  subroutine hf_energy(infos)
    use io_constants, only: iw
    use basis_tools, only: basis_set
    use messages, only: show_message
    use scf, only: scf_driver
    use dft, only: dft_initialize, dftclean
    use types, only: information
    use oqp_tagarray_driver
    use strings, only: Cstring, fstring
    use mod_dft_molgrid, only: dft_grid_t
    use printing, only: print_module_info

    implicit none

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

    type(information), target, intent(inout) :: infos

    integer :: nbf2, nbf, nsh2
    logical :: urohf, dft
    type(basis_set), pointer :: basis
    type(dft_grid_t) :: molGrid

    urohf = infos%control%scftype == 2 .or. infos%control%scftype == 3
    dft = infos%control%hamilton == 20

!   3. LOG: Write: Main output file
    open (unit=iw, file=infos%log_filename, position="append")
!
    call print_module_info('HF_DFT_Energy','Computing HF/DFT SCF Energy')

!   Load basis set
    basis => infos%basis
    basis%atoms => infos%atoms

!   Allocate memory
    nbf = basis%nbf
    nbf2 = nbf*(nbf+1)/2
    nsh2 = (basis%nshell**2+basis%nshell)/2

    ! clean data
    call infos%dat%remove_records((/ character(len=80) :: OQP_FOCK_A, OQP_FOCK_B /))

    call infos%dat%reserve_data(OQP_FOCK_A, TA_TYPE_REAL64, nbf2, comment=OQP_FOCK_A_comment)
    call check_status(infos%dat%get_status(), module_name, subroutine_name, OQP_FOCK_A)

    if (urohf) then
      call infos%dat%reserve_data(OQP_FOCK_B, TA_TYPE_REAL64, nbf2, comment=OQP_FOCK_B_comment)
      call check_status(infos%dat%get_status(), module_name, subroutine_name, OQP_FOCK_B)
    end if

!   Prepare dft grid
    if (dft) call dft_initialize(infos, basis, molGrid, verbose=.true.)

!   Run HF/DFT calculation
    call scf_driver(basis, infos, molGrid)

!   Cleanup
    if (dft) call dftclean(infos)

    close(iw)

  end subroutine hf_energy

end module hf_energy_mod