Actual source code: sregis.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/matimpl.h>     /*I       "petscmat.h"   I*/

  4: PETSC_EXTERN PetscErrorCode MatGetOrdering_Natural(Mat,MatOrderingType,IS*,IS*);
  5: PETSC_EXTERN PetscErrorCode MatGetOrdering_ND(Mat,MatOrderingType,IS*,IS*);
  6: PETSC_EXTERN PetscErrorCode MatGetOrdering_1WD(Mat,MatOrderingType,IS*,IS*);
  7: PETSC_EXTERN PetscErrorCode MatGetOrdering_QMD(Mat,MatOrderingType,IS*,IS*);
  8: PETSC_EXTERN PetscErrorCode MatGetOrdering_RCM(Mat,MatOrderingType,IS*,IS*);
  9: PETSC_EXTERN PetscErrorCode MatGetOrdering_RowLength(Mat,MatOrderingType,IS*,IS*);
 10: PETSC_EXTERN PetscErrorCode MatGetOrdering_DSC(Mat,MatOrderingType,IS*,IS*);
 11: PETSC_EXTERN PetscErrorCode MatGetOrdering_WBM(Mat,MatOrderingType,IS*,IS*);
 12: PETSC_EXTERN PetscErrorCode MatGetOrdering_Spectral(Mat,MatOrderingType,IS*,IS*);
 13: #if defined(PETSC_HAVE_SUITESPARSE)
 14: PETSC_EXTERN PetscErrorCode MatGetOrdering_AMD(Mat,MatOrderingType,IS*,IS*);
 15: #endif

 19: /*@C
 20:   MatOrderingRegisterAll - Registers all of the matrix
 21:   reordering routines in PETSc.

 23:   Not Collective

 25:   Level: developer

 27:   Adding new methods:
 28:   To add a new method to the registry. Copy this routine and
 29:   modify it to incorporate a call to MatReorderRegister() for
 30:   the new method, after the current list.

 32:   Restricting the choices: To prevent all of the methods from being
 33:   registered and thus save memory, copy this routine and comment out
 34:   those orderigs you do not wish to include.  Make sure that the
 35:   replacement routine is linked before libpetscmat.a.

 37: .keywords: matrix, reordering, register, all

 39: .seealso: MatOrderingRegister(), MatOrderingRegisterDestroy()
 40: @*/
 41: PetscErrorCode  MatOrderingRegisterAll(void)
 42: {

 46:   MatOrderingRegisterAllCalled = PETSC_TRUE;

 48:   MatOrderingRegister(MATORDERINGNATURAL,  MatGetOrdering_Natural);
 49:   MatOrderingRegister(MATORDERINGND,       MatGetOrdering_ND);
 50:   MatOrderingRegister(MATORDERING1WD,      MatGetOrdering_1WD);
 51:   MatOrderingRegister(MATORDERINGRCM,      MatGetOrdering_RCM);
 52:   MatOrderingRegister(MATORDERINGQMD,      MatGetOrdering_QMD);
 53:   MatOrderingRegister(MATORDERINGROWLENGTH,MatGetOrdering_RowLength);
 54: #if defined(PETSC_HAVE_SUPERLU_DIST)
 55:   MatOrderingRegister(MATORDERINGWBM,      MatGetOrdering_WBM);
 56: #endif
 57:   MatOrderingRegister(MATORDERINGSPECTRAL, MatGetOrdering_Spectral);
 58: #if defined(PETSC_HAVE_SUITESPARSE)
 59:   MatOrderingRegister(MATORDERINGAMD,      MatGetOrdering_AMD);
 60: #endif
 61:   return(0);
 62: }