evalDiff_forPoint Subroutine

private recursive subroutine evalDiff_forPoint(fun, varSys, point, time, tree, nPnts, res)

Same as evalDiff_forElement except it evaluate diff from points

The interface has to comply to the abstract interface tem_varSys_module#tem_varSys_proc_point

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(:,:)

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

Point in time at which to evaluate the variable.

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

global treelm mesh info

integer, intent(in) :: nPnts

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


Contents

Source Code


Source Code

  recursive subroutine evalDiff_forPoint( fun, varsys, point, time, tree, &
    &                                     nPnts, 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

    !> 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(:,:)

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

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

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

    !> 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 :: iDep, posDepVar
    integer :: nTotal
    ! first dimension:nElems,
    ! second dimension:size of difference variable in input_varname
    real(kind=rk), allocatable :: input_varRes(:,:)
    ! ---------------------------------------------------------------------- !

    nTotal = nPnts*fun%nComponents

    ! nInputs must be two
    allocate(input_varRes(nPnts*fun%nComponents, fun%nInputs))

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

      ! derive dependent variable
      call varSys%method%val(posDepVar)%get_point(                         &
        &                                   varSys = varSys,               &
        &                                   point  = point,                &
        &                                   time   = time,                 &
        &                                   tree   = tree,                 &
        &                                   nPnts  = nPnts,                &
        &                                   res    = input_varRes(:, iDep) )

    end do

    res(:nTotal) = input_varRes(:, 1) - input_varRes(:, 2)

    deallocate(input_varRes)

  end subroutine evalDiff_forPoint