Octopus
|
This module is intended to contain "only mathematical" functions and procedures. More...
This module is intended to contain "only mathematical" functions and procedures.
Data Types | |
interface | diagonal_matrix |
interface | interpolate |
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of the classical formula of Lagrange). More... | |
interface | is_close |
interface | log2 |
interface | lower_triangular_to_hermitian |
interface | pad |
interface | upper_triangular_to_hermitian |
Functions/Subroutines | |
elemental logical function | dis_close_scalar (x, y, rtol, atol) |
Are \(x\) and \(y\) equal within a tolerance. More... | |
elemental logical function | zis_close_scalar (x, y, rtol, atol) |
Same as dis_close_scalar for complex numbers. More... | |
pure integer function, dimension(dim, dim) | idiagonal_matrix (dim, diag) |
Currently only returns a matrix whose diagonal elements are all the same. Note that the real and complex versions are in math_inc.F90. More... | |
recursive real(real64) function, public | hermite (n, x) |
recursive integer function, public | factorial (n) |
subroutine, public | ylmr_cmplx (xx, li, mi, ylm) |
Computes spherical harmonics ylm at position (x, y, z) More... | |
subroutine, public | ylmr_real (xx, li, mi, ylm) |
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x, y, z): ylm = c * plm( cos(theta) ) * sin(m*phi) for m < 0 ylm = c * plm( cos(theta) ) * cos(m*phi) for m >= 0 with (theta,phi) the polar angles of r, c a positive normalization. More... | |
subroutine, public | weights (N, M, cc, side) |
Compute the weights for finite-difference calculations: More... | |
real(real64) pure function, public | ddelta (i, j) |
subroutine, public | make_idx_set (n, out, length, in) |
Construct out(1:length) = (/1, ..., n/) if in is not present, out(1:length) = in otherwise. More... | |
logical function, public | member (n, a) |
Considers a(1:ubound(a, 1)) as an integer set and checks if n is a member of it. More... | |
subroutine, public | interpolation_coefficients (nn, xa, xx, cc) |
logical pure function, public | even (n) |
Returns if n is even. More... | |
logical pure function, public | odd (n) |
Returns if n is odd. More... | |
subroutine, public | cartesian2hyperspherical (x, u) |
Performs a transformation of an n-dimensional vector from Cartesian coordinates to hyperspherical coordinates. More... | |
subroutine, public | hyperspherical2cartesian (u, x) |
Performs the inverse transformation of cartesian2hyperspherical. More... | |
subroutine, public | hypersphere_grad_matrix (grad_matrix, r, x) |
Gives the hyperspherical gradient matrix, which contains the derivatives of the Cartesian coordinates with respect to the hyperspherical angles. More... | |
integer(int64) pure function | pad88 (size, blk) |
integer(int64) pure function | pad48 (size, blk) |
integer(int64) pure function | pad8 (size, blk) |
integer pure function | pad4 (size, blk) |
integer pure function, public | pad_pow2 (size) |
create array size, which is padded to powers of 2 More... | |
real(real64) pure function | dlog2 (xx) |
integer pure function | ilog2 (xx) |
integer(int64) pure function | llog2 (xx) |
complex(real64) pure function, public | exponential (z) |
Wrapper for exponential. More... | |
complex(real64) pure function, public | phi1 (z) |
Compute phi1(z) = (exp(z)-1)/z. More... | |
complex(real64) pure function, public | phi2 (z) |
Compute phi2(z) = (phi1(z)-1)/z = (exp(z) - z - 1)/z^2. More... | |
real(real64) pure function, public | square_root (x) |
Wrapper for sqrt. More... | |
logical function, public | is_prime (n) |
subroutine, public | generate_rotation_matrix (R, ff, tt) |
Generates a rotation matrix R to rotate a vector f to t. More... | |
subroutine, public | numder_ridders (x, h, res, err, f) |
Numerical derivative (Ridder`s algorithm). More... | |
pure complex(real64) function, dimension(1:3), public | dzcross_product (a, b) |
pure complex(real64) function, dimension(1:3), public | zdcross_product (a, b) |
subroutine, public | generalized_laguerre_polynomial (np, nn, mm, xx, cx) |
subroutine | dupper_triangular_to_hermitian (nn, aa) |
subroutine | zupper_triangular_to_hermitian (nn, aa) |
subroutine | dlower_triangular_to_hermitian (nn, aa) |
subroutine | zlower_triangular_to_hermitian (nn, aa) |
subroutine, public | dsymmetrize_matrix (nn, aa) |
subroutine, public | dzero_small_elements_matrix (nn, aa, tol) |
pure complex(real64) function, dimension(dim, dim) | zdiagonal_matrix (dim, diag) |
Currently only returns a matrix whose diagonal elements are all the same. Note that the integer version is in math.F90. More... | |
subroutine | zinterpolate_2 (xa, ya, x, y) |
subroutine | zinterpolate_1 (xa, ya, x, y) |
subroutine | zinterpolate_0 (xa, ya, x, y) |
pure complex(real64) function, dimension(1:3), public | zcross_product (a, b) |
pure real(real64) function, dimension(dim, dim) | ddiagonal_matrix (dim, diag) |
Currently only returns a matrix whose diagonal elements are all the same. Note that the integer version is in math.F90. More... | |
subroutine | dinterpolate_2 (xa, ya, x, y) |
subroutine | dinterpolate_1 (xa, ya, x, y) |
subroutine | dinterpolate_0 (xa, ya, x, y) |
pure real(real64) function, dimension(1:3), public | dcross_product (a, b) |
|
private |
Are \(x\) and \(y\) equal within a tolerance.
The function evaluates the expression:
\[ |x - y| \leq (atol * rtol) * |y| \]
The tolerance values are positive, typically very small numbers. The relative difference \((rtol * |y|)\) and the absolute difference \( atol \) are added together to compare against the absolute difference between \(x\) and \(y\). Default tolerances are based on numpy''s [implementation](https:
[in] | x | Scalar. |
[in] | y | Scalar. |
[in] | rtol | Optional, relative tolerance. |
[in] | atol | Optional, absolute tolerance. |
|
private |
|
private |
recursive real(real64) function, public math_oct_m::hermite | ( | integer, intent(in) | n, |
real(real64), intent(in) | x | ||
) |
recursive integer function, public math_oct_m::factorial | ( | integer, intent(in) | n | ) |
subroutine, public math_oct_m::ylmr_cmplx | ( | real(real64), dimension(3), intent(in) | xx, |
integer, intent(in) | li, | ||
integer, intent(in) | mi, | ||
complex(real64), intent(out) | ylm | ||
) |
subroutine, public math_oct_m::ylmr_real | ( | real(real64), dimension(3), intent(in) | xx, |
integer, intent(in) | li, | ||
integer, intent(in) | mi, | ||
real(real64), intent(out) | ylm | ||
) |
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x, y, z): ylm = c * plm( cos(theta) ) * sin(m*phi) for m < 0 ylm = c * plm( cos(theta) ) * cos(m*phi) for m >= 0 with (theta,phi) the polar angles of r, c a positive normalization.
subroutine, public math_oct_m::weights | ( | integer, intent(in) | N, |
integer, intent(in) | M, | ||
real(real64), dimension(0:,0:,0:), intent(out) | cc, | ||
integer, intent(in), optional | side | ||
) |
Compute the weights for finite-difference calculations:
N -> highest order of the derivative to be approximated M -> number of grid points to be used in the approximation.
c(j,k,i) -> ith order derivative at kth-order approximation j=0,k: the coefficients acting of each point
side -> -1 left-sided, +1 right-sided, 0 centered (default)
[out] | cc | (0:M, 0:M, 0:N) |
real(real64) pure function, public math_oct_m::ddelta | ( | integer, intent(in) | i, |
integer, intent(in) | j | ||
) |
subroutine, public math_oct_m::make_idx_set | ( | integer, intent(in) | n, |
integer, dimension(:), intent(out), allocatable | out, | ||
integer, intent(out) | length, | ||
integer, dimension(:), intent(in), optional | in | ||
) |
logical function, public math_oct_m::member | ( | integer, intent(in) | n, |
integer, dimension(:), intent(in) | a | ||
) |
subroutine, public math_oct_m::interpolation_coefficients | ( | integer, intent(in) | nn, |
real(real64), dimension(:), intent(in) | xa, | ||
real(real64), intent(in) | xx, | ||
real(real64), dimension(:), intent(out) | cc | ||
) |
logical pure function, public math_oct_m::even | ( | integer, intent(in) | n | ) |
logical pure function, public math_oct_m::odd | ( | integer, intent(in) | n | ) |
subroutine, public math_oct_m::cartesian2hyperspherical | ( | real(real64), dimension(:), intent(in) | x, |
real(real64), dimension(:), intent(out) | u | ||
) |
subroutine, public math_oct_m::hyperspherical2cartesian | ( | real(real64), dimension(:), intent(in) | u, |
real(real64), dimension(:), intent(out) | x | ||
) |
subroutine, public math_oct_m::hypersphere_grad_matrix | ( | real(real64), dimension(:,:), intent(out) | grad_matrix, |
real(real64), intent(in) | r, | ||
real(real64), dimension(:), intent(in) | x | ||
) |
|
private |
|
private |
|
private |
|
private |
integer pure function, public math_oct_m::pad_pow2 | ( | integer, intent(in) | size | ) |
|
private |
|
private |
|
private |
complex(real64) pure function, public math_oct_m::exponential | ( | complex(real64), intent(in) | z | ) |
complex(real64) pure function, public math_oct_m::phi1 | ( | complex(real64), intent(in) | z | ) |
complex(real64) pure function, public math_oct_m::phi2 | ( | complex(real64), intent(in) | z | ) |
Compute phi2(z) = (phi1(z)-1)/z = (exp(z) - z - 1)/z^2.
Use the power series \( \sum^{4}_{k=0} z^k/(k+2)! \) for small arguments to avoid round-off errors. The expression and the cut value are similar to those from [GSL](https: (there is a difference by a factor of 2 between the two definitions).
real(real64) pure function, public math_oct_m::square_root | ( | real(real64), intent(in) | x | ) |
logical function, public math_oct_m::is_prime | ( | integer, intent(in) | n | ) |
subroutine, public math_oct_m::generate_rotation_matrix | ( | real(real64), dimension(:,:), intent(out) | R, |
real(real64), dimension(:), intent(in) | ff, | ||
real(real64), dimension(:), intent(in) | tt | ||
) |
subroutine, public math_oct_m::numder_ridders | ( | real(real64), intent(in) | x, |
real(real64), intent(in) | h, | ||
real(real64), intent(out) | res, | ||
real(real64), intent(out) | err, | ||
f | |||
) |
Numerical derivative (Ridder`s algorithm).
This is an alternative to "loct_numerical_derivative" (which is just an interface to the GSL numerical derivative). This version is an implementation of Ridders algorithm [C. J. F. Ridders, Adv. Eng. Software 4, 75 (1982); also described in Numerical Recipes]. It is more precise, but also typically more expensive, than the simpler 4-point algorithm implemented in the GSL library.
pure complex(real64) function, dimension(1:3), public math_oct_m::dzcross_product | ( | real(real64), dimension(:), intent(in) | a, |
complex(real64), dimension(:), intent(in) | b | ||
) |
pure complex(real64) function, dimension(1:3), public math_oct_m::zdcross_product | ( | complex(real64), dimension(:), intent(in) | a, |
real(real64), dimension(:), intent(in) | b | ||
) |
subroutine, public math_oct_m::generalized_laguerre_polynomial | ( | integer, intent(in) | np, |
integer, intent(in) | nn, | ||
integer, intent(in) | mm, | ||
real(real64), dimension(np), intent(in) | xx, | ||
real(real64), dimension(np), intent(out) | cx | ||
) |
|
private |
|
private |
|
private |
|
private |
subroutine, public math_oct_m::dsymmetrize_matrix | ( | integer, intent(in) | nn, |
real(real64), dimension(:, :), intent(inout) | aa | ||
) |
subroutine, public math_oct_m::dzero_small_elements_matrix | ( | integer, intent(in) | nn, |
real(real64), dimension(:, :), intent(inout) | aa, | ||
real(real64) | tol | ||
) |
|
private |
|
private |
|
private |
|
private |
pure complex(real64) function, dimension(1:3), public math_oct_m::zcross_product | ( | complex(real64), dimension(:), intent(in) | a, |
complex(real64), dimension(:), intent(in) | b | ||
) |
|
private |
|
private |
|
private |
|
private |