append_da_real Subroutine

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

appending a value to the dynamic array

with this subroutine, a given value can be added to the dynamic array. the actual position of this value 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_realarray_type) :: me
real(kind=rk), intent(in) :: val
integer, intent(in), optional :: length

optional length to expand the array

integer, intent(out), optional :: pos

position in the array, if the value is found

logical, intent(out), optional :: wasadded

flag to indicate, if val was newly added


Calls

proc~~append_da_real~~CallsGraph proc~append_da_real append_da_real interface~expand~27 expand proc~append_da_real->interface~expand~27 interface~sortedposofval~5 sortedposofval proc~append_da_real->interface~sortedposofval~5 proc~expand_da_label expand_da_label interface~expand~27->proc~expand_da_label proc~sortposofval_label sortposofval_label interface~sortedposofval~5->proc~sortposofval_label

Called by

proc~~append_da_real~~CalledByGraph proc~append_da_real append_da_real interface~append~28 append interface~append~28->proc~append_da_real

Source Code

  subroutine append_da_real(me, val, length, pos, wasadded )
    !------------------------------------------------------------------------
    type(dyn_realarray_type) :: me   !< array to append the value to
    real(kind=rk), intent(in)     :: val  !< value to append
    !> optional length to expand the array
    integer, intent(in), optional :: length
    !> position in the array, if the value is found
    integer, intent(out), optional :: pos
    !> flag to indicate, if val was newly added
    logical, intent(out), optional :: wasadded
    !------------------------------------------------------------------------
    integer :: foundpos
    integer :: i
    !------------------------------------------------------------------------

    ! do a binary search on existing entries (returns closest entry next to
    ! it if not found).
    foundpos = sortedposofval(me, val, .true.)
    if( present( wasadded ) ) wasadded = .false.

    ! if it found the value, the position is smaller than nvals
    if (foundpos <= me%nvals) then

      ! the returned position might actually be the right entry already or
      ! not, check for this here.
      if ( me%val(me%sorted(foundpos)) == val ) then

        ! found the value in a list of unique values,
        ! nothing to do, just return its position.
        if( present( pos ) ) pos = me%sorted(foundpos)

      else

        ! need to append a new value!

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

        if( present( wasadded ) ) wasadded = .true.
        if (me%nvals == me%containersize) then

          ! container is full, need to expand it
          call expand(me = me, length = length)
        end if
        me%nvals = me%nvals + 1

        ! put the new value into the last position in the
        ! array.
        me%val(me%nvals) = val
        do while( foundpos < me%nvals )
          if(me%val(me%sorted(foundpos)) /= val) then
            exit
          end if
          ! in case of multiple entries with the same value
          ! move on to the first differing entry.
          foundpos = foundpos + 1
        end do
        ! shift the sorted list of indices, to create a
        ! whole for the value to be inserted, at position
        ! foundpos.
        do i=me%nvals-1,foundpos,-1
          me%sorted(i+1) = me%sorted(i)
        end do
        ! put the index of the new value into the
        ! sorted list at the now freed position.
        me%sorted(foundpos) = me%nvals
        if( present( pos ) ) pos = me%nvals

      end if

    else

      ! value to append is larger than all existing ones,
      ! just put it to the end of the list, this captures
      ! also the case of empty lists.
      ! in this case foundpos = me%nvals + 1 holds.
      if( present( wasadded ) ) wasadded = .true.
      if (foundpos > me%containersize) then
        ! expand the array, if its boundary is reached
        call expand(me = me, length = length)
      end if
      me%nvals = foundpos
      me%val(foundpos) = val
      me%sorted(foundpos) = foundpos
      if( present( pos ) ) pos = foundpos

    end if

  end subroutine append_da_real