adi_criterium Function

public function adi_criterium(ilut, nI, sgn, ex, run) result(staticInit)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:NIfTot)
integer, intent(in) :: nI(nel)
real(kind=dp), intent(in) :: sgn(lenof_sign)
integer, intent(in) :: ex
integer, intent(in) :: run

Return Value logical


Contents

Source Code


Source Code

    function adi_criterium(ilut, nI, sgn, ex, run) result(staticInit)
        ! This is the adi-initiator criterium expansion
        ! I expect it to grow further
        use adi_references, only: initialize_c_caches, update_coherence_check, &
                                  eval_coherence
        use adi_data, only: tWeakCoherentDoubles, tAvCoherentDoubles, &
                            tUseCaches, exLvlRef
        use DetBitOps, only: FindBitExcitLevel
        implicit none
        integer(n_int), intent(in) :: ilut(0:NIfTot)
        real(dp), intent(in) :: sgn(lenof_sign)
        integer, intent(in) :: run, nI(nel), ex
        integer :: exLevel, i
        logical :: staticInit, tCCache, tCouplingPossible, tSI
        ! cache for the weak coherence check
        HElement_t(dp) :: signedCache
        real(dp) :: unsignedCache
        integer :: connections

        tSI = .false.
        staticInit = .false.
        tCCache = tWeakCoherentDoubles .or. tAvCoherentDoubles

        exLevel = 0
        if (tCCache) call initialize_c_caches(signedCache, unsignedCache, connections)
        ! Important : Only compare to the already initialized reference
        do i = 1, nRefs
            ! First, check if the excitation level differs by more than 2
            tCouplingPossible = .true.
            if (tUseCaches) then
                if (abs(ex - exLvlRef(i)) > max_ex_level) tCouplingPossible = .false.
            end if
            nExChecks = nExChecks + 1
            if (tCouplingPossible) then
                exLevel = FindBitExcitLevel(ilutRefAdi(:, i), ilut)
                ! We only need to do this if the excitation level is below 3
                if (exLevel < max_ex_level + 1) then
                    ! Check if we are sign-coherent if this is desired
                    if (unset_incoherent_initiator(exLevel, ilut, nI, sgn, i, staticInit, run)) &
                        ! If we find the determinant to be incoherent, we do not apply
                        ! the ADI rules and instead
                        return

                    if (tCCache) then
                        call update_coherence_check(ilut, nI, i, &
                                                    signedCache, unsignedCache, connections)
                        if (exLevel == 0) tSI = .true.
                    end if

                    ! Set the doubles to initiators
                    call set_si_initiator(exLevel, staticInit)
                    ! We have to keep looping over the SIs, even if we already have a double,
                    ! if we want to compute Xi. Else, we can exit if we found a valid double
                    if (staticInit .and. .not. tCCache) exit
                end if
            else
                nExCheckFails = nExCheckFails + 1
            end if
        end do

        ! superinitiators themselves are always static initiators -> no coherence check then
        if (tCCache .and. staticInit .and. .not. tSI) &
            call eval_coherence(signedCache, unsignedCache, sgn(run), connections, staticInit)

    end function adi_criterium