tem_reduction_spatial_append Subroutine

public subroutine tem_reduction_spatial_append(me, chunk, nElems, varSys, varPos, tree, treeID, nDofs)

Local chunk-wise reduction

Chunk-wise local reduction, which has to be performed for each chunk until all elements have been treated within the reduction. After all chunks are treated (=_append has been finished), the global reduction has to be performed by tem_reduction_spatial_close. NOTE: Reduction is applied only on 1st DOF

Arguments

Type IntentOptional Attributes Name
type(tem_reduction_spatial_type), intent(inout) :: me(:)

The reduction type to work on. All definitions should be present in here

real(kind=rk), intent(in) :: chunk(:)

chunk of results to reduce

integer, intent(in) :: nElems

number of elements the chunk has

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

global variable system defined in solver

integer, intent(in) :: varPos(:)

position of variable to reduce in the global varSys

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

the global tree

integer(kind=long_k), intent(in), optional :: treeID(nElems)

The list of treeIDs of the current chunk

integer, intent(in), optional :: nDofs

Number of degrees of freedom.


Calls

proc~~tem_reduction_spatial_append~~CallsGraph proc~tem_reduction_spatial_append tem_reduction_spatial_append proc~tem_elemsize tem_ElemSize proc~tem_reduction_spatial_append->proc~tem_elemsize proc~tem_levelof tem_LevelOf proc~tem_reduction_spatial_append->proc~tem_levelof proc~tem_elemsize->proc~tem_levelof proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_elemsize->proc~tem_elemsizelevel

Called by

proc~~tem_reduction_spatial_append~~CalledByGraph proc~tem_reduction_spatial_append tem_reduction_spatial_append proc~hvs_ascii_dump_elem_data hvs_ascii_dump_elem_data proc~hvs_ascii_dump_elem_data->proc~tem_reduction_spatial_append proc~hvs_ascii_dump_point_data hvs_ascii_dump_point_data proc~hvs_ascii_dump_point_data->proc~tem_reduction_spatial_append proc~tem_convergence_check_element tem_convergence_check_element proc~tem_convergence_check_element->proc~tem_reduction_spatial_append proc~tem_convergence_check_point tem_convergence_check_point proc~tem_convergence_check_point->proc~tem_reduction_spatial_append proc~hvs_output_write hvs_output_write proc~hvs_output_write->proc~hvs_ascii_dump_elem_data proc~hvs_output_write->proc~hvs_ascii_dump_point_data proc~tem_convergence_check tem_convergence_check proc~tem_convergence_check->proc~tem_convergence_check_element proc~tem_convergence_check->proc~tem_convergence_check_point proc~tem_tracker tem_tracker proc~tem_tracker->proc~hvs_output_write

