extract_fromIndex Subroutine

private recursive subroutine extract_fromIndex(fun, varSys, time, iLevel, idx, idxLen, nVals, res)

Same as extract_from Point except it extract from points via indices which are setup before

The interface has to comply to the abstract interface tem_varSys_module#tem_varSys_proc_getvalofIndex.

Arguments

Type IntentOptional Attributes Name
class(tem_varSys_op_type), intent(in) :: fun

Description of the method to obtain the variables, here some preset values might be stored, like the space time function to use or the required variables.

type(tem_varSys_type), intent(in) :: varSys

The variable system to obtain the variable from.

type(tem_time_type), intent(in) :: time

Point in time at which to evaluate the variable.

integer, intent(in) :: iLevel

Level on which values are requested

integer, intent(in) :: idx(:)

Index of points in the growing array and variable val array to return. Size: nVals

integer, intent(in), optional :: idxLen(:)

With idx as start index in contiguous memory, idxLength defines length of each contiguous memory Size: nVals

integer, intent(in) :: nVals

Number of values to obtain for this variable (vectorized access).

real(kind=rk), intent(out) :: res(:)

Resulting values for the requested variable.

Dimension: n requested entries x nComponents of this variable Access: (iElem-1)*fun%nComponents + iComp


Calls

proc~~extract_fromindex~~CallsGraph proc~extract_fromindex extract_fromIndex proc~tem_varsys_check_inargs tem_varSys_check_inArgs proc~extract_fromindex->proc~tem_varsys_check_inargs

Source Code

  recursive subroutine extract_fromIndex( fun, varSys, time, iLevel, idx, &
    &                                     idxLen, nVals, res              )
    ! ---------------------------------------------------------------------- !
    !> Description of the method to obtain the variables, here some preset
    !! values might be stored, like the space time function to use or the
    !! required variables.
    class(tem_varSys_op_type), intent(in) :: fun

    !> The variable system to obtain the variable from.
    type(tem_varSys_type), intent(in) :: varSys

    !> Point in time at which to evaluate the variable.
    type(tem_time_type), intent(in) :: time

    !> Level on which values are requested
    integer, intent(in) :: iLevel

    !> Index of points in the growing array and variable val array to
    !! return.
    !! Size: nVals
    integer, intent(in) :: idx(:)

    !> With idx as start index in contiguous memory,
    !! idxLength defines length of each contiguous memory
    !! Size: nVals
    integer, optional, intent(in) :: idxLen(:)

    !> Number of values to obtain for this variable (vectorized access).
    integer, intent(in) :: nVals

    !> Resulting values for the requested variable.
    !!
    !! Dimension: n requested entries x nComponents of this variable
    !! Access: (iElem-1)*fun%nComponents + iComp
    real(kind=rk), intent(out) :: res(:)
    ! ---------------------------------------------------------------------- !
    integer :: iVal, depVar_pos, depVar_nComps, iComp, comp_index
    real(kind=rk), allocatable :: input_varRes(:)
    ! ---------------------------------------------------------------------- !
    type(tem_varSys_op_data_type), pointer :: opData
    ! ---------------------------------------------------------------------- !

    call C_F_POINTER( fun%method_Data, opData )

    call tem_varSys_check_inArgs( fun, varSys, time, iLevel, idx, idxLen, &
      &    nVals, label = 'extract_fromIndex'                             )

    ! get the position of dependent variable, assuming only one
    depVar_pos = fun%input_varPos(1)

    ! get the nComp of dependent variable
    depVar_nComps = varSys%method%val( depVar_pos )%nComponents

    ! allocate array to save results
    allocate(input_varRes( nVals * depVar_nComps ))
    ! derive dependent variable
    call varSys%method%val(depVar_pos)%get_valOfIndex( &
        & varSys  = varSys,                            &
        & time    = time,                              &
        & iLevel  = iLevel,                            &
        & idx     = opData%input_pntIndex(1)           &
        &           %indexLvl(iLevel)%val( idx(:) ),   &
        & nVals   = nVals,                             &
        & res     = input_varRes                       )

    ! Extract component index from input_varRes
    do iVal = 1, nVals
      do iComp = 1, fun%nComponents
        comp_index = fun%input_varIndex(iComp)
        res((( ival-1)* fun%ncomponents+icomp) )                  &
          & = input_varRes( (( ival-1)* depvar_ncomps+comp_index) )
      end do
    end do

    deallocate(input_varRes)

  end subroutine extract_fromIndex