adjust_search_direction Subroutine

private subroutine adjust_search_direction(me, x, p, pnew, modified)

if necessary, adjust the search direction vector p so that bounds are not violated.

Reference

  • https://openmdao.org/newdocs/versions/latest/features/building_blocks/solvers/bounds_enforce.html

Type Bound

nlesolver_type

Arguments

Type IntentOptional Attributes Name
class(nlesolver_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(me%n) :: x

initial x. it is assumed that this does not violate any bounds.

real(kind=wp), intent(in), dimension(me%n) :: p

search direction p = xnew - x

real(kind=wp), intent(out), dimension(me%n) :: pnew

new search direction

logical, intent(out), dimension(me%n) :: modified

indicates the elements of p that were modified


Calls

proc~~adjust_search_direction~~CallsGraph proc~adjust_search_direction nlesolver_type%adjust_search_direction proc~int2str int2str proc~adjust_search_direction->proc~int2str proc~set_status nlesolver_type%set_status proc~adjust_search_direction->proc~set_status proc~set_status->proc~int2str proc~real2str real2str proc~set_status->proc~real2str

Called by

proc~~adjust_search_direction~~CalledByGraph proc~adjust_search_direction nlesolver_type%adjust_search_direction proc~backtracking_linesearch backtracking_linesearch proc~backtracking_linesearch->proc~adjust_search_direction proc~exact_linesearch exact_linesearch proc~exact_linesearch->proc~adjust_search_direction proc~fixed_point_linesearch fixed_point_linesearch proc~fixed_point_linesearch->proc~adjust_search_direction proc~simple_step simple_step proc~simple_step->proc~adjust_search_direction

Source Code

    subroutine adjust_search_direction(me,x,p,pnew,modified)

    implicit none

    class(nlesolver_type),intent(inout) :: me
    real(wp),dimension(me%n),intent(in) :: x  !! initial `x`. it is assumed that this does not violate any bounds.
    real(wp),dimension(me%n),intent(in) :: p  !! search direction `p = xnew - x`
    real(wp),dimension(me%n),intent(out) :: pnew  !! new search direction
    logical,dimension(me%n),intent(out) :: modified  !! indicates the elements of `p` that were modified

    integer :: i  !! counter
    real(wp),dimension(:),allocatable :: xnew  !! `x + pnew`
    logical :: search_direction_modifed !! if `p` has been modified
    real(wp) :: t !! indep var in the step equation: `xnew = x + t*p`

    if (me%bounds_mode==NLESOLVER_IGNORE_BOUNDS) then ! no change
        modified = .false.
        pnew = p
        return
    end if

    allocate(xnew(me%n))
    xnew = x + p  ! this is the full Newton step
    search_direction_modifed = .false.  ! will be set to true if any bounds are violated
    if (me%bounds_mode==NLESOLVER_VECTOR_BOUNDS) t = 1.0_wp ! for mode 2, start with full newton step
                                                            ! note: have to check each variable and
                                                            ! choose the smallest t. don't need to
                                                            ! keep track of xnew.
    modified = .false.

    do i = 1, me%n
        if (xnew(i)<me%xlow(i)) then
            search_direction_modifed = .true.
            modified(i) = .true.
            if (me%verbose) write(me%iunit, '(A)') 'x('//int2str(i)//') < xlow(i) : adjusting to lower bound'
            select case (me%bounds_mode)
            case(NLESOLVER_SCALAR_BOUNDS,NLESOLVER_WALL_BOUNDS); xnew(i) = me%xlow(i)
            case(NLESOLVER_VECTOR_BOUNDS);                       t = min(t,(me%xlow(i)-x(i))/p(i))
            end select
        else if (xnew(i)>me%xupp(i)) then
            search_direction_modifed = .true.
            modified(i) = .true.
            if (me%verbose) write(me%iunit, '(A)') 'x('//int2str(i)//') > xupp(i) : adjusting to upper bound'
            select case (me%bounds_mode)
            case(NLESOLVER_SCALAR_BOUNDS,NLESOLVER_WALL_BOUNDS); xnew(i) = me%xupp(i)
            case(NLESOLVER_VECTOR_BOUNDS);                       t = min(t,(me%xupp(i)-x(i))/p(i))
            end select
        end if
    end do

    if (search_direction_modifed) then
        ! adjust the search direction:
        select case (me%bounds_mode)
        case (NLESOLVER_SCALAR_BOUNDS,NLESOLVER_WALL_BOUNDS)
            pnew = xnew - x  ! here we have changed the search direction vector
            if (all(pnew==0.0_wp)) call me%set_status(istat = -17, string = 'Error adjusting line search direction for bounds')
        case (NLESOLVER_VECTOR_BOUNDS)
            ! here was are staying on the original search direction vector, just walking back
            if (t <= 0.0_wp) then ! something wrong
                call me%set_status(istat = -17, string = 'Error adjusting line search direction for bounds')
                pnew = p
            else
                pnew = p / t
            end if
        end select
        if (me%verbose) write(me%iunit, '(A)') 'Search direction modified to be within bounds'
    else
        pnew = p
    end if

    end subroutine adjust_search_direction