xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
17d0a6c19SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
34c1414c8SBarry Smith 
47087cfbeSBarry Smith extern PetscErrorCode MatInodeAdjustForInodes_SeqAIJ_Inode(Mat, IS *, IS *);
57087cfbeSBarry Smith extern PetscErrorCode MatInodeGetInodeSizes_SeqAIJ_Inode(Mat, PetscInt *, PetscInt *[], PetscInt *);
64c1414c8SBarry Smith 
79371c9d4SSatish Balay PetscErrorCode MatView_SeqAIJ_Inode(Mat A, PetscViewer viewer) {
84c1414c8SBarry Smith   Mat_SeqAIJ       *a = (Mat_SeqAIJ *)A->data;
9ace3abfcSBarry Smith   PetscBool         iascii;
104c1414c8SBarry Smith   PetscViewerFormat format;
114c1414c8SBarry Smith 
124c1414c8SBarry Smith   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
144c1414c8SBarry Smith   if (iascii) {
159566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
164c1414c8SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
174c1414c8SBarry Smith       if (a->inode.size) {
189566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "using I-node routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n", a->inode.node_count, a->inode.limit));
194c1414c8SBarry Smith       } else {
209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "not using I-node routines\n"));
214c1414c8SBarry Smith       }
224c1414c8SBarry Smith     }
234c1414c8SBarry Smith   }
244c1414c8SBarry Smith   PetscFunctionReturn(0);
254c1414c8SBarry Smith }
264c1414c8SBarry Smith 
279371c9d4SSatish Balay PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode) {
2871f1c65dSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
294c1414c8SBarry Smith 
304c1414c8SBarry Smith   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJCheckInode(A));
3271f1c65dSBarry Smith   a->inode.ibdiagvalid = PETSC_FALSE;
334c1414c8SBarry Smith   PetscFunctionReturn(0);
344c1414c8SBarry Smith }
354c1414c8SBarry Smith 
369371c9d4SSatish Balay PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A) {
374c1414c8SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
384c1414c8SBarry Smith 
394c1414c8SBarry Smith   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->inode.size));
419566063dSJacob Faibussowitsch   PetscCall(PetscFree3(a->inode.ibdiag, a->inode.bdiag, a->inode.ssor_work));
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatInodeAdjustForInodes_C", NULL));
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatInodeGetInodeSizes_C", NULL));
444c1414c8SBarry Smith   PetscFunctionReturn(0);
454c1414c8SBarry Smith }
464c1414c8SBarry Smith 
474108e4d5SBarry Smith /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
484c1414c8SBarry Smith /* It is also not registered as a type for use within MatSetType.                             */
494c1414c8SBarry Smith /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
504c1414c8SBarry Smith /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
514c1414c8SBarry Smith /* Maybe this is a bad idea. (?) */
529371c9d4SSatish Balay PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B) {
534c1414c8SBarry Smith   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
54ace3abfcSBarry Smith   PetscBool   no_inode, no_unroll;
554c1414c8SBarry Smith 
564c1414c8SBarry Smith   PetscFunctionBegin;
57da9f3051SSatish Balay   no_inode             = PETSC_FALSE;
58da9f3051SSatish Balay   no_unroll            = PETSC_FALSE;
59ec710b6aSStefano Zampini   b->inode.checked     = PETSC_FALSE;
604c1414c8SBarry Smith   b->inode.node_count  = 0;
61f4259b30SLisandro Dalcin   b->inode.size        = NULL;
624c1414c8SBarry Smith   b->inode.limit       = 5;
634c1414c8SBarry Smith   b->inode.max_limit   = 5;
6471f1c65dSBarry Smith   b->inode.ibdiagvalid = PETSC_FALSE;
65f4259b30SLisandro Dalcin   b->inode.ibdiag      = NULL;
66f4259b30SLisandro Dalcin   b->inode.bdiag       = NULL;
674c1414c8SBarry Smith 
68d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "Options for SEQAIJ matrix", "Mat");
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_unroll", "Do not optimize for inodes (slower)", NULL, no_unroll, &no_unroll, NULL));
70*48a46eb9SPierre Jolivet   if (no_unroll) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_unroll\n"));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_inode", "Do not optimize for inodes -slower-", NULL, no_inode, &no_inode, NULL));
72*48a46eb9SPierre Jolivet   if (no_inode) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_inode\n"));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_inode_limit", "Do not use inodes larger then this value", NULL, b->inode.limit, &b->inode.limit, NULL));
74d0609cedSBarry Smith   PetscOptionsEnd();
7526fbe8dcSKarl Rupp 
76ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
77d74fe821SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
78d74fe821SBarry Smith 
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatInodeAdjustForInodes_C", MatInodeAdjustForInodes_SeqAIJ_Inode));
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatInodeGetInodeSizes_C", MatInodeGetInodeSizes_SeqAIJ_Inode));
814c1414c8SBarry Smith   PetscFunctionReturn(0);
824c1414c8SBarry Smith }
834c1414c8SBarry Smith 
849371c9d4SSatish Balay PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A, MatOption op, PetscBool flg) {
854c1414c8SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
86d74fe821SBarry Smith 
874c1414c8SBarry Smith   PetscFunctionBegin;
884c1414c8SBarry Smith   switch (op) {
899371c9d4SSatish Balay   case MAT_USE_INODES: a->inode.use = flg; break;
909371c9d4SSatish Balay   default: break;
914c1414c8SBarry Smith   }
924c1414c8SBarry Smith   PetscFunctionReturn(0);
934c1414c8SBarry Smith }
94