tem_convergence_check_element Subroutine

private subroutine tem_convergence_check_element(me, time, status, varSys, tree, res)

This routine runs over convergence check using get_element interface

Arguments

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

convergence descriptions

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

current simulation time

type(tem_status_type), intent(inout) :: status

Status bits

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

global variable system

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

global tree

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

Output data size: io_buffer_size


Calls

proc~~tem_convergence_check_element~~CallsGraph proc~tem_convergence_check_element tem_convergence_check_element proc~tem_convergence_evaluate tem_convergence_evaluate proc~tem_convergence_check_element->proc~tem_convergence_evaluate proc~tem_get_element_chunk tem_get_element_chunk proc~tem_convergence_check_element->proc~tem_get_element_chunk proc~tem_reduction_spatial_append tem_reduction_spatial_append proc~tem_convergence_check_element->proc~tem_reduction_spatial_append proc~tem_reduction_spatial_close tem_reduction_spatial_close proc~tem_convergence_check_element->proc~tem_reduction_spatial_close proc~tem_reduction_spatial_open tem_reduction_spatial_open proc~tem_convergence_check_element->proc~tem_reduction_spatial_open proc~evaluate_residual evaluate_residual proc~tem_convergence_evaluate->proc~evaluate_residual proc~tem_comparator tem_comparator proc~tem_convergence_evaluate->proc~tem_comparator 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 mpi_reduce mpi_reduce proc~tem_reduction_spatial_close->mpi_reduce proc~tem_elemsize->proc~tem_levelof proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_elemsize->proc~tem_elemsizelevel

Called by

proc~~tem_convergence_check_element~~CalledByGraph proc~tem_convergence_check_element tem_convergence_check_element proc~tem_convergence_check tem_convergence_check proc~tem_convergence_check->proc~tem_convergence_check_element

Source Code

  subroutine tem_convergence_check_element( me, time, status, varSys, tree, &
    &                                       res )
    ! -------------------------------------------------------------------- !
    !> convergence descriptions
    type(tem_convergence_type), target, intent(inout)   :: me
    !> current simulation time
    type(tem_time_type), intent(in)     :: time
    !> Status bits
    type(tem_status_type), intent(inout) :: status
    !> global variable system
    type(tem_varSys_type), intent(in)         :: varSys
    !> global tree
    type(treelmesh_type ), intent(in)         :: tree
    !> Output data
    !! size: io_buffer_size
    real(kind=rk), intent(out) :: res(:)
    ! -------------------------------------------------------------------- !
    integer :: nVars, nElems, nScalars, elemOff, nChunkElems
    integer :: iElem, iChunk
    integer :: buf_start, buf_end
    integer, allocatable :: elemPos(:)
    logical :: isConverged
    ! -------------------------------------------------------------------- !
    ! number of variables in varMap
    nVars = me%varMap%varPos%nVals

    ! Number of scalars in current output
    nScalars = me%varMap%nScalars

    if (me%subTree%useGlobalMesh) then
      nElems = tree%nElems
    else
      nElems = me%subTree%nElems
    end if

    ! open spatial reduction
    call tem_reduction_spatial_open( me     = me%redSpatial,               &
      &                              varSys = varSys,                      &
      &                              varPos = me%varMap%varPos%val(:nVars) )

    ! allocate elemPos to size of chunkSize
    allocate(elemPos(me%chunkSize))

    ! Process all chunks to derive the quantities defined in the tracking object
    do iChunk = 1, me%nChunks
      ! Number of elements read so far in previous chunks.
      elemOff = ((iChunk-1)*me%chunkSize)

      ! number of elements written to THIS chunk
      nChunkElems = min(me%chunkSize, nElems-elemOff)

      ! Compute the element lower and upper bound for the current chunk
      buf_start = elemOff + 1
      buf_end = elemOff + nChunkElems

      if (me%subTree%useGlobalMesh) then
        elemPos(1:nChunkElems) = (/ (iElem, iElem=buf_start, buf_end) /)
      else
        elemPos(1:nChunkElems) = me%subTree                     &
          &                        %map2Global(buf_start:buf_end)
      end if

      ! evaluate all variables on current chunk
      call tem_get_element_chunk(varSys  = varSys,                 &
        &                        varPos  = me%varMap%varPos        &
        &                                           %val(:nVars),  &
        &                        elemPos = elemPos(1:nChunkElems), &
        &                        time    = time,                   &
        &                        tree    = tree,                   &
        &                        nElems  = nChunkElems,            &
        &                        nDofs   = me%nDofs,               &
        &                        res     = res                     )

      ! preform spatial reduction
      call tem_reduction_spatial_append( me     = me%redSpatial,            &
        &                                chunk  = res,                      &
        &                                nElems = nChunkElems,              &
        &                                treeID = tree%treeID(              &
        &                                         elemPos(1:nChunkElems) ), &
        &                                tree   = tree,                     &
        &                                varSys = varSys,                   &
        &                                nDofs  = me%nDofs,                 &
        &                                varPos = me%varMap%varPos          &
        &                                                  %val(:nVars)     )

    end do ! iChunk

    deallocate(elemPos)

    ! Close reduction for current convergence
    call tem_reduction_spatial_close( me   = me%redSpatial, &
      &                               proc = me%proc        )

    ! evaluate residual and check convergence for each scalar
    ! and set steady state flag only in root process of convergence
    ! communicator
    if (me%proc%rank == 0) then
      ! Now check for convergence for each reduced scalars
      call tem_convergence_evaluate( me       = me,         &
        &                            achieved = isConverged )

      ! set steady state if any convergence object is converged
      status%bits(tem_stat_steady_state) = status%bits(tem_stat_steady_state)&
        &                                  .or. isConverged
      if (isConverged) then
        write(logUnit(5),*) 'Reached steady state ', time%iter, isConverged
      end if
    end if

  end subroutine tem_convergence_check_element