cmerge Subroutine

private subroutine cmerge(n, a, i, m, h)

given the upper convex hulls of two consecutive sets of pairs (j,a(j)), compute the upper convex hull of their union

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

length of the vector a

real(kind=wp), intent(in) :: a(n)

vector defining the points (j,a(j))

integer, intent(in) :: i

abscissa of the common vertex of the two sets

integer, intent(in) :: m

the number of elements of each set is m+1

logical, intent(out) :: h(n)

vector defining the vertices of the convex hull, i.e., h(j) is .true. if (j,a(j)) is a vertex of the convex hull this vector is used also as output.


Calls

proc~~cmerge~~CallsGraph proc~cmerge polyroots_module::cmerge proc~ctest polyroots_module::ctest proc~cmerge->proc~ctest proc~left polyroots_module::left proc~cmerge->proc~left proc~right polyroots_module::right proc~cmerge->proc~right

Called by

proc~~cmerge~~CalledByGraph proc~cmerge polyroots_module::cmerge proc~cnvex polyroots_module::cnvex proc~cnvex->proc~cmerge proc~start polyroots_module::start proc~start->proc~cnvex proc~polzeros polyroots_module::polzeros proc~polzeros->proc~start

Source Code

    subroutine cmerge(n, a, i, m, h)

        !! given the upper convex hulls of two consecutive sets of pairs
        !! (j,a(j)), compute the upper convex hull of their union

        implicit none

        integer,intent(in) :: n !! length of the vector a
        real(wp),intent(in) :: a(n) !! vector defining the points (j,a(j))
        integer,intent(in) :: i !! abscissa of the common vertex of the two sets
        integer,intent(in) :: m !! the number of elements of each set is m+1
        logical,intent(out) :: h(n) !! vector defining the vertices of the convex hull, i.e.,
                                    !! h(j) is .true. if (j,a(j)) is a vertex of the convex hull
                                    !! this vector is used also as output.

        integer :: ir, il, irr, ill
        logical :: tstl, tstr

        ! at the left and the right of the common vertex (i,a(i)) determine
        ! the abscissae il,ir, of the closest vertices of the upper convex
        ! hull of the left and right sets, respectively
        call left(n, h, i, il)
        call right(n, h, i, ir)
        ! check the convexity of the angle formed by il,i,ir
        if (ctest(n, a, il, i, ir)) then
            return
        else
            ! continue the search of a pair of vertices in the left and right
            ! sets which yield the upper convex hull
            h(i) = .false.
            do
                if (il == (i - m)) then
                    tstl = .true.
                else
                    call left(n, h, il, ill)
                    tstl = ctest(n, a, ill, il, ir)
                end if
                if (ir == min(n, i + m)) then
                    tstr = .true.
                else
                    call right(n, h, ir, irr)
                    tstr = ctest(n, a, il, ir, irr)
                end if
                h(il) = tstl
                h(ir) = tstr
                if (tstl .and. tstr) return
                if (.not. tstl) il = ill
                if (.not. tstr) ir = irr
            end do
        end if

    end subroutine cmerge