Actual source code: aijperm.c
1: /*
2: Defines basic operations for the MATSEQAIJPERM matrix class.
3: This class is derived from the MATSEQAIJ class and retains the
4: compressed row storage (aka Yale sparse matrix format) but augments
5: it with some permutation information that enables some operations
6: to be more vectorizable. A physically rearranged copy of the matrix
7: may be stored if the user desires.
9: Eventually a variety of permutations may be supported.
10: */
12: #include <../src/mat/impls/aij/seq/aij.h>
14: #if defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
15: #include <immintrin.h>
17: #if !defined(_MM_SCALE_8)
18: #define _MM_SCALE_8 8
19: #endif
20: #if !defined(_MM_SCALE_4)
21: #define _MM_SCALE_4 4
22: #endif
23: #endif
25: #define NDIM 512
26: /* NDIM specifies how many rows at a time we should work with when
27: * performing the vectorized mat-vec. This depends on various factors
28: * such as vector register length, etc., and I really need to add a
29: * way for the user (or the library) to tune this. I'm setting it to
30: * 512 for now since that is what Ed D'Azevedo was using in his Fortran
31: * routines. */
33: typedef struct {
34: PetscObjectState nonzerostate; /* used to determine if the nonzero structure has changed and hence the permutations need updating */
36: PetscInt ngroup;
37: PetscInt *xgroup;
38: /* Denotes where groups of rows with same number of nonzeros
39: * begin and end, i.e., xgroup[i] gives us the position in iperm[]
40: * where the ith group begins. */
42: PetscInt *nzgroup; /* how many nonzeros each row that is a member of group i has. */
43: PetscInt *iperm; /* The permutation vector. */
45: /* Some of this stuff is for Ed's recursive triangular solve.
46: * I'm not sure what I need yet. */
47: PetscInt blocksize;
48: PetscInt nstep;
49: PetscInt *jstart_list;
50: PetscInt *jend_list;
51: PetscInt *action_list;
52: PetscInt *ngroup_list;
53: PetscInt **ipointer_list;
54: PetscInt **xgroup_list;
55: PetscInt **nzgroup_list;
56: PetscInt **iperm_list;
57: } Mat_SeqAIJPERM;
59: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJPERM_SeqAIJ(Mat A, MatType type, MatReuse reuse, Mat *newmat)
60: {
61: /* This routine is only called to convert a MATAIJPERM to its base PETSc type, */
62: /* so we will ignore 'MatType type'. */
63: Mat B = *newmat;
64: Mat_SeqAIJPERM *aijperm = (Mat_SeqAIJPERM *)A->spptr;
66: PetscFunctionBegin;
67: if (reuse == MAT_INITIAL_MATRIX) {
68: PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
69: aijperm = (Mat_SeqAIJPERM *)B->spptr;
70: }
72: /* Reset the original function pointers. */
73: B->ops->assemblyend = MatAssemblyEnd_SeqAIJ;
74: B->ops->destroy = MatDestroy_SeqAIJ;
75: B->ops->duplicate = MatDuplicate_SeqAIJ;
76: B->ops->mult = MatMult_SeqAIJ;
77: B->ops->multadd = MatMultAdd_SeqAIJ;
79: PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaijperm_seqaij_C", NULL));
81: /* Free everything in the Mat_SeqAIJPERM data structure.*/
82: PetscCall(PetscFree(aijperm->xgroup));
83: PetscCall(PetscFree(aijperm->nzgroup));
84: PetscCall(PetscFree(aijperm->iperm));
85: PetscCall(PetscFree(B->spptr));
87: /* Change the type of B to MATSEQAIJ. */
88: PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJ));
90: *newmat = B;
91: PetscFunctionReturn(PETSC_SUCCESS);
92: }
94: static PetscErrorCode MatDestroy_SeqAIJPERM(Mat A)
95: {
96: Mat_SeqAIJPERM *aijperm = (Mat_SeqAIJPERM *)A->spptr;
98: PetscFunctionBegin;
99: if (aijperm) {
100: /* If MatHeaderMerge() was used then this SeqAIJPERM matrix will not have a spprt. */
101: PetscCall(PetscFree(aijperm->xgroup));
102: PetscCall(PetscFree(aijperm->nzgroup));
103: PetscCall(PetscFree(aijperm->iperm));
104: PetscCall(PetscFree(A->spptr));
105: }
106: /* Change the type of A back to SEQAIJ and use MatDestroy_SeqAIJ()
107: * to destroy everything that remains. */
108: PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATSEQAIJ));
109: /* Note that I don't call MatSetType(). I believe this is because that
110: * is only to be called when *building* a matrix. I could be wrong, but
111: * that is how things work for the SuperLU matrix class. */
112: PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqaijperm_seqaij_C", NULL));
113: PetscCall(MatDestroy_SeqAIJ(A));
114: PetscFunctionReturn(PETSC_SUCCESS);
115: }
117: static PetscErrorCode MatDuplicate_SeqAIJPERM(Mat A, MatDuplicateOption op, Mat *M)
118: {
119: Mat_SeqAIJPERM *aijperm = (Mat_SeqAIJPERM *)A->spptr;
120: Mat_SeqAIJPERM *aijperm_dest;
121: PetscBool perm;
123: PetscFunctionBegin;
124: PetscCall(MatDuplicate_SeqAIJ(A, op, M));
125: PetscCall(PetscObjectTypeCompare((PetscObject)*M, MATSEQAIJPERM, &perm));
126: if (perm) {
127: aijperm_dest = (Mat_SeqAIJPERM *)(*M)->spptr;
128: PetscCall(PetscFree(aijperm_dest->xgroup));
129: PetscCall(PetscFree(aijperm_dest->nzgroup));
130: PetscCall(PetscFree(aijperm_dest->iperm));
131: } else {
132: PetscCall(PetscNew(&aijperm_dest));
133: (*M)->spptr = (void *)aijperm_dest;
134: PetscCall(PetscObjectChangeTypeName((PetscObject)*M, MATSEQAIJPERM));
135: PetscCall(PetscObjectComposeFunction((PetscObject)*M, "MatConvert_seqaijperm_seqaij_C", MatConvert_SeqAIJPERM_SeqAIJ));
136: }
137: PetscCall(PetscArraycpy(aijperm_dest, aijperm, 1));
138: /* Allocate space for, and copy the grouping and permutation info.
139: * I note that when the groups are initially determined in
140: * MatSeqAIJPERM_create_perm, xgroup and nzgroup may be sized larger than
141: * necessary. But at this point, we know how large they need to be, and
142: * allocate only the necessary amount of memory. So the duplicated matrix
143: * may actually use slightly less storage than the original! */
144: PetscCall(PetscMalloc1(A->rmap->n, &aijperm_dest->iperm));
145: PetscCall(PetscMalloc1(aijperm->ngroup + 1, &aijperm_dest->xgroup));
146: PetscCall(PetscMalloc1(aijperm->ngroup, &aijperm_dest->nzgroup));
147: PetscCall(PetscArraycpy(aijperm_dest->iperm, aijperm->iperm, A->rmap->n));
148: PetscCall(PetscArraycpy(aijperm_dest->xgroup, aijperm->xgroup, aijperm->ngroup + 1));
149: PetscCall(PetscArraycpy(aijperm_dest->nzgroup, aijperm->nzgroup, aijperm->ngroup));
150: PetscFunctionReturn(PETSC_SUCCESS);
151: }
153: static PetscErrorCode MatSeqAIJPERM_create_perm(Mat A)
154: {
155: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
156: Mat_SeqAIJPERM *aijperm = (Mat_SeqAIJPERM *)A->spptr;
157: PetscInt m; /* Number of rows in the matrix. */
158: PetscInt *ia; /* From the CSR representation; points to the beginning of each row. */
159: PetscInt maxnz; /* Maximum number of nonzeros in any row. */
160: PetscInt *rows_in_bucket;
161: /* To construct the permutation, we sort each row into one of maxnz
162: * buckets based on how many nonzeros are in the row. */
163: PetscInt nz;
164: PetscInt *nz_in_row; /* the number of nonzero elements in row k. */
165: PetscInt *ipnz;
166: /* When constructing the iperm permutation vector,
167: * ipnz[nz] is used to point to the next place in the permutation vector
168: * that a row with nz nonzero elements should be placed.*/
169: PetscInt i, ngroup, istart, ipos;
171: PetscFunctionBegin;
172: if (aijperm->nonzerostate == A->nonzerostate) PetscFunctionReturn(PETSC_SUCCESS); /* permutation exists and matches current nonzero structure */
173: aijperm->nonzerostate = A->nonzerostate;
174: /* Free anything previously put in the Mat_SeqAIJPERM data structure. */
175: PetscCall(PetscFree(aijperm->xgroup));
176: PetscCall(PetscFree(aijperm->nzgroup));
177: PetscCall(PetscFree(aijperm->iperm));
179: m = A->rmap->n;
180: ia = a->i;
182: /* Allocate the arrays that will hold the permutation vector. */
183: PetscCall(PetscMalloc1(m, &aijperm->iperm));
185: /* Allocate some temporary work arrays that will be used in
186: * calculating the permutation vector and groupings. */
187: PetscCall(PetscMalloc1(m, &nz_in_row));
189: /* Now actually figure out the permutation and grouping. */
191: /* First pass: Determine number of nonzeros in each row, maximum
192: * number of nonzeros in any row, and how many rows fall into each
193: * "bucket" of rows with same number of nonzeros. */
194: maxnz = 0;
195: for (i = 0; i < m; i++) {
196: nz_in_row[i] = ia[i + 1] - ia[i];
197: if (nz_in_row[i] > maxnz) maxnz = nz_in_row[i];
198: }
199: PetscCall(PetscMalloc1(PetscMax(maxnz, m) + 1, &rows_in_bucket));
200: PetscCall(PetscMalloc1(PetscMax(maxnz, m) + 1, &ipnz));
202: for (i = 0; i <= maxnz; i++) rows_in_bucket[i] = 0;
203: for (i = 0; i < m; i++) {
204: nz = nz_in_row[i];
205: rows_in_bucket[nz]++;
206: }
208: /* Allocate space for the grouping info. There will be at most (maxnz + 1)
209: * groups. (It is maxnz + 1 instead of simply maxnz because there may be
210: * rows with no nonzero elements.) If there are (maxnz + 1) groups,
211: * then xgroup[] must consist of (maxnz + 2) elements, since the last
212: * element of xgroup will tell us where the (maxnz + 1)th group ends.
213: * We allocate space for the maximum number of groups;
214: * that is potentially a little wasteful, but not too much so.
215: * Perhaps I should fix it later. */
216: PetscCall(PetscMalloc1(maxnz + 2, &aijperm->xgroup));
217: PetscCall(PetscMalloc1(maxnz + 1, &aijperm->nzgroup));
219: /* Second pass. Look at what is in the buckets and create the groupings.
220: * Note that it is OK to have a group of rows with no non-zero values. */
221: ngroup = 0;
222: istart = 0;
223: for (i = 0; i <= maxnz; i++) {
224: if (rows_in_bucket[i] > 0) {
225: aijperm->nzgroup[ngroup] = i;
226: aijperm->xgroup[ngroup] = istart;
227: ngroup++;
228: istart += rows_in_bucket[i];
229: }
230: }
232: aijperm->xgroup[ngroup] = istart;
233: aijperm->ngroup = ngroup;
235: /* Now fill in the permutation vector iperm. */
236: ipnz[0] = 0;
237: for (i = 0; i < maxnz; i++) ipnz[i + 1] = ipnz[i] + rows_in_bucket[i];
239: for (i = 0; i < m; i++) {
240: nz = nz_in_row[i];
241: ipos = ipnz[nz];
242: aijperm->iperm[ipos] = i;
243: ipnz[nz]++;
244: }
246: /* Clean up temporary work arrays. */
247: PetscCall(PetscFree(rows_in_bucket));
248: PetscCall(PetscFree(ipnz));
249: PetscCall(PetscFree(nz_in_row));
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: static PetscErrorCode MatAssemblyEnd_SeqAIJPERM(Mat A, MatAssemblyType mode)
254: {
255: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
257: PetscFunctionBegin;
258: if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
260: /* Since a MATSEQAIJPERM matrix is really just a MATSEQAIJ with some
261: * extra information, call the AssemblyEnd routine for a MATSEQAIJ.
262: * I'm not sure if this is the best way to do this, but it avoids
263: * a lot of code duplication.
264: * I also note that currently MATSEQAIJPERM doesn't know anything about
265: * the Mat_CompressedRow data structure that SeqAIJ now uses when there
266: * are many zero rows. If the SeqAIJ assembly end routine decides to use
267: * this, this may break things. (Don't know... haven't looked at it.) */
268: a->inode.use = PETSC_FALSE;
269: PetscCall(MatAssemblyEnd_SeqAIJ(A, mode));
271: /* Now calculate the permutation and grouping information. */
272: PetscCall(MatSeqAIJPERM_create_perm(A));
273: PetscFunctionReturn(PETSC_SUCCESS);
274: }
276: static PetscErrorCode MatMult_SeqAIJPERM(Mat A, Vec xx, Vec yy)
277: {
278: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
279: const PetscScalar *x;
280: PetscScalar *y;
281: const MatScalar *aa;
282: const PetscInt *aj, *ai;
283: PetscInt i, j;
284: #if defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
285: __m512d vec_x, vec_y, vec_vals;
286: __m256i vec_idx, vec_ipos, vec_j;
287: __mmask8 mask;
288: #endif
290: /* Variables that don't appear in MatMult_SeqAIJ. */
291: Mat_SeqAIJPERM *aijperm = (Mat_SeqAIJPERM *)A->spptr;
292: PetscInt *iperm; /* Points to the permutation vector. */
293: PetscInt *xgroup;
294: /* Denotes where groups of rows with same number of nonzeros
295: * begin and end in iperm. */
296: PetscInt *nzgroup;
297: PetscInt ngroup;
298: PetscInt jstart, jend;
299: /* jstart is used in loops to denote the position in iperm where a
300: * group starts; jend denotes the position where it ends.
301: * (jend + 1 is where the next group starts.) */
302: PetscInt iold, nz;
303: PetscInt istart, iend, isize;
304: PetscInt ipos;
305: PetscScalar yp[NDIM];
306: PetscInt ip[NDIM]; /* yp[] and ip[] are treated as vector "registers" for performing the mat-vec. */
308: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
309: #pragma disjoint(*x, *y, *aa)
310: #endif
312: PetscFunctionBegin;
313: PetscCall(VecGetArrayRead(xx, &x));
314: PetscCall(VecGetArray(yy, &y));
315: aj = a->j; /* aj[k] gives column index for element aa[k]. */
316: aa = a->a; /* Nonzero elements stored row-by-row. */
317: ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */
319: /* Get the info we need about the permutations and groupings. */
320: iperm = aijperm->iperm;
321: ngroup = aijperm->ngroup;
322: xgroup = aijperm->xgroup;
323: nzgroup = aijperm->nzgroup;
325: for (PetscInt igroup = 0; igroup < ngroup; igroup++) {
326: jstart = xgroup[igroup];
327: jend = xgroup[igroup + 1] - 1;
328: nz = nzgroup[igroup];
330: /* Handle the special cases where the number of nonzeros per row
331: * in the group is either 0 or 1. */
332: if (nz == 0) {
333: for (i = jstart; i <= jend; i++) y[iperm[i]] = 0.0;
334: } else if (nz == 1) {
335: for (i = jstart; i <= jend; i++) {
336: iold = iperm[i];
337: ipos = ai[iold];
338: y[iold] = aa[ipos] * x[aj[ipos]];
339: }
340: } else {
341: /* We work our way through the current group in chunks of NDIM rows
342: * at a time. */
344: for (istart = jstart; istart <= jend; istart += NDIM) {
345: /* Figure out where the chunk of 'isize' rows ends in iperm.
346: * 'isize may of course be less than NDIM for the last chunk. */
347: iend = istart + (NDIM - 1);
349: if (iend > jend) iend = jend;
351: isize = iend - istart + 1;
353: /* Initialize the yp[] array that will be used to hold part of
354: * the permuted results vector, and figure out where in aa each
355: * row of the chunk will begin. */
356: for (i = 0; i < isize; i++) {
357: iold = iperm[istart + i];
358: /* iold is a row number from the matrix A *before* reordering. */
359: ip[i] = ai[iold];
360: /* ip[i] tells us where the ith row of the chunk begins in aa. */
361: yp[i] = (PetscScalar)0.0;
362: }
364: /* If the number of zeros per row exceeds the number of rows in
365: * the chunk, we should vectorize along nz, that is, perform the
366: * mat-vec one row at a time as in the usual CSR case. */
367: if (nz > isize) {
368: #if defined(PETSC_HAVE_CRAY_VECTOR)
369: #pragma _CRI preferstream
370: #endif
371: for (i = 0; i < isize; i++) {
372: #if defined(PETSC_HAVE_CRAY_VECTOR)
373: #pragma _CRI prefervector
374: #endif
376: #if defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
377: vec_y = _mm512_setzero_pd();
378: ipos = ip[i];
379: for (j = 0; j < (nz >> 3); j++) {
380: vec_idx = _mm256_loadu_si256((__m256i const *)&aj[ipos]);
381: vec_vals = _mm512_loadu_pd(&aa[ipos]);
382: vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
383: vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
384: ipos += 8;
385: }
386: if ((nz & 0x07) > 2) {
387: mask = (__mmask8)(0xff >> (8 - (nz & 0x07)));
388: vec_idx = _mm256_loadu_si256((__m256i const *)&aj[ipos]);
389: vec_vals = _mm512_loadu_pd(&aa[ipos]);
390: vec_x = _mm512_mask_i32gather_pd(vec_x, mask, vec_idx, x, _MM_SCALE_8);
391: vec_y = _mm512_mask3_fmadd_pd(vec_x, vec_vals, vec_y, mask);
392: } else if ((nz & 0x07) == 2) {
393: yp[i] += aa[ipos] * x[aj[ipos]];
394: yp[i] += aa[ipos + 1] * x[aj[ipos + 1]];
395: } else if ((nz & 0x07) == 1) {
396: yp[i] += aa[ipos] * x[aj[ipos]];
397: }
398: yp[i] += _mm512_reduce_add_pd(vec_y);
399: #else
400: for (j = 0; j < nz; j++) {
401: ipos = ip[i] + j;
402: yp[i] += aa[ipos] * x[aj[ipos]];
403: }
404: #endif
405: }
406: } else {
407: /* Otherwise, there are enough rows in the chunk to make it
408: * worthwhile to vectorize across the rows, that is, to do the
409: * matvec by operating with "columns" of the chunk. */
410: for (j = 0; j < nz; j++) {
411: #if defined(PETSC_USE_AVX512_KERNELS) && defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
412: vec_j = _mm256_set1_epi32(j);
413: for (i = 0; i < ((isize >> 3) << 3); i += 8) {
414: vec_y = _mm512_loadu_pd(&yp[i]);
415: vec_ipos = _mm256_loadu_si256((__m256i const *)&ip[i]);
416: vec_ipos = _mm256_add_epi32(vec_ipos, vec_j);
417: vec_idx = _mm256_i32gather_epi32(aj, vec_ipos, _MM_SCALE_4);
418: vec_vals = _mm512_i32gather_pd(vec_ipos, aa, _MM_SCALE_8);
419: vec_x = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8);
420: vec_y = _mm512_fmadd_pd(vec_x, vec_vals, vec_y);
421: _mm512_storeu_pd(&yp[i], vec_y);
422: }
423: for (i = isize - (isize & 0x07); i < isize; i++) {
424: ipos = ip[i] + j;
425: yp[i] += aa[ipos] * x[aj[ipos]];
426: }
427: #else
428: for (i = 0; i < isize; i++) {
429: ipos = ip[i] + j;
430: yp[i] += aa[ipos] * x[aj[ipos]];
431: }
432: #endif
433: }
434: }
436: #if defined(PETSC_HAVE_CRAY_VECTOR)
437: #pragma _CRI ivdep
438: #endif
439: /* Put results from yp[] into non-permuted result vector y. */
440: for (i = 0; i < isize; i++) y[iperm[istart + i]] = yp[i];
441: } /* End processing chunk of isize rows of a group. */
442: } /* End handling matvec for chunk with nz > 1. */
443: } /* End loop over igroup. */
444: PetscCall(PetscLogFlops(PetscMax(2.0 * a->nz - A->rmap->n, 0)));
445: PetscCall(VecRestoreArrayRead(xx, &x));
446: PetscCall(VecRestoreArray(yy, &y));
447: PetscFunctionReturn(PETSC_SUCCESS);
448: }
450: /* MatMultAdd_SeqAIJPERM() calculates yy = ww + A * xx.
451: * Note that the names I used to designate the vectors differs from that
452: * used in MatMultAdd_SeqAIJ(). I did this to keep my notation consistent
453: * with the MatMult_SeqAIJPERM() routine, which is very similar to this one. */
454: /*
455: I hate having virtually identical code for the mult and the multadd!!!
456: */
457: static PetscErrorCode MatMultAdd_SeqAIJPERM(Mat A, Vec xx, Vec ww, Vec yy)
458: {
459: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
460: const PetscScalar *x;
461: PetscScalar *y, *w;
462: const MatScalar *aa;
463: const PetscInt *aj, *ai;
464: PetscInt i, j;
466: /* Variables that don't appear in MatMultAdd_SeqAIJ. */
467: Mat_SeqAIJPERM *aijperm;
468: PetscInt *iperm; /* Points to the permutation vector. */
469: PetscInt *xgroup;
470: /* Denotes where groups of rows with same number of nonzeros
471: * begin and end in iperm. */
472: PetscInt *nzgroup;
473: PetscInt ngroup;
474: PetscInt jstart, jend;
475: /* jstart is used in loops to denote the position in iperm where a
476: * group starts; jend denotes the position where it ends.
477: * (jend + 1 is where the next group starts.) */
478: PetscInt iold, nz;
479: PetscInt istart, iend, isize;
480: PetscInt ipos;
481: PetscScalar yp[NDIM];
482: PetscInt ip[NDIM];
483: /* yp[] and ip[] are treated as vector "registers" for performing
484: * the mat-vec. */
486: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
487: #pragma disjoint(*x, *y, *aa)
488: #endif
490: PetscFunctionBegin;
491: PetscCall(VecGetArrayRead(xx, &x));
492: PetscCall(VecGetArrayPair(yy, ww, &y, &w));
494: aj = a->j; /* aj[k] gives column index for element aa[k]. */
495: aa = a->a; /* Nonzero elements stored row-by-row. */
496: ai = a->i; /* ai[k] is the position in aa and aj where row k starts. */
498: /* Get the info we need about the permutations and groupings. */
499: aijperm = (Mat_SeqAIJPERM *)A->spptr;
500: iperm = aijperm->iperm;
501: ngroup = aijperm->ngroup;
502: xgroup = aijperm->xgroup;
503: nzgroup = aijperm->nzgroup;
505: for (PetscInt igroup = 0; igroup < ngroup; igroup++) {
506: jstart = xgroup[igroup];
507: jend = xgroup[igroup + 1] - 1;
509: nz = nzgroup[igroup];
511: /* Handle the special cases where the number of nonzeros per row
512: * in the group is either 0 or 1. */
513: if (nz == 0) {
514: for (i = jstart; i <= jend; i++) {
515: iold = iperm[i];
516: y[iold] = w[iold];
517: }
518: } else if (nz == 1) {
519: for (i = jstart; i <= jend; i++) {
520: iold = iperm[i];
521: ipos = ai[iold];
522: y[iold] = w[iold] + aa[ipos] * x[aj[ipos]];
523: }
524: }
525: /* For the general case: */
526: else {
527: /* We work our way through the current group in chunks of NDIM rows
528: * at a time. */
530: for (istart = jstart; istart <= jend; istart += NDIM) {
531: /* Figure out where the chunk of 'isize' rows ends in iperm.
532: * 'isize may of course be less than NDIM for the last chunk. */
533: iend = istart + (NDIM - 1);
534: if (iend > jend) iend = jend;
535: isize = iend - istart + 1;
537: /* Initialize the yp[] array that will be used to hold part of
538: * the permuted results vector, and figure out where in aa each
539: * row of the chunk will begin. */
540: for (i = 0; i < isize; i++) {
541: iold = iperm[istart + i];
542: /* iold is a row number from the matrix A *before* reordering. */
543: ip[i] = ai[iold];
544: /* ip[i] tells us where the ith row of the chunk begins in aa. */
545: yp[i] = w[iold];
546: }
548: /* If the number of zeros per row exceeds the number of rows in
549: * the chunk, we should vectorize along nz, that is, perform the
550: * mat-vec one row at a time as in the usual CSR case. */
551: if (nz > isize) {
552: #if defined(PETSC_HAVE_CRAY_VECTOR)
553: #pragma _CRI preferstream
554: #endif
555: for (i = 0; i < isize; i++) {
556: #if defined(PETSC_HAVE_CRAY_VECTOR)
557: #pragma _CRI prefervector
558: #endif
559: for (j = 0; j < nz; j++) {
560: ipos = ip[i] + j;
561: yp[i] += aa[ipos] * x[aj[ipos]];
562: }
563: }
564: }
565: /* Otherwise, there are enough rows in the chunk to make it
566: * worthwhile to vectorize across the rows, that is, to do the
567: * matvec by operating with "columns" of the chunk. */
568: else {
569: for (j = 0; j < nz; j++) {
570: for (i = 0; i < isize; i++) {
571: ipos = ip[i] + j;
572: yp[i] += aa[ipos] * x[aj[ipos]];
573: }
574: }
575: }
577: #if defined(PETSC_HAVE_CRAY_VECTOR)
578: #pragma _CRI ivdep
579: #endif
580: /* Put results from yp[] into non-permuted result vector y. */
581: for (i = 0; i < isize; i++) y[iperm[istart + i]] = yp[i];
582: } /* End processing chunk of isize rows of a group. */
584: } /* End handling matvec for chunk with nz > 1. */
585: } /* End loop over igroup. */
587: PetscCall(PetscLogFlops(2.0 * a->nz));
588: PetscCall(VecRestoreArrayRead(xx, &x));
589: PetscCall(VecRestoreArrayPair(yy, ww, &y, &w));
590: PetscFunctionReturn(PETSC_SUCCESS);
591: }
593: /* MatConvert_SeqAIJ_SeqAIJPERM converts a SeqAIJ matrix into a
594: * SeqAIJPERM matrix. This routine is called by the MatCreate_SeqAIJPERM()
595: * routine, but can also be used to convert an assembled SeqAIJ matrix
596: * into a SeqAIJPERM one. */
597: PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJPERM(Mat A, MatType type, MatReuse reuse, Mat *newmat)
598: {
599: Mat B = *newmat;
600: Mat_SeqAIJPERM *aijperm;
601: PetscBool sametype;
603: PetscFunctionBegin;
604: if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
605: PetscCall(PetscObjectTypeCompare((PetscObject)A, type, &sametype));
606: if (sametype) PetscFunctionReturn(PETSC_SUCCESS);
608: PetscCall(PetscNew(&aijperm));
609: B->spptr = (void *)aijperm;
611: /* Set function pointers for methods that we inherit from AIJ but override. */
612: B->ops->duplicate = MatDuplicate_SeqAIJPERM;
613: B->ops->assemblyend = MatAssemblyEnd_SeqAIJPERM;
614: B->ops->destroy = MatDestroy_SeqAIJPERM;
615: B->ops->mult = MatMult_SeqAIJPERM;
616: B->ops->multadd = MatMultAdd_SeqAIJPERM;
618: aijperm->nonzerostate = -1; /* this will trigger the generation of the permutation information the first time through MatAssembly()*/
619: /* If A has already been assembled, compute the permutation. */
620: if (A->assembled) PetscCall(MatSeqAIJPERM_create_perm(B));
622: PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqaijperm_seqaij_C", MatConvert_SeqAIJPERM_SeqAIJ));
624: PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQAIJPERM));
625: *newmat = B;
626: PetscFunctionReturn(PETSC_SUCCESS);
627: }
629: /*@C
630: MatCreateSeqAIJPERM - Creates a sparse matrix of type `MATSEQAIJPERM`.
632: Collective
634: Input Parameters:
635: + comm - MPI communicator, set to `PETSC_COMM_SELF`
636: . m - number of rows
637: . n - number of columns
638: . nz - number of nonzeros per row (same for all rows), ignored if `nnz` is given
639: - nnz - array containing the number of nonzeros in the various rows (possibly different for each row) or `NULL`
641: Output Parameter:
642: . A - the matrix
644: Level: intermediate
646: Notes:
647: This type inherits from `MATSEQAIJ`, but calculates some additional permutation information
648: that is used to allow better vectorization of some operations. At the cost of increased
649: storage, the `MATSEQAIJ` formatted matrix can be copied to a format in which pieces of the
650: matrix are stored in ELLPACK format, allowing the vectorized matrix multiply routine to use
651: stride-1 memory accesses.
653: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateMPIAIJPERM()`, `MatSetValues()`
654: @*/
655: PetscErrorCode MatCreateSeqAIJPERM(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt nz, const PetscInt nnz[], Mat *A)
656: {
657: PetscFunctionBegin;
658: PetscCall(MatCreate(comm, A));
659: PetscCall(MatSetSizes(*A, m, n, m, n));
660: PetscCall(MatSetType(*A, MATSEQAIJPERM));
661: PetscCall(MatSeqAIJSetPreallocation_SeqAIJ(*A, nz, nnz));
662: PetscFunctionReturn(PETSC_SUCCESS);
663: }
665: PETSC_EXTERN PetscErrorCode MatCreate_SeqAIJPERM(Mat A)
666: {
667: PetscFunctionBegin;
668: PetscCall(MatSetType(A, MATSEQAIJ));
669: PetscCall(MatConvert_SeqAIJ_SeqAIJPERM(A, MATSEQAIJPERM, MAT_INPLACE_MATRIX, &A));
670: PetscFunctionReturn(PETSC_SUCCESS);
671: }