make_hermitian_rdm Subroutine

public subroutine make_hermitian_rdm(rdm, nrdms_standard, spawn, rdm_recv)

Arguments

Type IntentOptional Attributes Name
type(rdm_list_t), intent(in) :: rdm
integer, intent(in) :: nrdms_standard
type(rdm_spawn_t), intent(inout) :: spawn
type(rdm_list_t), intent(inout) :: rdm_recv

Contents

Source Code


Source Code

    subroutine make_hermitian_rdm(rdm, nrdms_standard, spawn, rdm_recv)

        ! Take the RDM in the rdm object, and output a new RDM object which is
        ! the same but with hermiticy applied to the non-transition RDMs, i.e.,
        ! the elements above and below the diagonal are averaged appropriately.
        ! Transition RDM signs are set to zero, so transition RDMs are not
        ! included in the output RDM.

        use rdm_data_utils, only: annihilate_rdm_list

        type(rdm_list_t), intent(in) :: rdm
        integer, intent(in) :: nrdms_standard
        type(rdm_spawn_t), intent(inout) :: spawn
        type(rdm_list_t), intent(inout) :: rdm_recv

        integer(int_rdm) :: ijkl
        integer :: ielem, ij, kl, i, j, k, l
        real(dp) :: rdm_sign(rdm%sign_length)
        logical :: nearly_full, finished, all_finished
#ifdef DEBUG_
        character(*), parameter :: this_routine = "make_hermitian_rdm"
#endif
        integer(int_rdm) :: ij_, kl_

        ! If we're about to fill up the spawn list, perform a communication.
        nearly_full = .false.
        ! Have we finished adding RDM elements to the spawned list?
        finished = .false.
        rdm_recv%nelements = 0

        do ielem = 1, rdm%nelements
            ! If the spawned list is nearly full, perform a communication.
            if (nearly_full) then
                call communicate_rdm_spawn_t_wrapper(spawn, rdm_recv, finished, all_finished)
                nearly_full = .false.
            end if

            ijkl = rdm%elements(0, ielem)

            call extract_sign_rdm(rdm%elements(:, ielem), rdm_sign)
            ! Set sign for transition RDMs to 0.
            rdm_sign(nrdms_standard + 1:) = 0.0_dp

            if (tGUGA) then
                call extract_2_rdm_ind(ijkl, i, j, k, l, ij_, kl_)
                ij = int(ij_)
                kl = int(kl_)

                rdm_sign = rdm_sign / 2.0_dp
                call add_to_rdm_spawn_t(spawn, k, l, i, j, rdm_sign, .true., nearly_full)
                call add_to_rdm_spawn_t(spawn, i, j, k, l, rdm_sign, .true., nearly_full)

            else
                ! Obtain spin orbital labels and the RDM element.
                call calc_separate_rdm_labels(ijkl, ij, kl, i, j, k, l)

                ! Factor of a half to account for prevent double-counting, and
                ! instead average elements from above and below the diagonal.
                if (ij /= kl) rdm_sign = 0.5_dp * rdm_sign

                ! If in the lower half of the RDM, reflect to the upper half.
                if (ij > kl) then
                    call add_to_rdm_spawn_t(spawn, k, l, i, j, rdm_sign, .false., nearly_full)
                else
                    call add_to_rdm_spawn_t(spawn, i, j, k, l, rdm_sign, .false., nearly_full)
                end if
            end if

        end do

        finished = .true.
        ! Keep performing communications until all RDM spawnings on every
        ! processor have been communicated.
        do
            call communicate_rdm_spawn_t_wrapper(spawn, rdm_recv, finished, all_finished)
            if (all_finished) exit
        end do
        call annihilate_rdm_list(rdm_recv)

    end subroutine make_hermitian_rdm