xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision e87b5d96567786df49b4af65cff3120eb2cdf5c2)
1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
24c1414c8SBarry Smith 
37087cfbeSBarry Smith extern PetscErrorCode MatInodeAdjustForInodes_SeqAIJ_Inode(Mat, IS *, IS *);
47087cfbeSBarry Smith extern PetscErrorCode MatInodeGetInodeSizes_SeqAIJ_Inode(Mat, PetscInt *, PetscInt *[], PetscInt *);
54c1414c8SBarry Smith 
6d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqAIJ_Inode(Mat A, PetscViewer viewer)
7d71ae5a4SJacob Faibussowitsch {
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) {
174d12350bSJunchao Zhang       if (a->inode.size_csr) {
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   }
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254c1414c8SBarry Smith }
264c1414c8SBarry Smith 
27d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
28d71ae5a4SJacob Faibussowitsch {
2971f1c65dSBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
304c1414c8SBarry Smith 
314c1414c8SBarry Smith   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJCheckInode(A));
3371f1c65dSBarry Smith   a->inode.ibdiagvalid = PETSC_FALSE;
343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354c1414c8SBarry Smith }
364c1414c8SBarry Smith 
37d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
38d71ae5a4SJacob Faibussowitsch {
394c1414c8SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
404c1414c8SBarry Smith 
414c1414c8SBarry Smith   PetscFunctionBegin;
424d12350bSJunchao Zhang   PetscCall(PetscFree(a->inode.size_csr));
439566063dSJacob Faibussowitsch   PetscCall(PetscFree3(a->inode.ibdiag, a->inode.bdiag, a->inode.ssor_work));
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatInodeAdjustForInodes_C", NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatInodeGetInodeSizes_C", NULL));
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
474c1414c8SBarry Smith }
484c1414c8SBarry Smith 
494108e4d5SBarry Smith /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
504c1414c8SBarry Smith /* It is also not registered as a type for use within MatSetType.                             */
514c1414c8SBarry Smith /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
524c1414c8SBarry Smith /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
534c1414c8SBarry Smith /* Maybe this is a bad idea. (?) */
54d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
55d71ae5a4SJacob Faibussowitsch {
564c1414c8SBarry Smith   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
57ace3abfcSBarry Smith   PetscBool   no_inode, no_unroll;
584c1414c8SBarry Smith 
594c1414c8SBarry Smith   PetscFunctionBegin;
60da9f3051SSatish Balay   no_inode             = PETSC_FALSE;
61da9f3051SSatish Balay   no_unroll            = PETSC_FALSE;
62ec710b6aSStefano Zampini   b->inode.checked     = PETSC_FALSE;
634c1414c8SBarry Smith   b->inode.node_count  = 0;
644d12350bSJunchao Zhang   b->inode.size_csr    = NULL;
654c1414c8SBarry Smith   b->inode.limit       = 5;
664c1414c8SBarry Smith   b->inode.max_limit   = 5;
6771f1c65dSBarry Smith   b->inode.ibdiagvalid = PETSC_FALSE;
68f4259b30SLisandro Dalcin   b->inode.ibdiag      = NULL;
69f4259b30SLisandro Dalcin   b->inode.bdiag       = NULL;
704c1414c8SBarry Smith 
71d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B), ((PetscObject)B)->prefix, "Options for SEQAIJ matrix", "Mat");
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_unroll", "Do not optimize for inodes (slower)", NULL, no_unroll, &no_unroll, NULL));
7348a46eb9SPierre Jolivet   if (no_unroll) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_unroll\n"));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_inode", "Do not optimize for inodes -slower-", NULL, no_inode, &no_inode, NULL));
7548a46eb9SPierre Jolivet   if (no_inode) PetscCall(PetscInfo(B, "Not using Inode routines due to -mat_no_inode\n"));
76*e87b5d96SPierre Jolivet   PetscCall(PetscOptionsInt("-mat_inode_limit", "Do not use inodes larger than this value", NULL, b->inode.limit, &b->inode.limit, NULL));
77d0609cedSBarry Smith   PetscOptionsEnd();
7826fbe8dcSKarl Rupp 
79ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
80d74fe821SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
81d74fe821SBarry Smith 
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatInodeAdjustForInodes_C", MatInodeAdjustForInodes_SeqAIJ_Inode));
839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatInodeGetInodeSizes_C", MatInodeGetInodeSizes_SeqAIJ_Inode));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
854c1414c8SBarry Smith }
864c1414c8SBarry Smith 
87d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A, MatOption op, PetscBool flg)
88d71ae5a4SJacob Faibussowitsch {
894c1414c8SBarry Smith   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
90d74fe821SBarry Smith 
914c1414c8SBarry Smith   PetscFunctionBegin;
924c1414c8SBarry Smith   switch (op) {
93d71ae5a4SJacob Faibussowitsch   case MAT_USE_INODES:
94d71ae5a4SJacob Faibussowitsch     a->inode.use = flg;
95d71ae5a4SJacob Faibussowitsch     break;
96d71ae5a4SJacob Faibussowitsch   default:
97d71ae5a4SJacob Faibussowitsch     break;
984c1414c8SBarry Smith   }
993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1004c1414c8SBarry Smith }
101