tem_get_point_chunk Subroutine

public subroutine tem_get_point_chunk(varsys, varPos, point, time, tree, nPnts, res)

A routine to evaluate chunk of points for given list of variables

If subTree is present, it will use map2Global from subTree else map2Global is created for current chunk of global mesh

Arguments

Type IntentOptional Attributes Name
type(tem_varSys_type), intent(in) :: varsys

Variable system describing available data.

integer, intent(in) :: varPos(:)

Position of variables to evaluate in varSys

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

Three-dimensional coordinates at which the variable should be evaluated. Only useful for variables provided as space-time functions.

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

Time information for the current data.

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

Mesh definition of the input data.

integer, intent(in) :: nPnts

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

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

Output data size: io_buffer_size


Called by

proc~~tem_get_point_chunk~~CalledByGraph proc~tem_get_point_chunk tem_get_point_chunk proc~combine_forpoint combine_forPoint proc~combine_forpoint->proc~tem_get_point_chunk proc~hvs_asciispatial_dump_point_data hvs_asciiSpatial_dump_point_data proc~hvs_asciispatial_dump_point_data->proc~tem_get_point_chunk proc~tem_convergence_check_point tem_convergence_check_point proc~tem_convergence_check_point->proc~tem_get_point_chunk proc~hvs_ascii_dump_point_data hvs_ascii_dump_point_data proc~hvs_ascii_dump_point_data->proc~tem_get_point_chunk proc~hvs_output_write hvs_output_write proc~hvs_output_write->proc~hvs_asciispatial_dump_point_data proc~hvs_output_write->proc~hvs_ascii_dump_point_data proc~tem_convergence_check tem_convergence_check proc~tem_convergence_check->proc~tem_convergence_check_point proc~tem_tracker tem_tracker proc~tem_tracker->proc~hvs_output_write

Contents

Source Code


Source Code

  subroutine tem_get_point_chunk(varSys, varPos, point, time, tree, nPnts, &
    &                            res)
    ! --------------------------------------------------------------------------!
    !> Variable system describing available data.
    type(tem_varsys_type), intent(in) :: varsys

    !> Position of variables to evaluate in varSys
    integer, intent(in) :: varPos(:)

    !> Mesh definition of the input data.
    type(treelmesh_type), intent(in) :: tree

    !> Time information for the current data.
    type(tem_time_type), intent(in) :: time

    !> Three-dimensional coordinates at which the variable should be
    !! evaluated. Only useful for variables provided as space-time functions.
    real(kind=rk), intent(in) :: point(:,:)

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

    !> Output data
    !! size: io_buffer_size
    real(kind=rk), intent(out) :: res(:)
    ! --------------------------------------------------------------------------!
    integer :: maxComponents, nScalars, nComponents, compOff
    integer :: iPnt, iVar, nVars, res_size
    integer :: e_start, t_start
    real(kind=rk), allocatable :: tmpdat(:)
    ! --------------------------------------------------------------------------!
    ! number of variables to evaluate
    nVars = size(varPos)

    ! Number of scalars in current output
    nScalars = sum(varSys%method%val(varPos(:))%nComponents)

    ! Need to obtain the data variable for variable, and store it in an
    ! intermediate array, because all components should be put together in the
    ! res array.
    ! The temporary array therefore needs to be sufficiently large to store the
    ! maximal number of components.
    maxComponents = maxval(varSys%method%val(varPos(:))%nComponents)

    ! Using a temporary array to store the variables and transfer them to res
    ! in the correct ordering afterwards.
    allocate(tmpdat(nPnts*maxComponents))

    compOff = 0
    vars: do iVar=1, nVars

      ! get the number of components for variable iVar
      nComponents = varSys%method%val(varPos(iVar))%nComponents

      ! get the size of the needed part of the res array
      res_size = nPnts * nComponents

      ! derive the quantities for all the elements in the current chunk
      call varSys%method%val(varpos(iVar))%get_point(               &
        &                                varSys = varSys,           &
        &                                point  = point,            &
        &                                time   = time,             &
        &                                tree   = tree,             &
        &                                nPnts  = nPnts,            &
        &                                res    = tmpdat(:res_size) )

      ! copy the information to the right positions in the result array
      ! res contains results for all variables,
      ! tmpdat is only for one variable
      do iPnt=1,nPnts
        e_start = (iPnt-1)*nScalars + compOff
        t_start = (iPnt-1)*nComponents
        res( (e_start+1) : (e_start+nComponents) )         &
          &  = tmpdat( t_start + 1 : t_start + nComponents )
      end do
      ! Increase the component offset for the next variables.
      compOff = compOff + nComponents
    end do vars

  end subroutine tem_get_point_chunk