Actual source code: inode2.c

  1: #define PETSCMAT_DLL
 2:  #include src/mat/impls/aij/seq/aij.h

  4: EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth);
  6: EXTERN PetscErrorCode  MatInodeAdjustForInodes_Inode(Mat,IS*,IS*);
  7: EXTERN PetscErrorCode  MatInodeGetInodeSizes_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);

 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;

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

 37: PetscErrorCode MatAssemblyEnd_Inode(Mat A, MatAssemblyType mode)
 38: {
 40:   PetscTruth     samestructure;

 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:   Mat_CheckInode(A,samestructure);
 47:   return(0);
 48: }

 52: PetscErrorCode MatDestroy_Inode(Mat A)
 53: {
 55:   Mat_SeqAIJ      *a=(Mat_SeqAIJ*)A->data;

 58:   PetscFree(a->inode.size);
 59:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeAdjustForInodes_C","",PETSC_NULL);
 60:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeGetInodeSizes_C","",PETSC_NULL);
 61:   return(0);
 62: }

 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. (?) */
 71: PetscErrorCode MatCreate_Inode(Mat B)
 72: {
 73:   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
 75:   PetscTruth no_inode,no_unroll;

 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;

 85:   PetscOptionsBegin(B->comm,B->prefix,"Options for SEQAIJ matrix","Mat");
 86:     PetscOptionsTruth("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);
 87:     if (no_unroll) {PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");}
 88:     PetscOptionsTruth("-mat_no_inode","Do not optimize for inodes (slower)",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);
 89:     if (no_inode) {PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");}
 90:     PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",PETSC_NULL,b->inode.limit,&b->inode.limit,PETSC_NULL);
 91:   PetscOptionsEnd();
 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;

 95:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeAdjustForInodes_C",
 96:                                      "MatInodeAdjustForInodes_Inode",
 97:                                       MatInodeAdjustForInodes_Inode);
 98:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeGetInodeSizes_C",
 99:                                      "MatInodeGetInodeSizes_Inode",
100:                                       MatInodeGetInodeSizes_Inode);
101:   return(0);
102: }

106: PetscErrorCode MatSetOption_Inode(Mat A,MatOption op)
107: {
108:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;

112:   switch(op) {
113:     case MAT_USE_INODES:
114:       a->inode.use         = PETSC_TRUE;
115:       break;
116:     case MAT_DO_NOT_USE_INODES:
117:       a->inode.use         = PETSC_FALSE;
118:       PetscInfo(A,"Not using Inode routines due to MatSetOption(MAT_DO_NOT_USE_INODES\n");
119:       break;
120:     case MAT_INODE_LIMIT_1:
121:       a->inode.limit  = 1;
122:       break;
123:     case MAT_INODE_LIMIT_2:
124:       a->inode.limit  = 2;
125:       break;
126:     case MAT_INODE_LIMIT_3:
127:       a->inode.limit  = 3;
128:       break;
129:     case MAT_INODE_LIMIT_4:
130:       a->inode.limit  = 4;
131:       break;
132:     case MAT_INODE_LIMIT_5:
133:       a->inode.limit  = 5;
134:       break;
135:     default:
136:       break;
137:   }
138:   return(0);
139: }

143: PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *C)
144: {
145:   Mat            B=*C;
146:   Mat_SeqAIJ      *c=(Mat_SeqAIJ*)B->data,*a=(Mat_SeqAIJ*)A->data;
148:   PetscInt       m=A->rmap.n;


152:   c->inode.use          = a->inode.use;
153:   c->inode.limit        = a->inode.limit;
154:   c->inode.max_limit    = a->inode.max_limit;
155:   if (a->inode.size){
156:     PetscMalloc((m+1)*sizeof(PetscInt),&c->inode.size);
157:     c->inode.node_count = a->inode.node_count;
158:     PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(PetscInt));
159:   } else {
160:     c->inode.size       = 0;
161:     c->inode.node_count = 0;
162:   }
163:   return(0);
164: }

168: PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
169: {

173:     /* check for identical nodes. If found, use inode functions */
174:   Mat_CheckInode(*fact,PETSC_FALSE);
175:   return(0);
176: }

180: PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
181: {

185:     /* check for identical nodes. If found, use inode functions */
186:   Mat_CheckInode(*fact,PETSC_FALSE);
187:   return(0);
188: }

192: PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact)
193: {

197:     /* check for identical nodes. If found, use inode functions */
198:   Mat_CheckInode(*fact,PETSC_FALSE);
199:   return(0);
200: }