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