2: #include <../src/mat/impls/baij/seq/baij.h>
3: #include <../src/mat/impls/sbaij/seq/sbaij.h>
5: /*
6: This function is used before applying a
7: symmetric reordering to matrix A that is
8: in SBAIJ format.
10: The permutation is assumed to be symmetric, i.e.,
11: P = P^T (= inv(P)),
12: so the permuted matrix P*A*inv(P)=P*A*P^T is ensured to be symmetric.
13: - a wrong assumption! This code needs rework! -- Hong
15: The function is modified from sro.f of YSMP. The description from YSMP:
16: C THE NONZERO ENTRIES OF THE MATRIX M ARE ASSUMED TO BE STORED 17: C SYMMETRICALLY IN (IA,JA,A) FORMAT (I.E., NOT BOTH M(I,J) AND M(J,I)
18: C ARE STORED IF I NE J).
19: C
20: C SRO DOES NOT REARRANGE THE ORDER OF THE ROWS, BUT DOES MOVE
21: C NONZEROES FROM ONE ROW TO ANOTHER TO ENSURE THAT IF M(I,J) WILL BE 22: C IN THE UPPER TRIANGLE OF M WITH RESPECT TO THE NEW ORDERING, THEN 23: C M(I,J) IS STORED IN ROW I (AND THUS M(J,I) IS NOT STORED); WHEREAS
24: C IF M(I,J) WILL BE IN THE STRICT LOWER TRIANGLE OF M, THEN M(J,I) IS 25: C STORED IN ROW J (AND THUS M(I,J) IS NOT STORED).
28: -- output: new index set (inew, jnew) for A and a map a2anew that maps
29: values a to anew, such that all
30: nonzero A_(perm(i),iperm(k)) will be stored in the upper triangle.
31: Note: matrix A is not permuted by this function!
32: */
35: PetscErrorCode MatReorderingSeqSBAIJ(Mat A,IS perm) 36: {
37: Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
38: const PetscInt mbs=a->mbs;
41: if (!mbs) return(0);
42: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
43: #if 0
45: const PetscInt *rip,*riip;
46: PetscInt *ai,*aj,*r;
47: PetscInt *nzr,nz,jmin,jmax,j,k,ajk,i;
48: IS iperm; /* inverse of perm */
49: ISGetIndices(perm,&rip);
51: ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
52: ISGetIndices(iperm,&riip);
54: for (i=0; i<mbs; i++) {
55: if (rip[i] != riip[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, use symmetric permutation for symmetric matrices");
56: }
57: ISRestoreIndices(iperm,&riip);
58: ISDestroy(&iperm);
60: if (!a->inew) {
61: PetscMalloc2(mbs+1,&ai, 2*a->i[mbs],&aj);
62: } else {
63: ai = a->inew; aj = a->jnew;
64: }
65: PetscMemcpy(ai,a->i,(mbs+1)*sizeof(PetscInt));
66: PetscMemcpy(aj,a->j,(a->i[mbs])*sizeof(PetscInt));
68: /*
69: Phase 1: Find row index r in which to store each nonzero.
70: Initialize count of nonzeros to be stored in each row (nzr).
71: At the end of this phase, a nonzero a(*,*)=a(r(),aj())
72: s.t. a(perm(r),perm(aj)) will fall into upper triangle part.
73: */
75: PetscMalloc1(mbs,&nzr);
76: PetscMalloc1(ai[mbs],&r);
77: for (i=0; i<mbs; i++) nzr[i] = 0;
78: for (i=0; i<ai[mbs]; i++) r[i] = 0;
80: /* for each nonzero element */
81: for (i=0; i<mbs; i++) {
82: nz = ai[i+1] - ai[i];
83: j = ai[i];
84: /* printf("nz = %d, j=%d\n",nz,j); */
85: while (nz--) {
86: /* --- find row (=r[j]) and column (=aj[j]) in which to store a[j] ...*/
87: k = aj[j]; /* col. index */
88: /* printf("nz = %d, k=%d\n", nz,k); */
89: /* for entry that will be permuted into lower triangle, swap row and col. index */
90: if (rip[k] < rip[i]) aj[j] = i;
91: else k = i;
93: r[j] = k; j++;
94: nzr[k]++; /* increment count of nonzeros in that row */
95: }
96: }
98: /* Phase 2: Find new ai and permutation to apply to (aj,a).
99: Determine pointers (r) to delimit rows in permuted (aj,a).
100: Note: r is different from r used in phase 1.
101: At the end of this phase, (aj[j],a[j]) will be stored in
102: (aj[r(j)],a[r(j)]).
103: */
104: for (i=0; i<mbs; i++) {
105: ai[i+1] = ai[i] + nzr[i];
106: nzr[i] = ai[i+1];
107: }
109: /* determine where each (aj[j], a[j]) is stored in new (aj,a)
110: for each nonzero element (in reverse order) */
111: jmin = ai[0]; jmax = ai[mbs];
112: nz = jmax - jmin;
113: j = jmax-1;
114: while (nz--) {
115: i = r[j]; /* row value */
116: if (aj[j] == i) r[j] = ai[i]; /* put diagonal nonzero at beginning of row */
117: else { /* put off-diagonal nonzero in last unused location in row */
118: nzr[i]--; r[j] = nzr[i];
119: }
120: j--;
121: }
123: a->a2anew = aj + ai[mbs];
124: PetscMemcpy(a->a2anew,r,ai[mbs]*sizeof(PetscInt));
126: /* Phase 3: permute (aj,a) to upper triangular form (wrt new ordering) */
127: for (j=jmin; j<jmax; j++) {
128: while (r[j] != j) {
129: k = r[j]; r[j] = r[k]; r[k] = k;
130: ajk = aj[k]; aj[k] = aj[j]; aj[j] = ajk;
131: /* ak = aa[k]; aa[k] = aa[j]; aa[j] = ak; */
132: }
133: }
134: ierr= ISRestoreIndices(perm,&rip);
136: a->inew = ai;
137: a->jnew = aj;
139: ISDestroy(&a->row);
140: ISDestroy(&a->icol);
141: PetscObjectReference((PetscObject)perm);
142: ISDestroy(&a->row);
143: a->row = perm;
144: PetscObjectReference((PetscObject)perm);
145: ISDestroy(&a->icol);
146: a->icol = perm;
148: PetscFree(nzr);
149: PetscFree(r);
150: return(0);
151: #endif
152: }