1 #define PETSCMAT_DLL 2 #include "src/mat/impls/aij/seq/aij.h" 3 4 EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth); 5 EXTERN_C_BEGIN 6 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes_Inode(Mat,IS*,IS*); 7 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*); 8 EXTERN_C_END 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "MatView_Inode" 12 PetscErrorCode MatView_Inode(Mat A,PetscViewer viewer) 13 { 14 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 15 PetscErrorCode ierr; 16 PetscTruth iascii; 17 PetscViewerFormat format; 18 19 PetscFunctionBegin; 20 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 21 if (iascii) { 22 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 23 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) { 24 if (a->inode.size) { 25 ierr = PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n", 26 a->inode.node_count,a->inode.limit);CHKERRQ(ierr); 27 } else { 28 ierr = PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");CHKERRQ(ierr); 29 } 30 } 31 } 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "MatAssemblyEnd_Inode" 37 PetscErrorCode MatAssemblyEnd_Inode(Mat A, MatAssemblyType mode) 38 { 39 PetscErrorCode ierr; 40 PetscTruth samestructure; 41 42 PetscFunctionBegin; 43 /* info.nz_unneeded of zero denotes no structural change was made to the matrix during Assembly */ 44 samestructure = (PetscTruth)(!A->info.nz_unneeded); 45 /* check for identical nodes. If found, use inode functions */ 46 ierr = Mat_CheckInode(A,samestructure);CHKERRQ(ierr); 47 PetscFunctionReturn(0); 48 } 49 50 #undef __FUNCT__ 51 #define __FUNCT__ "MatDestroy_Inode" 52 PetscErrorCode MatDestroy_Inode(Mat A) 53 { 54 PetscErrorCode ierr; 55 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 56 57 PetscFunctionBegin; 58 ierr = PetscFree(a->inode.size);CHKERRQ(ierr); 59 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeAdjustForInodes_C","",PETSC_NULL);CHKERRQ(ierr); 60 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeGetInodeSizes_C","",PETSC_NULL);CHKERRQ(ierr); 61 PetscFunctionReturn(0); 62 } 63 64 /* MatCreate_Inode is not DLLEXPORTed because it is not a constructor for a complete type. */ 65 /* It is also not registered as a type for use within MatSetType. */ 66 /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should */ 67 /* inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */ 68 /* Maybe this is a bad idea. (?) */ 69 #undef __FUNCT__ 70 #define __FUNCT__ "MatCreate_Inode" 71 PetscErrorCode MatCreate_Inode(Mat B) 72 { 73 Mat_SeqAIJ *b=(Mat_SeqAIJ*)B->data; 74 PetscErrorCode ierr; 75 PetscTruth no_inode,no_unroll; 76 77 PetscFunctionBegin; 78 no_inode = PETSC_FALSE; 79 no_unroll = PETSC_FALSE; 80 b->inode.node_count = 0; 81 b->inode.size = 0; 82 b->inode.limit = 5; 83 b->inode.max_limit = 5; 84 85 ierr = PetscOptionsBegin(B->comm,B->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr); 86 ierr = PetscOptionsTruth("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);CHKERRQ(ierr); 87 if (no_unroll) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);} 88 ierr = PetscOptionsTruth("-mat_no_inode","Do not optimize for inodes (slower)",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);CHKERRQ(ierr); 89 if (no_inode) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);} 90 ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",PETSC_NULL,b->inode.limit,&b->inode.limit,PETSC_NULL);CHKERRQ(ierr); 91 ierr = PetscOptionsEnd();CHKERRQ(ierr); 92 b->inode.use = (PetscTruth)(!(no_unroll || no_inode)); 93 if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit; 94 95 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeAdjustForInodes_C", 96 "MatInodeAdjustForInodes_Inode", 97 MatInodeAdjustForInodes_Inode);CHKERRQ(ierr); 98 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeGetInodeSizes_C", 99 "MatInodeGetInodeSizes_Inode", 100 MatInodeGetInodeSizes_Inode);CHKERRQ(ierr); 101 PetscFunctionReturn(0); 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "MatSetOption_Inode" 106 PetscErrorCode MatSetOption_Inode(Mat A,MatOption op,PetscTruth flg) 107 { 108 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 109 110 PetscFunctionBegin; 111 switch(op) { 112 case MAT_USE_INODES: 113 a->inode.use = flg; 114 break; 115 case MAT_INODE_LIMIT_1: 116 a->inode.limit = 1; 117 break; 118 case MAT_INODE_LIMIT_2: 119 a->inode.limit = 2; 120 break; 121 case MAT_INODE_LIMIT_3: 122 a->inode.limit = 3; 123 break; 124 case MAT_INODE_LIMIT_4: 125 a->inode.limit = 4; 126 break; 127 case MAT_INODE_LIMIT_5: 128 a->inode.limit = 5; 129 break; 130 default: 131 break; 132 } 133 PetscFunctionReturn(0); 134 } 135 136 #undef __FUNCT__ 137 #define __FUNCT__ "MatDuplicate_Inode" 138 PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *C) 139 { 140 Mat B=*C; 141 Mat_SeqAIJ *c=(Mat_SeqAIJ*)B->data,*a=(Mat_SeqAIJ*)A->data; 142 PetscErrorCode ierr; 143 PetscInt m=A->rmap.n; 144 145 PetscFunctionBegin; 146 147 c->inode.use = a->inode.use; 148 c->inode.limit = a->inode.limit; 149 c->inode.max_limit = a->inode.max_limit; 150 if (a->inode.size){ 151 ierr = PetscMalloc((m+1)*sizeof(PetscInt),&c->inode.size);CHKERRQ(ierr); 152 c->inode.node_count = a->inode.node_count; 153 ierr = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 154 } else { 155 c->inode.size = 0; 156 c->inode.node_count = 0; 157 } 158 PetscFunctionReturn(0); 159 } 160 161 #undef __FUNCT__ 162 #define __FUNCT__ "MatILUDTFactor_Inode" 163 PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 164 { 165 PetscErrorCode ierr; 166 167 PetscFunctionBegin; 168 /* check for identical nodes. If found, use inode functions */ 169 ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 170 PetscFunctionReturn(0); 171 } 172 173 #undef __FUNCT__ 174 #define __FUNCT__ "MatLUFactorSymbolic_Inode" 175 PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 176 { 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 /* check for identical nodes. If found, use inode functions */ 181 ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 182 PetscFunctionReturn(0); 183 } 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "MatILUFactorSymbolic_Inode" 187 PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 188 { 189 PetscErrorCode ierr; 190 191 PetscFunctionBegin; 192 /* check for identical nodes. If found, use inode functions */ 193 ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 198