if necessary, adjust the search direction vector p so that bounds are not violated.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(nlesolver_type), | intent(inout) | :: | me | |||
| real(kind=wp), | intent(in), | dimension(me%n) | :: | x |
initial |
|
| real(kind=wp), | intent(in), | dimension(me%n) | :: | p |
search direction |
|
| real(kind=wp), | intent(out), | dimension(me%n) | :: | pnew |
new search direction |
|
| logical, | intent(out), | dimension(me%n) | :: | modified |
indicates the elements of |
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