nonlin 1.5.2
A library that provides routines to compute the solutions to systems of nonlinear equations.
Loading...
Searching...
No Matches
nonlin_core::fcnnvar_helper Type Reference

Defines a type capable of encapsulating an equation of N variables. More...

Public Member Functions

procedure, public fcn fnh_fcn
 Executes the routine containing the function to evaluate.
 
procedure, public is_fcn_defined fnh_is_fcn_defined
 Tests if the pointer to the function has been assigned.
 
procedure, public set_fcn fnh_set_fcn
 Establishes a pointer to the routine containing the function.
 
procedure, public get_variable_count fnh_get_nvar
 Gets the number of variables in this system.
 
procedure, public set_gradient_fcn fnh_set_grad
 Establishes a pointer to the routine containing the gradient vector of the function.
 
procedure, public is_gradient_defined fnh_is_grad_defined
 Tests if the pointer to the routine containing the gradient has been assigned.
 
procedure, public gradient fnh_grad_fcn
 Computes the gradient of the function.
 

Public Attributes

integer(int32) m_nvar = 0
 The number of variables in m_fcn.
 

Static Public Attributes

procedure(gradientfcn), pointer, nopass m_grad => null()
 A pointer to the gradient routine.
 

Static Private Attributes

procedure(fcnnvar), pointer, nopass m_fcn => null()
 A pointer to the target fcnnvar routine.
 

Detailed Description

Defines a type capable of encapsulating an equation of N variables.

Example
The following example illustrates how to use this type to find the minimum of a function of multiple variables. In this instance the Nelder-Mead simplex method is utilized to minimize the Rosenbrock function in two variables.
program example
use iso_fortran_env
implicit none
! Local Variables
type(nelder_mead) :: solver
type(fcnnvar_helper) :: obj
procedure(fcnnvar), pointer :: fcn
real(real64) :: x(2), f
! Initialization
fcn => rosenbrock
call obj%set_fcn(fcn, 2)
! Define an initial guess - the solution is (1, 1)
call random_number(x)
! Call the solver
call solver%solve(obj, x, f)
! Display the output
print '(AF7.5AF7.5A)', "Minimum: (", x(1), ", ", x(2), ")"
print '(AE9.3)', "Function Value: ", f
contains
! Rosenbrock's Function
function rosenbrock(x) result(f)
real(real64), intent(in), dimension(:) :: x
real(real64) :: f
f = 1.0d2 * (x(2) - x(1)**2)**2 + (x(1) - 1.0d0)**2
end function
end program
nonlin_core
nonlin_optimize
The above program produces the following output.
Minimum: (1.00000, 1.00000)
Function Value: 0.213E-12

Definition at line 896 of file nonlin_core.f90.

Member Function/Subroutine Documentation

◆ fcn()

procedure, public nonlin_core::fcnnvar_helper::fcn

Executes the routine containing the function to evaluate.

Syntax
real(real64) function fcn(class(fcnnvar_helper) this, real(real64) x(:))
Parameters
[in]thisThe fcnnvar_helper object.
[in]xThe value of the independent variable at which the function should be evaluated.
Returns
The value of the function at x.

Definition at line 916 of file nonlin_core.f90.

◆ get_variable_count()

procedure, public nonlin_core::fcnnvar_helper::get_variable_count

Gets the number of variables in this system.

Syntax
integer(int32) function get_variable_count(class(fcnnvar_helper) this)
Parameters
[in]thisThe fcnnvar_helper object.
Returns
The number of variables.

Definition at line 993 of file nonlin_core.f90.

◆ gradient()

procedure, public nonlin_core::fcnnvar_helper::gradient

Computes the gradient of the function.

Syntax
subroutine gradient(class(fcnnvar_helper) this, real(real64) x(:), &
real(real64) g(:), optional real(real64) fv(:), &
optional integer(int32) err)
Parameters
[in]thisThe fcnnvar_helper object.
[in,out]xAn N-element array containing the independent variables defining the point about which the derivatives will be calculated. This array is restored upon output.
[out]gAn N-element array where the gradient will be written upon output.
[in]fvAn optional input providing the function value at x.
[out]errAn optional integer output that can be used to determine error status. If not used, and an error is encountered, the routine simply returns silently. If used, the following error codes identify error status:
  • 0: No error has occurred.
  • n: A positive integer denoting the index of an invalid input.

Definition at line 1104 of file nonlin_core.f90.

◆ is_fcn_defined()

procedure, public nonlin_core::fcnnvar_helper::is_fcn_defined

Tests if the pointer to the function has been assigned.

Syntax
logical function is_fcn_defined(class(fcnnvar_helper) this)
Parameters
[in]thisThe fcnnvar_helper object.
Returns
Returns true if the pointer has been assigned; else, false.

Definition at line 926 of file nonlin_core.f90.

◆ is_gradient_defined()

procedure, public nonlin_core::fcnnvar_helper::is_gradient_defined

Tests if the pointer to the routine containing the gradient has been assigned.

Syntax
logical function is_gradient_defined(class(fcnnvar_helper) this)
Parameters
[in]thisThe fcnnvar_helper object.
Returns
Returns true if the pointer has been assigned; else, false.

Definition at line 1081 of file nonlin_core.f90.

