append_da_veclabel Subroutine

public subroutine append_da_veclabel(me, val, length, pos, wasadded)

appending a sorted list of values to the dynamic array

with this subroutine, a given list of sorted values can be added to the dynamic array. the actual positions of these values in the dynamic array will be returned, so it can be found again easily later. with the wasadded flag, it is indicated,\n wasadded = true, if this entry had to be added,\n wasadded = false, if this was already found in the array.

Arguments

Type IntentOptional Attributes Name
type(dyn_labelarray_type) :: me
character(len=*), intent(in) :: val(:)
integer, intent(in), optional :: length

optional length to expand the array

integer, intent(out), optional :: pos(:)

position in the array, the values are found at.

logical, intent(out), optional :: wasadded(:)

flag to indicate, if val was newly added


Calls

proc~~append_da_veclabel~~CallsGraph proc~append_da_veclabel append_da_veclabel interface~expand~27 expand proc~append_da_veclabel->interface~expand~27 proc~expand_da_label expand_da_label interface~expand~27->proc~expand_da_label

Called by

proc~~append_da_veclabel~~CalledByGraph proc~append_da_veclabel append_da_veclabel interface~append~29 append interface~append~29->proc~append_da_veclabel proc~append_pointdata append_pointData proc~append_pointdata->interface~append~29 proc~append_singlepnt2grwpoints append_singlePnt2grwPoints proc~append_singlepnt2grwpoints->interface~append~29 proc~append_vectorpnt2grwpoints append_vectorPnt2grwPoints proc~append_vectorpnt2grwpoints->interface~append~29 proc~setup_indices_spacetime setup_indices_spacetime proc~setup_indices_spacetime->interface~append~29 proc~sorttruncate_da_int sorttruncate_da_int proc~sorttruncate_da_int->interface~append~29 proc~sorttruncate_da_label sorttruncate_da_label proc~sorttruncate_da_label->interface~append~29 proc~sorttruncate_da_long sorttruncate_da_long proc~sorttruncate_da_long->interface~append~29 proc~sorttruncate_da_real sorttruncate_da_real proc~sorttruncate_da_real->interface~append~29 proc~tem_addtimer tem_addTimer proc~tem_addtimer->interface~append~29 proc~tem_appendtimers tem_appendTimers proc~tem_appendtimers->interface~append~29 proc~tem_cano_checkneigh tem_cano_checkNeigh proc~tem_cano_checkneigh->interface~append~29 proc~tem_cano_initsubtree tem_cano_initSubTree proc~tem_cano_initsubtree->interface~append~29 proc~tem_cano_storepntsinsubtree tem_cano_storePntsInSubTree proc~tem_cano_storepntsinsubtree->interface~append~29 proc~tem_comm_createbuffer tem_comm_createBuffer proc~tem_comm_createbuffer->interface~append~29 proc~tem_commbuf_int_gatherindexed tem_commbuf_int_gatherindexed proc~tem_commbuf_int_gatherindexed->interface~append~29 proc~tem_commbuf_long_gatherindexed tem_commbuf_long_gatherindexed proc~tem_commbuf_long_gatherindexed->interface~append~29 proc~tem_commbuf_real_gatherindexed tem_commbuf_real_gatherindexed proc~tem_commbuf_real_gatherindexed->interface~append~29 proc~tem_shape_findelembybclabels tem_shape_findElemByBCLabels proc~tem_shape_findelembybclabels->interface~append~29 proc~tem_shape_initbylevels tem_shape_initByLevels proc~tem_shape_initbylevels->interface~append~29 proc~tem_shape_initpropelements tem_shape_initPropElements proc~tem_shape_initpropelements->interface~append~29 proc~tem_shape_subtreefromgeominters tem_shape_subTreeFromGeomInters proc~tem_shape_subtreefromgeominters->interface~append~29 proc~tem_timer_loadconfig tem_timer_loadconfig proc~tem_timer_loadconfig->interface~append~29 proc~tem_varsys_append_auxfieldvar tem_varSys_append_auxFieldVar proc~tem_varsys_append_auxfieldvar->interface~append~29 proc~tem_varsys_append_dervar tem_varSys_append_derVar proc~tem_varsys_append_dervar->interface~append~29 proc~tem_varsys_append_statevar tem_varSys_append_stateVar proc~tem_varsys_append_statevar->interface~append~29 proc~tem_varsys_append_stfun_raw tem_varSys_append_stFun_raw proc~tem_varsys_append_stfun_raw->interface~append~29 proc~tem_varsys_append_stfunvar tem_varSys_append_stFunVar proc~tem_varsys_append_stfunvar->interface~append~29 proc~tem_varsys_load_single tem_varSys_load_single proc~tem_varsys_load_single->interface~append~29

