|
FGSL
Fortran interface for the GNU scientific library
|

Functions/Subroutines | |
| type(fgsl_vector) function | fgsl_vector_init (type) |
| Initialize a GSL vector object. This is invoked via the generic fgsl_vector_init. More... | |
| integer(fgsl_int) function | fgsl_vector_align (array, len, fvec, size, offset, stride) |
| Wrap a rank 1 Fortran array slice inside a double precision real GSL vector object. This is invoked via the generic fgsl_vector_align. More... | |
| integer(fgsl_int) function | fgsl_vector_pointer_align (ptr, fvec) |
Associate a Fortran pointer with the data stored inside a GSL vector object. This is invoked via the generic fgsl_vector_align. Objects of type gsl_vector which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist. More... | |
| subroutine | fgsl_vector_to_array (result, source) |
| The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a GSL vector into a Fortran array. More... | |
| subroutine | fgsl_vector_free (fvec) |
| Free the resources inside a GSL vector object previously established by a call to fgsl_vector_init(). This is invoked via the generic fgsl_vector_free. More... | |
| subroutine | fgsl_vector_c_ptr (res, src) |
| logical function | fgsl_vector_status (vector) |
| integer(fgsl_size_t) function | fgsl_sizeof_vector (w) |
| Inquire the size of a double precision real GSL vector object. More... | |
| type(fgsl_vector_complex) function | fgsl_vector_complex_init (type) |
| Initialize a complex GSL vector object. This is invoked via the generic fgsl_vector_init. More... | |
| integer(fgsl_int) function | fgsl_vector_complex_align (array, len, fvec, size, offset, stride) |
| Wrap a rank 1 Fortran array slice inside a double precision complex real GSL vector object. This is invoked via the generic fgsl_vector_align. More... | |
| integer(fgsl_int) function | fgsl_vector_complex_pointer_align (ptr, fvec) |
Associate a Fortran pointer with the data stored inside a GSL vector object. This is invoked via the generic fgsl_vector_align. Objects of type gsl_vector_complex which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist. More... | |
| subroutine | fgsl_vector_complex_to_array (result, source) |
| The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a complex GSL vector into a Fortran array. More... | |
| subroutine | fgsl_vector_complex_free (fvec) |
| Free the resources inside a complex GSL vector object previously established by a call to fgsl_vector_complex_init(). This is invoked via the generic fgsl_vector_free. More... | |
| subroutine | fgsl_vector_complex_c_ptr (res, src) |
| logical function | fgsl_vector_complex_status (vector_complex) |
| integer(fgsl_size_t) function | fgsl_sizeof_vector_complex (w) |
| Inquire the size of a double precision complex GSL vector object. More... | |
| type(fgsl_matrix) function | fgsl_matrix_init (type) |
| Initialize a GSL matrix object. This is invoked via the generic fgsl_matrix_init. More... | |
| integer(fgsl_int) function | fgsl_matrix_align (array, lda, n, m, fmat) |
| Wrap a rank 2 Fortran array inside a double precision real GSL matrix object. This is invoked via the generic fgsl_matrix_align. More... | |
| integer(fgsl_int) function | fgsl_matrix_pointer_align (ptr, fmat) |
Associate a Fortran pointer with the data stored inside a GSL matrix object. This is invoked via the generic fgsl_matrix_align. Objects of type gsl_matrix which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist. More... | |
| subroutine | fgsl_matrix_to_array (result, source) |
| The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a GSL matrix into a rank 2 Fortran array. More... | |
| subroutine | fgsl_matrix_free (fvec) |
| Free the resources inside a GSL matrix object previously established by a call to fgsl_matrix_init(). This is invoked via the generic fgsl_matrix_free. More... | |
| subroutine | fgsl_matrix_c_ptr (res, src) |
| logical function | fgsl_matrix_status (matrix) |
| integer(fgsl_size_t) function | fgsl_sizeof_matrix (w) |
| Inquire the number of elements in a double precision real GSL matrix object. More... | |
| type(fgsl_matrix_complex) function | fgsl_matrix_complex_init (type) |
| Initialize a GSL matrix object. This is invoked via the generic fgsl_matrix_init. More... | |
| integer(fgsl_int) function | fgsl_matrix_complex_align (array, lda, n, m, fmat) |
| Wrap a rank 2 Fortran array inside a double precision complex GSL matrix object. This is invoked via the generic fgsl_matrix_align. More... | |
| integer(fgsl_int) function | fgsl_matrix_complex_pointer_align (ptr, fmat) |
Associate a Fortran pointer with the data stored inside a complex GSL matrix object. This is invoked via the generic fgsl_matrix_align. Objects of type gsl_matrix_complex which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist. More... | |
| subroutine | fgsl_matrix_complex_to_array (result, source) |
| The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a complex GSL matrix into a rank 2 Fortran array. More... | |
| subroutine | fgsl_matrix_complex_free (fvec) |
| Free the resources inside a complex GSL matrix object previously established by a call to fgsl_matrix_complex_init(). This is invoked via the generic fgsl_matrix_free. More... | |
| subroutine | fgsl_matrix_complex_c_ptr (res, src) |
| logical function | fgsl_matrix_complex_status (matrix_complex) |
| integer(fgsl_size_t) function | fgsl_sizeof_matrix_complex (w) |
| Inquire the number of elements in a double precision complex GSL matrix object. More... | |
| integer(fgsl_size_t) function | fgsl_vector_get_size (vec) |
| integer(fgsl_size_t) function | fgsl_vector_get_stride (vec) |
| integer(fgsl_size_t) function | fgsl_matrix_get_size1 (matr) |
| integer(fgsl_size_t) function | fgsl_matrix_get_size2 (matr) |
| integer(fgsl_size_t) function | fgsl_matrix_get_tda (matr) |
| integer(fgsl_int) function fgsl_matrix_align | ( | real(fgsl_double), dimension(lda, m), intent(in), target | array, |
| integer(fgsl_size_t), intent(in) | lda, | ||
| integer(fgsl_size_t), intent(in) | n, | ||
| integer(fgsl_size_t), intent(in) | m, | ||
| type(fgsl_matrix), intent(inout) | fmat | ||
| ) |
Wrap a rank 2 Fortran array inside a double precision real GSL matrix object. This is invoked via the generic fgsl_matrix_align.
| array | - requires the actual argument to have the TARGET attribute. Otherwise being passed by reference is not guaranteed by the Fortran standard. |
| lda | - leading dimension of the rank 2 array |
| n | - number of rows in array |
| m | - number of columns in array |
| fmat | - previously initialized double precision GSL matrix object |
| subroutine fgsl_matrix_c_ptr | ( | type(fgsl_matrix), intent(out) | res, |
| type(c_ptr), intent(in) | src | ||
| ) |
| integer(fgsl_int) function fgsl_matrix_complex_align | ( | complex(fgsl_double_complex), dimension(lda, m), intent(in), target | array, |
| integer(fgsl_size_t), intent(in) | lda, | ||
| integer(fgsl_size_t), intent(in) | n, | ||
| integer(fgsl_size_t), intent(in) | m, | ||
| type(fgsl_matrix_complex), intent(inout) | fmat | ||
| ) |
Wrap a rank 2 Fortran array inside a double precision complex GSL matrix object. This is invoked via the generic fgsl_matrix_align.
| array | - requires the actual argument to have the TARGET attribute. Otherwise being passed by reference is not guaranteed by the Fortran standard. |
| lda | - leading dimension of the rank 2 array |
| n | - number of rows in array |
| m | - number of columns in array |
| fmat | - previously initialized double precision complex GSL matrix object |
| subroutine fgsl_matrix_complex_c_ptr | ( | type(fgsl_matrix_complex), intent(out) | res, |
| type(c_ptr), intent(in) | src | ||
| ) |
| subroutine fgsl_matrix_complex_free | ( | type(fgsl_matrix_complex), intent(inout) | fvec | ) |
Free the resources inside a complex GSL matrix object previously established by a call to fgsl_matrix_complex_init(). This is invoked via the generic fgsl_matrix_free.
| type(fgsl_matrix_complex) function fgsl_matrix_complex_init | ( | complex(fgsl_double_complex), intent(in) | type | ) |
Initialize a GSL matrix object. This is invoked via the generic fgsl_matrix_init.
| type | - determine intrinsic type of vector object |
| integer(fgsl_int) function fgsl_matrix_complex_pointer_align | ( | complex(fgsl_double_complex), dimension(:,:), intent(out), pointer | ptr, |
| type(fgsl_matrix_complex), intent(in) | fmat | ||
| ) |
Associate a Fortran pointer with the data stored inside a complex GSL matrix object. This is invoked via the generic fgsl_matrix_align. Objects of type gsl_matrix_complex which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist.
| ptr | - rank 2 Fortran pointer |
| fmat | - double precision complex GSL matrix |
| logical function fgsl_matrix_complex_status | ( | type(fgsl_matrix_complex), intent(in) | matrix_complex | ) |
| subroutine fgsl_matrix_complex_to_array | ( | complex(fgsl_double_complex), dimension(:,:), intent(inout) | result, |
| type(fgsl_matrix_complex), intent(in) | source | ||
| ) |
The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a complex GSL matrix into a rank 2 Fortran array.
| subroutine fgsl_matrix_free | ( | type(fgsl_matrix), intent(inout) | fvec | ) |
Free the resources inside a GSL matrix object previously established by a call to fgsl_matrix_init(). This is invoked via the generic fgsl_matrix_free.
| integer(fgsl_size_t) function fgsl_matrix_get_size1 | ( | type(fgsl_matrix), intent(in) | matr | ) |
| integer(fgsl_size_t) function fgsl_matrix_get_size2 | ( | type(fgsl_matrix), intent(in) | matr | ) |
| integer(fgsl_size_t) function fgsl_matrix_get_tda | ( | type(fgsl_matrix), intent(in) | matr | ) |
| type(fgsl_matrix) function fgsl_matrix_init | ( | real(fgsl_double), intent(in) | type | ) |
Initialize a GSL matrix object. This is invoked via the generic fgsl_matrix_init.
| type | - determine intrinsic type of vector object |
| integer(fgsl_int) function fgsl_matrix_pointer_align | ( | real(fgsl_double), dimension(:,:), intent(out), pointer | ptr, |
| type(fgsl_matrix), intent(in) | fmat | ||
| ) |
Associate a Fortran pointer with the data stored inside a GSL matrix object. This is invoked via the generic fgsl_matrix_align. Objects of type gsl_matrix which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist.
| ptr | - rank 2 Fortran pointer |
| fmat | - double precision real GSL matrix |
| logical function fgsl_matrix_status | ( | type(fgsl_matrix), intent(in) | matrix | ) |
| subroutine fgsl_matrix_to_array | ( | real(fgsl_double), dimension(:,:), intent(inout) | result, |
| type(fgsl_matrix), intent(in) | source | ||
| ) |
The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a GSL matrix into a rank 2 Fortran array.
| integer(fgsl_size_t) function fgsl_sizeof_matrix | ( | type(fgsl_matrix), intent(in) | w | ) |
Inquire the number of elements in a double precision real GSL matrix object.
| integer(fgsl_size_t) function fgsl_sizeof_matrix_complex | ( | type(fgsl_matrix_complex), intent(in) | w | ) |
Inquire the number of elements in a double precision complex GSL matrix object.
| integer(fgsl_size_t) function fgsl_sizeof_vector | ( | type(fgsl_vector), intent(in) | w | ) |
Inquire the size of a double precision real GSL vector object.
| integer(fgsl_size_t) function fgsl_sizeof_vector_complex | ( | type(fgsl_vector_complex), intent(in) | w | ) |
Inquire the size of a double precision complex GSL vector object.
| integer(fgsl_int) function fgsl_vector_align | ( | real(fgsl_double), dimension(len), intent(in), target | array, |
| integer(fgsl_size_t), intent(in) | len, | ||
| type(fgsl_vector), intent(inout) | fvec, | ||
| integer(fgsl_size_t), intent(in) | size, | ||
| integer(fgsl_size_t), intent(in) | offset, | ||
| integer(fgsl_size_t), intent(in) | stride | ||
| ) |
Wrap a rank 1 Fortran array slice inside a double precision real GSL vector object. This is invoked via the generic fgsl_vector_align.
| array | - requires the actual argument to have the TARGET attribute. Otherwise being passed by reference is not guaranteed by the Fortran standard. |
| len | - number of elements of the rank 1 array |
| fvec | - previously initialized GSL vector object |
| size | - number of elements from array wrapped inside fvec |
| offset | - index of first element of array to be mapped to fvec |
| stride | - stride in array for successive elements of fvec |
| subroutine fgsl_vector_c_ptr | ( | type(fgsl_vector), intent(out) | res, |
| type(c_ptr), intent(in) | src | ||
| ) |
| integer(fgsl_int) function fgsl_vector_complex_align | ( | complex(fgsl_double_complex), dimension(len), intent(in), target | array, |
| integer(fgsl_size_t), intent(in) | len, | ||
| type(fgsl_vector_complex), intent(inout) | fvec, | ||
| integer(fgsl_size_t), intent(in) | size, | ||
| integer(fgsl_size_t), intent(in) | offset, | ||
| integer(fgsl_size_t), intent(in) | stride | ||
| ) |
Wrap a rank 1 Fortran array slice inside a double precision complex real GSL vector object. This is invoked via the generic fgsl_vector_align.
| array | - requires the actual argument to have the TARGET attribute. Otherwise being passed by reference is not guaranteed by the Fortran standard. |
| len | - number of elements of the rank 1 array |
| fvec | - previously initialized complex GSL vector object |
| size | - number of elements from array wrapped inside fvec |
| offset | - index of first element of array to be mapped to fvec |
| stride | - stride in array for successive elements of fvec |
| subroutine fgsl_vector_complex_c_ptr | ( | type(fgsl_vector_complex), intent(out) | res, |
| type(c_ptr), intent(in) | src | ||
| ) |
| subroutine fgsl_vector_complex_free | ( | type(fgsl_vector_complex), intent(inout) | fvec | ) |
Free the resources inside a complex GSL vector object previously established by a call to fgsl_vector_complex_init(). This is invoked via the generic fgsl_vector_free.
| type(fgsl_vector_complex) function fgsl_vector_complex_init | ( | complex(fgsl_double_complex), intent(in) | type | ) |
Initialize a complex GSL vector object. This is invoked via the generic fgsl_vector_init.
| type | - determine intrinsic type of vector object |
| integer(fgsl_int) function fgsl_vector_complex_pointer_align | ( | complex(fgsl_double_complex), dimension(:), intent(out), pointer | ptr, |
| type(fgsl_vector_complex), intent(in) | fvec | ||
| ) |
Associate a Fortran pointer with the data stored inside a GSL vector object. This is invoked via the generic fgsl_vector_align. Objects of type gsl_vector_complex which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist.
| ptr | - rank 1 Fortran pointer |
| fvec | - double precision complex GSL vector |
| logical function fgsl_vector_complex_status | ( | type(fgsl_vector_complex), intent(in) | vector_complex | ) |
| subroutine fgsl_vector_complex_to_array | ( | complex(fgsl_double_complex), dimension(:), intent(inout) | result, |
| type(fgsl_vector_complex), intent(in) | source | ||
| ) |
The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a complex GSL vector into a Fortran array.
| subroutine fgsl_vector_free | ( | type(fgsl_vector), intent(inout) | fvec | ) |
Free the resources inside a GSL vector object previously established by a call to fgsl_vector_init(). This is invoked via the generic fgsl_vector_free.
| integer(fgsl_size_t) function fgsl_vector_get_size | ( | type(fgsl_vector), intent(in) | vec | ) |
| integer(fgsl_size_t) function fgsl_vector_get_stride | ( | type(fgsl_vector), intent(in) | vec | ) |
| type(fgsl_vector) function fgsl_vector_init | ( | real(fgsl_double), intent(in) | type | ) |
Initialize a GSL vector object. This is invoked via the generic fgsl_vector_init.
| type | - determine intrinsic type of vector object |
| integer(fgsl_int) function fgsl_vector_pointer_align | ( | real(fgsl_double), dimension(:), intent(out), pointer | ptr, |
| type(fgsl_vector), intent(in) | fvec | ||
| ) |
Associate a Fortran pointer with the data stored inside a GSL vector object. This is invoked via the generic fgsl_vector_align. Objects of type gsl_vector which are returned by GSL routines often are persistent subobjects of other GSL objects. A Fortran pointer aligned with a subobject hence will remain up-to-date throughout the lifetime of the object; it may become undefined once the object ceases to exist.
| ptr | - rank 1 Fortran pointer |
| fvec | - double precision real GSL vector |
| logical function fgsl_vector_status | ( | type(fgsl_vector), intent(in) | vector | ) |
| subroutine fgsl_vector_to_array | ( | real(fgsl_double), dimension(:), intent(inout) | result, |
| type(fgsl_vector), intent(in) | source | ||
| ) |
The assignment operator (see interface/generics.finc) is overloaded to enable copying of the content of a GSL vector into a Fortran array.
1.8.9.1