Sparse Matrix Representation and Operations. More...
#include "cs_defs.h"
#include <chrono>
#include <stdarg.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <assert.h>
#include <math.h>
#include "bft_mem.h"
#include "bft_error.h"
#include "bft_printf.h"
#include "cs_base.h"
#include "cs_blas.h"
#include "cs_dispatch.h"
#include "cs_halo.h"
#include "cs_halo_perio.h"
#include "cs_log.h"
#include "cs_numbering.h"
#include "cs_prototypes.h"
#include "cs_sort.h"
#include "cs_timer.h"
#include "cs_matrix.h"
#include "cs_matrix_priv.h"
#include "cs_matrix_spmv.h"
Functions | |
cs_matrix_structure_t * | cs_matrix_structure_create (cs_matrix_type_t type, cs_lnum_t n_rows, cs_lnum_t n_cols_ext, cs_lnum_t n_edges, const cs_lnum_2_t *restrict edges, const cs_halo_t *halo, const cs_numbering_t *numbering) |
Create a matrix structure. More... | |
cs_matrix_structure_t * | cs_matrix_structure_create_msr (cs_matrix_type_t type, bool transfer, bool have_diag, bool ordered, cs_lnum_t n_rows, cs_lnum_t n_cols_ext, cs_lnum_t **row_index, cs_lnum_t **col_id, const cs_halo_t *halo, const cs_numbering_t *numbering) |
Create a matrix structure based on a MSR connectivity definition. More... | |
cs_matrix_structure_t * | cs_matrix_structure_create_msr_shared (bool direct_assembly, cs_lnum_t n_rows, cs_lnum_t n_cols_ext, const cs_lnum_t *row_index, const cs_lnum_t *col_id, const cs_halo_t *halo, const cs_numbering_t *numbering) |
Create an MSR matrix structure sharing an existing connectivity definition. More... | |
cs_matrix_structure_t * | cs_matrix_structure_create_from_assembler (cs_matrix_type_t type, cs_matrix_assembler_t *ma) |
Create a matrix structure using a matrix assembler. More... | |
void | cs_matrix_structure_destroy (cs_matrix_structure_t **ms) |
Destroy a matrix structure. More... | |
cs_matrix_t * | cs_matrix_create (const cs_matrix_structure_t *ms) |
Create a matrix container using a given structure. More... | |
cs_matrix_t * | cs_matrix_create_from_assembler (cs_matrix_type_t type, cs_matrix_assembler_t *ma) |
Create a matrix directly from assembler. More... | |
cs_matrix_t * | cs_matrix_create_by_copy (cs_matrix_t *src) |
Create a matrix container by copying another. More... | |
cs_matrix_t * | cs_matrix_create_by_local_restrict (const cs_matrix_t *src) |
Create a matrix based on the local restriction of a base matrix. More... | |
void | cs_matrix_destroy (cs_matrix_t **matrix) |
cs_matrix_type_t | cs_matrix_get_type (const cs_matrix_t *matrix) |
Return matrix type. More... | |
const char * | cs_matrix_get_type_name (const cs_matrix_t *matrix) |
const char * | cs_matrix_get_type_fullname (const cs_matrix_t *matrix) |
cs_lnum_t | cs_matrix_get_n_columns (const cs_matrix_t *matrix) |
Return number of columns in a matrix. More... | |
cs_lnum_t | cs_matrix_get_n_rows (const cs_matrix_t *matrix) |
Return the number of rows in matrix. More... | |
cs_lnum_t | cs_matrix_get_n_entries (const cs_matrix_t *matrix) |
Return the number of entries in matrix. More... | |
cs_lnum_t | cs_matrix_get_diag_block_size (const cs_matrix_t *matrix) |
Return the size of the diagonal block for the given matrix. More... | |
cs_lnum_t | cs_matrix_get_extra_diag_block_size (const cs_matrix_t *matrix) |
Return the size of the extra-diagonal block for the given matrix. More... | |
const cs_matrix_assembler_t * | cs_matrix_get_assembler (const cs_matrix_t *matrix) |
Return matrix assembler if present. More... | |
const cs_halo_t * | cs_matrix_get_halo (const cs_matrix_t *matrix) |
Return the pointer to the halo structure for the given matrix. More... | |
const cs_gnum_t * | cs_matrix_get_l_range (const cs_matrix_t *matrix) |
Return a pointer to local global row range associated with a matrix, if available. More... | |
cs_alloc_mode_t | cs_matrix_get_alloc_mode (const cs_matrix_t *matrix) |
Query matrix allocation mode. More... | |
void | cs_matrix_set_alloc_mode (cs_matrix_t *matrix, cs_alloc_mode_t alloc_mode) |
Set matrix allocation mode. More... | |
bool | cs_matrix_get_need_xa (const cs_matrix_t *matrix) |
Query matrix allocation mode. More... | |
void | cs_matrix_set_need_xa (cs_matrix_t *matrix, bool need_xa) |
Indicate whether matrix will need xa coefficients. More... | |
cs_matrix_fill_type_t | cs_matrix_get_fill_type (bool symmetric, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size) |
Get matrix fill type, depending on block sizes. More... | |
void | cs_matrix_set_coefficients (cs_matrix_t *matrix, bool symmetric, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size, const cs_lnum_t n_edges, const cs_lnum_2_t edges[], const cs_real_t *da, const cs_real_t *xa) |
Set matrix coefficients defined relative to a "native" edge graph, sharing arrays with the caller when possible. More... | |
void | cs_matrix_transfer_coefficients (cs_matrix_t *matrix, bool symmetric, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size, const cs_lnum_t n_edges, const cs_lnum_2_t edges[], cs_real_t **d_val, cs_real_t **x_val) |
void | cs_matrix_transfer_coefficients_msr (cs_matrix_t *matrix, bool symmetric, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size, const cs_lnum_t row_index[], const cs_lnum_t col_id[], cs_real_t **d_val, cs_real_t **x_val) |
Set matrix coefficients in an MSR format, transferring the property of those arrays to the matrix. More... | |
void | cs_matrix_get_coefficients_msr_w (cs_matrix_t *matrix, bool symmetric, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size, cs_real_t **d_val, cs_real_t **e_val) |
Directly access matrix coefficients in an MSR format for writing. More... | |
void | cs_matrix_release_coefficients (cs_matrix_t *matrix) |
Release shared matrix coefficients. More... | |
cs_matrix_assembler_values_t * | cs_matrix_assembler_values_init (cs_matrix_t *matrix, cs_lnum_t diag_block_size, cs_lnum_t extra_diag_block_size) |
Create and initialize a CSR matrix assembler values structure. More... | |
void | cs_matrix_copy_diagonal (const cs_matrix_t *matrix, cs_real_t *restrict da) |
Copy matrix diagonal values. More... | |
bool | cs_matrix_is_symmetric (const cs_matrix_t *matrix) |
Query matrix coefficients symmetry. More... | |
bool | cs_matrix_is_mapped_from_native (const cs_matrix_t *matrix) |
Indicate whether coefficients were mapped from native face-based arrays. More... | |
const cs_real_t * | cs_matrix_get_diagonal (const cs_matrix_t *matrix) |
Get matrix diagonal values. More... | |
const cs_real_t * | cs_matrix_get_extra_diagonal (const cs_matrix_t *matrix) |
Get pointer to matrix extra-diagonal values in "native" format. More... | |
void | cs_matrix_row_init (cs_matrix_row_info_t *r) |
Initialize row info for a given matrix. More... | |
void | cs_matrix_row_finalize (cs_matrix_row_info_t *r) |
Finalize row info for a given matrix. More... | |
void | cs_matrix_get_row (const cs_matrix_t *matrix, const cs_lnum_t row_id, cs_matrix_row_info_t *r) |
Get row values for a given matrix. More... | |
void | cs_matrix_get_native_arrays (const cs_matrix_t *matrix, bool *symmetric, cs_lnum_t *n_edges, const cs_lnum_2_t **edges, const cs_real_t **d_val, const cs_real_t **x_val) |
Get arrays describing a matrix in native format. More... | |
void | cs_matrix_get_csr_arrays (const cs_matrix_t *matrix, const cs_lnum_t **row_index, const cs_lnum_t **col_id, const cs_real_t **val) |
Get arrays describing a matrix in CSR format. More... | |
void | cs_matrix_get_msr_arrays (const cs_matrix_t *matrix, const cs_lnum_t **row_index, const cs_lnum_t **col_id, const cs_real_t **d_val, const cs_real_t **x_val) |
Get arrays describing a matrix in MSR format. More... | |
void | cs_matrix_set_mesh_association (cs_matrix_t *matrix, const cs_lnum_t *c2f_idx, const cs_lnum_t *c2f, const short int *c2f_sgn, const cs_real_3_t *cell_cen, const cs_real_t *cell_vol, const cs_real_3_t *face_normal) |
Associate mesh information with a matrix. More... | |
void | cs_matrix_get_mesh_association (const cs_matrix_t *matrix, const cs_lnum_t **c2f_idx, const cs_lnum_t **c2f, const short int **c2f_sgn, const cs_real_3_t **cell_cen, const cs_real_t **cell_vol, const cs_real_3_t **face_normal) |
Query mesh information that me be associated with a matrix. More... | |
void | cs_matrix_vector_multiply (const cs_matrix_t *matrix, cs_real_t *restrict x, cs_real_t *restrict y) |
Matrix.vector product y = A.x. More... | |
void | cs_matrix_vector_multiply_nosync (const cs_matrix_t *matrix, cs_real_t *restrict x, cs_real_t *restrict y) |
Matrix.vector product y = A.x with no prior halo update of x. More... | |
void | cs_matrix_vector_multiply_partial (const cs_matrix_t *matrix, cs_matrix_spmv_type_t op_type, cs_real_t *restrict x, cs_real_t *restrict y) |
Partial matrix.vector product. More... | |
cs_matrix_variant_t * | cs_matrix_variant_create (cs_matrix_t *m) |
Build matrix variant. More... | |
void | cs_matrix_variant_build_list (const cs_matrix_t *m, int *n_variants, cs_matrix_variant_t **m_variant) |
Build list of variants for tuning or testing. More... | |
void | cs_matrix_variant_destroy (cs_matrix_variant_t **mv) |
Destroy a matrix variant structure. More... | |
void | cs_matrix_variant_apply (cs_matrix_t *m, cs_matrix_variant_t *mv) |
Apply a variant to a given matrix. More... | |
void | cs_matrix_variant_apply_tuned (cs_matrix_t *m, cs_matrix_variant_t *mv) |
Apply variants defined by tuning to a given matrix. More... | |
void | cs_matrix_variant_set_func (cs_matrix_variant_t *mv, cs_matrix_fill_type_t fill_type, cs_matrix_spmv_type_t spmv_type, const cs_numbering_t *numbering, const char *func_name) |
Select the sparse matrix-vector product function to be used by a matrix variant for a given fill type. More... | |
cs_matrix_type_t | cs_matrix_variant_type (const cs_matrix_variant_t *mv) |
Get the type associated with a matrix variant. More... | |
Sparse Matrix Representation and Operations.
Notes:
The aim of these structures and associated functions is multiple:
Choice of a given matrix structure may be strongly related to the choice of linear solver or smoother:
Though many external libraries assume a classical CSR or block-CSR structure, separate storage of diagonal elements (a.k.a. MSR format) is very useful for Jacobi and Gauss-Seidel solvers and Jacobi or polynomial preconditioners, in which the diagonal elements need to be accessed and/or multiplied separately, so we select that approach.
In distributed parallel mode, some matrix entries reference distant rows (and ghost values). For a given local row i and column j, if j > n_local_rows, a_ij references a distant value. For the native (graph edge-based) format, symmetric storage is used, but for other (matrix row-based) formats, only local rows are stored (i.e.* a_ij is not stored if i is not in the local row range), so the matrix storage is locally rectangular, even for a logically square matrix.
An additional optimization (also used in some external libraries) is to store the matrix entries referencing distant rows (and ghost values) separately, to simplify handling of computation-communication overlap. The "distributed" format is thus a variation of the MSR format with separate distant elements.
The specific access requirements of Gauss-Seidel solvers and smoothers lead us to only consider the MSR format for their implementation. When requesting a Gauss-Seidel solver or smoother for another storage format, a Jacobi solver or smoother may be substituted.
Please refer to the matrix section of the theory guide for more informations.
cs_matrix_assembler_values_t * cs_matrix_assembler_values_init | ( | cs_matrix_t * | matrix, |
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size | ||
) |
Create and initialize a CSR matrix assembler values structure.
The associated matrix's structure must have been created using cs_matrix_structure_create_from_assembler.
[in,out] | matrix | pointer to matrix structure |
[in] | diag_block_size | block sizes for diagonal |
[in] | extra_diag_block_size | block sizes for extra diagonal |
void cs_matrix_copy_diagonal | ( | const cs_matrix_t * | matrix, |
cs_real_t *restrict | da | ||
) |
Copy matrix diagonal values.
In case of matrixes with block diagonal coefficients, only the true diagonal values are copied.
[in] | matrix | pointer to matrix structure |
[out] | da | diagonal (pre-allocated, size: n_rows*block_size) |
cs_matrix_t * cs_matrix_create | ( | const cs_matrix_structure_t * | ms | ) |
Create a matrix container using a given structure.
Note that the matrix container maps to the assigned structure, so it must be destroyed before that structure.
[in] | ms | associated matrix structure |
cs_matrix_t * cs_matrix_create_by_copy | ( | cs_matrix_t * | src | ) |
Create a matrix container by copying another.
Note that the matrix containers share the same assigned structure, so they must be both destroyed before that structure.
If assigned, coefficients are not copied.
[in] | src | reference matrix structure |
cs_matrix_t * cs_matrix_create_by_local_restrict | ( | const cs_matrix_t * | src | ) |
Create a matrix based on the local restriction of a base matrix.
Coefficients are copied. Some coefficients may be shared with the parent matrix, so the base matrix must not be destroyed before the restriction matrix.
[in] | src | reference matrix structure |
cs_matrix_t * cs_matrix_create_from_assembler | ( | cs_matrix_type_t | type, |
cs_matrix_assembler_t * | ma | ||
) |
Create a matrix directly from assembler.
Only CSR and MSR formats are handled.
[in] | type | type of matrix considered |
[in] | ma | pointer to matrix assembler structure |
void cs_matrix_destroy | ( | cs_matrix_t ** | matrix | ) |
cs_alloc_mode_t cs_matrix_get_alloc_mode | ( | const cs_matrix_t * | matrix | ) |
Query matrix allocation mode.
[in] | matrix | pointer to matrix structure |
const cs_matrix_assembler_t * cs_matrix_get_assembler | ( | const cs_matrix_t * | matrix | ) |
Return matrix assembler if present.
[in] | matrix | pointer to matrix structure |
void cs_matrix_get_coefficients_msr_w | ( | cs_matrix_t * | matrix, |
bool | symmetric, | ||
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size, | ||
cs_real_t ** | d_val, | ||
cs_real_t ** | e_val | ||
) |
Directly access matrix coefficients in an MSR format for writing.
The matrix's fil type is also set by this function.
The associated matrix must be in MSR format, and the associated row index and column ids known.
[in,out] | matrix | pointer to matrix structure |
[in] | symmetric | indicates if matrix coefficients are symmetric |
[in] | diag_block_size | block sizes for diagonal |
[in] | extra_diag_block_size | block sizes for extra diagonal |
[in,out] | d_val | diagonal values (nullptr if zero) |
[in,out] | e_val | extradiagonal values (nullptr if zero) |
void cs_matrix_get_csr_arrays | ( | const cs_matrix_t * | matrix, |
const cs_lnum_t ** | row_index, | ||
const cs_lnum_t ** | col_id, | ||
const cs_real_t ** | val | ||
) |
Get arrays describing a matrix in CSR format.
This function only works for an CSR matrix (i.e. there is no automatic conversion from another matrix type).
Matrix block sizes can be obtained by cs_matrix_get_diag_block_size() and cs_matrix_get_extra_diag_block_size().
[in] | matrix | pointer to matrix structure |
[out] | row_index | CSR row index |
[out] | col_id | CSR column id |
[out] | val | values |
cs_lnum_t cs_matrix_get_diag_block_size | ( | const cs_matrix_t * | matrix | ) |
Return the size of the diagonal block for the given matrix.
[in] | matrix | pointer to matrix structure |
const cs_real_t * cs_matrix_get_diagonal | ( | const cs_matrix_t * | matrix | ) |
Get matrix diagonal values.
In case of matrixes with block diagonal coefficients, a pointer to the complete block diagonal is returned.
[in] | matrix | pointer to matrix structure |
cs_lnum_t cs_matrix_get_extra_diag_block_size | ( | const cs_matrix_t * | matrix | ) |
Return the size of the extra-diagonal block for the given matrix.
[in] | matrix | pointer to matrix structure |
const cs_real_t * cs_matrix_get_extra_diagonal | ( | const cs_matrix_t * | matrix | ) |
Get pointer to matrix extra-diagonal values in "native" format.
This function only functions if the coefficients were mapped from native coefficients using cs_matrix_set_coefficients(), in which case the pointer returned is the same as the one passed to that function.
It is used in the current multigrid code, but should be removed as soon as the dependency to the native format is removed.
[in] | matrix | pointer to matrix structure |
cs_matrix_fill_type_t cs_matrix_get_fill_type | ( | bool | symmetric, |
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size | ||
) |
Get matrix fill type, depending on block sizes.
[in] | symmetric | indicates if matrix coefficients are symmetric |
[in] | diag_block_size | block sizes for diagonal |
[in] | extra_diag_block_size | block sizes for extra diagonal |
const cs_halo_t * cs_matrix_get_halo | ( | const cs_matrix_t * | matrix | ) |
Return the pointer to the halo structure for the given matrix.
[in] | matrix | pointer to matrix structure |
const cs_gnum_t * cs_matrix_get_l_range | ( | const cs_matrix_t * | matrix | ) |
Return a pointer to local global row range associated with a matrix, if available.
[in] | matrix | pointer to matrix structure |
void cs_matrix_get_mesh_association | ( | const cs_matrix_t * | matrix, |
const cs_lnum_t ** | c2f_idx, | ||
const cs_lnum_t ** | c2f, | ||
const short int ** | c2f_sgn, | ||
const cs_real_3_t ** | cell_cen, | ||
const cs_real_t ** | cell_vol, | ||
const cs_real_3_t ** | face_normal | ||
) |
Query mesh information that me be associated with a matrix.
This may be useful for multigrid smoothing.
[in] | matrix | pointer to matrix structure |
[out] | c2f_idx | cell to faces index, or nullptr |
[out] | c2f | cell to faces adjacency, or nullptr |
[out] | c2f_sgn | cell to faces adjacency sign, or nullptr |
[out] | cell_cen | cell center coordinates, or nullptr |
[out] | cell_vol | cell volumes, or nullptr |
[out] | face_normal | face normas, or nullptr |
void cs_matrix_get_msr_arrays | ( | const cs_matrix_t * | matrix, |
const cs_lnum_t ** | row_index, | ||
const cs_lnum_t ** | col_id, | ||
const cs_real_t ** | d_val, | ||
const cs_real_t ** | x_val | ||
) |
Get arrays describing a matrix in MSR format.
This function only works for an MSR matrix (i.e. there is no automatic conversion from another matrix type).
Matrix block sizes can be obtained by cs_matrix_get_diag_block_size() and cs_matrix_get_extra_diag_block_size().
[in] | matrix | pointer to matrix structure |
[out] | row_index | MSR row index |
[out] | col_id | MSR column id |
[out] | d_val | diagonal values |
[out] | x_val | extra-diagonal values |
cs_lnum_t cs_matrix_get_n_columns | ( | const cs_matrix_t * | matrix | ) |
Return number of columns in a matrix.
[in] | matrix | pointer to matrix structure |
cs_lnum_t cs_matrix_get_n_entries | ( | const cs_matrix_t * | matrix | ) |
Return the number of entries in matrix.
When the block size is > 1, the number reported is the number of entry blocks, not individual entries.
[in] | matrix | pointer to matrix structure |
cs_lnum_t cs_matrix_get_n_rows | ( | const cs_matrix_t * | matrix | ) |
Return the number of rows in matrix.
[in] | matrix | pointer to matrix structure |
void cs_matrix_get_native_arrays | ( | const cs_matrix_t * | matrix, |
bool * | symmetric, | ||
cs_lnum_t * | n_edges, | ||
const cs_lnum_2_t ** | edges, | ||
const cs_real_t ** | d_val, | ||
const cs_real_t ** | x_val | ||
) |
Get arrays describing a matrix in native format.
This function works for matrix in native format.
Matrix block sizes can be obtained by cs_matrix_get_diag_block_size() and cs_matrix_get_extra_diag_block_size().
[in] | matrix | pointer to matrix structure |
[out] | symmetric | true if symmetric |
[out] | n_edges | number of associated faces |
[out] | edges | edges (symmetric row <-> column) connectivity |
[out] | d_val | diagonal values |
[out] | x_val | extra-diagonal values |
bool cs_matrix_get_need_xa | ( | const cs_matrix_t * | matrix | ) |
Query matrix allocation mode.
[in] | matrix | pointer to matrix structure |
void cs_matrix_get_row | ( | const cs_matrix_t * | matrix, |
const cs_lnum_t | row_id, | ||
cs_matrix_row_info_t * | r | ||
) |
Get row values for a given matrix.
This function may not work for all matrix types.
In the case of blocked matrices, the true (non-blocked) values are returned.
The row information structure must have been previously initialized using cs_matrix_row_init, and should be finalized using using cs_matrix_row_finalize, so as to free buffers it may have built for certain matrix formats.
[in] | matrix | pointer to matrix structure |
[in] | row_id | id of row to query |
[in,out] | r | row info structure |
cs_matrix_type_t cs_matrix_get_type | ( | const cs_matrix_t * | matrix | ) |
Return matrix type.
[in] | matrix | pointer to matrix structure |
const char * cs_matrix_get_type_fullname | ( | const cs_matrix_t * | matrix | ) |
const char * cs_matrix_get_type_name | ( | const cs_matrix_t * | matrix | ) |
bool cs_matrix_is_mapped_from_native | ( | const cs_matrix_t * | matrix | ) |
Indicate whether coefficients were mapped from native face-based arrays.
It is used in the current multgrid code, but should be removed as soon as the dependency to the native format is removed.
[in] | matrix | pointer to matrix structure |
bool cs_matrix_is_symmetric | ( | const cs_matrix_t * | matrix | ) |
Query matrix coefficients symmetry.
[in] | matrix | pointer to matrix structure |
void cs_matrix_release_coefficients | ( | cs_matrix_t * | matrix | ) |
Release shared matrix coefficients.
Pointers to mapped coefficients are set to nullptr, while coefficient copies owned by the matrix are not modified.
This simply ensures the matrix does not maintain pointers to nonexistant data.
[in,out] | matrix | pointer to matrix structure |
void cs_matrix_row_finalize | ( | cs_matrix_row_info_t * | r | ) |
Finalize row info for a given matrix.
[in,out] | r | row info structure |
void cs_matrix_row_init | ( | cs_matrix_row_info_t * | r | ) |
Initialize row info for a given matrix.
[out] | r | row info structure |
void cs_matrix_set_alloc_mode | ( | cs_matrix_t * | matrix, |
cs_alloc_mode_t | alloc_mode | ||
) |
Set matrix allocation mode.
[in,out] | matrix | pointer to matrix structure |
[in] | alloc_mode | host/device allocation mode |
void cs_matrix_set_coefficients | ( | cs_matrix_t * | matrix, |
bool | symmetric, | ||
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size, | ||
const cs_lnum_t | n_edges, | ||
const cs_lnum_2_t | edges[], | ||
const cs_real_t * | da, | ||
const cs_real_t * | xa | ||
) |
Set matrix coefficients defined relative to a "native" edge graph, sharing arrays with the caller when possible.
With shared arrays, the matrix becomes unusable if the arrays passed as arguments are not be modified (its coefficients should be unset first to mark this).
Depending on current options and initialization, values will be copied or simply mapped.
[in,out] | matrix | pointer to matrix structure |
[in] | symmetric | indicates if matrix coefficients are symmetric |
[in] | diag_block_size | block sizes for diagonal |
[in] | extra_diag_block_size | block sizes for extra diagonal |
[in] | n_edges | local number of graph edges |
[in] | edges | edges (row <-> column) connectivity |
[in] | da | diagonal values (nullptr if zero) |
[in] | xa | extradiagonal values (nullptr if zero) casts as: xa[n_edges] if symmetric, xa[n_edges][2] if non symmetric |
void cs_matrix_set_mesh_association | ( | cs_matrix_t * | matrix, |
const cs_lnum_t * | c2f_idx, | ||
const cs_lnum_t * | c2f, | ||
const short int * | c2f_sgn, | ||
const cs_real_3_t * | cell_cen, | ||
const cs_real_t * | cell_vol, | ||
const cs_real_3_t * | face_normal | ||
) |
Associate mesh information with a matrix.
This may be useful for multigrid smoothing.
At least cell centers and volumes are needed for relaxation, and face adjacency and normals are needed for the "classical" option.
Note that cells and faces here do not need to be primary mesh elements, but could be dual mesh elements of some sort.
The arrays passed to the matrix are shared, so should have a lifetime at least as long as the matrix.
[in,out] | matrix | pointer to matrix structure |
[in] | c2f_idx | cell to faces index, or nullptr |
[in] | c2f | cell to faces adjacency, or nullptr |
[in] | c2f_sgn | cell to faces adjacency sign, or nullptr |
[in] | cell_cen | cell center coordinates |
[in] | cell_vol | cell volumes |
[in] | face_normal | face normal, or nullptr |
void cs_matrix_set_need_xa | ( | cs_matrix_t * | matrix, |
bool | need_xa | ||
) |
Indicate whether matrix will need xa coefficients.
[in,out] | matrix | pointer to matrix structure |
[in] | need_xa | is thr face-based xa array needed ? |
cs_matrix_structure_t * cs_matrix_structure_create | ( | cs_matrix_type_t | type, |
cs_lnum_t | n_rows, | ||
cs_lnum_t | n_cols_ext, | ||
cs_lnum_t | n_edges, | ||
const cs_lnum_2_t *restrict | edges, | ||
const cs_halo_t * | halo, | ||
const cs_numbering_t * | numbering | ||
) |
Create a matrix structure.
Note that the structure created usually maps to the given existing cell global number, face -> cell connectivity arrays, and cell halo structure, so it must be destroyed before they are freed (usually along with the code's main face -> cell structure).
Note that the resulting matrix structure will contain a full main main diagonal, and that the extra-diagonal structure is always symmetric (though the coefficients my not be, and we may choose a matrix format that does not exploit this symmetry). If the edges connectivity argument is nullptr, the matrix will be purely diagonal.
[in] | type | type of matrix considered |
[in] | n_rows | local number of rows |
[in] | n_cols_ext | number of local + ghost columns |
[in] | n_edges | local number of (undirected) graph edges |
[in] | edges | edges (symmetric row <-> column) connectivity |
[in] | halo | halo structure associated with cells, or nullptr |
[in] | numbering | vectorization or thread-related numbering info, or nullptr |
cs_matrix_structure_t * cs_matrix_structure_create_from_assembler | ( | cs_matrix_type_t | type, |
cs_matrix_assembler_t * | ma | ||
) |
Create a matrix structure using a matrix assembler.
Only CSR and MSR formats are handled.
[in] | type | type of matrix considered |
[in] | ma | pointer to matrix assembler structure |
cs_matrix_structure_t * cs_matrix_structure_create_msr | ( | cs_matrix_type_t | type, |
bool | transfer, | ||
bool | have_diag, | ||
bool | ordered, | ||
cs_lnum_t | n_rows, | ||
cs_lnum_t | n_cols_ext, | ||
cs_lnum_t ** | row_index, | ||
cs_lnum_t ** | col_id, | ||
const cs_halo_t * | halo, | ||
const cs_numbering_t * | numbering | ||
) |
Create a matrix structure based on a MSR connectivity definition.
Only CSR and MSR formats are handled.
col_id is sorted row by row during the creation of this structure.
In case the property of the row index and col_id arrays are transferred to the structure, the arrays pointers passed as arguments are set to nullptr, to help ensure the caller does not use the original arrays directly after this call.
[in] | type | type of matrix considered |
[in] | transfer | transfer property of row_index and col_id if true, map them otherwise |
[in] | have_diag | indicates if the structure includes the diagonal (should be the same for all rows) |
[in] | ordered | indicates if column ids are already ordered |
[in] | n_rows | local number of rows |
[in] | n_cols_ext | local number of columns + ghosts |
[in] | row_index | pointer to index on rows |
[in,out] | col_id | pointer to array of column ids related to the row index |
[in] | halo | halo structure for synchronization, or nullptr |
[in] | numbering | vectorization or thread-related numbering info, or nullptr |
cs_matrix_structure_t * cs_matrix_structure_create_msr_shared | ( | bool | direct_assembly, |
cs_lnum_t | n_rows, | ||
cs_lnum_t | n_cols_ext, | ||
const cs_lnum_t * | row_index, | ||
const cs_lnum_t * | col_id, | ||
const cs_halo_t * | halo, | ||
const cs_numbering_t * | numbering | ||
) |
Create an MSR matrix structure sharing an existing connectivity definition.
Note that as the structure created maps to the given existing cell global number, face -> cell connectivity arrays, and cell halo structure, it must be destroyed before they are freed (usually along with the code's main face -> cell structure).
[in] | direct_assembly | true if each value corresponds to a unique face |
[in] | n_rows | local number of rows |
[in] | n_cols_ext | local number of columns + ghosts |
[in] | row_index | index on rows |
[in] | col_id | array of column ids related to the row index |
[in] | halo | halo structure for synchronization, or nullptr |
[in] | numbering | vectorization or thread-related numbering info, or nullptr |
void cs_matrix_structure_destroy | ( | cs_matrix_structure_t ** | ms | ) |
Destroy a matrix structure.
[in,out] | ms | pointer to matrix structure pointer |
void cs_matrix_transfer_coefficients | ( | cs_matrix_t * | matrix, |
bool | symmetric, | ||
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size, | ||
const cs_lnum_t | n_edges, | ||
const cs_lnum_2_t | edges[], | ||
cs_real_t ** | d_val, | ||
cs_real_t ** | x_val | ||
) |
void cs_matrix_transfer_coefficients_msr | ( | cs_matrix_t * | matrix, |
bool | symmetric, | ||
cs_lnum_t | diag_block_size, | ||
cs_lnum_t | extra_diag_block_size, | ||
const cs_lnum_t | row_index[], | ||
const cs_lnum_t | col_id[], | ||
cs_real_t ** | d_val, | ||
cs_real_t ** | x_val | ||
) |
Set matrix coefficients in an MSR format, transferring the property of those arrays to the matrix.
If the matrix is also in MSR format, this avoids an extra copy. If it is in a different format, values are copied to the structure, and the original arrays freed. In any case, the arrays pointers passed as arguments are set to nullptr, to help ensure the caller does not use the original arrays directly after this call.
[in,out] | matrix | pointer to matrix structure |
[in] | symmetric | indicates if matrix coefficients are symmetric |
[in] | diag_block_size | block sizes for diagonal |
[in] | extra_diag_block_size | block sizes for extra diagonal |
[in] | row_index | MSR row index (0 to n-1) |
[in] | col_id | MSR column id (0 to n-1) |
[in,out] | d_val | diagonal values (nullptr if zero) |
[in,out] | x_val | extradiagonal values (nullptr if zero) |
void cs_matrix_variant_apply | ( | cs_matrix_t * | m, |
cs_matrix_variant_t * | mv | ||
) |
Apply a variant to a given matrix.
[in,out] | m | pointer to matrix |
[in] | mv | pointer to matrix variant pointer |
void cs_matrix_variant_apply_tuned | ( | cs_matrix_t * | m, |
cs_matrix_variant_t * | mv | ||
) |
Apply variants defined by tuning to a given matrix.
[in,out] | m | pointer to matrix |
[in] | mv | pointer to matrix variant pointer |
void cs_matrix_variant_build_list | ( | const cs_matrix_t * | m, |
int * | n_variants, | ||
cs_matrix_variant_t ** | m_variant | ||
) |
Build list of variants for tuning or testing.
The matrix coefficients should be assigned, so the fill type can be determined.
[in] | m | associated matrix |
[out] | n_variants | number of variants |
[out] | m_variant | array of matrix variants |
cs_matrix_variant_t * cs_matrix_variant_create | ( | cs_matrix_t * | m | ) |
Build matrix variant.
The variant will initially use default matrix-vector functions, which can be later modified using cs_matrix_variant_set_func().
[in] | m | pointer to matrix |
void cs_matrix_variant_destroy | ( | cs_matrix_variant_t ** | mv | ) |
Destroy a matrix variant structure.
[in,out] | mv | pointer to matrix variant pointer |
void cs_matrix_variant_set_func | ( | cs_matrix_variant_t * | mv, |
cs_matrix_fill_type_t | fill_type, | ||
cs_matrix_spmv_type_t | spmv_type, | ||
const cs_numbering_t * | numbering, | ||
const char * | func_name | ||
) |
Select the sparse matrix-vector product function to be used by a matrix variant for a given fill type.
Currently, possible variant functions are:
CS_MATRIX_NATIVE (all fill types) default baseline omp (for OpenMP with compatible numbering) omp_atomic (for OpenMP with atomics) vector (For vector machine with compatible numbering)
CS_MATRIX_CSR (for CS_MATRIX_SCALAR or CS_MATRIX_SCALAR_SYM) default mkl (with MKL)
CS_MATRIX_MSR (all fill types) default mkl (with MKL, for CS_MATRIX_SCALAR or CS_MATRIX_SCALAR_SYM) omp_sched (For OpenMP with scheduling)
parameters: mv <-> pointer to matrix variant numbering <– mesh numbering info, or nullptr fill type <– matrix fill type to merge from spmv_type <– SpMV operation type (full or sub-matrix) (all types if CS_MATRIX_SPMV_N_TYPES) func_name <– function type name
cs_matrix_type_t cs_matrix_variant_type | ( | const cs_matrix_variant_t * | mv | ) |
Get the type associated with a matrix variant.
[in] | mv | pointer to matrix variant structure |
void cs_matrix_vector_multiply | ( | const cs_matrix_t * | matrix, |
cs_real_t *restrict | x, | ||
cs_real_t *restrict | y | ||
) |
Matrix.vector product y = A.x.
This function includes a halo update of x prior to multiplication by A.
[in] | matrix | pointer to matrix structure |
[in,out] | x | multiplying vector values (ghost values updated) |
[out] | y | resulting vector |
void cs_matrix_vector_multiply_nosync | ( | const cs_matrix_t * | matrix, |
cs_real_t *restrict | x, | ||
cs_real_t *restrict | y | ||
) |
Matrix.vector product y = A.x with no prior halo update of x.
This function does not include a halo update of x prior to multiplication by A, so it should be called only when the halo of x is known to already be up to date (in which case we avoid the performance penalty of a redundant update by using this variant of the matrix.vector product).
[in] | matrix | pointer to matrix structure |
[in] | x | multiplying vector values |
[out] | y | resulting vector |
void cs_matrix_vector_multiply_partial | ( | const cs_matrix_t * | matrix, |
cs_matrix_spmv_type_t | op_type, | ||
cs_real_t *restrict | x, | ||
cs_real_t *restrict | y | ||
) |
Partial matrix.vector product.
This function includes a halo update of x prior to multiplication, except for the CS_MATRIX_SPMV_L operation type, which does not require it, as halo adjacencies are only present and useful in the upper-diagonal part..
[in] | matrix | pointer to matrix structure |
[in] | op_type | SpMV operation type |
[in,out] | x | multiplying vector values (ghost values updated) |
[out] | y | resulting vector |