get_from_cache Subroutine

private subroutine get_from_cache(me, x, ifs, i, f, xfound, ffound)

Check if the x vector is in the cache, if so return f. Note that only some of the elements may be present, so it will return the ones there are there, and indicate which ones were found.

Type Bound

function_cache

Arguments

Type IntentOptional Attributes Name
class(function_cache), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

independant variable vector

integer, intent(in), dimension(:) :: ifs

elements of f needed

integer, intent(out) :: i

index in the hash table

real(kind=wp), intent(out), dimension(:) :: f

f(x) from the cache (if it was found)

logical, intent(out) :: xfound

if x was found in the cache

logical, intent(out), dimension(size(ifs)) :: ffound

which ifs were found in the cache


Calls

proc~~get_from_cache~~CallsGraph proc~get_from_cache numdiff_cache_module::function_cache%get_from_cache proc~vector_djb_hash numdiff_cache_module::vector_djb_hash proc~get_from_cache->proc~vector_djb_hash

Called by

proc~~get_from_cache~~CalledByGraph proc~get_from_cache numdiff_cache_module::function_cache%get_from_cache proc~compute_function_with_cache numerical_differentiation_module::compute_function_with_cache proc~compute_function_with_cache->proc~get_from_cache

Source Code

    subroutine get_from_cache(me,x,ifs,i,f,xfound,ffound)

    implicit none

    class(function_cache),intent(inout)      :: me
    real(wp),dimension(:),intent(in)         :: x      !! independant variable vector
    integer,dimension(:),intent(in)          :: ifs    !! elements of `f` needed
    integer,intent(out)                      :: i      !! index in the hash table
    real(wp),dimension(:),intent(out)        :: f      !! `f(x)` from the cache (if it was found)
    logical,intent(out)                      :: xfound !! if `x` was found in the cache
    logical,dimension(size(ifs)),intent(out) :: ffound !! which `ifs` were found in the cache

    integer :: j !! counter

    ! initialize:
    xfound = .false.
    ffound = .false.

    if (allocated(me%c)) then

        ! get index in the hash table:
        i = mod( abs(vector_djb_hash(x)), int(size(me%c),ip) )

        ! check the table:
        if (allocated(me%c(i)%x) .and. allocated(me%c(i)%f)) then
            if (size(me%c(i)%x)==size(x) .and. &
                size(me%c(i)%f)==size(f)) then
                if (all(me%c(i)%x==x)) then
                    xfound = .true.
                    ! return the elements that were found in the cache:
                    f(me%c(i)%ifs) = me%c(i)%f(me%c(i)%ifs)
                    ! what indices are in the cache?
                    do j = 1, size(ifs)
                        ffound(j) = any(ifs(j)==me%c(i)%ifs)
                    end do
                end if
            end if
        end if

    else
        error stop 'Error: the cache has not been initialized.'
    end if

    end subroutine get_from_cache