calc_determin_hamil_full Subroutine

public subroutine calc_determin_hamil_full(hamil, rep)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(out), allocatable :: hamil(:,:)
type(core_space_t) :: rep

Contents


Source Code

    subroutine calc_determin_hamil_full(hamil, rep)
        use guga_data, only: ExcitationInformation_t
        use guga_matrixElements, only: calc_guga_matrix_element
        type(core_space_t) :: rep
        type(CSF_Info_t) :: csf_i, csf_j
        type(ExcitationInformation_t) :: excitInfo

        HElement_t(dp), allocatable, intent(out) :: hamil(:, :)
        integer :: i, j, nI(nel), nJ(nel)

        allocate(hamil(rep%determ_space_size, rep%determ_space_size))

        hamil = h_cast(0.0_dp)

        do i = 1, rep%determ_space_size
            call decode_bit_det(nI, rep%core_space(:, i))

            if (tGUGA) csf_i = CSF_Info_t(rep%core_space(0:nifd, i))

            if (tHPHF) then
                hamil(i, i) = hphf_diag_helement(nI, rep%core_space(:, i))
            else
                hamil(i, i) = get_helement(nI, nI, 0)
            end if

            do j = 1, rep%determ_space_size

                if (i == j) cycle

                call decode_bit_det(nJ, rep%core_space(:, j))

                if (tGUGA) csf_j = CSF_Info_t(rep%core_space(:, j))

                if (tHPHF) then
                    hamil(i, j) = hphf_off_diag_helement(nI, nJ, &
                                     rep%core_space(:, i), rep%core_space(:, j))
                else if (tGUGA) then
                    call calc_guga_matrix_element(&
                        rep%core_space(:, i), csf_i, &
                        rep%core_space(:, j), csf_j, excitInfo, hamil(i, j), .true.)
                else
                    hamil(i, j) = get_helement(nI, nJ, rep%core_space(:, i), &
                        rep%core_space(:, j))
                end if

            end do
        end do

    end subroutine calc_determin_hamil_full