tem_evalMag_forElement Subroutine

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

Evaluate magnitude of any vectorial variable. In lua file, first define new variable with varType operation kind as "magnitude" and provide name of the variable from which magnitude to be derived in input_varname. If input_varname variable is not part of predefined solver variables then add also that variable via variable table.

\verbatim -- in lua file, one can define as following: variable = {{ name = 'velMag', ncomponents = 1, vartype = "operation", operation = {kind='magnitude',input_varname={'velocity'}} }, } tracking = { variable = {'velMag'}, folder = 'tracking/', shape = {kind = 'canoND', object = {origin = {3.0,3.1,3.0} } }, format = 'ascii', time = {min = 0, max = tmax, interval = 1}, } \endverbatim

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


Source Code

  recursive subroutine tem_evalMag_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 :: iElem, depVar_pos, depVar_nComps, offset, iDof
    real(kind=rk), allocatable :: input_varRes(:)
    ! ---------------------------------------------------------------------- !

    ! get the position of dependent variable, assuming only one
    depVar_pos = fun%input_varPos(1)

    ! get the nComp of dependent variable
    depVar_nComps = varSys%method%val( depVar_pos )%nComponents

    ! allocate array to save results
    allocate(input_varRes( nElems * nDofs * depVar_nComps ))

    ! derive dependent variable
    call varSys%method%val(depVar_pos)%get_element(                 &
      &                                   varSys  = varSys,         &
      &                                   elemPos = elemPos,        &
      &                                   time    = time,           &
      &                                   tree    = tree,           &
      &                                   nElems  = nElems,         &
      &                                   nDofs   = nDofs,          &
      &                                   res     = input_varRes(:) )

    ! compute magnitude
    ! assuming fun%nComponents = 1
    input_varRes = input_varRes * input_varRes
    do iElem = 1, nElems
      do iDof = 1, nDofs
        offset = (( ielem-1)* depvar_ncomps* ndofs+( idof-1)* depvar_ncomps+0)
        res((( ielem-1)* 1* ndofs+( idof-1)* 1+1) )                    &
          & = sqrt(sum(input_varRes(offset+1 : offset+depVar_nComps)))
      end do
    end do

    deallocate(input_varRes)

  end subroutine tem_evalMag_forElement