evalLogicalLess_forElement Subroutine

public subroutine evalLogicalLess_forElement(fun, varSys, elempos, time, tree, nElems, nDofs, res)

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~~evallogicalless_forelement~~CallsGraph proc~evallogicalless_forelement evalLogicalLess_forElement proc~logicaltoreal logicalToReal proc~evallogicalless_forelement->proc~logicaltoreal

Source Code

  subroutine evalLogicalLess_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, iComp, iElem, lIdx, uIdx
    ! first dimension:nElems,
    ! second dimension:size of difference variable in input_varname
    real(kind=rk), allocatable :: input_varRes(:,:)
    !--------------------------------------------------------------------------!


    ! nInputs must be two
    allocate(input_varRes(nElems*nDofs*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_element(                        &
        &                                   varSys  = varSys,               &
        &                                   elemPos = elemPos,              &
        &                                   time    = time,                 &
        &                                   tree    = tree,                 &
        &                                   nElems  = nElems,               &
        &                                   nDofs   = nDofs,                &
        &                                   res     = input_varRes(:, iDep) )

    end do

    do iComp = 1, fun%nComponents
      do iElem = 1, nElems
        ! We only take the first mode into consideration, but store the result
        ! value for all modes.
        lIdx = (( ielem-1)* fun%ncomponents* ndofs+( 1-1)* fun%ncomponents+icomp)
        uIdx = (( ielem-1)* fun%ncomponents* ndofs+( ndofs-1)* fun%ncomponents+icomp)

        res( lIdx : uIdx ) = logicalToReal(                       &
          & input_varRes( lIdx, 1 ) .flt. input_varRes( lIdx, 2 ) )
      end do
    end do

    deallocate(input_varRes)

  end subroutine evalLogicalLess_forElement