Source Code

  subroutine append_da_veclabel(me, val, length, pos, wasadded )
    !------------------------------------------------------------------------
    type(dyn_labelarray_type) :: me   !< array to append the value to
    character(len=*), intent(in)     :: val(:)  !< values to append
    !> optional length to expand the array
    integer, intent(in), optional :: length
    !> position in the array, the values are found at.
    integer, intent(out), optional :: pos(:)
    !> flag to indicate, if val was newly added
    logical, intent(out), optional :: wasadded(:)
    !------------------------------------------------------------------------
    character(len=labellen) :: lastval
    logical :: addedval(size(val))
    integer :: i
    integer :: veclen
    integer :: maxlen
    integer :: nappend
    integer :: rem_app
    integer :: curval, ival, iold, iadd
    integer, allocatable :: newsorted(:)
    !------------------------------------------------------------------------

    if (size(val) == 0) return

    veclen = size(val)
    maxlen = veclen + me%nvals

    allocate(newsorted(maxlen))

    addedval = .false.

    iold = 1
    iadd = 1

    nappend = 0
    curval = 0

    ! select the first entry before the loop unconditionally without checks
    ! for uniqueness (nothing to check against yet).
    if ( me%val(me%sorted(iold)) <= val(iadd) ) then
      curval = curval + 1
      newsorted(curval) = me%sorted(iold)
      lastval = me%val(me%sorted(iold))
      iold = iold + 1
    else
      curval = curval + 1
      nappend = nappend + 1
      newsorted(curval) = me%nvals + nappend
      lastval = val(iadd)
      if (present(pos))  pos(iadd) = newsorted(curval)
      addedval(iadd) = .true.
      iadd = iadd + 1
    end if

    do ival=2,maxlen

      if ( (iadd <= veclen) .and. (iold <= me%nvals) ) then

        if ( me%val(me%sorted(iold)) <= val(iadd) ) then

          ! the original list's values are appended to newsorted before
          ! the additional list is appended.
          curval = curval + 1
          newsorted(curval) = me%sorted(iold)
          lastval = me%val(me%sorted(iold))
          iold = iold + 1

        else

          ! only append the value to unique lists, if it is not yet in the list.
          ! (if it is already in the list, it has to be the previous (curval-1)
          !  entry.)
          if ( lastval < val(iadd) ) then
            nappend = nappend + 1
            curval = curval + 1
            newsorted(curval) = me%nvals + nappend
            lastval = val(iadd)
            addedval(iadd) = .true.
          end if
          if (present(pos)) pos(iadd) = newsorted(curval)
          iadd = iadd + 1

        end if

      else

        ! reached the end of one or both of the sorted lists.
        exit

      end if

    end do

    if (iold <= me%nvals) then
      ! still some values from the original list left.
      newsorted(curval+1:me%nvals+nappend) = me%sorted(iold:me%nvals)
    end if

    if (iadd <= veclen) then
      ! still some values from the list to append left.
      rem_app = iadd
      do i = rem_app,veclen
        if ( lastval < val(iadd) ) then
          nappend = nappend + 1
          curval = curval + 1
          newsorted(curval) = me%nvals + nappend
          lastval = val(iadd)
          addedval(iadd) = .true.
        end if
        if (present(pos)) pos(iadd) = newsorted(curval)
        iadd = iadd + 1
      end do
    end if

    if (me%nvals > huge(me%nvals)-nappend) then
       write(*,*) "reached end of integer range for dynamic array!"
       write(*,*) "aborting!!"
       stop
    end if

    if (me%nvals + nappend > me%containersize) then
      call expand( me        = me,      &
        &          increment = nappend, &
        &          length    = length   )
    end if
    me%sorted(:me%nvals+nappend) = newsorted(:me%nvals+nappend)
    curval = me%nvals
    do iadd=1,veclen
      if (addedval(iadd)) then
        curval = curval + 1
        me%val(curval) = val(iadd)
      end if
    end do
    me%nvals = me%nvals + nappend

    if( present( wasadded ) ) wasadded = addedval

  end subroutine append_da_veclabel