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