An exact linesearch that uses a derivative-free minimizer to
find the minimum value of f(x) between
x = xold + p * alpha_min and
x = xold + p * alpha_max.
Usually this is overkill and not necessary, but is here as an option for testing.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(nlesolver_type), | intent(inout) | :: | me | |||
| real(kind=wp), | intent(in), | dimension(me%n) | :: | xold |
previous value of |
|
| real(kind=wp), | intent(in), | dimension(me%n) | :: | p |
search direction |
|
| real(kind=wp), | intent(out), | dimension(me%n) | :: | x |
new |
|
| real(kind=wp), | intent(inout) | :: | f |
magnitude of |
||
| real(kind=wp), | intent(inout), | dimension(me%m) | :: | fvec |
function vector |
|
| real(kind=wp), | intent(in), | optional, | dimension(:,:) | :: | fjac |
jacobian matrix [dense] |
| real(kind=wp), | intent(in), | optional, | dimension(:) | :: | fjac_sparse |
jacobian matrix [sparse] |
subroutine exact_linesearch(me,xold,p,x,f,fvec,fjac,fjac_sparse) implicit none class(nlesolver_type),intent(inout) :: me real(wp),dimension(me%n),intent(in) :: xold !! previous value of `x` real(wp),dimension(me%n),intent(in) :: p !! search direction real(wp),dimension(me%n),intent(out) :: x !! new `x` real(wp),intent(inout) :: f !! magnitude of `fvec` real(wp),dimension(me%m),intent(inout) :: fvec !! function vector real(wp),dimension(:,:),intent(in),optional :: fjac !! jacobian matrix [dense] real(wp),dimension(:),intent(in),optional :: fjac_sparse !! jacobian matrix [sparse] real(wp),dimension(:),allocatable :: xnew !! used in [[func_for_fmin]] real(wp) :: alpha_min real(wp),dimension(:),allocatable :: search_direction !! search direction to use (may be modified from `p` if bounds are violated) logical,dimension(:),allocatable :: modified !! indicates the elements of p that were modified allocate(xnew(me%n)) allocate(search_direction(me%n)) allocate(modified(me%n)) ! find the minimum value of f in the range of alphas: alpha_min = fmin(func_for_fmin,me%alpha_min,me%alpha_max,me%fmin_tol) if (me%verbose) write(me%iunit,'(1P,*(A,1X,E16.6))') ' alpha_min = ', alpha_min call me%adjust_search_direction(xold,p,search_direction,modified) call me%compute_next_step(xold, search_direction, alpha_min, modified, x) if (all(x==xnew)) then ! already computed in the func else call me%func(x,fvec) f = me%norm(fvec) end if contains real(wp) function func_for_fmin(alpha) !! function for [[fmin]] implicit none real(wp),intent(in) :: alpha !! indep variable call me%compute_next_step(xold, search_direction, alpha, modified, xnew) call me%func(xnew,fvec) func_for_fmin = me%norm(fvec) ! return result f = func_for_fmin ! just in case this is the solution end function func_for_fmin end subroutine exact_linesearch