SpM Handbook 1.2.4
Loading...
Searching...
No Matches
s_spm_sort.c
Go to the documentation of this file.
1/**
2 *
3 * @file s_spm_sort.c
4 *
5 * SParse Matrix package precision dependent sort routines.
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 Mathieu Faverge
12 * @author Tony Delarue
13 * @author Alycia Lisito
14 * @date 2024-06-25
15 *
16 * @generated from /builds/2mk6rsew/0/fpruvost/spm/src/z_spm_sort.c, normal z -> s, Fri Nov 29 11:34:30 2024
17 *
18 **/
19#include "common.h"
20
21/**
22 *******************************************************************************
23 *
24 * @ingroup spm_dev_check
25 *
26 * @brief This routine sorts the single dof spm matrix.
27 *
28 * For the CSC and CSR formats, the subarray of edges for each vertex are sorted.
29 * For the IJV format, the edges are storted first by column indexes, and then
30 * by row indexes. To perform a sort first by row, second by column, please swap
31 * the colptr and rowptr of the structure before calling the subroutine.
32 *
33 * @warning This function should NOT be called if dof is greater than 1.
34 *
35 *******************************************************************************
36 *
37 * @param[inout] spm
38 * On entry, the pointer to the sparse matrix structure.
39 * On exit, the same sparse matrix with subarrays of edges sorted by
40 * ascending order.
41 *
42 *******************************************************************************/
43static inline void
45{
46 spm_int_t *colptr = spm->colptr;
47 spm_int_t *rowptr = spm->rowptr;
48 float *values = spm->values;
49 void *sortptr[3];
50 spm_int_t n = spm->n;
51 spm_int_t i, size;
52 (void)sortptr;
53
54 /* Sort in place each subset */
55 if ( spm->fmttype == SpmCSC ) {
56 for (i=0; i<n; i++, colptr++)
57 {
58 size = colptr[1] - colptr[0];
59
60#if defined(PRECISION_p)
61 spmIntSort1Asc1( rowptr, size );
62#else
63 sortptr[0] = rowptr;
64 sortptr[1] = values;
65 s_spmIntFltSortAsc( sortptr, size );
66#endif
67 rowptr += size;
68 values += size;
69 }
70 }
71 else if ( spm->fmttype == SpmCSR ) {
72 for (i=0; i<n; i++, rowptr++)
73 {
74 size = rowptr[1] - rowptr[0];
75
76#if defined(PRECISION_p)
77 spmIntSort1Asc1( colptr, size );
78#else
79 sortptr[0] = colptr;
80 sortptr[1] = values;
81 s_spmIntFltSortAsc( sortptr, size );
82#endif
83 colptr += size;
84 values += size;
85 }
86 }
87 else if ( spm->fmttype == SpmIJV ) {
88 size = spm->nnz;
89
90 sortptr[0] = colptr;
91 sortptr[1] = rowptr;
92
93#if defined(PRECISION_p)
94 spmIntMSortIntAsc( sortptr, size );
95#else
96 sortptr[2] = values;
97 s_spmIntIntFltSortAsc( sortptr, size );
98#endif
99 }
100 (void) values;
101}
102
103/**
104 *******************************************************************************
105 *
106 * @brief Apply the permutations on the values array for a CSX spm
107 * The values array holds the permutation indexes of the new values array.
108 *
109 *******************************************************************************
110 *
111 * @param[in] spm
112 * The pointer to the spm.
113 *
114 * @param[in] values
115 * The original values array.
116 *
117 * @param[inout] newval
118 * The new values array with the correct permutations.
119 * Must be allocated before this routine.
120 *
121 *******************************************************************************/
122static inline void
124 const float *values,
125 float *newval )
126{
127 spm_int_t i, j, ig, jg, index;
128 spm_int_t size, baseval, dof;
129 spm_int_t dofi, dofj, dof2;
130
131 const spm_int_t *colptr = (spm->fmttype == SpmCSC) ? spm->colptr: spm->rowptr;
132 const spm_int_t *rowptr = (spm->fmttype == SpmCSC) ? spm->rowptr: spm->colptr;
133 const spm_int_t *indexes = spm->values;
134 const spm_int_t *dofs = spm->dofs;
135 const spm_int_t *loc2glob = spm->loc2glob;
136
137 float *valtmp = newval;
138
139 size = spm->n;
140 baseval = spm->baseval;
141 dof = spm->dof;
142 for ( j = 0; j < size; j++, colptr++, loc2glob++ )
143 {
144 jg = spm->replicated ? j : *loc2glob - baseval;
145 dofj = (dof > 0) ? dof : dofs[jg+1] - dofs[jg];
146
147 for ( i = colptr[0]; i < colptr[1]; i++, rowptr++, indexes++ )
148 {
149 ig = *rowptr - baseval;
150 dofi = (dof > 0) ? dof : dofs[ig+1] - dofs[ig];
151 dof2 = dofi * dofj;
152
153 index = *indexes;
154 memcpy( valtmp, values + index, dof2 * sizeof(float) );
155 valtmp += dof2;
156 }
157 }
158}
159
160/**
161 *******************************************************************************
162 *
163 * @brief Apply the permutations on the values array fon an IJV spm.
164 * The values array holds the permutation indexes of the new values array.
165 *
166 *******************************************************************************
167 *
168 * @param[in] spm
169 * The pointer to the spm.
170 *
171 * @param[in] values
172 * The original values array.
173 *
174 * @param[inout] newval
175 * The new values array with the correct permutations.
176 * Must be allocated before this routine.
177 *
178 *******************************************************************************/
179static inline void
181 const float *values,
182 float *newval )
183{
184 spm_int_t i, ig, jg, index;
185 spm_int_t size, baseval, dof;
186 spm_int_t dofi, dofj, dof2;
187
188 const spm_int_t *colptr = spm->colptr;
189 const spm_int_t *rowptr = spm->rowptr;
190 const spm_int_t *indexes = spm->values;
191 const spm_int_t *dofs;
192
193 float *valtmp = newval;
194
195 size = spm->nnz;
196 baseval = spm->baseval;
197 dof = spm->dof;
198 dofs = spm->dofs - baseval;
199 for ( i = 0; i < size; i++, colptr++, rowptr++, indexes++ )
200 {
201 jg = *colptr;
202 dofj = (dof > 0) ? dof : dofs[jg+1] - dofs[jg];
203 ig = *rowptr;
204 dofi = (dof > 0) ? dof : dofs[ig+1] - dofs[ig];
205 dof2 = dofi * dofj;
206
207 index = *indexes;
208
209 memcpy( valtmp, values + index, dof2 * sizeof(float) );
210 valtmp += dof2;
211 }
212}
213
214/**
215 *******************************************************************************
216 *
217 * @ingroup spm_dev_check
218 *
219 * @brief This routine sorts the multiple dof spm matrix.
220 *
221 * For the CSC and CSR formats, the subarray of edges for each vertex are sorted.
222 * For the IJV format, the edges are sorted first by column indexes, and then
223 * by row indexes. To perform a sort first by row, second by column, please swap
224 * the colptr and rowptr of the structure before calling the subroutine.
225 * This routine is used for multidof matrices. It's way less efficient than the
226 * single dof one.
227 *
228 *******************************************************************************
229 *
230 * @param[inout] spm
231 * On entry, the pointer to the sparse matrix structure.
232 * On exit, the same sparse matrix with subarrays of edges sorted by
233 * ascending order.
234 *
235 *******************************************************************************/
236static inline void
238{
239 spm_int_t dof;
240 spm_coeftype_t flttype;
241 float *values = spm->values;
242 /* This array will store the solution */
243 float *newval = malloc( spm->nnzexp * sizeof(float) );
244
245 /* Create a tmp array composed by the multidof indexes of the valptr */
246 spm_int_t *indexes = spm_get_value_idx_by_elt( spm );
247
248 /*
249 * Sort the spm as a single dof matrix.
250 * The value array will represent the permutations.
251 */
252 dof = spm->dof;
253 flttype = spm->flttype;
254
255 spm->values = indexes;
256 spm->dof = 1;
257
258 if ( sizeof(spm_int_t) == sizeof(spm_fixdbl_t) ) {
259 spm->flttype = 3; /* SpmFloat */
260 }
261 else {
262 spm->flttype = 2; /* SpmFloat */
263 }
264 spmSort( spm );
265
266 spm->dof = dof;
267 spm->flttype = flttype;
268
269 /* Apply the permutations and copy datas in the newval */
270 if ( spm->fmttype != SpmIJV ) {
271 s_spm_sort_multidof_csx_values( spm, values, newval );
272 }
273 else {
274 s_spm_sort_multidof_ijv_values( spm, values, newval );
275 }
276 free(indexes);
277 free(values);
278
279 spm->values = newval;
280}
281
282/**
283 *******************************************************************************
284 *
285 * @ingroup spm_dev_check
286 *
287 * @brief This routine sorts the spm matrix.
288 *
289 * For the CSC and CSR formats, the subarray of edges for each vertex are sorted.
290 * For the IJV format, the edges are storted first by column indexes, and then
291 * by row indexes. To perform a sort first by row, second by column, please swap
292 * the colptr and rowptr of the structure before calling the subroutine.
293 *
294 *******************************************************************************
295 *
296 * @param[inout] spm
297 * On entry, the pointer to the sparse matrix structure.
298 * On exit, the same sparse matrix with subarrays of edges sorted by
299 * ascending order.
300 *
301 *******************************************************************************/
302void
304{
305 int swapped = 0;
306
307 if ( spm->fmttype == SpmIJV ) {
308 int distribution;
309
311 distribution = spm_get_distribution( spm );
312
313 if ( distribution == SpmDistByRow ) {
314 spm_int_t *tmp = spm->colptr;
315 spm->colptr = spm->rowptr;
316 spm->rowptr = tmp;
317 swapped = 1;
318 }
319 }
320
321 if ( (spm->dof != 1) && (spm->flttype != SpmPattern) ) {
322 s_spmSortMultidof( spm );
323 }
324 else {
325 s_spmSortNoDof( spm );
326 }
327
328 if ( swapped ) {
329 spm_int_t *tmp = spm->colptr;
330 spm->colptr = spm->rowptr;
331 spm->rowptr = tmp;
332 }
333}
double spm_fixdbl_t
Double datatype that is not converted through precision generator functions.
Definition datatypes.h:150
#define SpmDistByRow
Definition const.h:43
enum spm_coeftype_e spm_coeftype_t
Arithmetic types.
@ SpmPattern
Definition const.h:62
@ SpmCSC
Definition const.h:73
@ SpmCSR
Definition const.h:74
@ SpmIJV
Definition const.h:75
static void s_spmSortNoDof(spmatrix_t *spm)
This routine sorts the single dof spm matrix.
Definition s_spm_sort.c:44
spm_int_t * spm_get_value_idx_by_elt(const spmatrix_t *spm)
Create an array that represents the shift for each sub-element of the original multidof value array.
Definition spm.c:2211
void s_spmSort(spmatrix_t *spm)
This routine sorts the spm matrix.
Definition s_spm_sort.c:303
static void s_spmSortMultidof(spmatrix_t *spm)
This routine sorts the multiple dof spm matrix.
Definition s_spm_sort.c:237
void s_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....
void s_spmIntFltSortAsc(void **const pbase, const spm_int_t n)
Integer routines.
spm_coeftype_t flttype
Definition spm.h:67
spm_int_t * dofs
Definition spm.h:85
void * values
Definition spm.h:92
spm_int_t nnz
Definition spm.h:75
spm_int_t nnzexp
Definition spm.h:80
spm_int_t * rowptr
Definition spm.h:90
spm_fmttype_t fmttype
Definition spm.h:69
spm_int_t dof
Definition spm.h:82
spm_int_t n
Definition spm.h:73
spm_int_t * loc2glob
Definition spm.h:91
spm_int_t baseval
Definition spm.h:71
spm_int_t * colptr
Definition spm.h:89
int replicated
Definition spm.h:98
int spm_int_t
The main integer datatype used in spm arrays.
Definition datatypes.h:70
int spmSort(spmatrix_t *spm)
Sort the subarray of edges of each vertex.
Definition spm.c:746
The sparse matrix data structure.
Definition spm.h:64
static void s_spm_sort_multidof_csx_values(const spmatrix_t *spm, const float *values, float *newval)
Apply the permutations on the values array for a CSX spm The values array holds the permutation index...
Definition s_spm_sort.c:123
static void s_spm_sort_multidof_ijv_values(const spmatrix_t *spm, const float *values, float *newval)
Apply the permutations on the values array fon an IJV spm. The values array holds the permutation ind...
Definition s_spm_sort.c:180
spm_int_t * spm_getandset_glob2loc(spmatrix_t *spm)
Computes the glob2loc array if needed, and returns it.
Definition spm.c:2007
int spm_get_distribution(const spmatrix_t *spm)
Search the distribution pattern used in the spm structure.
Definition spm.c:2049