nlesolver_test_1 Program

Uses

  • program~~nlesolver_test_1~~UsesGraph program~nlesolver_test_1 nlesolver_test_1 module~nlesolver_module nlesolver_module program~nlesolver_test_1->module~nlesolver_module iso_fortran_env iso_fortran_env module~nlesolver_module->iso_fortran_env module~fmin_module fmin_module module~nlesolver_module->module~fmin_module module~fmin_module->iso_fortran_env

Test of a small, square (n=m) problem.


Calls

program~~nlesolver_test_1~~CallsGraph program~nlesolver_test_1 nlesolver_test_1 proc~get_status nlesolver_module::nlesolver_type%get_status program~nlesolver_test_1->proc~get_status proc~initialize_nlesolver_variables nlesolver_module::nlesolver_type%initialize_nlesolver_variables program~nlesolver_test_1->proc~initialize_nlesolver_variables proc~nlesolver_solver nlesolver_module::nlesolver_type%nlesolver_solver program~nlesolver_test_1->proc~nlesolver_solver proc~set_status nlesolver_module::nlesolver_type%set_status proc~initialize_nlesolver_variables->proc~set_status proc~initialize_ez lsqr_module::lsqr_solver_ez%initialize_ez proc~nlesolver_solver->proc~initialize_ez proc~linear_solver nlesolver_module::linear_solver proc~nlesolver_solver->proc~linear_solver proc~nlesolver_solver->proc~set_status proc~solve_ez lsqr_module::lsqr_solver_ez%solve_ez proc~nlesolver_solver->proc~solve_ez proc~lsqr lsqr_module::lsqr_solver%LSQR proc~solve_ez->proc~lsqr aprod aprod proc~lsqr->aprod proc~d2norm lsqr_module::d2norm proc~lsqr->proc~d2norm proc~dcopy lsqpblas_module::dcopy proc~lsqr->proc~dcopy proc~dnrm2 lsqpblas_module::dnrm2 proc~lsqr->proc~dnrm2 proc~dscal lsqpblas_module::dscal proc~lsqr->proc~dscal

Variables

Type Attributes Name Initial
integer, parameter :: n = 2
integer, parameter :: m = 2
integer, parameter :: max_iter = 100
real(kind=wp), parameter :: tol = 1.0e-8_wp
logical, parameter :: verbose = .false.
type(nlesolver_type) :: solver
real(kind=wp) :: alpha
logical :: use_broyden
integer :: step_mode
integer :: n_intervals
integer :: istat

Integer status code.

character(len=:), allocatable :: message

Text status message

real(kind=wp), dimension(n) :: x
integer :: f_evals
integer :: i
character(len=:), allocatable :: description
real(kind=wp) :: fmin_tol

Subroutines

subroutine func(me, x, f)

compute the function

Arguments

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

subroutine grad(me, x, g)

compute the gradient of the function (Jacobian):

Arguments

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

subroutine export(me, x, f, iter)

export an iteration:

Arguments

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

iteration number

Source Code

    program nlesolver_test_1

    use nlesolver_module, wp => nlesolver_rk

    implicit none

    integer,parameter :: n = 2
    integer,parameter :: m = 2
    integer,parameter :: max_iter = 100
    real(wp),parameter :: tol = 1.0e-8_wp
    logical,parameter :: verbose = .false.

    type(nlesolver_type) :: solver
    real(wp) :: alpha
    logical :: use_broyden
    integer :: step_mode
    integer :: n_intervals
    integer :: istat !! Integer status code.
    character(len=:),allocatable :: message  !! Text status message
    real(wp),dimension(n) :: x
    integer :: f_evals
    integer :: i
    character(len=:),allocatable :: description
    real(wp) :: fmin_tol

    fmin_tol = 1.0e-2_wp ! don't need a tight tol for this
    n_intervals = 2
    alpha = 1.0_wp

    write(*,*) ''
    write(*,*) '***********************'
    write(*,*) '* nlesolver_test_1    *'
    write(*,*) '***********************'
    write(*,*) ''
    do i = 1, 8

        select case (i)
        case(1)
            step_mode = 1
            use_broyden = .false.
            f_evals = 0
            description = 'Constant alpha'
        case(2)
            step_mode = 1
            use_broyden = .true.
            f_evals = 0
            description = 'Constant alpha + broyden'
        case(3)
            step_mode = 2
            use_broyden = .false.
            f_evals = 0
            description = 'Backtracking line search'
        case(4)
            step_mode = 2
            use_broyden = .true.
            f_evals = 0
            description = 'Backtracking line search + broyden'
        case(5)
            step_mode = 3
            use_broyden = .false.
            f_evals = 0
            description = 'Exact line search'
        case(6)
            step_mode = 3
            use_broyden = .true.
            f_evals = 0
            description = 'Exact line search + broyden'
        case(7)
            step_mode = 4
            use_broyden = .false.
            f_evals = 0
            description = 'Fixed point search'
        case(8)
            step_mode = 4
            use_broyden = .true.
            f_evals = 0
            description = 'Fixed point search + broyden'
        case default
            error stop 'invalid case'
        end select

        write(*,*) '-------------------------------------------------------'
        write(*,'(A,I3,A,A)') 'Case ', i, ' : ', description
        write(*,*) ''

        call solver%initialize( n = n, &
                                m = m, &
                                max_iter = max_iter, &
                                tol = tol, &
                                func = func, &
                                grad = grad, &
                                step_mode = step_mode,&
                                use_broyden = use_broyden,&
                                export_iteration = export,&
                                n_intervals = n_intervals, &
                                fmin_tol = fmin_tol, &
                                verbose = verbose)
        call solver%status(istat, message)
        write(*,'(I3,1X,A)') istat, message
        if (istat /= 0) error stop

        x = [1.0_wp, 2.0_wp]
        call solver%solve(x)

        call solver%status(istat, message)
        write(*,'(I3,1X,A)') istat, message
        write(*,*) ''

    end do

    contains

    subroutine func(me,x,f)
        !! compute the function
        implicit none
        class(nlesolver_type),intent(inout) :: me
        real(wp),dimension(:),intent(in)    :: x
        real(wp),dimension(:),intent(out)   :: f

        f_evals = f_evals + 1

        f(1) = x(1)**2 + x(2) - 0.1_wp
        f(2) = x(2) + 0.2_wp

        ! root is 5.477226E-01  -2.000000E-01

    end subroutine func

    subroutine grad(me,x,g)
        !! compute the gradient of the function (Jacobian):
        implicit none
        class(nlesolver_type),intent(inout) :: me
        real(wp),dimension(:),intent(in)    :: x
        real(wp),dimension(:,:),intent(out) :: g

        f_evals = f_evals + 2   ! to approximate forward diff derivatives

        g(1,1) = 2.0_wp * x(1)  !df(1)/dx
        g(2,1) = 0.0_wp         !df(2)/dx

        g(1,2) = 1.0_wp         !df(1)/dy
        g(2,2) = 1.0_wp         !df(2)/dy

    end subroutine grad

    subroutine export(me,x,f,iter)
        !! export an iteration:
        implicit none
        class(nlesolver_type),intent(inout) :: me
        real(wp),dimension(:),intent(in)    :: x
        real(wp),dimension(:),intent(in)    :: f
        integer,intent(in)                  :: iter !! iteration number

        write(*,'(1P,I3,1X,A,I3,A,*(E15.6))') iter, '(',f_evals,')', x, norm2(f)

    end subroutine export

!******************************************************************************************************
    end program nlesolver_test_1