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