root_scalar_by_name Subroutine

private subroutine root_scalar_by_name(method, fun, ax, bx, xzero, fzero, iflag, ftol, rtol, atol, maxiter, fax, fbx, bisect_on_failure)

Non-object-oriented wrapper.

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: method

the method to use

procedure(func2) :: fun

user function to find the root of

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, -999=invalid method)

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

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.


Calls

proc~~root_scalar_by_name~~CallsGraph proc~root_scalar_by_name root_module::root_scalar_by_name interface~root_scalar root_module::root_scalar proc~root_scalar_by_name->interface~root_scalar proc~root_name_to_method root_module::root_name_to_method proc~root_scalar_by_name->proc~root_name_to_method interface~root_scalar->proc~root_scalar_by_name proc~root_scalar_by_type root_module::root_scalar_by_type interface~root_scalar->proc~root_scalar_by_type proc~lowercase root_module::lowercase proc~root_name_to_method->proc~lowercase proc~initialize_root_solver root_module::root_solver%initialize_root_solver proc~root_scalar_by_type->proc~initialize_root_solver proc~solve root_module::root_solver%solve proc~root_scalar_by_type->proc~solve proc~solve->interface~root_scalar find_root find_root proc~solve->find_root proc~get_fa_fb root_module::root_solver%get_fa_fb proc~solve->proc~get_fa_fb proc~solution root_module::root_solver%solution proc~solve->proc~solution

Called by

proc~~root_scalar_by_name~~CalledByGraph proc~root_scalar_by_name root_module::root_scalar_by_name interface~root_scalar root_module::root_scalar proc~root_scalar_by_name->interface~root_scalar interface~root_scalar->proc~root_scalar_by_name proc~root_scalar_by_type root_module::root_scalar_by_type interface~root_scalar->proc~root_scalar_by_type proc~solve root_module::root_solver%solve proc~solve->interface~root_scalar proc~root_scalar_by_type->proc~solve

Source Code

    subroutine root_scalar_by_name(method,fun,ax,bx,xzero,fzero,iflag,&
                                   ftol,rtol,atol,maxiter,fax,fbx,&
                                   bisect_on_failure)

    implicit none

    character(len=*),intent(in)   :: method   !! the method to use
    procedure(func2)              :: fun      !! user function to find the root of
    real(wp),intent(in)           :: ax       !! left endpoint of initial interval
    real(wp),intent(in)           :: bx       !! right endpoint of initial interval
    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 (`-1`=error, `0`=root found, `-999`=invalid method)
    real(wp),intent(in),optional  :: ftol     !! absolute tolerance for `f=0`
    real(wp),intent(in),optional  :: rtol     !! relative tol for x
    real(wp),intent(in),optional  :: atol     !! absolute tol for x
    integer,intent(in),optional   :: maxiter  !! maximum number of iterations
    real(wp),intent(in),optional  :: fax      !! if `f(ax)` is already known, it can be input here
    real(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.

    type(root_method) :: r

    r = root_name_to_method(method)

    if (r%id /= 0) then
        call root_scalar(r,fun,ax,bx,xzero,fzero,iflag,&
                         ftol,rtol,atol,maxiter,fax,fbx,&
                         bisect_on_failure)
    else
        iflag = -999    ! invalid method
        return
    end if

    end subroutine root_scalar_by_name