Source Code

  subroutine tem_reduction_spatial_append(me, chunk, nElems, varSys, varPos, &
    &                                     tree, treeID, nDofs)
    ! ---------------------------------------------------------------------------
    !> The reduction type to work on. All definitions should be
    !! present in here
    type( tem_reduction_spatial_type ), intent(inout) :: me(:)
    !> number of elements the chunk has
    integer, intent(in) :: nElems
    !> global variable system defined in solver
    type( tem_varSys_type ), intent(in)          :: varSys
    !> position of variable to reduce in the global varSys
    integer, intent(in)                          :: varPos(:)
    !> chunk of results to reduce
    real(kind=rk), intent(in) :: chunk(:)
    !> Number of degrees of freedom.
    integer, optional, intent(in) :: nDofs
    !> the global tree
    type(treelmesh_type), intent(in) :: tree
    !> The list of treeIDs of the current chunk
    integer(kind=long_k), optional, intent(in) :: treeID( nElems )
    ! ---------------------------------------------------------------------------
    real(kind=rk) :: temp
    real(kind=rk) :: dx, V, Vloc
    integer :: iVar, iComp, iPos, chunkVarPos, iElem, nDofs_L
    integer :: nPerElem
    integer :: minLevel, myLevel, volumeFac, nScalars
    ! ---------------------------------------------------------------------------
    minLevel = tree%global%minLevel
    nScalars = sum(varSys%method%val(varPos(:))%nComponents)

    ! spatial reduction is applied only on 1st Degrees of freedom
    if (present(nDofs)) then
      nDofs_L = nDofs
    else
      nDofs_L = 1
    end if
    nPerElem = nDofs_L*nScalars

    chunkVarPos = 0
    ! Run over all variables in the chunk
    do iVar = 1, size( me )
      me( iVar )%nElems = me( iVar )%nElems + nElems
      ! count up the me%nElems here as well
      ! update the me%val
      do iComp = 1, me( iVar )%nComponents
        chunkVarPos = chunkVarPos + 1
        temp = 0._rk
        select case( trim(me( iVar )%reduceType) )
          ! Upon extending the reductions, make sure to add them in the
          ! tem_load_spatial_reduction  routine
          case('sum', 'average')
            ! Take the values corresponding to the current component of the
            ! variable with step (nScalars)
            temp = sum( chunk( chunkVarPos : nElems*nPerElem &
              &                            : nPerElem )      )
            me( iVar )%val( iComp ) = me( iVar )%val( iComp ) + temp

          case('l2normalized')
            Vloc = 0.0_rk
            ! Sum up the squares of elemets of chunk into a temp variable
            ! Total error L2 on the domain normalize by Volum or number of points
            if (present(treeID)) then
              do iPos = chunkVarPos, nElems*nPerElem, nPerElem
                iElem  = ((iPos-1)/nPerElem) + 1
                dx = tem_Elemsize(tree, treeID(iElem))
                V = dx**3
                Vloc = Vloc + V

                temp = temp + chunk(iPos)*chunk(iPos) * V ! real(volumeFac,kind=rk)
              end do
            else
              do iPos = chunkVarPos, nElems*nPerElem, nPerElem
                temp = temp + chunk(iPos)*chunk(iPos)
                Vloc = Vloc + 1
              end do
            end if
            me( iVar )%val( iComp ) = me( iVar)%val( iComp ) + temp
            me( iVar )%Vloc = Vloc

          case('l2norm')
            ! Sum up the squares of elemets of chunk into a temp variable
            ! Total error L2 on the domain Omega is with the error function u
            ! $L2(\Omega) =  \sqrt(\int_\Omega |u|^2 ))
            !             =  \sqrt(\sum_i \int_\Omega_i |u|^2 ))
            !             =  \sqrt(\sum_i \Omega_i \cdot |u|^2 ))
            ! Value of elements on level other than minLevel is scale by a
            ! factor of 8^{myLevel-minLevel} to obtain a consistant results for
            ! multi-level.
            if (present(treeID)) then
              do iPos = chunkVarPos, nElems*nPerElem, nPerElem
                iElem  = ((iPos-1)/nPerElem) + 1
                ! get level, then get scale factor
                !myLevel = tem_levelOf( treeID(iElem) )
                !volumeFac = 8**( myLevel - minLevel )
                dx = tem_Elemsize(tree, treeID(iElem))
                V = dx**3

                temp = temp + chunk(iPos)*chunk(iPos) * V ! real(volumeFac,kind=rk)
              end do
            else
              do iPos = chunkVarPos, nElems*nPerElem, nPerElem
                temp = temp + chunk(iPos)*chunk(iPos)
              end do
            end if
            me( iVar )%val( iComp ) = me( iVar)%val( iComp ) + temp

          case('linfnorm')
            ! L_inf norm evaluates the maximum of the absolute occurring value
            me( iVar )%val( iComp ) = max( me( iVar )%val( iComp ),   &
              & maxval( abs(chunk( chunkVarPos : nElems*nPerElem      &
              &                                : nPerElem         ) )))

          case('max')
            ! maximum between previous value and current
            me( iVar )%val( iComp ) = max(me( iVar )%val( iComp ), &
              & maxval( chunk( chunkVarPos : nElems*nPerElem       &
              &                            : nScalars*nDofs_L ) )  )

          case('min')
            ! minimum between previous value and current
            me( iVar )%val( iComp ) = min(me(iVar)%val(iComp),  &
              & minval( chunk( chunkVarPos : nElems*nPerElem    &
              &                            : nPerElem         )))

          case('weighted_sum')
            ! Sum up scaled results in chunk into a temp variable
            ! Value of elements on level other than minLevel is scale by a
            ! factor of 8^{myLevel-minLevel} to obtain a consistant results for
            ! multi-level.
            if (present(treeID)) then
              do iPos = chunkVarPos, nElems*nPerElem, nPerElem
                iElem  = ((iPos-1)/nPerElem) + 1
                ! get level, then get scale factor
                myLevel = tem_levelOf( treeID(iElem) )
                volumeFac = 8**( myLevel - minLevel )
                temp = temp + chunk(iPos) / real(volumeFac,rk)
              end do
            else
              temp = sum( chunk( chunkVarPos : nElems*nPerElem &
                &                            : nPerElem )      )
            end if
            me( iVar )%val( iComp ) = me( iVar)%val( iComp ) + temp

          case('none')
            me( iVar )%val( iComp ) = 0._rk

        end select
      enddo
    enddo

  end subroutine tem_reduction_spatial_append