SpM Handbook 1.2.4
Loading...
Searching...
No Matches
d_spm.h
Go to the documentation of this file.
1/**
2 *
3 * @file d_spm.h
4 *
5 * SParse Matrix package precision dependent header.
6 *
7 * @copyright 2016-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8 * Univ. Bordeaux. All rights reserved.
9 *
10 * @version 1.2.4
11 * @author Pierre Ramet
12 * @author Mathieu Faverge
13 * @author Alban Bellot
14 * @author Tony Delarue
15 * @author Alycia Lisito
16 * @date 2024-05-29
17 *
18 * @generated from /builds/2mk6rsew/0/fpruvost/spm/include/spm/z_spm.h, normal z -> d, Fri Nov 29 11:34:27 2024
19 *
20 * @addtogroup spm_dev_integer
21 * @{
22 *
23 **/
24#ifndef _d_spm_h_
25#define _d_spm_h_
26
27/**
28 * @brief Integer routines
29 */
30void d_spmIntFltSortAsc( void **const pbase,
31 const spm_int_t n );
32void d_spmIntIntFltSortAsc( void ** const pbase,
33 const spm_int_t n );
34
35/**
36 * @}
37 * @addtogroup spm_dev_convert
38 * @{
39 *
40 * @brief Conversion routines
41 */
48
49void d_spm2dense( const spmatrix_t *spm, double *A );
50
51/**
52 * @}
53 * @addtogroup spm_dev_matvec
54 * @{
55 *
56 * @brief Matrix-Vector and matrix-matrix product routines
57 */
58int spm_dspmv( spm_trans_t trans,
59 double alpha,
60 const spmatrix_t *A,
61 const double *x,
62 spm_int_t incx,
63 double beta,
64 double *y,
65 spm_int_t incy );
66int spm_dspmm( spm_side_t side,
67 spm_trans_t transA,
68 spm_trans_t transB,
69 spm_int_t K,
70 double alpha,
71 const spmatrix_t *A,
72 const double *B,
73 spm_int_t ldb,
74 double beta,
75 double *C,
76 spm_int_t ldc );
77
78/**
79 * @}
80 * @addtogroup spm_dev_norm
81 * @{
82 *
83 * @brief Norm computation routines
84 */
85double d_spmNorm( spm_normtype_t ntype,
86 const spmatrix_t *spm );
87double d_spmNormMat( spm_normtype_t ntype,
88 const spmatrix_t *spm,
89 spm_int_t n,
90 const double *A,
91 spm_int_t lda );
92
93/**
94 * @}
95 * @addtogroup spm_dev_check
96 * @{
97 *
98 * @brief Extra routines
99 */
100void d_spmSort( spmatrix_t *spm );
102
103/**
104 * @}
105 * @addtogroup spm_dev_rhs
106 * @{
107 *
108 * @brief RHS routines
109 */
110int d_spmGenRHS( spm_rhstype_t type,
111 int nrhs,
112 const spmatrix_t *spm,
113 void *x,
114 int ldx,
115 void *b,
116 int ldb );
117int d_spmGenMat( spm_rhstype_t type,
118 int nrhs,
119 const spmatrix_t *spm,
120 void *alpha,
121 unsigned long long int seed,
122 void *A,
123 int lda );
125 int nrhs,
126 const spmatrix_t *spm,
127 void *x0,
128 int ldx0,
129 void *b,
130 int ldb,
131 const void *x,
132 int ldx );
133
134void d_spmGatherRHS( int nrhs,
135 const spmatrix_t *spm,
136 const double *x,
137 spm_int_t ldx,
138 int root,
139 double *gx,
140 spm_int_t ldgx );
141
142void d_spmReduceRHS( int nrhs,
143 const spmatrix_t *spm,
144 double *bglob,
145 spm_int_t ldbg,
146 double *bloc,
147 spm_int_t ldbl );
148void d_spmExtractLocalRHS( int nrhs,
149 const spmatrix_t *spm,
150 const double *bglob,
151 spm_int_t ldbg,
152 double *bloc,
153 spm_int_t ldbl );
154
155/**
156 * @}
157 * @addtogroup spm_dev_print
158 * @{
159 *
160 * @brief Output routines
161 */
162void d_spmPrint( FILE *f,
163 const spmatrix_t *spm );
164void d_spmPrintRHS( FILE *f,
165 const spmatrix_t *spm,
166 int nrhs,
167 const void *x,
168 spm_int_t ldx );
169void d_spmDensePrint( FILE *f,
170 spm_int_t m,
171 spm_int_t n,
172 const double *A,
173 spm_int_t lda );
174
175/**
176 * @}
177 * @addtogroup spm_dev_dof
178 * @{
179 *
180 * @brief DOF routines
181 */
182void d_spmExpand( const spmatrix_t *spm_in,
183 spmatrix_t *spm_out );
184void d_spmDofExtend( spmatrix_t *spm );
185
186/**
187 * @}
188 * @addtogroup spm_dev_scal
189 * @{
190 *
191 * @brief Scaling routines
192 */
193void d_spmScal( const double alpha,
194 spmatrix_t *spm );
195
196/**
197 * @}
198 * @addtogroup spm_dev_rhs
199 * @{
200 *
201 * @brief Subroutines for random vector generation to be used in testings
202 */
203int d_spmRhsGenRndShm( const spmatrix_t *spm,
204 double scale,
205 spm_int_t n,
206 double *A,
207 spm_int_t lda,
208 int shift,
209 unsigned long long int seed );
210int d_spmRhsGenRndDist( const spmatrix_t *spm,
211 double scale,
212 spm_int_t n,
213 double *A,
214 spm_int_t lda,
215 int shift,
216 unsigned long long int seed );
217
218/**
219 * @}
220 */
221#endif /* _d_spm_h_ */
double spm_fixdbl_t
Double datatype that is not converted through precision generator functions.
Definition datatypes.h:150
enum spm_rhstype_e spm_rhstype_t
How to generate RHS.
enum spm_normtype_e spm_normtype_t
Norms.
enum spm_trans_e spm_trans_t
Transpostion.
enum spm_side_e spm_side_t
Side of the operation.
void d_spmSort(spmatrix_t *spm)
This routine sorts the spm matrix.
Definition d_spm_sort.c:303
spm_int_t d_spmMergeDuplicate(spmatrix_t *spm)
This routine merge the multiple entries in a sparse matrix by summing their values together.
int d_spmConvertIJV2CSC(spmatrix_t *spm)
convert a matrix in IJV format to a matrix in CSC format.
int d_spmConvertCSR2IJV(spmatrix_t *spm)
convert a matrix in CSR format to a matrix in IJV format.
int d_spmConvertCSR2CSC(spmatrix_t *spm)
convert a matrix in CSR format to a matrix in CSC format.
int d_spmConvertCSC2IJV(spmatrix_t *spm)
Convert a matrix in CSC format to a matrix in IJV format.
int d_spmConvertCSC2CSR(spmatrix_t *spm)
convert a matrix in CSC format to a matrix in CSR format.
int d_spmConvertIJV2CSR(spmatrix_t *spm)
convert a matrix in IJV format to a matrix in CSR format.
void d_spm2dense(const spmatrix_t *spm, double *A)
Convert a sparse matrix into a dense matrix.
void d_spmExpand(const spmatrix_t *spm_in, spmatrix_t *spm_out)
Expand a single dof sparse matrix to a multi-dofs sparse matrix.
void d_spmDofExtend(spmatrix_t *spm)
Extend a single dof sparse matrix to a multi-dof sparse matrix.
void d_spmIntFltSortAsc(void **const pbase, const spm_int_t n)
Integer routines.
void d_spmIntIntFltSortAsc(void **const pbase, const spm_int_t n)
Sort 2 arrays simultaneously, the first array is an array of spm_int_t and used as key for sorting....
int spm_dspmm(spm_side_t side, spm_trans_t transA, spm_trans_t transB, spm_int_t K, double alpha, const spmatrix_t *A, const double *B, spm_int_t ldb, double beta, double *C, spm_int_t ldc)
Compute a matrix-matrix product.
int spm_dspmv(spm_trans_t trans, double alpha, const spmatrix_t *A, const double *x, spm_int_t incx, double beta, double *y, spm_int_t incy)
compute the matrix-vector product:
double d_spmNorm(spm_normtype_t ntype, const spmatrix_t *spm)
Compute the norm of an spm matrix.
double d_spmNormMat(spm_normtype_t ntype, const spmatrix_t *spm, spm_int_t n, const double *A, spm_int_t lda)
Compute the norm of a dense matrix that follows the distribution of an spm matrix.
void d_spmDensePrint(FILE *f, spm_int_t m, spm_int_t n, const double *A, spm_int_t lda)
Print a dense matrix to the given file.
void d_spmPrintRHS(FILE *f, const spmatrix_t *spm, int nrhs, const void *x, spm_int_t ldx)
Write into a file the vectors associated to a spm.
void d_spmPrint(FILE *f, const spmatrix_t *spm)
Write a spm matrix in a file.
int d_spmRhsGenRndShm(const spmatrix_t *spm, double scale, spm_int_t n, double *A, spm_int_t lda, int shift, unsigned long long int seed)
Generate a set of vectors of random values in shared memory.
int d_spmCheckAxb(spm_fixdbl_t eps, int nrhs, const spmatrix_t *spm, void *x0, int ldx0, void *b, int ldb, const void *x, int ldx)
Check the backward error, and the forward error if x0 is provided.
int d_spmGenRHS(spm_rhstype_t type, int nrhs, const spmatrix_t *spm, void *x, int ldx, void *b, int ldb)
Generate nrhs right hand side vectors associated to a given matrix to test a problem with a solver.
int d_spmGenMat(spm_rhstype_t type, int nrhs, const spmatrix_t *spm, void *alpha, unsigned long long int seed, void *A, int lda)
Generate nrhs right hand side vectors associated to a given matrix to test a problem with a solver.
void d_spmReduceRHS(int nrhs, const spmatrix_t *spm, double *bglob, spm_int_t ldbg, double *bloc, spm_int_t ldbl)
Reduce all the global coefficients of a rhs and store the local ones.
Definition d_spm_rhs.c:109
int d_spmRhsGenRndDist(const spmatrix_t *spm, double scale, spm_int_t n, double *A, spm_int_t lda, int shift, unsigned long long int seed)
Generate a set of vectors of random values in distributed memory.
void d_spmGatherRHS(int nrhs, const spmatrix_t *spm, const double *x, spm_int_t ldx, int root, double *gx, spm_int_t ldgx)
Gather all the global C coefficients and store the good ones in local.
Definition d_spm_rhs.c:165
void d_spmExtractLocalRHS(int nrhs, const spmatrix_t *spm, const double *bglob, spm_int_t ldbg, double *bloc, spm_int_t ldbl)
Stores the local values of a global RHS in the local one.
Definition d_spm_rhs.c:51
void d_spmScal(const double alpha, spmatrix_t *spm)
Scal the spm: A = alpha * A.
Definition d_spm_scal.c:39
int spm_int_t
The main integer datatype used in spm arrays.
Definition datatypes.h:70
The sparse matrix data structure.
Definition spm.h:64