Actual source code: inode2.c
petsc-3.7.3 2016-08-01
2: #include <../src/mat/impls/aij/seq/aij.h>
4: extern PetscErrorCode MatSeqAIJCheckInode(Mat);
5: extern PetscErrorCode MatInodeAdjustForInodes_SeqAIJ_Inode(Mat,IS*,IS*);
6: extern PetscErrorCode MatInodeGetInodeSizes_SeqAIJ_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);
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;
18: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
19: if (iascii) {
20: PetscViewerGetFormat(viewer,&format);
21: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
22: if (a->inode.size) {
23: PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n",
24: a->inode.node_count,a->inode.limit);
25: } else {
26: PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");
27: }
28: }
29: }
30: return(0);
31: }
35: PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
36: {
37: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
41: MatSeqAIJCheckInode(A);
42: a->inode.ibdiagvalid = PETSC_FALSE;
43: return(0);
44: }
48: PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
49: {
51: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
54: PetscFree(a->inode.size);
55: PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);
56: PetscObjectComposeFunction((PetscObject)A,"MatInodeAdjustForInodes_C",NULL);
57: PetscObjectComposeFunction((PetscObject)A,"MatInodeGetInodeSizes_C",NULL);
58: return(0);
59: }
61: /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type. */
62: /* It is also not registered as a type for use within MatSetType. */
63: /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should */
64: /* inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
65: /* Maybe this is a bad idea. (?) */
68: PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
69: {
70: Mat_SeqAIJ *b=(Mat_SeqAIJ*)B->data;
72: PetscBool no_inode,no_unroll;
75: no_inode = PETSC_FALSE;
76: no_unroll = PETSC_FALSE;
77: b->inode.node_count = 0;
78: b->inode.size = 0;
79: b->inode.limit = 5;
80: b->inode.max_limit = 5;
81: b->inode.ibdiagvalid = PETSC_FALSE;
82: b->inode.ibdiag = 0;
83: b->inode.bdiag = 0;
85: PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");
86: PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);
87: if (no_unroll) {
88: PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");
89: }
90: PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",NULL,no_inode,&no_inode,NULL);
91: if (no_inode) {
92: PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");
93: }
94: PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);
95: PetscOptionsEnd();
97: b->inode.use = (PetscBool)(!(no_unroll || no_inode));
98: if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
100: PetscObjectComposeFunction((PetscObject)B,"MatInodeAdjustForInodes_C",MatInodeAdjustForInodes_SeqAIJ_Inode);
101: PetscObjectComposeFunction((PetscObject)B,"MatInodeGetInodeSizes_C",MatInodeGetInodeSizes_SeqAIJ_Inode);
102: return(0);
103: }
107: PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool flg)
108: {
109: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
112: switch (op) {
113: case MAT_USE_INODES:
114: a->inode.use = flg;
115: break;
116: default:
117: break;
118: }
119: return(0);
120: }