Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Sparse matrix linear solver utilities.
This part of selalib is derived from SPM library developed by Ahmed Ratnani (http://ratnani.org/spm_doc/html/)
Michel Mehrenberger did the implementation of UMFPACK option
Derived types and interfaces | |
type | sll_t_csr_matrix |
type for CSR format More... | |
interface | sll_f_new_csr_matrix |
interface | sll_o_delete |
Functions/Subroutines | |
subroutine, public | sll_s_free_csr_matrix (mat) |
type(sll_t_csr_matrix) function, pointer | new_csr_matrix_with_dof (num_rows, num_cols, num_elts, local_to_global_row, num_local_dof_row, local_to_global_col, num_local_dof_col, solver) |
allocates the memory space for a new CSR matrix, More... | |
subroutine | sll_s_csr_matrix_init (mat, num_rows, num_cols, num_elts, local_to_global_row, num_local_dof_row, local_to_global_col, num_local_dof_col) |
initialization of CSR matrix type thanks to the global index of each local dof of each element More... | |
subroutine, public | sll_s_csr_matrix_with_constraint_init (mat, mat_a) |
type(sll_t_csr_matrix) function, pointer, public | sll_f_new_csr_matrix_with_constraint (mat_a) |
subroutine, public | sll_s_csr_add_one_constraint (ia_in, ja_in, a_in, num_rows_in, num_nz_in, constraint_vec, ia_out, ja_out, a_out) |
subroutine, public | sll_s_factorize_csr_matrix (mat) |
subroutine, public | sll_s_mult_csr_matrix_vector (mat, input, output) |
subroutine, public | sll_s_add_to_csr_matrix (mat, val, row, col) |
subroutine, public | sll_s_solve_csr_matrix (mat, b, u) |
subroutine, public | sll_s_solve_csr_matrix_perper (mat, b, u, masse_tot) |
subroutine, public | sll_s_csr_todense (mat, dense_matrix) |
subroutine | solve_with_cg (mat, b, u) |
subroutine | print_solver_type (mat) |
subroutine, public | sll_s_delete_csr_matrix (csr_mat) |
type(sll_t_csr_matrix) function, pointer, public | sll_f_new_csr_matrix (num_rows, num_cols, num_elts, local_to_global_row, num_local_dof_row, local_to_global_col, num_local_dof_col) |
allocates the memory space for a new CSR type on the heap, initializes it with the given arguments and returns a pointer to the object. param[in] num_rows : number of rows param[in] num_cols : number of columns param[in] num_element : number of elements param[in] local_to_global_row : local_to_global_row(\ell,i) gives the global row index of the matrix, for the element i and local degree of freedom \ell param[in] num_local_dof_row : number of local degrees of freedom for the rows param[in] local_to_global_col : local_to_global_col(\ell,i) gives the global column index of the matrix, for the element i and local degree of freedom \ell param[in] num_local_dof_col : number of local degrees of freedom for the columns return a pointer to the newly allocated object. More... | |
subroutine, public | sll_s_csr_matrix_init_with_constraint (mat, mat_a) |
subroutine, public | sll_s_csr_add_one_constraint (ia_in, ja_in, a_in, num_rows_in, num_nz_in, constraint_vec, ia_out, ja_out, a_out) |
subroutine, public | sll_s_mult_csr_matrix_vector (mat, input, output) |
subroutine | set_values_csr_matrix (mat, val) |
subroutine, public | sll_s_solve_csr_matrix (mat, apr_B, apr_U) |
integer(kind=i32) function | sll_count_non_zero_elts (ai_nR, ai_nC, ai_nel, api_LM_1, ai_nen_1, api_LM_2, ai_nen_2, api_columns, api_occ) |
subroutine | init_sparsematrix (ai_nel, api_LM_1, ai_nen_1, api_LM_2, ai_nen_2, api_columns, api_occ, row_ptr, col_ind, num_rows) |
Variables | |
integer(kind=i32), parameter, public | sll_p_umfpack = 1 |
integer(kind=i32), parameter, public | sll_p_pastix = 2 |
integer(kind=i32), parameter, public | sll_p_mumps = 3 |
|
private |
Definition at line 484 of file sll_m_sparse_matrix_pastix.F90.
|
private |
allocates the memory space for a new CSR matrix,
initializes it with the given arguments and returns a pointer to the object.
[in] | num_rows | : number of rows |
[in] | num_cols | : number of columns |
[in] | num_element | : number of elements |
[in] | local_to_global_row | : local_to_global_row(\ell,i) gives the global row index of the matrix, for the element i and local degree of freedom \ell |
[in] | num_local_dof_row | : number of local degrees of freedom for the rows |
[in] | local_to_global_col | : local_to_global_col(\ell,i) gives the global column index of the matrix, for the element i and local degree of freedom \ell |
[in] | num_local_dof_col | : number of local degrees of freedom for the columns |
Definition at line 123 of file sll_m_sparse_matrix.F90.
|
private |
|
private |
Definition at line 368 of file sll_m_sparse_matrix_pastix.F90.
|
private |
Definition at line 399 of file sll_m_sparse_matrix_pastix.F90.
type(sll_t_csr_matrix) function, pointer, public sll_m_sparse_matrix::sll_f_new_csr_matrix | ( | integer(kind=i32), intent(in) | num_rows, |
integer(kind=i32), intent(in) | num_cols, | ||
integer(kind=i32), intent(in) | num_elts, | ||
integer(kind=i32), dimension(:, :), intent(in) | local_to_global_row, | ||
integer(kind=i32), intent(in) | num_local_dof_row, | ||
integer(kind=i32), dimension(:, :), intent(in) | local_to_global_col, | ||
integer(kind=i32), intent(in) | num_local_dof_col | ||
) |
allocates the memory space for a new CSR type on the heap, initializes it with the given arguments and returns a pointer to the object. param[in] num_rows : number of rows param[in] num_cols : number of columns param[in] num_element : number of elements param[in] local_to_global_row : local_to_global_row(\ell,i) gives the global row index of the matrix, for the element i and local degree of freedom \ell param[in] num_local_dof_row : number of local degrees of freedom for the rows param[in] local_to_global_col : local_to_global_col(\ell,i) gives the global column index of the matrix, for the element i and local degree of freedom \ell param[in] num_local_dof_col : number of local degrees of freedom for the columns return a pointer to the newly allocated object.
Definition at line 91 of file sll_m_sparse_matrix_pastix.F90.
type(sll_t_csr_matrix) function, pointer, public sll_m_sparse_matrix::sll_f_new_csr_matrix_with_constraint | ( | type(sll_t_csr_matrix) | mat_a | ) |
subroutine public sll_s_add_to_csr_matrix | ( | type(sll_t_csr_matrix), intent(inout) | mat, |
real(kind=f64), intent(in) | val, | ||
integer(kind=i32), intent(in) | row, | ||
integer(kind=i32), intent(in) | col | ||
) |
Definition at line 573 of file sll_m_sparse_matrix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_csr_add_one_constraint | ( | integer, dimension(:), intent(in) | ia_in, |
integer, dimension(:), intent(in) | ja_in, | ||
real(8), dimension(:), intent(in) | a_in, | ||
integer, intent(in) | num_rows_in, | ||
integer, intent(in) | num_nz_in, | ||
real(8), dimension(:), intent(in) | constraint_vec, | ||
integer, dimension(:), intent(out) | ia_out, | ||
integer, dimension(:), intent(out) | ja_out, | ||
real(8), dimension(:), intent(out) | a_out | ||
) |
Definition at line 378 of file sll_m_sparse_matrix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_csr_add_one_constraint | ( | integer(kind=i32), dimension(:), intent(in) | ia_in, |
integer(kind=i32), dimension(:), intent(in) | ja_in, | ||
real(8), dimension(:), intent(in) | a_in, | ||
integer(kind=i32), intent(in) | num_rows_in, | ||
integer(kind=i32), intent(in) | num_nz_in, | ||
real(8), dimension(:), intent(in) | constraint_vec, | ||
integer(kind=i32), dimension(:), intent(out) | ia_out, | ||
integer(kind=i32), dimension(:), intent(out) | ja_out, | ||
real(8), dimension(:), intent(out) | a_out | ||
) |
Definition at line 233 of file sll_m_sparse_matrix_pastix.F90.
|
private |
initialization of CSR matrix type thanks to the global index of each local dof of each element
initialization of CSR matrix type thanks to the global index of each local dof of each element param[inout] mat : CSR matrix structure param[in] num_rows : number of rows param[in] num_cols : number of columns param[in] num_element : number of elements param[in] local_to_global_row : local_to_global_row(\ell,i) gives the global row index of the matrix, for the element i and local degree of freedom \ell param[in] num_local_dof_row : number of local degrees of freedom for the rows param[in] local_to_global_col : local_to_global_col(\ell,i) gives the global column index of the matrix, for the element i and local degree of freedom \ell param[in] num_local_dof_col : number of local degrees of freedom for the columns
[in,out] | mat | : CSR matrix structure |
[in] | num_rows | : number of rows |
[in] | num_cols | : number of columns |
[in] | num_element | : number of elements |
[in] | local_to_global_row | : local_to_global_row(\ell,i) gives the global row index of the matrix, for the element i and local degree of freedom \ell |
[in] | num_local_dof_row | : number of local degrees of freedom for the rows |
[in] | local_to_global_col | : local_to_global_col(\ell,i) gives the global column index of the matrix, for the element i and local degree of freedom \ell |
[in] | num_local_dof_col | : number of local degrees of freedom for the columns |
Definition at line 185 of file sll_m_sparse_matrix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_csr_matrix_init_with_constraint | ( | type(sll_t_csr_matrix), intent(inout) | mat, |
type(sll_t_csr_matrix), intent(in) | mat_a | ||
) |
Definition at line 203 of file sll_m_sparse_matrix_pastix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_csr_matrix_with_constraint_init | ( | type(sll_t_csr_matrix), intent(inout) | mat, |
type(sll_t_csr_matrix), intent(in) | mat_a | ||
) |
subroutine, public sll_m_sparse_matrix::sll_s_csr_todense | ( | type(sll_t_csr_matrix) | mat, |
real(kind=f64), dimension(:, :) | dense_matrix | ||
) |
Definition at line 697 of file sll_m_sparse_matrix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_delete_csr_matrix | ( | type(sll_t_csr_matrix), pointer | csr_mat | ) |
Definition at line 69 of file sll_m_sparse_matrix_pastix.F90.
subroutine public sll_s_factorize_csr_matrix | ( | type(sll_t_csr_matrix), intent(inout) | mat | ) |
subroutine, public sll_m_sparse_matrix::sll_s_free_csr_matrix | ( | type(sll_t_csr_matrix), pointer | mat | ) |
Definition at line 104 of file sll_m_sparse_matrix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_mult_csr_matrix_vector | ( | type(sll_t_csr_matrix), intent(in) | mat, |
real(kind=f64), dimension(:), intent(in) | input, | ||
real(kind=f64), dimension(:), intent(out) | output | ||
) |
subroutine, public sll_m_sparse_matrix::sll_s_mult_csr_matrix_vector | ( | type(sll_t_csr_matrix) | mat, |
real(kind=f64), dimension(:), intent(in) | input, | ||
real(kind=f64), dimension(:), intent(out) | output | ||
) |
Definition at line 325 of file sll_m_sparse_matrix_pastix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_solve_csr_matrix | ( | type(sll_t_csr_matrix) | mat, |
real(kind=f64), dimension(:) | apr_B, | ||
real(kind=f64), dimension(:) | apr_U | ||
) |
Definition at line 386 of file sll_m_sparse_matrix_pastix.F90.
subroutine, public sll_m_sparse_matrix::sll_s_solve_csr_matrix | ( | type(sll_t_csr_matrix), intent(in) | mat, |
real(kind=f64), dimension(:), intent(inout) | b, | ||
real(kind=f64), dimension(:), intent(out) | u | ||
) |
subroutine public sll_s_solve_csr_matrix_perper | ( | type(sll_t_csr_matrix) | mat, |
real(kind=f64), dimension(:) | b, | ||
real(kind=f64), dimension(:) | u, | ||
real(kind=f64), dimension(:), pointer | masse_tot | ||
) |
|
private |
Definition at line 725 of file sll_m_sparse_matrix.F90.
integer(kind=i32), parameter, public sll_p_mumps = 3 |
Definition at line 94 of file sll_m_sparse_matrix.F90.
integer(kind=i32), parameter, public sll_p_pastix = 2 |
Definition at line 93 of file sll_m_sparse_matrix.F90.
integer(kind=i32), parameter, public sll_p_umfpack = 1 |
Definition at line 92 of file sll_m_sparse_matrix.F90.