xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
74108e4d5SBarry Smith PetscErrorCode MatView_SeqAIJ_Inode(Mat A,PetscViewer viewer)
84c1414c8SBarry Smith {
94c1414c8SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
10ace3abfcSBarry Smith   PetscBool         iascii;
114c1414c8SBarry Smith   PetscViewerFormat format;
124c1414c8SBarry Smith 
134c1414c8SBarry Smith   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
154c1414c8SBarry Smith   if (iascii) {
169566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
174c1414c8SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
184c1414c8SBarry Smith       if (a->inode.size) {
199566063dSJacob 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));
204c1414c8SBarry Smith       } else {
219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"not using I-node routines\n"));
224c1414c8SBarry Smith       }
234c1414c8SBarry Smith     }
244c1414c8SBarry Smith   }
254c1414c8SBarry Smith   PetscFunctionReturn(0);
264c1414c8SBarry Smith }
274c1414c8SBarry Smith 
284108e4d5SBarry Smith PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
294c1414c8SBarry Smith {
3071f1c65dSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
314c1414c8SBarry Smith 
324c1414c8SBarry Smith   PetscFunctionBegin;
339566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJCheckInode(A));
3471f1c65dSBarry Smith   a->inode.ibdiagvalid = PETSC_FALSE;
354c1414c8SBarry Smith   PetscFunctionReturn(0);
364c1414c8SBarry Smith }
374c1414c8SBarry Smith 
384108e4d5SBarry Smith PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
394c1414c8SBarry Smith {
404c1414c8SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
414c1414c8SBarry Smith 
424c1414c8SBarry Smith   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->inode.size));
449566063dSJacob Faibussowitsch   PetscCall(PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work));
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatInodeAdjustForInodes_C",NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatInodeGetInodeSizes_C",NULL));
474c1414c8SBarry Smith   PetscFunctionReturn(0);
484c1414c8SBarry Smith }
494c1414c8SBarry Smith 
504108e4d5SBarry Smith /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
514c1414c8SBarry Smith /* It is also not registered as a type for use within MatSetType.                             */
524c1414c8SBarry Smith /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
534c1414c8SBarry Smith /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
544c1414c8SBarry Smith /* Maybe this is a bad idea. (?) */
554108e4d5SBarry Smith PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
564c1414c8SBarry Smith {
574c1414c8SBarry Smith   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
58ace3abfcSBarry Smith   PetscBool      no_inode,no_unroll;
594c1414c8SBarry Smith 
604c1414c8SBarry Smith   PetscFunctionBegin;
61da9f3051SSatish Balay   no_inode             = PETSC_FALSE;
62da9f3051SSatish Balay   no_unroll            = PETSC_FALSE;
63ec710b6aSStefano Zampini   b->inode.checked     = PETSC_FALSE;
644c1414c8SBarry Smith   b->inode.node_count  = 0;
65f4259b30SLisandro Dalcin   b->inode.size        = NULL;
664c1414c8SBarry Smith   b->inode.limit       = 5;
674c1414c8SBarry Smith   b->inode.max_limit   = 5;
6871f1c65dSBarry Smith   b->inode.ibdiagvalid = PETSC_FALSE;
69f4259b30SLisandro Dalcin   b->inode.ibdiag      = NULL;
70f4259b30SLisandro Dalcin   b->inode.bdiag       = NULL;
714c1414c8SBarry Smith 
72*d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL));
742205254eSKarl Rupp   if (no_unroll) {
759566063dSJacob Faibussowitsch     PetscCall(PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n"));
762205254eSKarl Rupp   }
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",NULL,no_inode,&no_inode,NULL));
782205254eSKarl Rupp   if (no_inode) {
799566063dSJacob Faibussowitsch     PetscCall(PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n"));
802205254eSKarl Rupp   }
819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL));
82*d0609cedSBarry Smith   PetscOptionsEnd();
8326fbe8dcSKarl Rupp 
84ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
85d74fe821SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
86d74fe821SBarry Smith 
879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatInodeAdjustForInodes_C",MatInodeAdjustForInodes_SeqAIJ_Inode));
889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatInodeGetInodeSizes_C",MatInodeGetInodeSizes_SeqAIJ_Inode));
894c1414c8SBarry Smith   PetscFunctionReturn(0);
904c1414c8SBarry Smith }
914c1414c8SBarry Smith 
92ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool flg)
934c1414c8SBarry Smith {
944c1414c8SBarry Smith   Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
95d74fe821SBarry Smith 
964c1414c8SBarry Smith   PetscFunctionBegin;
974c1414c8SBarry Smith   switch (op) {
984c1414c8SBarry Smith   case MAT_USE_INODES:
994e0d8c25SBarry Smith     a->inode.use = flg;
1004c1414c8SBarry Smith     break;
1014c1414c8SBarry Smith   default:
1024c1414c8SBarry Smith     break;
1034c1414c8SBarry Smith   }
1044c1414c8SBarry Smith   PetscFunctionReturn(0);
1054c1414c8SBarry Smith }
106