FMUTIL
0.1
Fortran Miscellaneous UTILities
|
RootFinding Module. More...
Data Types | |
interface | func |
Interface of the function whose root needs to be computed. More... | |
Functions/Subroutines | |
subroutine | polyroots (c, r, error, BalanceOn) |
Computes the complex roots of the given polynomial. More... | |
pure subroutine | root (a, b, ya, yb, tol, f, x, fx, niter, error) |
Subroutine for computing the root using Brent algorithm. More... | |
Variables | |
integer, parameter | max_iterations = 50 |
Max iterations for all the root-finding. More... | |
RootFinding Module.
This module provides functions for root-finding of polynomials by solving associated eigenvalue problems and nonlinear equations using algorithms like Brent's Algorithm.
subroutine rootfinding::polyroots | ( | complex(wp), dimension(1:), intent(in) | c, |
complex(wp), dimension(:), intent(out), allocatable | r, | ||
integer, intent(out) | error, | ||
logical, intent(in), optional | BalanceOn | ||
) |
Computes the complex roots of the given polynomial.
PolyRoots finds the roots of the polynomial defined by complex coefficients given in array "c" (c(1) is the highest degree coeff.). It forms the companion matrix and then solves its eigenvalues using LAPACk routines from Intel MKL library. The algorithm is same as that used by MATLAB's "roots" command. References:
[in] | c | Complex coefficients of the polynomial defined by p(x) = c(1)*x^n + ... + c(n)*x + c(n+1) |
[out] | r | Complex roots array |
[out] | error | Error status returned to the user. "0" : Success, valid results in "r" "-1": Error, input coefficients array has length less than two "-2": Error, all the coefficients are zero "-3": Error, companion matrix balancing failed "-4": Error, LAPACK routine failed to compute eigenvalues |
[in] | balanceon | If True, then the companion matrix is balanced first before finding its eigenvalues, otherwise not. Default is True. It may produce inaccurate results in case of very poorly-scaled matrices. |
Definition at line 66 of file rootfinding.f90.
References fmutilbase::inf.
pure subroutine rootfinding::root | ( | real(wp), intent(in) | a, |
real(wp), intent(in) | b, | ||
real(wp), intent(in) | ya, | ||
real(wp), intent(in) | yb, | ||
real(wp), intent(in) | tol, | ||
procedure(func) | f, | ||
real(wp), intent(out) | x, | ||
real(wp), intent(out) | fx, | ||
integer, intent(out) | niter, | ||
integer, intent(out) | error | ||
) |
Subroutine for computing the root using Brent algorithm.
For details on the algorithm, see Numerical receipes in C book at https://www2.units.it/ipl/students_area/imm2/files/Numerical_Recipes.pdf
[in] | a | Lower bound on the root |
[in] | b | Upper bound on the root |
[in] | ya | f(a) |
[in] | yb | f(b) |
[in] | tol | tolerance |
f | function | |
[out] | x | root |
[out] | fx | value at root |
[out] | niter | no. of iterations |
[out] | error | Error status returned to the user. "0" means root has been found. "1" means provided bounds do not contain the root, i.e., root is not bracketed. "2" means maximum number of interations reached. |
Definition at line 208 of file rootfinding.f90.
References max_iterations.
integer, parameter rootfinding::max_iterations = 50 |
Max iterations for all the root-finding.
Definition at line 42 of file rootfinding.f90.
Referenced by root().