Actual source code: sorder.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:      Provides the code that allows PETSc users to register their own
  4:   sequential matrix Ordering routines.
  5: */
  6:  #include <petsc/private/matimpl.h>
  7:  #include <petscmat.h>

  9: PetscFunctionList MatOrderingList              = 0;
 10: PetscBool         MatOrderingRegisterAllCalled = PETSC_FALSE;

 12: extern PetscErrorCode MatGetOrdering_Flow_SeqAIJ(Mat,MatOrderingType,IS*,IS*);

 14: PetscErrorCode MatGetOrdering_Flow(Mat mat,MatOrderingType type,IS *irow,IS *icol)
 15: {
 17:   SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot do default flow ordering for matrix type");
 18: #if !defined(PETSC_USE_DEBUG)
 19:   return(0);
 20: #endif
 21: }

 23: PETSC_INTERN PetscErrorCode MatGetOrdering_Natural(Mat mat,MatOrderingType type,IS *irow,IS *icol)
 24: {
 26:   PetscInt       n,i,*ii;
 27:   PetscBool      done;
 28:   MPI_Comm       comm;

 31:   PetscObjectGetComm((PetscObject)mat,&comm);
 32:   MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done);
 33:   MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,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:     PetscMalloc1(n,&ii);
 42:     for (i=0; i<n; i++) ii[i] = i;
 43:     ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow);
 44:     ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol);
 45:   } else {
 46:     PetscInt start,end;

 48:     MatGetOwnershipRange(mat,&start,&end);
 49:     ISCreateStride(comm,end-start,start,1,irow);
 50:     ISCreateStride(comm,end-start,start,1,icol);
 51:   }
 52:   ISSetIdentity(*irow);
 53:   ISSetIdentity(*icol);
 54:   return(0);
 55: }

 57: /*
 58:      Orders the rows (and columns) by the lengths of the rows.
 59:    This produces a symmetric Ordering but does not require a
 60:    matrix with symmetric non-zero structure.
 61: */
 62: PETSC_INTERN PetscErrorCode MatGetOrdering_RowLength(Mat mat,MatOrderingType type,IS *irow,IS *icol)
 63: {
 65:   PetscInt       n,*permr,*lens,i;
 66:   const PetscInt *ia,*ja;
 67:   PetscBool      done;

 70:   MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,&ia,&ja,&done);
 71:   if (!done) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get rows for matrix");

 73:   PetscMalloc2(n,&lens,n,&permr);
 74:   for (i=0; i<n; i++) {
 75:     lens[i]  = ia[i+1] - ia[i];
 76:     permr[i] = i;
 77:   }
 78:   MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,&ia,&ja,&done);

 80:   PetscSortIntWithPermutation(n,lens,permr);

 82:   ISCreateGeneral(PETSC_COMM_SELF,n,permr,PETSC_COPY_VALUES,irow);
 83:   ISCreateGeneral(PETSC_COMM_SELF,n,permr,PETSC_COPY_VALUES,icol);
 84:   PetscFree2(lens,permr);
 85:   return(0);
 86: }

 88: /*@C
 89:    MatOrderingRegister - Adds a new sparse matrix ordering to the matrix package.

 91:    Not Collective

 93:    Input Parameters:
 94: +  sname - name of ordering (for example MATORDERINGND)
 95: -  function - function pointer that creates the ordering

 97:    Level: developer

 99:    Sample usage:
100: .vb
101:    MatOrderingRegister("my_order", MyOrder);
102: .ve

104:    Then, your partitioner can be chosen with the procedural interface via
105: $     MatOrderingSetType(part,"my_order)
106:    or at runtime via the option
107: $     -pc_factor_mat_ordering_type my_order

109: .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
110: @*/
111: PetscErrorCode  MatOrderingRegister(const char sname[],PetscErrorCode (*function)(Mat,MatOrderingType,IS*,IS*))
112: {

116:   MatInitializePackage();
117:   PetscFunctionListAdd(&MatOrderingList,sname,function);
118:   return(0);
119: }

121:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
122: /*@C
123:    MatGetOrdering - Gets a reordering for a matrix to reduce fill or to
124:    improve numerical stability of LU factorization.

126:    Collective on Mat

128:    Input Parameters:
129: +  mat - the matrix
130: -  type - type of reordering, one of the following:
131: $      MATORDERINGNATURAL_OR_ND - Nested dissection unless matrix is SBAIJ then it is natural
132: $      MATORDERINGNATURAL - Natural
133: $      MATORDERINGND - Nested Dissection
134: $      MATORDERING1WD - One-way Dissection
135: $      MATORDERINGRCM - Reverse Cuthill-McKee
136: $      MATORDERINGQMD - Quotient Minimum Degree

138:    Output Parameters:
139: +  rperm - row permutation indices
140: -  cperm - column permutation indices

142:    Options Database Key:
143: + -mat_view_ordering draw - plots matrix nonzero structure in new ordering
144: - -pc_factor_mat_ordering_type <nd,natural,..> - ordering to use with PCs based on factorization, LU, ILU, Cholesky, ICC

146:    Level: intermediate

148:    Notes:
149:       This DOES NOT actually reorder the matrix; it merely returns two index sets
150:    that define a reordering. This is usually not used directly, rather use the
151:    options PCFactorSetMatOrderingType()

153:    The user can define additional orderings; see MatOrderingRegister().

155:    These are generally only implemented for sequential sparse matrices.

157:    The external packages that PETSc can use for direct factorization such as SuperLU do not accept orderings provided by
158:    this call.


161:            fill, reordering, natural, Nested Dissection,
162:            One-way Dissection, Cholesky, Reverse Cuthill-McKee,
163:            Quotient Minimum Degree

165: .seealso:   MatOrderingRegister(), PCFactorSetMatOrderingType(), MatColoring, MatColoringCreate()
166: @*/
167: PetscErrorCode  MatGetOrdering(Mat mat,MatOrderingType type,IS *rperm,IS *cperm)
168: {
170:   PetscInt       mmat,nmat,mis,m;
171:   PetscErrorCode (*r)(Mat,MatOrderingType,IS*,IS*);
172:   PetscBool      flg = PETSC_FALSE,isseqdense,ismpidense,ismpiaij,ismpibaij,ismpisbaij,ismpiaijcusparse,iselemental,flg1;

178:   if (!mat->assembled) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
179:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

181:   PetscStrcmp(type,MATORDERINGNATURAL_OR_ND,&flg1);
182:   if (flg1) {
183:     PetscBool isseqsbaij;
184:     PetscObjectTypeCompareAny((PetscObject)mat,&isseqsbaij,MATSEQSBAIJ,MATSEQBAIJ,NULL);
185:     if (isseqsbaij) {
186:       type = MATORDERINGNATURAL;
187:     } else {
188:       type = MATORDERINGND;
189:     }
190:   }

192:   /* This code is terrible. MatGetOrdering() multiple dispatch should use matrix and this code should move to impls/aij/mpi. */
193:   PetscObjectTypeCompare((PetscObject)mat,MATMPIAIJ,&ismpiaij);
194:   if (ismpiaij) {               /* Reorder using diagonal block */
195:     Mat            Ad,Ao;
196:     const PetscInt *colmap;
197:     IS             lrowperm,lcolperm;
198:     PetscInt       i,rstart,rend,*idx;
199:     const PetscInt *lidx;

201:     MatMPIAIJGetSeqAIJ(mat,&Ad,&Ao,&colmap);
202:     MatGetOrdering(Ad,type,&lrowperm,&lcolperm);
203:     MatGetOwnershipRange(mat,&rstart,&rend);
204:     /* Remap row index set to global space */
205:     ISGetIndices(lrowperm,&lidx);
206:     PetscMalloc1(rend-rstart,&idx);
207:     for (i=0; i+rstart<rend; i++) idx[i] = rstart + lidx[i];
208:     ISRestoreIndices(lrowperm,&lidx);
209:     ISDestroy(&lrowperm);
210:     ISCreateGeneral(PetscObjectComm((PetscObject)mat),rend-rstart,idx,PETSC_OWN_POINTER,rperm);
211:     ISSetPermutation(*rperm);
212:     /* Remap column index set to global space */
213:     ISGetIndices(lcolperm,&lidx);
214:     PetscMalloc1(rend-rstart,&idx);
215:     for (i=0; i+rstart<rend; i++) idx[i] = rstart + lidx[i];
216:     ISRestoreIndices(lcolperm,&lidx);
217:     ISDestroy(&lcolperm);
218:     ISCreateGeneral(PetscObjectComm((PetscObject)mat),rend-rstart,idx,PETSC_OWN_POINTER,cperm);
219:     ISSetPermutation(*cperm);
220:     return(0);
221:   }

223:   /* this chunk of code is REALLY bad, should maybe get the ordering from the factor matrix,
224:      then those that don't support orderings will handle their cases themselves. */
225:   PetscObjectTypeCompare((PetscObject)mat,MATSEQDENSE,&isseqdense);
226:   PetscObjectTypeCompare((PetscObject)mat,MATMPIDENSE,&ismpidense);
227:   PetscObjectTypeCompare((PetscObject)mat,MATMPIAIJCUSPARSE,&ismpiaijcusparse);
228:   PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&ismpibaij);
229:   PetscObjectTypeCompare((PetscObject)mat,MATMPISBAIJ,&ismpisbaij);
230:   PetscObjectTypeCompare((PetscObject)mat,MATELEMENTAL,&iselemental);
231:   if (isseqdense || ismpidense || ismpibaij || ismpisbaij || ismpiaijcusparse || iselemental) {
232:     MatGetLocalSize(mat,&m,NULL);
233:     /*
234:        These matrices only give natural ordering
235:     */
236:     ISCreateStride(PETSC_COMM_SELF,m,0,1,cperm);
237:     ISCreateStride(PETSC_COMM_SELF,m,0,1,rperm);
238:     ISSetIdentity(*cperm);
239:     ISSetIdentity(*rperm);
240:     return(0);
241:   }

