brent_module Module

Brent algorithms for minimization and root solving without derivatives.

See also

  1. R. Brent, "Algorithms for Minimization Without Derivatives", Prentice-Hall, Inc., 1973.

Uses

  • module~~brent_module~~UsesGraph module~brent_module brent_module module~kind_module kind_module module~brent_module->module~kind_module module~numbers_module numbers_module module~brent_module->module~numbers_module iso_fortran_env iso_fortran_env module~kind_module->iso_fortran_env module~numbers_module->module~kind_module

Used by

  • module~~brent_module~~UsedByGraph module~brent_module brent_module module~fortran_astrodynamics_toolkit fortran_astrodynamics_toolkit module~fortran_astrodynamics_toolkit->module~brent_module proc~integrate_to_event rk_module_variable_step::rk_variable_step_class%integrate_to_event proc~integrate_to_event->module~brent_module proc~integrate_to_event~2 rk_module::rk_class%integrate_to_event proc~integrate_to_event~2->module~brent_module

Abstract Interfaces

abstract interface

  • private function func(me, x) result(f)

    Interface to the function to be minimized. It should evaluate f(x) for any x in the interval (ax,bx)

    Arguments

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

    Return Value real(kind=wp)


Derived Types

type, public ::  brent_class

the main class

Components

Type Visibility Attributes Name Initial
procedure(func), public, pointer :: f => null()

function to be minimized

Type-Bound Procedures

procedure, public :: set_function
procedure, public :: minimize => fmin
procedure, public :: find_zero => zeroin

Functions

private function fmin(me, ax, bx, tol) result(xmin)

An approximation x to the point where f attains a minimum on the interval (ax,bx) is determined.

Read more…

Arguments

Type IntentOptional Attributes Name
class(brent_class), 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) :: tol

desired length of the interval of uncertainty of the final result (>=0)

Return Value real(kind=wp)

abcissa approximating the point where f attains a minimum


Subroutines

private subroutine set_function(me, f)

Author
Jacob Williams
Date
7/19/2014

Set the function to be minimized.

Arguments

Type IntentOptional Attributes Name
class(brent_class), intent(inout) :: me
procedure(func) :: f

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

Find a zero of the function in the given interval to within a tolerance , where is the relative machine precision defined as the smallest representable number such that .

Read more…

Arguments

Type IntentOptional Attributes Name
class(brent_class), 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) :: tol

desired length of the interval of uncertainty of the final result (>=0)

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)

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(ax) is already known, it can be input here

public subroutine brent_test()

Author
Jacob Williams
Date
7/16/2014

Test of the fmin and zeroin functions.

Arguments

None