Actual source code: inode2.c

petsc-3.4.2 2013-07-02
  2: #include <../src/mat/impls/aij/seq/aij.h>

  4: extern PetscErrorCode Mat_CheckInode(Mat,PetscBool);
  5: extern PetscErrorCode  MatInodeAdjustForInodes_SeqAIJ_Inode(Mat,IS*,IS*);
  6: extern PetscErrorCode  MatInodeGetInodeSizes_SeqAIJ_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);

 10: PetscErrorCode MatView_SeqAIJ_Inode(Mat A,PetscViewer viewer)
 11: {
 12:   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
 13:   PetscErrorCode    ierr;
 14:   PetscBool         iascii;
 15:   PetscViewerFormat format;

 18:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
 19:   if (iascii) {
 20:     PetscViewerGetFormat(viewer,&format);
 21:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
 22:       if (a->inode.size) {
 23:         PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n",
 24:                                       a->inode.node_count,a->inode.limit);
 25:       } else {
 26:         PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");
 27:       }
 28:     }
 29:   }
 30:   return(0);
 31: }

 35: PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
 36: {
 37:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
 39:   PetscBool      samestructure;

 42:   /* info.nz_unneeded of zero denotes no structural change was made to the matrix during Assembly */
 43:   samestructure = (PetscBool)(!A->info.nz_unneeded);
 44:   /* check for identical nodes. If found, use inode functions */
 45:   Mat_CheckInode(A,samestructure);

 47:   a->inode.ibdiagvalid = PETSC_FALSE;
 48:   return(0);
 49: }

 53: PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
 54: {
 56:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;

 59:   PetscFree(a->inode.size);
 60:   PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);
 61:   PetscObjectComposeFunction((PetscObject)A,"MatInodeAdjustForInodes_C",NULL);
 62:   PetscObjectComposeFunction((PetscObject)A,"MatInodeGetInodeSizes_C",NULL);
 63:   return(0);
 64: }

 66: /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
 67: /* It is also not registered as a type for use within MatSetType.                             */
 68: /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
 69: /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
 70: /* Maybe this is a bad idea. (?) */
 73: PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
 74: {
 75:   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
 77:   PetscBool      no_inode,no_unroll;

 80:   no_inode             = PETSC_FALSE;
 81:   no_unroll            = PETSC_FALSE;
 82:   b->inode.node_count  = 0;
 83:   b->inode.size        = 0;
 84:   b->inode.limit       = 5;
 85:   b->inode.max_limit   = 5;
 86:   b->inode.ibdiagvalid = PETSC_FALSE;
 87:   b->inode.ibdiag      = 0;
 88:   b->inode.bdiag       = 0;

 90:   PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");
 91:   PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
 92:   if (no_unroll) {
 93:     PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
 94:   }
 95:   PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",NULL,no_inode,&no_inode,NULL);
 96:   if (no_inode) {
 97:     PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
 98:   }
 99:   PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
100:   PetscOptionsEnd();

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

105:   PetscObjectComposeFunction((PetscObject)B,"MatInodeAdjustForInodes_C",MatInodeAdjustForInodes_SeqAIJ_Inode);
106:   PetscObjectComposeFunction((PetscObject)B,"MatInodeGetInodeSizes_C",MatInodeGetInodeSizes_SeqAIJ_Inode);
107:   return(0);
108: }

112: PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool flg)
113: {
114:   Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;

117:   switch (op) {
118:   case MAT_USE_INODES:
119:     a->inode.use = flg;
120:     break;
121:   default:
122:     break;
123:   }
124:   return(0);
125: }