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 76 PetscFunctionBegin; 77 b->inode.use = PETSC_TRUE; 78 b->inode.node_count = 0; 79 b->inode.size = 0; 80 b->inode.limit = 5; 81 b->inode.max_limit = 5; 82 83 ierr = PetscOptionsBegin(B->comm,B->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr); 84 ierr = PetscOptionsTruth("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,b->inode.use,&b->inode.use,PETSC_NULL);CHKERRQ(ierr); 85 if (!b->inode.use) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);} 86 ierr = PetscOptionsTruth("-mat_no_inode","Do not optimize for inodes (slower)",PETSC_NULL,b->inode.use,&b->inode.use,PETSC_NULL);CHKERRQ(ierr); 87 if (!b->inode.use) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);} 88 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); 89 ierr = PetscOptionsEnd();CHKERRQ(ierr); 90 if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit; 91 92 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeAdjustForInodes_C", 93 "MatInodeAdjustForInodes_Inode", 94 MatInodeAdjustForInodes_Inode);CHKERRQ(ierr); 95 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeGetInodeSizes_C", 96 "MatInodeGetInodeSizes_Inode", 97 MatInodeGetInodeSizes_Inode);CHKERRQ(ierr); 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "MatSetOption_Inode" 103 PetscErrorCode MatSetOption_Inode(Mat A,MatOption op) 104 { 105 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 106 PetscErrorCode ierr; 107 108 PetscFunctionBegin; 109 switch(op) { 110 case MAT_USE_INODES: 111 a->inode.use = PETSC_TRUE; 112 break; 113 case MAT_DO_NOT_USE_INODES: 114 a->inode.use = PETSC_FALSE; 115 ierr = PetscInfo(A,"Not using Inode routines due to MatSetOption(MAT_DO_NOT_USE_INODES\n");CHKERRQ(ierr); 116 break; 117 case MAT_INODE_LIMIT_1: 118 a->inode.limit = 1; 119 break; 120 case MAT_INODE_LIMIT_2: 121 a->inode.limit = 2; 122 break; 123 case MAT_INODE_LIMIT_3: 124 a->inode.limit = 3; 125 break; 126 case MAT_INODE_LIMIT_4: 127 a->inode.limit = 4; 128 break; 129 case MAT_INODE_LIMIT_5: 130 a->inode.limit = 5; 131 break; 132 default: 133 break; 134 } 135 PetscFunctionReturn(0); 136 } 137 138 #undef __FUNCT__ 139 #define __FUNCT__ "MatDuplicate_Inode" 140 PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *C) 141 { 142 Mat B=*C; 143 Mat_SeqAIJ *c=(Mat_SeqAIJ*)B->data,*a=(Mat_SeqAIJ*)A->data; 144 PetscErrorCode ierr; 145 PetscInt m=A->rmap.n; 146 147 PetscFunctionBegin; 148 149 c->inode.use = a->inode.use; 150 c->inode.limit = a->inode.limit; 151 c->inode.max_limit = a->inode.max_limit; 152 if (a->inode.size){ 153 ierr = PetscMalloc((m+1)*sizeof(PetscInt),&c->inode.size);CHKERRQ(ierr); 154 c->inode.node_count = a->inode.node_count; 155 ierr = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 156 } else { 157 c->inode.size = 0; 158 c->inode.node_count = 0; 159 } 160 PetscFunctionReturn(0); 161 } 162 163 #undef __FUNCT__ 164 #define __FUNCT__ "MatILUDTFactor_Inode" 165 PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 166 { 167 PetscErrorCode ierr; 168 169 PetscFunctionBegin; 170 /* check for identical nodes. If found, use inode functions */ 171 ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 172 PetscFunctionReturn(0); 173 } 174 175 #undef __FUNCT__ 176 #define __FUNCT__ "MatLUFactorSymbolic_Inode" 177 PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 178 { 179 PetscErrorCode ierr; 180 181 PetscFunctionBegin; 182 /* check for identical nodes. If found, use inode functions */ 183 ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 #undef __FUNCT__ 188 #define __FUNCT__ "MatILUFactorSymbolic_Inode" 189 PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact) 190 { 191 PetscErrorCode ierr; 192 193 PetscFunctionBegin; 194 /* check for identical nodes. If found, use inode functions */ 195 ierr = Mat_CheckInode(*fact,PETSC_FALSE);CHKERRQ(ierr); 196 PetscFunctionReturn(0); 197 } 198 199 200