◆ set_fcn()

procedure, public nonlin_core::fcnnvar_helper::set_fcn

Establishes a pointer to the routine containing the function.

Syntax
subroutine set_fcn(class(fcnnvar_helper) this, procedure(fcnnvar) pointer fcn, integer(int32) nvar)
Parameters
[in,out]thisThe fcnnvar_helper object.
[in]fcnThe function pointer.
[in]nvarThe number of variables in the function.
Example
The following example illustrates how to use this routine to inform the solver of the function to optimize. This particular example utilizes the Nelder-Mead simplex method to minimize the function.
program example
use iso_fortran_env
implicit none
! Local Variables
type(nelder_mead) :: solver
type(fcnnvar_helper) :: obj
procedure(fcnnvar), pointer :: fcn
real(real64) :: x(2), f
! Initialization
fcn => rosenbrock
call obj%set_fcn(fcn, 2)
! Define an initial guess - the solution is (1, 1)
call random_number(x)
! Call the solver
call solver%solve(obj, x, f)
! Display the output
print '(AF7.5AF7.5A)', "Minimum: (", x(1), ", ", x(2), ")"
print '(AE9.3)', "Function Value: ", f
contains
! Rosenbrock's Function
function rosenbrock(x) result(f)
real(real64), intent(in), dimension(:) :: x
real(real64) :: f
f = 1.0d2 * (x(2) - x(1)**2)**2 + (x(1) - 1.0d0)**2
end function
end program
The above program produces the following output.
Minimum: (1.00000, 1.00000)
Function Value: 0.213E-12

Definition at line 983 of file nonlin_core.f90.

◆ set_gradient_fcn()

procedure, public nonlin_core::fcnnvar_helper::set_gradient_fcn

Establishes a pointer to the routine containing the gradient vector of the function.

Syntax
subroutine set_gradient_fcn(class(fcnnvar_helper) this, procedure(gradientfcn) pointer fcn)
Parameters
[in,out]thisThe fcnnvar_helper object.
[in]fcnThe pointer to the gradient routine.
Example
The following example illustrates the use of a user-defined gradient function. This particular example utilizes the BFGS method to minimize Beale's function.
program example
use iso_fortran_env
implicit none
! Local Variables
type(bfgs) :: solver
type(fcnnvar_helper) :: obj
procedure(fcnnvar), pointer :: fcn
procedure(gradientfcn), pointer :: grad
real(real64) :: x(2), f
! Tell the solver where to find the function
fcn => beale
grad => bealegrad
call obj%set_fcn(fcn, 2)
call obj%set_gradient_fcn(grad)
! Define an initial guess
x = 1.0d0
! Compute the solution
call solver%solve(obj, x, f)
! Display the output
print '(AF7.5AF7.5A)', "Minimum: (", x(1), ", ", x(2), ")"
print '(AE9.3)', "Function Value: ", f
contains
! The Beale function:
! f(x) = (1.5 - x + xy)**2 + (2.25 - x + xy**2)**2 + (2.625 - x + xy**3)**2
! The minimum is at x = 3, y = 0.5, and f(3, 0.5) = 0
function beale(x) result(f)
real(real64), intent(in), dimension(:) :: x
real(real64) :: f
f = (1.5d0 - x(1) + x(1) * x(2))**2 + &
(2.25d0 - x(1) + x(1) * x(2)**2)**2 + &
(2.625d0 - x(1) + x(1) * x(2)**3)**2
end function
! The gradient
subroutine bealegrad(x, g)
real(real64), intent(in), dimension(:) :: x
real(real64), intent(out), dimension(:) :: g
g(1) = 2.0d0 * (x(2)**3 - 1.0d0) * (x(1) * x(2)**3 - x(1) + 2.625d0) + &
2.0d0 * (x(2)**2 - 1.0d0) * (x(1) * x(2)**2 - x(1) + 2.25d0) + &
2.0d0 * (x(2) - 1.0d0) * (x(1) * x(2) - x(1) + 1.5d0)
g(2) = 6.0d0 * x(1) * x(2)**2 * (x(1) * x(2)**3 - x(1) + 2.625d0) + &
4.0d0 * x(1) * x(2) * (x(1) * x(2)**2 - x(1) + 2.25d0) + &
2.0d0 * x(1) * (x(1) * x(2) - x(1) + 1.5d0)
end subroutine
end program
The above program produces the following output.
Minimum: (3.00000, 0.50000)
Function Value: 0.999E-28

Definition at line 1070 of file nonlin_core.f90.

Member Data Documentation

◆ m_fcn

procedure(fcnnvar), pointer, nopass nonlin_core::fcnnvar_helper::m_fcn => null()
staticprivate

A pointer to the target fcnnvar routine.

Definition at line 899 of file nonlin_core.f90.

◆ m_grad

procedure(gradientfcn), pointer, nopass nonlin_core::fcnnvar_helper::m_grad => null()
static

A pointer to the gradient routine.

Definition at line 901 of file nonlin_core.f90.

◆ m_nvar

integer(int32) nonlin_core::fcnnvar_helper::m_nvar = 0

The number of variables in m_fcn.

Definition at line 903 of file nonlin_core.f90.


The documentation for this type was generated from the following file: