Actual source code: inode2.c

  1: #include <../src/mat/impls/aij/seq/aij.h>

  3: extern PetscErrorCode MatInodeAdjustForInodes_SeqAIJ_Inode(Mat, IS *, IS *);
  4: extern PetscErrorCode MatInodeGetInodeSizes_SeqAIJ_Inode(Mat, PetscInt *, PetscInt *[], PetscInt *);

  6: PetscErrorCode MatView_SeqAIJ_Inode(Mat A, PetscViewer viewer)
  7: {
  8:   Mat_SeqAIJ       *a = (Mat_SeqAIJ *)A->data;
  9:   PetscBool         iascii;
 10:   PetscViewerFormat format;

 12:   PetscFunctionBegin;
 13:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
 14:   if (iascii) {
 15:     PetscCall(PetscViewerGetFormat(viewer, &format));
 16:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
 17:       if (a->inode.size) {
 18:         PetscCall(PetscViewerASCIIPrintf(viewer, "using I-node routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n", a->inode.node_count, a->inode.limit));
 19:       } else {
 20:         PetscCall(PetscViewerASCIIPrintf(viewer, "not using I-node routines\n"));
 21:       }
 22:     }
 23:   }
 24:   PetscFunctionReturn(PETSC_SUCCESS);
 25: }

 27: PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
 28: {
 29:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

 31:   PetscFunctionBegin;
 32:   PetscCall(MatSeqAIJCheckInode(A));
 33:   a->inode.ibdiagvalid = PETSC_FALSE;
 34:   PetscFunctionReturn(PETSC_SUCCESS);
 35: }

 37: PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
 38: {
 39:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

 41:   PetscFunctionBegin;
 42:   PetscCall(PetscFree(a->inode.size));
 43:   PetscCall(PetscFree3(a->inode.ibdiag, a->inode.bdiag, a->inode.ssor_work));
 44:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatInodeAdjustForInodes_C", NULL));
 45:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatInodeGetInodeSizes_C", NULL));
 46:   PetscFunctionReturn(PETSC_SUCCESS);
 47: }

 49: /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
 50: /* It is also not registered as a type for use within MatSetType.                             */
 51: /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
 52: /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
 53: /* Maybe this is a bad idea. (?) */
 54: PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
 55: {
 56:   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
 57:   PetscBool   no_inode, no_unroll;

 59:   PetscFunctionBegin;
 60:   no_inode             = PETSC_FALSE;
 61:   no_unroll            = PETSC_FALSE;
 62:   b->inode.checked     = PETSC_FALSE;
 63:   b->inode.node_count  = 0;
 64:   b->inode.size        = NULL;
 65:   b->inode.limit       = 5;
 66:   b->inode.max_limit   = 5;
 67:   b->inode.ibdiagvalid = PETSC_FALSE;
 68:   b->inode.ibdiag      = NULL;
 69:   b->inode.bdiag       = NULL;

 71:   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "Options for SEQAIJ matrix", "Mat");
 72:   PetscCall(PetscOptionsBool("-mat_no_unroll", "Do not optimize for inodes (slower)", NULL, no_unroll, &no_unroll, NULL));
 73:   if (no_unroll) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_unroll\n"));
 74:   PetscCall(PetscOptionsBool("-mat_no_inode", "Do not optimize for inodes -slower-", NULL, no_inode, &no_inode, NULL));
 75:   if (no_inode) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_inode\n"));
 76:   PetscCall(PetscOptionsInt("-mat_inode_limit", "Do not use inodes larger then this value", NULL, b->inode.limit, &b->inode.limit, NULL));
 77:   PetscOptionsEnd();

 79:   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
 80:   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;

 82:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatInodeAdjustForInodes_C", MatInodeAdjustForInodes_SeqAIJ_Inode));
 83:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatInodeGetInodeSizes_C", MatInodeGetInodeSizes_SeqAIJ_Inode));
 84:   PetscFunctionReturn(PETSC_SUCCESS);
 85: }

 87: PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A, MatOption op, PetscBool flg)
 88: {
 89:   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;

 91:   PetscFunctionBegin;
 92:   switch (op) {
 93:   case MAT_USE_INODES:
 94:     a->inode.use = flg;
 95:     break;
 96:   default:
 97:     break;
 98:   }
 99:   PetscFunctionReturn(PETSC_SUCCESS);
100: }