Actual source code: sorder.c
1: /*$Id: sorder.c,v 1.90 2001/06/21 21:17:30 bsmith Exp $*/
2: /*
3: Provides the code that allows PETSc users to register their own
4: sequential matrix Ordering routines.
5: */
6: #include src/mat/matimpl.h
7: #include petscsys.h
9: PetscFList MatOrderingList = 0;
10: PetscTruth MatOrderingRegisterAllCalled = PETSC_FALSE;
12: EXTERN int MatOrdering_Flow_SeqAIJ(Mat,const MatOrderingType,IS *,IS *);
16: int MatOrdering_Flow(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
17: {
19: SETERRQ(PETSC_ERR_SUP,"Cannot do default flow ordering for matrix type");
20: #if !defined(PETSC_USE_DEBUG)
21: return(0);
22: #endif
23: }
25: EXTERN_C_BEGIN
28: int MatOrdering_Natural(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
29: {
30: int n,ierr,i,*ii;
31: PetscTruth done;
32: MPI_Comm comm;
35: PetscObjectGetComm((PetscObject)mat,&comm);
36: MatGetRowIJ(mat,0,PETSC_FALSE,&n,PETSC_NULL,PETSC_NULL,&done);
37: MatRestoreRowIJ(mat,0,PETSC_FALSE,&n,PETSC_NULL,PETSC_NULL,&done);
38: if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */
39: /*
40: We actually create general index sets because this avoids mallocs to
41: to obtain the indices in the MatSolve() routines.
42: ISCreateStride(PETSC_COMM_SELF,n,0,1,irow);
43: ISCreateStride(PETSC_COMM_SELF,n,0,1,icol);
44: */
45: PetscMalloc(n*sizeof(int),&ii);
46: for (i=0; i<n; i++) ii[i] = i;
47: ISCreateGeneral(PETSC_COMM_SELF,n,ii,irow);
48: ISCreateGeneral(PETSC_COMM_SELF,n,ii,icol);
49: PetscFree(ii);
50: } else {
51: int start,end;
53: MatGetOwnershipRange(mat,&start,&end);
54: ISCreateStride(comm,end-start,start,1,irow);
55: ISCreateStride(comm,end-start,start,1,icol);
56: }
57: ISSetIdentity(*irow);
58: ISSetIdentity(*icol);
59: return(0);
60: }
61: EXTERN_C_END
63: EXTERN_C_BEGIN
64: /*
65: Orders the rows (and columns) by the lengths of the rows.
66: This produces a symmetric Ordering but does not require a
67: matrix with symmetric non-zero structure.
68: */
71: int MatOrdering_RowLength(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
72: {
73: int ierr,n,*ia,*ja,*permr,*lens,i;
74: PetscTruth done;
77: MatGetRowIJ(mat,0,PETSC_FALSE,&n,&ia,&ja,&done);
78: if (!done) SETERRQ(PETSC_ERR_SUP,"Cannot get rows for matrix");
80: PetscMalloc(2*n*sizeof(int),&lens);
81: permr = lens + n;
82: for (i=0; i<n; i++) {
83: lens[i] = ia[i+1] - ia[i];
84: permr[i] = i;
85: }
86: MatRestoreRowIJ(mat,0,PETSC_FALSE,&n,&ia,&ja,&done);
88: PetscSortIntWithPermutation(n,lens,permr);
90: ISCreateGeneral(PETSC_COMM_SELF,n,permr,irow);
91: ISCreateGeneral(PETSC_COMM_SELF,n,permr,icol);
92: PetscFree(lens);
93: return(0);
94: }
95: EXTERN_C_END
99: int MatOrderingRegister(const char sname[],const char path[],const char name[],int (*function)(Mat,const MatOrderingType,IS*,IS*))
100: {
101: int ierr;
102: char fullname[256];
105: PetscFListConcat(path,name,fullname);
106: PetscFListAdd(&MatOrderingList,sname,fullname,(void (*)(void))function);
107: return(0);
108: }
112: /*@C
113: MatOrderingRegisterDestroy - Frees the list of ordering routines.
115: Not collective
117: Level: developer
118:
119: .keywords: matrix, register, destroy
121: .seealso: MatOrderingRegisterDynamic(), MatOrderingRegisterAll()
122: @*/
123: int MatOrderingRegisterDestroy(void)
124: {
128: if (MatOrderingList) {
129: PetscFListDestroy(&MatOrderingList);
130: MatOrderingList = 0;
131: }
132: return(0);
133: }
135: EXTERN int MatAdjustForInodes(Mat,IS *,IS *);
137: #include src/mat/impls/aij/mpi/mpiaij.h
140: /*@C
141: MatGetOrdering - Gets a reordering for a matrix to reduce fill or to
142: improve numerical stability of LU factorization.
144: Collective on Mat
146: Input Parameters:
147: + mat - the matrix
148: - type - type of reordering, one of the following:
149: $ MATORDERING_NATURAL - Natural
150: $ MATORDERING_ND - Nested Dissection
151: $ MATORDERING_1WD - One-way Dissection
152: $ MATORDERING_RCM - Reverse Cuthill-McKee
153: $ MATORDERING_QMD - Quotient Minimum Degree
155: Output Parameters:
156: + rperm - row permutation indices
157: - cperm - column permutation indices
160: Options Database Key:
161: . -mat_view_ordering_draw - plots matrix nonzero structure in new ordering
163: Level: intermediate
164:
165: Notes:
166: This DOES NOT actually reorder the matrix; it merely returns two index sets
167: that define a reordering. This is usually not used directly, rather use the
168: options PCLUSetMatOrdering() or PCILUSetMatOrdering().
170: The user can define additional orderings; see MatOrderingRegisterDynamic().
172: .keywords: matrix, set, ordering, factorization, direct, ILU, LU,
173: fill, reordering, natural, Nested Dissection,
174: One-way Dissection, Cholesky, Reverse Cuthill-McKee,
175: Quotient Minimum Degree
177: .seealso: MatOrderingRegisterDynamic(), PCLUSetMatOrdering(), PCILUSetMatOrdering()
178: @*/
179: int MatGetOrdering(Mat mat,const MatOrderingType type,IS *rperm,IS *cperm)
180: {
181: int ierr,mmat,nmat,mis,m;
182: int (*r)(Mat,const MatOrderingType,IS*,IS*);
183: PetscTruth flg,isseqdense,ismpidense;
189: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
190: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
192: PetscTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);
193: PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&ismpidense);
194: if (isseqdense || ismpidense) {
195: MatGetLocalSize(mat,&m,PETSC_NULL);
196: /*
197: Dense matrices only give natural ordering
198: */
199: ISCreateStride(PETSC_COMM_SELF,0,m,1,cperm);
200: ISCreateStride(PETSC_COMM_SELF,0,m,1,rperm);
201: ISSetIdentity(*cperm);
202: ISSetIdentity(*rperm);
203: ISSetPermutation(*rperm);
204: ISSetPermutation(*cperm);
205: return(0);
206: }
208: if (!mat->M) { /* matrix has zero rows */
209: ISCreateStride(PETSC_COMM_SELF,0,0,1,cperm);
210: ISCreateStride(PETSC_COMM_SELF,0,0,1,rperm);
211: ISSetIdentity(*cperm);
212: ISSetIdentity(*rperm);
213: ISSetPermutation(*rperm);
214: ISSetPermutation(*cperm);
215: return(0);
216: }
218: if (!MatOrderingRegisterAllCalled) {
219: MatOrderingRegisterAll(PETSC_NULL);
220: }
222: PetscLogEventBegin(MAT_GetOrdering,mat,0,0,0);
223: PetscFListFind(mat->comm,MatOrderingList,type,(void (**)(void)) &r);
224: if (!r) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown or unregistered type: %s",type);}
226: (*r)(mat,type,rperm,cperm);
227: ISSetPermutation(*rperm);
228: ISSetPermutation(*cperm);
230: /*
231: Adjust for inode (reduced matrix ordering) only if row permutation
232: is smaller then matrix size
233: */
234: MatGetLocalSize(mat,&mmat,&nmat);
235: ISGetLocalSize(*rperm,&mis);
236: if (mmat > mis) {
237: MatAdjustForInodes(mat,rperm,cperm);
238: }
240: PetscLogEventEnd(MAT_GetOrdering,mat,0,0,0);
242: PetscOptionsHasName(PETSC_NULL,"-mat_view_ordering_draw",&flg);
243: if (flg) {
244: Mat tmat;
245: PetscOptionsHasName(PETSC_NULL,"-mat_view_contour",&flg);
246: if (flg) {
247: PetscViewerPushFormat(PETSC_VIEWER_DRAW_(mat->comm),PETSC_VIEWER_DRAW_CONTOUR);
248: }
249: MatPermute(mat,*rperm,*cperm,&tmat);
250: MatView(tmat,PETSC_VIEWER_DRAW_(mat->comm));
251: if (flg) {
252: PetscViewerPopFormat(PETSC_VIEWER_DRAW_(mat->comm));
253: }
254: MatDestroy(tmat);
255: }
257: return(0);
258: }
262: int MatGetOrderingList(PetscFList *list)
263: {
265: *list = MatOrderingList;
266: return(0);
267: }