A simple search that just evaluates the function at a specified number of points and picks the one with the minimum function value.
| 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 fixed_point_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] integer :: i !! counter integer :: n_points !! number of points to compute real(wp),dimension(:),allocatable :: alphas_to_try !! set of `alpha` values to try real(wp),dimension(:),allocatable :: x_tmp !! temp `x` real(wp),dimension(:),allocatable :: fvec_tmp !! temp `fvec` real(wp) :: f_tmp !! temp `f` real(wp) :: step_size !! step size for `alpha` integer :: n !! number of steps to divide the interval 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 ! 1 o-----------o ! 2 o-----o-----o ! 3 o---o---o---o n = me%n_intervals n_points = n + 1 allocate(alphas_to_try(n_points)) allocate(x_tmp(me%n)) allocate(fvec_tmp(me%m)) allocate(search_direction(me%n)) allocate(modified(me%n)) step_size = (me%alpha_max - me%alpha_min) / real(n,wp) ! compute the alphas: alphas_to_try(1) = me%alpha_min do i = 2, n alphas_to_try(i) = alphas_to_try(i-1) + step_size end do alphas_to_try(n_points) = me%alpha_max ! now compute the functions at these alphas: f = big call me%adjust_search_direction(xold,p,search_direction,modified) do i = 1, n_points call me%compute_next_step(xold, search_direction, alphas_to_try(i), modified, x_tmp) ! evaluate the function at this point: call me%func(x_tmp,fvec_tmp) f_tmp = me%norm(fvec_tmp) if (f_tmp<=f) then ! new best point x = x_tmp f = f_tmp fvec = fvec_tmp end if end do end subroutine fixed_point_linesearch