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,MatOrderingType,IS *,IS *);
14: int MatOrdering_Flow(Mat mat,MatOrderingType type,IS *irow,IS *icol)
15: {
17: SETERRQ(PETSC_ERR_SUP,"Cannot do default flow ordering for matrix type");
18: #if !defined(PETSC_USE_DEBUG)
19: return(0);
20: #endif
21: }
23: EXTERN_C_BEGIN
24: int MatOrdering_Natural(Mat mat,MatOrderingType type,IS *irow,IS *icol)
25: {
26: int n,ierr,i,*ii;
27: PetscTruth done;
28: MPI_Comm comm;
31: PetscObjectGetComm((PetscObject)mat,&comm);
32: MatGetRowIJ(mat,0,PETSC_FALSE,&n,PETSC_NULL,PETSC_NULL,&done);
33: MatRestoreRowIJ(mat,0,PETSC_FALSE,&n,PETSC_NULL,PETSC_NULL,&done);
34: if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */
35: /*
36: We actually create general index sets because this avoids mallocs to
37: to obtain the indices in the MatSolve() routines.
38: ISCreateStride(PETSC_COMM_SELF,n,0,1,irow);
39: ISCreateStride(PETSC_COMM_SELF,n,0,1,icol);
40: */
41: PetscMalloc(n*sizeof(int),&ii);
42: for (i=0; i<n; i++) ii[i] = i;
43: ISCreateGeneral(PETSC_COMM_SELF,n,ii,irow);
44: ISCreateGeneral(PETSC_COMM_SELF,n,ii,icol);
45: PetscFree(ii);
46: } else {
47: int start,end;
49: MatGetOwnershipRange(mat,&start,&end);
50: ISCreateStride(comm,end-start,start,1,irow);
51: ISCreateStride(comm,end-start,start,1,icol);
52: }
53: ISSetIdentity(*irow);
54: ISSetIdentity(*icol);
55: return(0);
56: }
57: EXTERN_C_END
59: EXTERN_C_BEGIN
60: /*
61: Orders the rows (and columns) by the lengths of the rows.
62: This produces a symmetric Ordering but does not require a
63: matrix with symmetric non-zero structure.
64: */
65: int MatOrdering_RowLength(Mat mat,MatOrderingType type,IS *irow,IS *icol)
66: {
67: int ierr,n,*ia,*ja,*permr,*lens,i;
68: PetscTruth done;
71: MatGetRowIJ(mat,0,PETSC_FALSE,&n,&ia,&ja,&done);
72: if (!done) SETERRQ(PETSC_ERR_SUP,"Cannot get rows for matrix");
74: ierr = PetscMalloc(2*n*sizeof(int),&lens);
75: permr = lens + n;
76: for (i=0; i<n; i++) {
77: lens[i] = ia[i+1] - ia[i];
78: permr[i] = i;
79: }
80: MatRestoreRowIJ(mat,0,PETSC_FALSE,&n,&ia,&ja,&done);
82: PetscSortIntWithPermutation(n,lens,permr);
84: ISCreateGeneral(PETSC_COMM_SELF,n,permr,irow);
85: ISCreateGeneral(PETSC_COMM_SELF,n,permr,icol);
86: PetscFree(lens);
87: return(0);
88: }
89: EXTERN_C_END
91: /*MC
92: MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
93: matrix package.
95: Synopsis:
96: int MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,int (*routine_create)(MatOrdering))
98: Not Collective
100: Input Parameters:
101: + sname - name of ordering (for example MATORDERING_ND)
102: . path - location of library where creation routine is
103: . name - name of function that creates the ordering type,a string
104: - function - function pointer that creates the ordering
106: Level: developer
108: If dynamic libraries are used, then the fourth input argument (function)
109: is ignored.
111: Sample usage:
112: .vb
113: MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
114: "MyOrder",MyOrder);
115: .ve
117: Then, your partitioner can be chosen with the procedural interface via
118: $ MatOrderingSetType(part,"my_order)
119: or at runtime via the option
120: $ -pc_ilu_mat_ordering_type my_order
121: $ -pc_lu_mat_ordering_type my_order
123: ${PETSC_ARCH} and ${BOPT} occuring in pathname will be replaced with appropriate values.
125: .keywords: matrix, ordering, register
127: .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
128: M*/
130: int MatOrderingRegister(char *sname,char *path,char *name,int (*function)(Mat,MatOrderingType,IS*,IS*))
131: {
132: int ierr;
133: char fullname[256];
136: PetscFListConcat(path,name,fullname);
137: PetscFListAdd(&MatOrderingList,sname,fullname,(void (*)(void))function);
138: return(0);
139: }
141: /*@C
142: MatOrderingRegisterDestroy - Frees the list of ordering routines.
144: Not collective
146: Level: developer
147:
148: .keywords: matrix, register, destroy
150: .seealso: MatOrderingRegisterDynamic(), MatOrderingRegisterAll()
151: @*/
152: int MatOrderingRegisterDestroy(void)
153: {
157: if (MatOrderingList) {
158: PetscFListDestroy(&MatOrderingList);
159: MatOrderingList = 0;
160: }
161: return(0);
162: }
164: EXTERN int MatAdjustForInodes(Mat,IS *,IS *);
166: #include src/mat/impls/aij/mpi/mpiaij.h
167: /*@C
168: MatGetOrdering - Gets a reordering for a matrix to reduce fill or to
169: improve numerical stability of LU factorization.
171: Collective on Mat
173: Input Parameters:
174: + mat - the matrix
175: - type - type of reordering, one of the following:
176: $ MATORDERING_NATURAL - Natural
177: $ MATORDERING_ND - Nested Dissection
178: $ MATORDERING_1WD - One-way Dissection
179: $ MATORDERING_RCM - Reverse Cuthill-McKee
180: $ MATORDERING_QMD - Quotient Minimum Degree
182: Output Parameters:
183: + rperm - row permutation indices
184: - cperm - column permutation indices
187: Options Database:
188: . -mat_view_ordering_draw - plots matrix nonzero structure in new ordering
190: Level: intermediate
191:
192: Notes:
193: This DOES NOT actually reorder the matrix; it merely returns two index sets
194: that define a reordering. This is usually not used directly, rather use the
195: options PCLUSetMatOrdering() or PCILUSetMatOrdering().
197: The user can define additional orderings; see MatOrderingRegisterDynamic().
199: .keywords: matrix, set, ordering, factorization, direct, ILU, LU,
200: fill, reordering, natural, Nested Dissection,
201: One-way Dissection, Cholesky, Reverse Cuthill-McKee,
202: Quotient Minimum Degree
204: .seealso: MatOrderingRegisterDynamic(), PCLUSetMatOrdering(), PCILUSetMatOrdering()
205: @*/
206: int MatGetOrdering(Mat mat,MatOrderingType type,IS *rperm,IS *cperm)
207: {
208: int ierr,mmat,nmat,mis,m;
209: int (*r)(Mat,MatOrderingType,IS*,IS*);
210: PetscTruth flg,isseqdense,ismpidense;
214: if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
215: if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
217: PetscTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);
218: PetscTypeCompare((PetscObject)mat,MATMPIDENSE,&ismpidense);
219: if (isseqdense || ismpidense) {
220: MatGetLocalSize(mat,&m,PETSC_NULL);
221: /*
222: Dense matrices only give natural ordering
223: */
224: ISCreateStride(PETSC_COMM_SELF,0,m,1,cperm);
225: ISCreateStride(PETSC_COMM_SELF,0,m,1,rperm);
226: ISSetIdentity(*cperm);
227: ISSetIdentity(*rperm);
228: ISSetPermutation(*rperm);
229: ISSetPermutation(*cperm);
230: return(0);
231: }
233: if (!mat->M) { /* matrix has zero rows */
234: ISCreateStride(PETSC_COMM_SELF,0,0,1,cperm);
235: ISCreateStride(PETSC_COMM_SELF,0,0,1,rperm);
236: ISSetIdentity(*cperm);
237: ISSetIdentity(*rperm);
238: ISSetPermutation(*rperm);
239: ISSetPermutation(*cperm);
240: return(0);
241: }
243: if (!MatOrderingRegisterAllCalled) {
244: MatOrderingRegisterAll(PETSC_NULL);
245: }
247: PetscLogEventBegin(MAT_GetOrdering,mat,0,0,0);
248: PetscFListFind(mat->comm,MatOrderingList,type,(void (**)(void)) &r);
249: if (!r) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Unknown or unregistered type: %s",type);}
251: (*r)(mat,type,rperm,cperm);
252: ISSetPermutation(*rperm);
253: ISSetPermutation(*cperm);
255: /*
256: Adjust for inode (reduced matrix ordering) only if row permutation
257: is smaller then matrix size
258: */
259: MatGetLocalSize(mat,&mmat,&nmat);
260: ISGetLocalSize(*rperm,&mis);
261: if (mmat > mis) {
262: MatAdjustForInodes(mat,rperm,cperm);
263: }
265: PetscLogEventEnd(MAT_GetOrdering,mat,0,0,0);
267: PetscOptionsHasName(PETSC_NULL,"-mat_view_ordering_draw",&flg);
268: if (flg) {
269: Mat tmat;
270: PetscOptionsHasName(PETSC_NULL,"-mat_view_contour",&flg);
271: if (flg) {
272: PetscViewerPushFormat(PETSC_VIEWER_DRAW_(mat->comm),PETSC_VIEWER_DRAW_CONTOUR);
273: }
274: MatPermute(mat,*rperm,*cperm,&tmat);
275: MatView(tmat,PETSC_VIEWER_DRAW_(mat->comm));
276: if (flg) {
277: PetscViewerPopFormat(PETSC_VIEWER_DRAW_(mat->comm));
278: }
279: MatDestroy(tmat);
280: }
282: return(0);
283: }
285: int MatGetOrderingList(PetscFList *list)
286: {
288: *list = MatOrderingList;
289: return(0);
290: }