muller_solver Derived Type

type, public, extends(root_solver) :: muller_solver

muller root solver


Inherits

type~~muller_solver~~InheritsGraph type~muller_solver muller_solver type~root_solver root_solver type~muller_solver->type~root_solver

Type-Bound Procedures

procedure, public :: initialize => initialize_root_solver

initialize the class [must be called first]

  • private subroutine initialize_root_solver(me, f, ftol, rtol, atol, maxiter)

    Initialize the root_solver class.

    Read more…

    Arguments

    Type IntentOptional Attributes Name
    class(root_solver), intent(out) :: me
    procedure(func) :: f

    user function f(x) to find the root of

    real(kind=wp), intent(in), optional :: ftol

    absolute tolerance for f=0

    real(kind=wp), intent(in), optional :: rtol

    relative tol for x

    real(kind=wp), intent(in), optional :: atol

    absolute tol for x

    integer, intent(in), optional :: maxiter

    maximum number of iterations

procedure, public :: solve

main routine for finding the root

  • private subroutine solve(me, ax, bx, xzero, fzero, iflag, fax, fbx, bisect_on_failure)

    Main wrapper routine for all the methods.

    Arguments

    Type IntentOptional Attributes Name
    class(root_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(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 (-1=error, 0=root found, -4=ax must be /= bx)

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

    if f(ax) is already known, it can be input here

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

    if f(bx) is already known, it can be input here

    logical, intent(in), optional :: bisect_on_failure

    if true, then if the specified method fails, it will be retried using the bisection method. (default is False). Note that this can use up to maxiter additional function evaluations.

procedure, public :: find_root => muller

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

    Improved Muller method (for real roots only). Will fall back to bisection if any step fails.

    Read more…

    Arguments

    Type IntentOptional Attributes Name
    class(muller_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)

Source Code

    type,extends(root_solver),public :: muller_solver
    !! muller root solver
    private
    contains
    private
    procedure,public :: find_root => muller
    end type muller_solver