Program listing for file numerics/src/tools/SparseBlockMatrix.h#
Return to documentation for this file
1#ifndef SparseBlockMatrix_H
2#define SparseBlockMatrix_H
3
4#include <stdio.h>
5#include "CSparseMatrix.h"
6#include "NumericsFwd.h"
7
8#include "SiconosConfig.h"
9#include "NumericsDataVersion.h"
10
11
12
13
14
15
16struct SparseBlockStructuredMatrix
17{
18
19 unsigned int nbblocks;
20 double **block;
21
22 unsigned int blocknumber0;
23
24 unsigned int blocknumber1;
25
26 unsigned int *blocksize0;
27
28 unsigned int *blocksize1;
29
30
31 size_t filled1;
32
33 size_t filled2;
34
35 size_t *index1_data;
36 size_t *index2_data;
37
38
39 unsigned int * diagonal_blocks;
40
41 NumericsDataVersion version;
42};
43
44struct SparseBlockCoordinateMatrix
45{
46
47 unsigned int nbblocks;
48
49
50 unsigned int blocknumber0;
51
52
53 unsigned int blocknumber1;
54
55
56 double **block;
57
58
59 unsigned int *blocksize0;
60
61
62
63 unsigned int *blocksize1;
64
65
66 unsigned int *row;
67
68
69 unsigned int *column;
70};
71
72struct SparseBlockStructuredMatrixPred
73{
74 int nbbldiag;
75 int **indic;
76 int **indicop;
77 double **submatlcp;
78 double **submatlcpop;
79 int **ipiv;
80 int *sizesublcp;
81 int *sizesublcpop;
82 double **subq;
83 double **bufz;
84 double **newz;
85 double **workspace;
86};
87
88struct SBM_index_by_column
89{
90
91 size_t filled3;
92
93 size_t filled4;
94
95 size_t * index3_data;
96 size_t * index4_data;
97 size_t * blockMap;
98};
99
100
101
102#define NUMERICS_SBM_FREE_BLOCK 4
103#define NUMERICS_SBM_FREE_SBM 8
104
105#if defined(__cplusplus) && !defined(BUILD_AS_CPP)
106extern "C"
107{
108#endif
109
110
111 SparseBlockStructuredMatrix* SBM_new(void);
112
113
114 void SBM_null(SparseBlockStructuredMatrix* sbm);
115
116
117 void SBM_gemv(unsigned int sizeX, unsigned int sizeY,
118 double alpha, const SparseBlockStructuredMatrix* const A,
119 const double* x, double beta, double* y);
120
121
122 void SBM_gemv_3x3(unsigned int sizeX, unsigned int sizeY,
123 const SparseBlockStructuredMatrix* const A,
124 double* const x, double* y);
125
126
127 void SBM_gemm_without_allocation(double alpha, const SparseBlockStructuredMatrix* const A,
128 const SparseBlockStructuredMatrix* const B,
129 double beta, SparseBlockStructuredMatrix* C);
130
131
132 SparseBlockStructuredMatrix* SBM_multiply(const SparseBlockStructuredMatrix* const A, const SparseBlockStructuredMatrix* const B);
133
134
135 SparseBlockStructuredMatrix* SBM_zero_matrix_for_multiply(const SparseBlockStructuredMatrix* const A, const SparseBlockStructuredMatrix* const B);
136
137
138 SparseBlockStructuredMatrix * SBM_add(SparseBlockStructuredMatrix * A, SparseBlockStructuredMatrix * B, double alpha, double beta);
139
140
141 void SBM_add_without_allocation(SparseBlockStructuredMatrix * A, SparseBlockStructuredMatrix * B,
142 double alpha, double beta,
143 SparseBlockStructuredMatrix * C,
144 double gamma);
145
146
147
148 void SBM_scal(double alpha, SparseBlockStructuredMatrix * A);
149
150
151
152 void SBM_row_prod(unsigned int sizeX, unsigned int sizeY, unsigned int currentRowNumber, const SparseBlockStructuredMatrix* const A, const double* const x, double* y, int init);
153
154
155 void SBM_row_prod_no_diag(unsigned int sizeX, unsigned int sizeY, unsigned int currentRowNumber, const SparseBlockStructuredMatrix* const A, const double* const x, double* y, int init);
156
157
158 void SBM_row_prod_no_diag_3x3(unsigned int sizeX, unsigned int sizeY, unsigned int currentRowNumber, const SparseBlockStructuredMatrix* const A, double* const x, double* y);
159 void SBM_row_prod_no_diag_2x2(unsigned int sizeX, unsigned int sizeY, unsigned int currentRowNumber, const SparseBlockStructuredMatrix* const A, double* const x, double* y);
160 void SBM_row_prod_no_diag_1x1(unsigned int sizeX, unsigned int sizeY, unsigned int currentRowNumber, const SparseBlockStructuredMatrix* const A, double* const x, double* y);
161
162
163 void SBM_extract_component_3x3(const SparseBlockStructuredMatrix* const A,
164 SparseBlockStructuredMatrix* B,
165 unsigned int *row_components, unsigned int row_components_size,
166 unsigned int *col_components, unsigned int col_components_size);
167
168
169 void SBM_clear(SparseBlockStructuredMatrix * blmat);
170
171
172 void SBMfree(SparseBlockStructuredMatrix* A, unsigned int level);
173
174
175 void SBM_print(const SparseBlockStructuredMatrix* const m);
176
177
178 void SBM_write_in_file(const SparseBlockStructuredMatrix* const m, FILE* file);
179
180
181 void SBM_read_in_file(SparseBlockStructuredMatrix* const M, FILE *file);
182
183
184 SparseBlockStructuredMatrix* SBM_new_from_file(FILE *file);
185
186
187 void SBM_write_in_fileForScilab(const SparseBlockStructuredMatrix* const M, FILE* file);
188
189
190
191 void SBM_write_in_filename(const SparseBlockStructuredMatrix* const M, const char *filename);
192
193
194 void SBM_read_in_filename(SparseBlockStructuredMatrix* const M, const char *filename);
195
196
197 void SBM_clear_pred(SparseBlockStructuredMatrixPred *blmatpred);
198
199
200 unsigned int * SBM_diagonal_block_indices(SparseBlockStructuredMatrix* const M);
201
202
203 unsigned int SBM_diagonal_block_index(SparseBlockStructuredMatrix* const M, unsigned int row);
204
205
206 int SBM_entry(SparseBlockStructuredMatrix* M, unsigned int row, unsigned int col, double val);
207
208
209 double SBM_get_value(const SparseBlockStructuredMatrix* const M, unsigned int row, unsigned int col);
210
211
212 int SBM_copy(const SparseBlockStructuredMatrix* const A, SparseBlockStructuredMatrix* B, unsigned int copyBlock);
213
214
215
216 int SBM_transpose(const SparseBlockStructuredMatrix* const A, SparseBlockStructuredMatrix* B);
217
218
219 int SBM_inverse_diagonal_block_matrix_in_place(const SparseBlockStructuredMatrix* M, int* ipiv);
220
221
222
223 void SBM_to_dense(const SparseBlockStructuredMatrix* const A, double *denseMat);
224
225
226 int SBM_to_sparse(const SparseBlockStructuredMatrix* const A, CSparseMatrix *outSparseMat);
227
228
229 int SBM_to_sparse_init_memory(const SparseBlockStructuredMatrix* const A, CSparseMatrix *sparseMat);
230
231
232 void SBM_row_to_dense(const SparseBlockStructuredMatrix* const A, int row, double *denseMat, int rowPos, int rowNb);
233
234
235 void SBM_row_permutation(unsigned int *rowIndex, SparseBlockStructuredMatrix* A, SparseBlockStructuredMatrix* C);
236
237
238 void SBM_column_permutation(unsigned int *colIndex, SparseBlockStructuredMatrix* A, SparseBlockStructuredMatrix* C);
239
240 void SBCM_null(SparseBlockCoordinateMatrix* MC);
241
242 SparseBlockCoordinateMatrix* SBCM_new(void);
243
244
245
246
247
248 SparseBlockCoordinateMatrix* SBCM_new_3x3(unsigned int m, unsigned int n,
249 unsigned int nbblocks,
250 unsigned int *row, unsigned int *column, double *block);
251
252
253
254 void SBCM_free_3x3(SparseBlockCoordinateMatrix *MC);
255
256
257 SparseBlockStructuredMatrix* SBCM_to_SBM(SparseBlockCoordinateMatrix* MC);
258
259
260
261 void SBM_free_from_SBCM(SparseBlockStructuredMatrix* M);
262
263
264
265 int SBM_from_csparse(int blocksize, const CSparseMatrix* const sparseMat, SparseBlockStructuredMatrix* outSBM);
266
267
268#if defined(__cplusplus) && !defined(BUILD_AS_CPP)
269}
270#endif
271
272#endif