1 2 #include <../src/mat/impls/aij/seq/aij.h> 3 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*); 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatView_SeqAIJ_Inode" 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; 16 17 PetscFunctionBegin; 18 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 19 if (iascii) { 20 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 21 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) { 22 if (a->inode.size) { 23 ierr = PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n", 24 a->inode.node_count,a->inode.limit);CHKERRQ(ierr); 25 } else { 26 ierr = PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");CHKERRQ(ierr); 27 } 28 } 29 } 30 PetscFunctionReturn(0); 31 } 32 33 #undef __FUNCT__ 34 #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_Inode" 35 PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode) 36 { 37 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 38 PetscErrorCode ierr; 39 PetscBool samestructure; 40 41 PetscFunctionBegin; 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 ierr = Mat_CheckInode(A,samestructure);CHKERRQ(ierr); 46 47 a->inode.ibdiagvalid = PETSC_FALSE; 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "MatDestroy_SeqAIJ_Inode" 53 PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A) 54 { 55 PetscErrorCode ierr; 56 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 57 58 PetscFunctionBegin; 59 ierr = PetscFree(a->inode.size);CHKERRQ(ierr); 60 ierr = PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);CHKERRQ(ierr); 61 ierr = PetscObjectComposeFunction((PetscObject)A,"MatInodeAdjustForInodes_C","",NULL);CHKERRQ(ierr); 62 ierr = PetscObjectComposeFunction((PetscObject)A,"MatInodeGetInodeSizes_C","",NULL);CHKERRQ(ierr); 63 PetscFunctionReturn(0); 64 } 65 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. (?) */ 71 #undef __FUNCT__ 72 #define __FUNCT__ "MatCreate_SeqAIJ_Inode" 73 PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B) 74 { 75 Mat_SeqAIJ *b=(Mat_SeqAIJ*)B->data; 76 PetscErrorCode ierr; 77 PetscBool no_inode,no_unroll; 78 79 PetscFunctionBegin; 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; 89 90 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr); 91 ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);CHKERRQ(ierr); 92 if (no_unroll) { 93 ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr); 94 } 95 ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",NULL,no_inode,&no_inode,NULL);CHKERRQ(ierr); 96 if (no_inode) { 97 ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr); 98 } 99 ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);CHKERRQ(ierr); 100 ierr = PetscOptionsEnd();CHKERRQ(ierr); 101 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; 104 105 ierr = PetscObjectComposeFunction((PetscObject)B,"MatInodeAdjustForInodes_C","MatInodeAdjustForInodes_SeqAIJ_Inode",MatInodeAdjustForInodes_SeqAIJ_Inode);CHKERRQ(ierr); 106 ierr = PetscObjectComposeFunction((PetscObject)B,"MatInodeGetInodeSizes_C","MatInodeGetInodeSizes_SeqAIJ_Inode",MatInodeGetInodeSizes_SeqAIJ_Inode);CHKERRQ(ierr); 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "MatSetOption_SeqAIJ_Inode" 112 PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool flg) 113 { 114 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 115 116 PetscFunctionBegin; 117 switch (op) { 118 case MAT_USE_INODES: 119 a->inode.use = flg; 120 break; 121 default: 122 break; 123 } 124 PetscFunctionReturn(0); 125 } 126 127 128 129 130 131