anderson_bjorck Subroutine

private subroutine anderson_bjorck(me, ax, bx, fax, fbx, xzero, fzero, iflag)

Compute the zero of the function f(x) in the interval ax,bx using the Anderson-Bjorck method.

See also

  • G.E. Mullges & F. Uhlig, "Numerical Algorithms with Fortran", Springer, 1996. Section 2.8.2, p 36.

Type Bound

anderson_bjorck_solver

Arguments

Type IntentOptional Attributes Name
class(anderson_bjorck_solver), intent(inout) :: me
real(kind=wp), intent(in) :: ax

left endpoint of initial interval

real(kind=wp), intent(in) :: bx

right endpoint of initial interval

real(kind=wp), intent(in) :: fax

f(ax)

real(kind=wp), intent(in) :: fbx

f(ax)

real(kind=wp), intent(out) :: xzero

abscissa approximating a zero of f in the interval ax,bx

real(kind=wp), intent(out) :: fzero

value of f at the root (f(xzero))

integer, intent(out) :: iflag

status flag (0=root found, -2=max iterations reached)


Calls

proc~~anderson_bjorck~~CallsGraph proc~anderson_bjorck root_module::anderson_bjorck_solver%anderson_bjorck f f proc~anderson_bjorck->f proc~choose_best root_module::choose_best proc~anderson_bjorck->proc~choose_best proc~converged root_module::root_solver%converged proc~anderson_bjorck->proc~converged proc~secant root_module::secant proc~anderson_bjorck->proc~secant proc~solution root_module::root_solver%solution proc~anderson_bjorck->proc~solution proc~bisect root_module::bisect proc~secant->proc~bisect

Source Code

    subroutine anderson_bjorck(me,ax,bx,fax,fbx,xzero,fzero,iflag)

    implicit none

    class(anderson_bjorck_solver),intent(inout) :: me
    real(wp),intent(in)    :: ax      !! left endpoint of initial interval
    real(wp),intent(in)    :: bx      !! right endpoint of initial interval
    real(wp),intent(in)    :: fax     !! `f(ax)`
    real(wp),intent(in)    :: fbx     !! `f(ax)`
    real(wp),intent(out)   :: xzero   !! abscissa approximating a zero of `f` in the interval `ax`,`bx`
    real(wp),intent(out)   :: fzero   !! value of `f` at the root (`f(xzero)`)
    integer,intent(out)    :: iflag   !! status flag (`0`=root found, `-2`=max iterations reached)

    integer :: i !! counter
    logical :: root_found !! convergence in x
    real(wp) :: x1,x2,x3,f1,f2,f3,g,f1tmp

    ! initialize:
    iflag = 0
    x1    = ax
    x2    = bx
    f1    = fax
    f2    = fbx

    ! main loop:
    do i = 1,me%maxiter

        x3 = secant(x1,x2,f1,f2,ax,bx)
        f3 = me%f(x3)
        if (me%solution(x3,f3,xzero,fzero)) return

        ! determine a new inclusion interval:
        if (f2*f3<0.0_wp) then
            ! zero lies between x2 and x3
            x1 = x2
            x2 = x3
            f1 = f2
            f2 = f3
            f1tmp = f1
        else
            ! zero lies between x1 and x3
            g = 1.0_wp - f3/f2
            if (g<=0.0_wp) g = 0.5_wp
            x2 = x3
            f1tmp = f1
            f1 = g*f1
            f2 = f3
        end if

        ! check for convergence:
        root_found = me%converged(x1,x2)
        if (root_found .or. i == me%maxiter) then
            call choose_best(x1,x2,f1tmp,f2,xzero,fzero)
            if (.not. root_found) iflag = -2  ! max iterations reached
            exit
        end if

    end do

    end subroutine anderson_bjorck