evalRelDiff_forElement Subroutine

private recursive subroutine evalRelDiff_forElement(fun, varSys, elempos, time, tree, nElems, nDofs, res)

Same as evalDiff but output of 2nd dependent variable is used to compute relative difference.

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

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.

integer, intent(in) :: elempos(:)

Position of the TreeID of the element to get the variable for in the global treeID list.

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) :: nElems

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

integer, intent(in) :: nDofs

Number of degrees of freedom within an element.

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

Resulting values for the requested variable.

Linearized array dimension: (n requested entries) x (nComponents of this variable) x (nDegrees of freedom) Access: (iElem-1)fun%nComponentsnDofs + (iDof-1)*fun%nComponents + iComp


Calls

proc~~evalreldiff_forelement~~CallsGraph proc~evalreldiff_forelement evalRelDiff_forElement proc~tem_abort tem_abort proc~evalreldiff_forelement->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Contents


Source Code

  recursive subroutine evalRelDiff_forElement( fun, varsys, elempos, time, &
    &                                          tree, nElems, nDofs, 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

    !> Position of the TreeID of the element to get the variable for in the
    !! global treeID list.
    integer, intent(in) :: elempos(:)

    !> 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) :: nElems

    !> Number of degrees of freedom within an element.
    integer, intent(in) :: nDofs

    !> Resulting values for the requested variable.
    !!
    !! Linearized array dimension:
    !! (n requested entries) x (nComponents of this variable)
    !! x (nDegrees of freedom)
    !! Access: (iElem-1)*fun%nComponents*nDofs +
    !!         (iDof-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 = nElems*nDofs*fun%nComponents
    ! nInputs must be two
    allocate(input_varRes(nTotal, fun%nInputs))

    if (nDofs > 1) then
      write(*,*) 'TODO: evalreldiff does not work for polynomial data yet'
      write(*,*) '      It makes use of a division, which can not directly'
      write(*,*) '      be done in modal space!'
      write(*,*) ''
      write(*,*) 'Need to replace this routine in Ateles!'
      write(*,*) 'Stopping now.'
      call tem_abort()
    end if

    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_element(                        &
        &                                   varSys  = varSys,               &
        &                                   elemPos = elemPos,              &
        &                                   time    = time,                 &
        &                                   tree    = tree,                 &
        &                                   nElems  = nElems,               &
        &                                   nDofs   = nDofs,                &
        &                                   res     = input_varRes(:, iDep) )

    end do

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

    deallocate(input_varRes)

  end subroutine evalRelDiff_forElement