Actual source code: sorder.c
1: /*
2: Provides the code that allows PETSc users to register their own
3: sequential matrix Ordering routines.
4: */
5: #include src/mat/matimpl.h
6: #include petscsys.h
8: PetscFList MatOrderingList = 0;
9: PetscTruth MatOrderingRegisterAllCalled = PETSC_FALSE;
11: EXTERN PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat,const MatOrderingType,IS *,IS *);
15: PetscErrorCode MatOrdering_Flow(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
16: {
18: SETERRQ(PETSC_ERR_SUP,"Cannot do default flow ordering for matrix type");
19: #if !defined(PETSC_USE_DEBUG)
20: return(0);
21: #endif
22: }
27: PetscErrorCode MatOrdering_Natural(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
28: {
30: PetscInt n,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(PetscInt),&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: PetscInt 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: }
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: PetscErrorCode MatOrdering_RowLength(Mat mat,const MatOrderingType type,IS *irow,IS *icol)
72: {
74: PetscInt n,*ia,*ja,*permr,*lens,i;
75: PetscTruth done;
78: MatGetRowIJ(mat,0,PETSC_FALSE,&n,&ia,&ja,&done);
79: if (!done) SETERRQ(PETSC_ERR_SUP,"Cannot get rows for matrix");
81: PetscMalloc(2*n*sizeof(PetscInt),&lens);
82: permr = lens + n;
83: for (i=0; i<n; i++) {
84: lens[i] = ia[i+1] - ia[i];
85: permr[i] = i;
86: }
87: MatRestoreRowIJ(mat,0,PETSC_FALSE,&n,&ia,&ja,&done);
89: PetscSortIntWithPermutation(n,lens,permr);
91: ISCreateGeneral(PETSC_COMM_SELF,n,permr,irow);
92: ISCreateGeneral(PETSC_COMM_SELF,n,permr,icol);
93: PetscFree(lens);
94: return(0);
95: }
100: PetscErrorCode MatOrderingRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(Mat,const MatOrderingType,IS*,IS*))
101: {
103: char fullname[PETSC_MAX_PATH_LEN];
106: PetscFListConcat(path,name,fullname);
107: PetscFListAdd(&MatOrderingList,sname,fullname,(void (*)(void))function);
108: return(0);
109: }
113: /*@C
114: MatOrderingRegisterDestroy - Frees the list of ordering routines.
116: Not collective
118: Level: developer
119:
120: .keywords: matrix, register, destroy
122: .seealso: MatOrderingRegisterDynamic(), MatOrderingRegisterAll()
123: @*/
124: PetscErrorCode MatOrderingRegisterDestroy(void)
125: {
129: if (MatOrderingList) {
130: PetscFListDestroy(&MatOrderingList);
131: MatOrderingList = 0;
132: }
133: return(0);
134: }
136: EXTERN PetscErrorCode MatAdjustForInodes(Mat,IS *,IS *);
138: #include src/mat/impls/aij/mpi/mpiaij.h
141: /*@C
142: MatGetOrdering - Gets a reordering for a matrix to reduce fill or to
143: improve numerical stability of LU factorization.
145: Collective on Mat
147: Input Parameters:
148: + mat - the matrix
149: - type - type of reordering, one of the following:
150: $ MATORDERING_NATURAL - Natural
151: $ MATORDERING_ND - Nested Dissection
152: $ MATORDERING_1WD - One-way Dissection
153: $ MATORDERING_RCM - Reverse Cuthill-McKee
154: $ MATORDERING_QMD - Quotient Minimum Degree
156: Output Parameters:
157: + rperm - row permutation indices
158: - cperm - column permutation indices
161: Options Database Key:
162: . -mat_view_ordering_draw - plots matrix nonzero structure in new ordering
164: Level: intermediate
165:
166: Notes:
167: This DOES NOT actually reorder the matrix; it merely returns two index sets
168: that define a reordering. This is usually not used directly, rather use the
169: options PCLUSetMatOrdering() or PCILUSetMatOrdering().
171: The user can define additional orderings; see MatOrderingRegisterDynamic().
173: .keywords: matrix, set, ordering, factorization, direct, ILU, LU,
174: fill, reordering, natural, Nested Dissection,
175: One-way Dissection, Cholesky, Reverse Cuthill-McKee,
176: Quotient Minimum Degree
178: .seealso: MatOrderingRegisterDynamic(), PCLUSetMatOrdering(), PCILUSetMatOrdering()
179: @*/
180: PetscErrorCode MatGetOrdering(Mat mat,const MatOrderingType type,IS *rperm,IS *cperm)
181: {
182: PetscErrorCode ierr;
183: PetscInt mmat,nmat,mis,m;
184: PetscErrorCode (*r)(Mat,const MatOrderingType,IS*,IS*);
185: PetscTruth flg,isseqdense,ismpidense;
191: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
192: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
194: PetscTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);
195: PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&ismpidense);
196: if (isseqdense || ismpidense) {
197: MatGetLocalSize(mat,&m,PETSC_NULL);
198: /*
199: Dense matrices only give natural ordering
200: */
201: ISCreateStride(PETSC_COMM_SELF,0,m,1,cperm);
202: ISCreateStride(PETSC_COMM_SELF,0,m,1,rperm);
203: ISSetIdentity(*cperm);
204: ISSetIdentity(*rperm);
205: ISSetPermutation(*rperm);
206: ISSetPermutation(*cperm);
207: return(0);
208: }
210: if (!mat->M) { /* matrix has zero rows */
211: ISCreateStride(PETSC_COMM_SELF,0,0,1,cperm);
212: ISCreateStride(PETSC_COMM_SELF,0,0,1,rperm);
213: ISSetIdentity(*cperm);
214: ISSetIdentity(*rperm);
215: ISSetPermutation(*rperm);
216: ISSetPermutation(*cperm);
217: return(0);
218: }
220: if (!MatOrderingRegisterAllCalled) {
221: MatOrderingRegisterAll(PETSC_NULL);
222: }
224: PetscLogEventBegin(MAT_GetOrdering,mat,0,0,0);
225: PetscFListFind(mat->comm,MatOrderingList,type,(void (**)(void)) &r);
226: if (!r) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown or unregistered type: %s",type);}
228: (*r)(mat,type,rperm,cperm);
229: ISSetPermutation(*rperm);
230: ISSetPermutation(*cperm);
232: /*
233: Adjust for inode (reduced matrix ordering) only if row permutation
234: is smaller then matrix size
235: */
236: MatGetLocalSize(mat,&mmat,&nmat);
237: ISGetLocalSize(*rperm,&mis);
238: if (mmat > mis) {
239: MatAdjustForInodes(mat,rperm,cperm);
240: }
242: PetscLogEventEnd(MAT_GetOrdering,mat,0,0,0);
244: PetscOptionsHasName(PETSC_NULL,"-mat_view_ordering_draw",&flg);
245: if (flg) {
246: Mat tmat;
247: PetscOptionsHasName(PETSC_NULL,"-mat_view_contour",&flg);
248: if (flg) {
249: PetscViewerPushFormat(PETSC_VIEWER_DRAW_(mat->comm),PETSC_VIEWER_DRAW_CONTOUR);
250: }
251: MatPermute(mat,*rperm,*cperm,&tmat);
252: MatView(tmat,PETSC_VIEWER_DRAW_(mat->comm));
253: if (flg) {
254: PetscViewerPopFormat(PETSC_VIEWER_DRAW_(mat->comm));
255: }
256: MatDestroy(tmat);
257: }
259: return(0);
260: }
264: PetscErrorCode MatGetOrderingList(PetscFList *list)
265: {
267: *list = MatOrderingList;
268: return(0);
269: }