243:   if (!mat->rmap->N) { /* matrix has zero rows */
244:     ISCreateStride(PETSC_COMM_SELF,0,0,1,cperm);
245:     ISCreateStride(PETSC_COMM_SELF,0,0,1,rperm);
246:     ISSetIdentity(*cperm);
247:     ISSetIdentity(*rperm);
248:     return(0);
249:   }

251:   MatGetLocalSize(mat,&mmat,&nmat);
252:   if (mmat != nmat) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",mmat,nmat);

254:   MatOrderingRegisterAll();
255:   PetscFunctionListFind(MatOrderingList,type,&r);
256:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown or unregistered type: %s",type);

258:   PetscLogEventBegin(MAT_GetOrdering,mat,0,0,0);
259:   (*r)(mat,type,rperm,cperm);
260:   ISSetPermutation(*rperm);
261:   ISSetPermutation(*cperm);
262:   /* Adjust for inode (reduced matrix ordering) only if row permutation is smaller the matrix size */
263:   ISGetLocalSize(*rperm,&mis);
264:   if (mmat > mis) {MatInodeAdjustForInodes(mat,rperm,cperm);}
265:   PetscLogEventEnd(MAT_GetOrdering,mat,0,0,0);


268:   PetscOptionsHasName(((PetscObject)mat)->options,((PetscObject)mat)->prefix,"-mat_view_ordering",&flg);
269:   if (flg) {
270:     Mat tmat;
271:     MatPermute(mat,*rperm,*cperm,&tmat);
272:     MatViewFromOptions(tmat,(PetscObject)mat,"-mat_view_ordering");
273:     MatDestroy(&tmat);
274:   }
275:   return(0);
276: }

278: PetscErrorCode MatGetOrderingList(PetscFunctionList *list)
279: {
281:   *list = MatOrderingList;
282:   return(0);
283: }