tem_convergence_check_point Subroutine

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

This routine runs over convergence check using get_point 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_point~~CallsGraph proc~tem_convergence_check_point tem_convergence_check_point proc~tem_baryofid tem_BaryOfId proc~tem_convergence_check_point->proc~tem_baryofid proc~tem_convergence_evaluate tem_convergence_evaluate proc~tem_convergence_check_point->proc~tem_convergence_evaluate proc~tem_get_point_chunk tem_get_point_chunk proc~tem_convergence_check_point->proc~tem_get_point_chunk proc~tem_reduction_spatial_append tem_reduction_spatial_append proc~tem_convergence_check_point->proc~tem_reduction_spatial_append proc~tem_reduction_spatial_close tem_reduction_spatial_close proc~tem_convergence_check_point->proc~tem_reduction_spatial_close proc~tem_reduction_spatial_open tem_reduction_spatial_open proc~tem_convergence_check_point->proc~tem_reduction_spatial_open proc~tem_coordofid tem_CoordOfId proc~tem_baryofid->proc~tem_coordofid proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_baryofid->proc~tem_elemsizelevel 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_coordofid->proc~tem_levelof proc~tem_elemsize->proc~tem_elemsizelevel proc~tem_elemsize->proc~tem_levelof

Called by

proc~~tem_convergence_check_point~~CalledByGraph proc~tem_convergence_check_point tem_convergence_check_point proc~tem_convergence_check tem_convergence_check proc~tem_convergence_check->proc~tem_convergence_check_point

Source Code

  subroutine tem_convergence_check_point( 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, nPoints, nScalars, pointsOff, nChunkPoints
    integer :: iPoint, iChunk, counter
    integer :: buf_start, buf_end
    real(kind=rk), allocatable :: points(:,:)
    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
      nPoints = tree%nElems
    else
      nPoints = me%subTree%nPoints
    end if

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

    ! allocate points to size of chunkSize
    allocate(points(me%chunkSize,3))

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

      ! number of points written to THIS chunk
      nChunkPoints = min(me%chunkSize, nPoints-pointsOff)

      ! Compute the points lower and upper bound for the current chunk
      buf_start = pointsOff + 1
      buf_end = pointsOff + nChunkPoints

      if (me%subTree%useGlobalMesh) then
        counter = 0
        do iPoint = buf_start, buf_end
          counter = counter + 1
          points(counter, :) = tem_BaryOfId( tree,               &
            &                                tree%treeID(iPoint) )
        end do
      else
        points(1:nChunkPoints,:) = me%subTree%points(buf_start:buf_end,:)
      end if

      ! evaluate all variables on current chunk
      call tem_get_point_chunk(varSys  = varSys,                   &
        &                      varPos  = me%varMap%varPos          &
        &                                         %val(:nVars),    &
        &                      point   = points(1:nChunkPoints,:), &
        &                      time    = time,                     &
        &                      tree    = tree,                     &
        &                      nPnts   = nChunkPoints,             &
        &                      res     = res                       )

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

    end do ! iChunk

    deallocate(points)

    ! 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(10),*) 'Reached steady state ', time%iter, isConverged
      end if
    end if

  end subroutine tem_convergence_check_point