tem_opVar_reduction_transient_update Subroutine

public subroutine tem_opVar_reduction_transient_update(redTransVarPos, varSys, tree, time)

Update all time reduction operation variables for entire domain

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: redTransVarPos(:)

Position of time reduction variables in varSys

type(tem_varSys_type), intent(in) :: varSys

Global Variable system

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

treelmesh_type

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

Current time


Calls

proc~~tem_opvar_reduction_transient_update~~CallsGraph proc~tem_opvar_reduction_transient_update tem_opVar_reduction_transient_update proc~tem_reduction_transient_update tem_reduction_transient_update proc~tem_opvar_reduction_transient_update->proc~tem_reduction_transient_update proc~tem_reduction_transient_reset tem_reduction_transient_reset proc~tem_reduction_transient_update->proc~tem_reduction_transient_reset

Contents


Source Code

  subroutine tem_opVar_reduction_transient_update(redTransVarPos, varSys, tree,&
    &                                             time)
    ! ---------------------------------------------------------------------- !
    !> Position of time reduction variables in varSys
    integer, intent(in) :: redTransVarPos(:)
    !> Global Variable system
    type(tem_varSys_type), intent(in) :: varSys
    !> treelmesh_type
    type(treelmesh_type), intent(in) :: tree
    !> Current time
    type(tem_time_type), intent(in) :: time
    ! ---------------------------------------------------------------------- !
    integer :: elemPos(tree%nElems)
    real(kind=rk), allocatable :: input_varRes(:)
    integer :: iVar, iElem, posDepVar, nCompMax, idxMax, nDofs
    type(tem_varSys_op_data_type), pointer :: opData
    ! ---------------------------------------------------------------------- !

    ! Only need to do anything here if at least one variable with reduction
    ! is to be computed.
    if (size(redTransVarPos) < 1) RETURN

    ! nDofs of solver are stored in opData%redTrans which is same for all
    ! variables so get from any reduction_transient variable
    call C_F_Pointer(varSys%method%val(redTransVarPos(1))%method_data, &
      &              opData                                            )
    nDofs = opData%redTrans%nDofs

    elemPos(1:tree%nElems) = (/ (iElem, iElem=1, tree%nElems) /)

    nCompMax = maxval(varSys%method%val(redTransVarPos(:))%nComponents)
    allocate(input_varRes(tree%nElems*nCompMax*nDofs))

    do iVar = 1, size(redTransVarPos)
      call C_F_Pointer(varSys%method%val(redTransVarPos(iVar))%method_data, &
        &              opData                                               )

      posDepVar = varSys%method%val(redTransVarPos(iVar))%input_varPos(1)
      call varSys%method%val(posDepVar)%get_element(                        &
        &                                   varSys  = varSys,               &
        &                                   elemPos = elemPos,              &
        &                                   time    = time,                 &
        &                                   tree    = tree,                 &
        &                                   nElems  = tree%nElems,          &
        &                                   nDofs   = nDofs,                &
        &                                   res     = input_varRes(:)       )

      idxMax = opData%redTrans%nEntries

      call tem_reduction_transient_update(me     = opData%redTrans,       &
        &                                 res    = input_varRes(1:idxMax) )
    end do
    deallocate(input_varRes)
  end subroutine tem_opVar_reduction_transient_update