tem_opVar_fill_inputIndex Subroutine

public recursive subroutine tem_opVar_fill_inputIndex(fun, varSys, point, offset_bit, iLevel, tree, nPnts, inputIndex)

subroutine to fill index for the setuo Index routine called for operation variables, it is also used by the solver

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.

real(kind=rk), intent(in) :: point(:,:)

List of space coordinate points to store as growing array in method_data

character, intent(in), optional :: offset_bit(:)

Offset bit encoded as character for every point.

Offset integer coord(3) is converted into a character with offset_bit = achar( (coord(1)+1) + (coord(2)+1)4 + (coord(3)+1)16 ) Backward transformation form character to 3 integer: coord(1) = mod(ichar(offset_bit),4) - 1 coord(2) = mod(ichar(offset_bit),16)/4 - 1 coord(3) = ichar(offset_bit)/16 - 1

If not present default is to center i.e offset_bit = achar(1+4+16)

integer, intent(in) :: iLevel

Level to which input points belong to

type(treelmesh_type), intent(in) :: tree

global treelm mesh info

integer, intent(in) :: nPnts

Number of points to add in method_data of this variable

type(grw_intarray_type), intent(out) :: inputIndex(:)

input index for dependent variables size: fun%nInputs


Calls

proc~~tem_opvar_fill_inputindex~~CallsGraph proc~tem_opvar_fill_inputindex tem_opVar_fill_inputIndex interface~append~23 append proc~tem_opvar_fill_inputindex->interface~append~23 interface~truncate~17 truncate proc~tem_opvar_fill_inputindex->interface~truncate~17 proc~append_singlega2d_real append_singlega2d_real interface~append~23->proc~append_singlega2d_real proc~append_arrayga2d_real append_arrayga2d_real interface~append~23->proc~append_arrayga2d_real proc~truncate_ga_char truncate_ga_char interface~truncate~17->proc~truncate_ga_char interface~expand~22 expand proc~append_singlega2d_real->interface~expand~22 proc~append_arrayga2d_real->interface~expand~22 proc~expand_ga2d_real expand_ga2d_real interface~expand~22->proc~expand_ga2d_real

Called by

proc~~tem_opvar_fill_inputindex~~CalledByGraph proc~tem_opvar_fill_inputindex tem_opVar_fill_inputIndex proc~tem_opvar_setupindices tem_opVar_setupIndices proc~tem_opvar_setupindices->proc~tem_opvar_fill_inputindex

Contents


Source Code

  recursive subroutine tem_opVar_fill_inputIndex( fun, varSys, point,       &
    &                                             offset_bit, iLevel, tree, &
    &                                             nPnts, inputIndex         )
    !---------------------------`----------------------------------------------!
    !> 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

    !> List of space coordinate points to store as growing array in
    !! method_data
    real(kind=rk), intent(in) :: point(:,:)

    !> Offset bit encoded as character for every point.
    !!
    !! Offset integer coord(3) is converted into a character with
    !! offset_bit = achar( (coord(1)+1) + (coord(2)+1)*4 + (coord(3)+1)*16 )
    !! Backward transformation form character to 3 integer:
    !! coord(1) = mod(ichar(offset_bit),4) - 1
    !! coord(2) = mod(ichar(offset_bit),16)/4 - 1
    !! coord(3) = ichar(offset_bit)/16 - 1
    !!
    !! If not present default is to center i.e offset_bit = achar(1+4+16)
    character, optional, intent(in) :: offset_bit(:)

    !> Level to which input points belong to
    integer, intent(in) :: iLevel

    !> global treelm mesh info
    type(treelmesh_type), intent(in) :: tree

    !> Number of points to add in method_data of this variable
    integer, intent(in) :: nPnts

    !> input index for dependent variables
    !! size: fun%nInputs
    type(grw_intArray_type), intent(out) :: inputIndex(:)
    !--------------------------------------------------------------------------!
    integer, allocatable :: idx_loc(:)
    integer :: iDep, posDepVar
    !--------------------------------------------------------------------------!
    ! allocate local index array, needed for setup_indices call
    allocate( idx_loc(nPnts) )

    do iDep = 1, fun%nInputs
      idx_loc = 0

      ! get position of dependent var in varSys
      posDepVar = fun%input_varPos(iDep)

      ! setup indices in dependent variable.
      ! idx output from any variables will be same so it just overwrites and
      ! returns the idx of last dependent variable
      call varSys%method%val(posDepVar)%setup_indices( &
        & varSys     = varSys,                         &
        & point      = point,                          &
        & offset_bit = offset_bit,                     &
        & iLevel     = iLevel,                         &
        & tree       = tree,                           &
        & nPnts      = nPnts,                          &
        & idx        = idx_loc                         )

      call append( me  = inputIndex(iDep), &
        &          val = idx_loc           )

      call truncate (inputIndex(iDep))
      write(logUnit(9),*) 'nIndex on level for input variable ',  &
        &                 trim(varSys%varname%val(posDepVar)),    &
        &                 ' on Lvl ', iLevel,                     &
        &                 ' are = ', inputIndex(iDep)%nVals
    end do
    deallocate(idx_loc)

  end subroutine tem_opVar_fill_inputIndex