GAS_singles_uniform_init Subroutine

private subroutine GAS_singles_uniform_init(this, GAS_spec, use_lookup, create_lookup)

Type Bound

GAS_singles_PC_uniform_ExcGenerator_t

Arguments

Type IntentOptional Attributes Name
class(GAS_singles_PC_uniform_ExcGenerator_t), intent(inout) :: this
class(GASSpec_t), intent(in) :: GAS_spec
logical, intent(in) :: use_lookup
logical, intent(in) :: create_lookup

Contents


Source Code

    subroutine GAS_singles_uniform_init(this, GAS_spec, use_lookup, create_lookup)
        class(GAS_singles_PC_uniform_ExcGenerator_t), intent(inout) :: this
        class(GASSpec_t), intent(in) :: GAS_spec
        logical, intent(in) :: use_lookup, create_lookup
        integer, allocatable :: supergroups(:, :)
        character(*), parameter :: this_routine = 'GAS_singles_uniform_init'

        integer :: i_sg, src, tgt

        this%GAS_spec = GAS_spec
        allocate(this%indexer, source=SuperGroupIndexer_t(GAS_spec, nEl))
        this%create_lookup = create_lookup
        if (create_lookup) then
            if (associated(lookup_supergroup_indexer)) then
                call stop_all(this_routine, 'Someone else is already managing the supergroup lookup.')
            else
                write(stdout, *) 'GAS singles is creating and managing the supergroup lookup'
                lookup_supergroup_indexer => this%indexer
            end if
        end if
        this%use_lookup = use_lookup
        if (use_lookup) write(stdout, *) 'GAS singles is using the supergroup lookup'
        ! possible supergroups
        supergroups = this%indexer%get_supergroups()

        allocate(this%allowed_holes(0 : nIfD, this%GAS_spec%n_spin_orbs(), size(supergroups, 2)), source=0_n_int)


        ! Find for each supergroup (i_sg)
        ! the allowed holes `tgt` for a given `src`.
        do i_sg = 1, size(supergroups, 2)
            ! Note that the loop cannot take a symmetric shape.
            ! It may be, that (1 -> 5) is allowed, but (5 -> 1) is not.
            do src = 1, this%GAS_spec%n_spin_orbs()
                do tgt = 1, this%GAS_spec%n_spin_orbs()
                    if (src == tgt) cycle
                    associate(exc => Excite_1_t(src, tgt))
                    if (spin_allowed(exc) &
                            .and. this%GAS_spec%is_allowed(exc, supergroups(:, i_sg)) &
                            .and. symmetry_allowed(exc)) then
                        call my_set_orb(this%allowed_holes(:, src, i_sg), tgt)
                    end if
                    end associate
                end do
            end do
        end do

    contains

            ! Ugly Fortran syntax rules forces us to not use the macro.
            ! Even associate would not help here
            ! https://stackoverflow.com/questions/65734764/non-one-indexed-array-from-associate
            pure subroutine my_set_orb(ilut, orb)
                integer(n_int), intent(inout) :: ilut(0 : nIfD)
                integer, intent(in) :: orb
                set_orb(ilut, orb)
            end subroutine

            ! For single excitations it is simple
            logical pure function symmetry_allowed(exc)
                use SymExcitDataMod, only: SpinOrbSymLabel
                type(Excite_1_t), intent(in) :: exc
                symmetry_allowed = SpinOrbSymLabel(exc%val(1, 1)) == SpinOrbSymLabel(exc%val(2, 1))
            end function
    end subroutine GAS_singles_uniform_init