xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision b9c12af58e66e36adbaa844832c74a52f523537f)
17d0a6c19SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
34c1414c8SBarry Smith 
4a02bda8eSBarry Smith extern PetscErrorCode MatSeqAIJCheckInode(Mat);
57087cfbeSBarry Smith extern PetscErrorCode MatInodeAdjustForInodes_SeqAIJ_Inode(Mat,IS*,IS*);
67087cfbeSBarry Smith extern PetscErrorCode MatInodeGetInodeSizes_SeqAIJ_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);
74c1414c8SBarry Smith 
84c1414c8SBarry Smith #undef __FUNCT__
94108e4d5SBarry Smith #define __FUNCT__ "MatView_SeqAIJ_Inode"
104108e4d5SBarry Smith PetscErrorCode MatView_SeqAIJ_Inode(Mat A,PetscViewer viewer)
114c1414c8SBarry Smith {
124c1414c8SBarry Smith   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
134c1414c8SBarry Smith   PetscErrorCode    ierr;
14ace3abfcSBarry Smith   PetscBool         iascii;
154c1414c8SBarry Smith   PetscViewerFormat format;
164c1414c8SBarry Smith 
174c1414c8SBarry Smith   PetscFunctionBegin;
18251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
194c1414c8SBarry Smith   if (iascii) {
204c1414c8SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
214c1414c8SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
224c1414c8SBarry Smith       if (a->inode.size) {
23*b9c12af5SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n",a->inode.node_count,a->inode.limit);CHKERRQ(ierr);
244c1414c8SBarry Smith       } else {
254c1414c8SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");CHKERRQ(ierr);
264c1414c8SBarry Smith       }
274c1414c8SBarry Smith     }
284c1414c8SBarry Smith   }
294c1414c8SBarry Smith   PetscFunctionReturn(0);
304c1414c8SBarry Smith }
314c1414c8SBarry Smith 
324c1414c8SBarry Smith #undef __FUNCT__
334108e4d5SBarry Smith #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_Inode"
344108e4d5SBarry Smith PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
354c1414c8SBarry Smith {
3671f1c65dSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
374c1414c8SBarry Smith   PetscErrorCode ierr;
384c1414c8SBarry Smith 
394c1414c8SBarry Smith   PetscFunctionBegin;
40a02bda8eSBarry Smith   ierr = MatSeqAIJCheckInode(A);CHKERRQ(ierr);
4171f1c65dSBarry Smith   a->inode.ibdiagvalid = PETSC_FALSE;
424c1414c8SBarry Smith   PetscFunctionReturn(0);
434c1414c8SBarry Smith }
444c1414c8SBarry Smith 
454c1414c8SBarry Smith #undef __FUNCT__
464108e4d5SBarry Smith #define __FUNCT__ "MatDestroy_SeqAIJ_Inode"
474108e4d5SBarry Smith PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
484c1414c8SBarry Smith {
494c1414c8SBarry Smith   PetscErrorCode ierr;
504c1414c8SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
514c1414c8SBarry Smith 
524c1414c8SBarry Smith   PetscFunctionBegin;
534c1414c8SBarry Smith   ierr = PetscFree(a->inode.size);CHKERRQ(ierr);
5489c6957cSBarry Smith   ierr = PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);CHKERRQ(ierr);
55bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatInodeAdjustForInodes_C",NULL);CHKERRQ(ierr);
56bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatInodeGetInodeSizes_C",NULL);CHKERRQ(ierr);
574c1414c8SBarry Smith   PetscFunctionReturn(0);
584c1414c8SBarry Smith }
594c1414c8SBarry Smith 
604108e4d5SBarry Smith /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
614c1414c8SBarry Smith /* It is also not registered as a type for use within MatSetType.                             */
624c1414c8SBarry Smith /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
634c1414c8SBarry Smith /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
644c1414c8SBarry Smith /* Maybe this is a bad idea. (?) */
654c1414c8SBarry Smith #undef __FUNCT__
664108e4d5SBarry Smith #define __FUNCT__ "MatCreate_SeqAIJ_Inode"
674108e4d5SBarry Smith PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
684c1414c8SBarry Smith {
694c1414c8SBarry Smith   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
704c1414c8SBarry Smith   PetscErrorCode ierr;
71ace3abfcSBarry Smith   PetscBool      no_inode,no_unroll;
724c1414c8SBarry Smith 
734c1414c8SBarry Smith   PetscFunctionBegin;
74da9f3051SSatish Balay   no_inode             = PETSC_FALSE;
75da9f3051SSatish Balay   no_unroll            = PETSC_FALSE;
764c1414c8SBarry Smith   b->inode.node_count  = 0;
774c1414c8SBarry Smith   b->inode.size        = 0;
784c1414c8SBarry Smith   b->inode.limit       = 5;
794c1414c8SBarry Smith   b->inode.max_limit   = 5;
8071f1c65dSBarry Smith   b->inode.ibdiagvalid = PETSC_FALSE;
8171f1c65dSBarry Smith   b->inode.ibdiag      = 0;
8271f1c65dSBarry Smith   b->inode.bdiag       = 0;
834c1414c8SBarry Smith 
84ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr);
850298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);CHKERRQ(ierr);
862205254eSKarl Rupp   if (no_unroll) {
872205254eSKarl Rupp     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);
882205254eSKarl Rupp   }
890298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",NULL,no_inode,&no_inode,NULL);CHKERRQ(ierr);
902205254eSKarl Rupp   if (no_inode) {
912205254eSKarl Rupp     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);
922205254eSKarl Rupp   }
930298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);CHKERRQ(ierr);
94d74fe821SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
9526fbe8dcSKarl Rupp 
96ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
97d74fe821SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
98d74fe821SBarry Smith 
99bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatInodeAdjustForInodes_C",MatInodeAdjustForInodes_SeqAIJ_Inode);CHKERRQ(ierr);
100bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatInodeGetInodeSizes_C",MatInodeGetInodeSizes_SeqAIJ_Inode);CHKERRQ(ierr);
1014c1414c8SBarry Smith   PetscFunctionReturn(0);
1024c1414c8SBarry Smith }
1034c1414c8SBarry Smith 
1044c1414c8SBarry Smith #undef __FUNCT__
1054108e4d5SBarry Smith #define __FUNCT__ "MatSetOption_SeqAIJ_Inode"
106ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool flg)
1074c1414c8SBarry Smith {
1084c1414c8SBarry Smith   Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
109d74fe821SBarry Smith 
1104c1414c8SBarry Smith   PetscFunctionBegin;
1114c1414c8SBarry Smith   switch (op) {
1124c1414c8SBarry Smith   case MAT_USE_INODES:
1134e0d8c25SBarry Smith     a->inode.use = flg;
1144c1414c8SBarry Smith     break;
1154c1414c8SBarry Smith   default:
1164c1414c8SBarry Smith     break;
1174c1414c8SBarry Smith   }
1184c1414c8SBarry Smith   PetscFunctionReturn(0);
1194c1414c8SBarry Smith }
1204c1414c8SBarry Smith 
1214c1414c8SBarry Smith 
1224c1414c8SBarry Smith 
1234c1414c8SBarry Smith 
1244c1414c8SBarry Smith 
125