xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision 4fc747eaadbeca11629f314a99edccbc2ed7b3d3) !
1 
2 #include <../src/mat/impls/aij/seq/aij.h>
3 
4 extern PetscErrorCode MatSeqAIJCheckInode(Mat);
5 extern PetscErrorCode MatInodeAdjustForInodes_SeqAIJ_Inode(Mat,IS*,IS*);
6 extern PetscErrorCode MatInodeGetInodeSizes_SeqAIJ_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "MatView_SeqAIJ_Inode"
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;
16 
17   PetscFunctionBegin;
18   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
19   if (iascii) {
20     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
21     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
22       if (a->inode.size) {
23         ierr = PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n",a->inode.node_count,a->inode.limit);CHKERRQ(ierr);
24       } else {
25         ierr = PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");CHKERRQ(ierr);
26       }
27     }
28   }
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_Inode"
34 PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
35 {
36   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
37   PetscErrorCode ierr;
38 
39   PetscFunctionBegin;
40   ierr = MatSeqAIJCheckInode(A);CHKERRQ(ierr);
41   a->inode.ibdiagvalid = PETSC_FALSE;
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "MatDestroy_SeqAIJ_Inode"
47 PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
48 {
49   PetscErrorCode ierr;
50   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
51 
52   PetscFunctionBegin;
53   ierr = PetscFree(a->inode.size);CHKERRQ(ierr);
54   ierr = PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);CHKERRQ(ierr);
55   ierr = PetscObjectComposeFunction((PetscObject)A,"MatInodeAdjustForInodes_C",NULL);CHKERRQ(ierr);
56   ierr = PetscObjectComposeFunction((PetscObject)A,"MatInodeGetInodeSizes_C",NULL);CHKERRQ(ierr);
57   PetscFunctionReturn(0);
58 }
59 
60 /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
61 /* It is also not registered as a type for use within MatSetType.                             */
62 /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
63 /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
64 /* Maybe this is a bad idea. (?) */
65 #undef __FUNCT__
66 #define __FUNCT__ "MatCreate_SeqAIJ_Inode"
67 PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
68 {
69   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
70   PetscErrorCode ierr;
71   PetscBool      no_inode,no_unroll;
72 
73   PetscFunctionBegin;
74   no_inode             = PETSC_FALSE;
75   no_unroll            = PETSC_FALSE;
76   b->inode.node_count  = 0;
77   b->inode.size        = 0;
78   b->inode.limit       = 5;
79   b->inode.max_limit   = 5;
80   b->inode.ibdiagvalid = PETSC_FALSE;
81   b->inode.ibdiag      = 0;
82   b->inode.bdiag       = 0;
83 
84   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr);
85   ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);CHKERRQ(ierr);
86   if (no_unroll) {
87     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);
88   }
89   ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",NULL,no_inode,&no_inode,NULL);CHKERRQ(ierr);
90   if (no_inode) {
91     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);
92   }
93   ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);CHKERRQ(ierr);
94   ierr = PetscOptionsEnd();CHKERRQ(ierr);
95 
96   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
97   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
98 
99   ierr = PetscObjectComposeFunction((PetscObject)B,"MatInodeAdjustForInodes_C",MatInodeAdjustForInodes_SeqAIJ_Inode);CHKERRQ(ierr);
100   ierr = PetscObjectComposeFunction((PetscObject)B,"MatInodeGetInodeSizes_C",MatInodeGetInodeSizes_SeqAIJ_Inode);CHKERRQ(ierr);
101   PetscFunctionReturn(0);
102 }
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "MatSetOption_SeqAIJ_Inode"
106 PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool 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   default:
116     break;
117   }
118   PetscFunctionReturn(0);
119 }
120 
121 
122 
123 
124 
125