xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision 7d0a6c19129e7069c8a40e210b34ed62989173db)
1 
2 #include "../src/mat/impls/aij/seq/aij.h"
3 
4 extern PetscErrorCode Mat_CheckInode(Mat,PetscBool );
5 EXTERN_C_BEGIN
6 extern PetscErrorCode  MatInodeAdjustForInodes_SeqAIJ_Inode(Mat,IS*,IS*);
7 extern PetscErrorCode  MatInodeGetInodeSizes_SeqAIJ_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);
8 EXTERN_C_END
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "MatView_SeqAIJ_Inode"
12 PetscErrorCode MatView_SeqAIJ_Inode(Mat A,PetscViewer viewer)
13 {
14   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data;
15   PetscErrorCode    ierr;
16   PetscBool         iascii;
17   PetscViewerFormat format;
18 
19   PetscFunctionBegin;
20   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&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_SeqAIJ_Inode"
37 PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
38 {
39   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
40   PetscErrorCode ierr;
41   PetscBool      samestructure;
42 
43   PetscFunctionBegin;
44   /* info.nz_unneeded of zero denotes no structural change was made to the matrix during Assembly */
45   samestructure = (PetscBool)(!A->info.nz_unneeded);
46   /* check for identical nodes. If found, use inode functions */
47   ierr = Mat_CheckInode(A,samestructure);CHKERRQ(ierr);
48   a->inode.ibdiagvalid = PETSC_FALSE;
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "MatDestroy_SeqAIJ_Inode"
54 PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
55 {
56   PetscErrorCode ierr;
57   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
58 
59   PetscFunctionBegin;
60   ierr = PetscFree(a->inode.size);CHKERRQ(ierr);
61   ierr = PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);CHKERRQ(ierr);
62   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeAdjustForInodes_C","",PETSC_NULL);CHKERRQ(ierr);
63   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeGetInodeSizes_C","",PETSC_NULL);CHKERRQ(ierr);
64   PetscFunctionReturn(0);
65 }
66 
67 /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
68 /* It is also not registered as a type for use within MatSetType.                             */
69 /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
70 /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
71 /* Maybe this is a bad idea. (?) */
72 #undef __FUNCT__
73 #define __FUNCT__ "MatCreate_SeqAIJ_Inode"
74 PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
75 {
76   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
77   PetscErrorCode ierr;
78   PetscBool      no_inode,no_unroll;
79 
80   PetscFunctionBegin;
81   no_inode             = PETSC_FALSE;
82   no_unroll            = PETSC_FALSE;
83   b->inode.node_count  = 0;
84   b->inode.size        = 0;
85   b->inode.limit       = 5;
86   b->inode.max_limit   = 5;
87   b->inode.ibdiagvalid = PETSC_FALSE;
88   b->inode.ibdiag      = 0;
89   b->inode.bdiag       = 0;
90 
91   ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr);
92     ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);CHKERRQ(ierr);
93     if (no_unroll) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);}
94     ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);CHKERRQ(ierr);
95     if (no_inode) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);}
96     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);
97   ierr = PetscOptionsEnd();CHKERRQ(ierr);
98   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
99   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
100 
101   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeAdjustForInodes_C",
102                                      "MatInodeAdjustForInodes_SeqAIJ_Inode",
103                                       MatInodeAdjustForInodes_SeqAIJ_Inode);CHKERRQ(ierr);
104   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeGetInodeSizes_C",
105                                      "MatInodeGetInodeSizes_SeqAIJ_Inode",
106                                       MatInodeGetInodeSizes_SeqAIJ_Inode);CHKERRQ(ierr);
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "MatSetOption_SeqAIJ_Inode"
112 PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool  flg)
113 {
114   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
115 
116   PetscFunctionBegin;
117   switch(op) {
118     case MAT_USE_INODES:
119       a->inode.use         = flg;
120       break;
121     default:
122       break;
123   }
124   PetscFunctionReturn(0);
125 }
126 
127 #undef __FUNCT__
128 #define __FUNCT__ "MatDuplicate_SeqAIJ_Inode"
129 PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat A,MatDuplicateOption cpvalues,Mat *C)
130 {
131   Mat            B=*C;
132   Mat_SeqAIJ     *c=(Mat_SeqAIJ*)B->data,*a=(Mat_SeqAIJ*)A->data;
133   PetscErrorCode ierr;
134   PetscInt       m=A->rmap->n;
135 
136   PetscFunctionBegin;
137 
138   c->inode.use          = a->inode.use;
139   c->inode.limit        = a->inode.limit;
140   c->inode.max_limit    = a->inode.max_limit;
141   if (a->inode.size){
142     ierr                = PetscMalloc((m+1)*sizeof(PetscInt),&c->inode.size);CHKERRQ(ierr);
143     c->inode.node_count = a->inode.node_count;
144     ierr                = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
145   } else {
146     c->inode.size       = 0;
147     c->inode.node_count = 0;
148   }
149   c->inode.ibdiagvalid = PETSC_FALSE;
150   c->inode.ibdiag      = 0;
151   c->inode.bdiag       = 0;
152   PetscFunctionReturn(0);
153 }
154 
155 
156 
157