xref: /petsc/src/mat/impls/aij/seq/inode2.c (revision 7d0a6c19129e7069c8a40e210b34ed62989173db)
1*7d0a6c19SBarry Smith 
27c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h"
34c1414c8SBarry Smith 
409573ac7SBarry Smith extern PetscErrorCode Mat_CheckInode(Mat,PetscBool );
54c1414c8SBarry Smith EXTERN_C_BEGIN
67087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes_SeqAIJ_Inode(Mat,IS*,IS*);
77087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes_SeqAIJ_Inode(Mat,PetscInt*,PetscInt*[],PetscInt*);
84c1414c8SBarry Smith EXTERN_C_END
94c1414c8SBarry Smith 
104c1414c8SBarry Smith #undef __FUNCT__
114108e4d5SBarry Smith #define __FUNCT__ "MatView_SeqAIJ_Inode"
124108e4d5SBarry Smith PetscErrorCode MatView_SeqAIJ_Inode(Mat A,PetscViewer viewer)
134c1414c8SBarry Smith {
144c1414c8SBarry Smith   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data;
154c1414c8SBarry Smith   PetscErrorCode    ierr;
16ace3abfcSBarry Smith   PetscBool         iascii;
174c1414c8SBarry Smith   PetscViewerFormat format;
184c1414c8SBarry Smith 
194c1414c8SBarry Smith   PetscFunctionBegin;
202692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
214c1414c8SBarry Smith   if (iascii) {
224c1414c8SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
234c1414c8SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL || format == PETSC_VIEWER_ASCII_INFO) {
244c1414c8SBarry Smith       if (a->inode.size) {
254c1414c8SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"using I-node routines: found %D nodes, limit used is %D\n",
264c1414c8SBarry Smith                                       a->inode.node_count,a->inode.limit);CHKERRQ(ierr);
274c1414c8SBarry Smith       } else {
284c1414c8SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"not using I-node routines\n");CHKERRQ(ierr);
294c1414c8SBarry Smith       }
304c1414c8SBarry Smith     }
314c1414c8SBarry Smith   }
324c1414c8SBarry Smith   PetscFunctionReturn(0);
334c1414c8SBarry Smith }
344c1414c8SBarry Smith 
354c1414c8SBarry Smith #undef __FUNCT__
364108e4d5SBarry Smith #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_Inode"
374108e4d5SBarry Smith PetscErrorCode MatAssemblyEnd_SeqAIJ_Inode(Mat A, MatAssemblyType mode)
384c1414c8SBarry Smith {
3971f1c65dSBarry Smith   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
404c1414c8SBarry Smith   PetscErrorCode ierr;
41ace3abfcSBarry Smith   PetscBool      samestructure;
424c1414c8SBarry Smith 
434c1414c8SBarry Smith   PetscFunctionBegin;
444c1414c8SBarry Smith   /* info.nz_unneeded of zero denotes no structural change was made to the matrix during Assembly */
45ace3abfcSBarry Smith   samestructure = (PetscBool)(!A->info.nz_unneeded);
464c1414c8SBarry Smith   /* check for identical nodes. If found, use inode functions */
474c1414c8SBarry Smith   ierr = Mat_CheckInode(A,samestructure);CHKERRQ(ierr);
4871f1c65dSBarry Smith   a->inode.ibdiagvalid = PETSC_FALSE;
494c1414c8SBarry Smith   PetscFunctionReturn(0);
504c1414c8SBarry Smith }
514c1414c8SBarry Smith 
524c1414c8SBarry Smith #undef __FUNCT__
534108e4d5SBarry Smith #define __FUNCT__ "MatDestroy_SeqAIJ_Inode"
544108e4d5SBarry Smith PetscErrorCode MatDestroy_SeqAIJ_Inode(Mat A)
554c1414c8SBarry Smith {
564c1414c8SBarry Smith   PetscErrorCode ierr;
574c1414c8SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
584c1414c8SBarry Smith 
594c1414c8SBarry Smith   PetscFunctionBegin;
604c1414c8SBarry Smith   ierr = PetscFree(a->inode.size);CHKERRQ(ierr);
6189c6957cSBarry Smith   ierr = PetscFree3(a->inode.ibdiag,a->inode.bdiag,a->inode.ssor_work);CHKERRQ(ierr);
624c1414c8SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeAdjustForInodes_C","",PETSC_NULL);CHKERRQ(ierr);
634c1414c8SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatInodeGetInodeSizes_C","",PETSC_NULL);CHKERRQ(ierr);
644c1414c8SBarry Smith   PetscFunctionReturn(0);
654c1414c8SBarry Smith }
664c1414c8SBarry Smith 
674108e4d5SBarry Smith /* MatCreate_SeqAIJ_Inode is not DLLEXPORTed because it is not a constructor for a complete type.    */
684c1414c8SBarry Smith /* It is also not registered as a type for use within MatSetType.                             */
694c1414c8SBarry Smith /* It is intended as a helper for the MATSEQAIJ class, so classes which desire Inodes should  */
704c1414c8SBarry Smith /*    inherit off of MATSEQAIJ instead by calling MatSetType(MATSEQAIJ) in their constructor. */
714c1414c8SBarry Smith /* Maybe this is a bad idea. (?) */
724c1414c8SBarry Smith #undef __FUNCT__
734108e4d5SBarry Smith #define __FUNCT__ "MatCreate_SeqAIJ_Inode"
744108e4d5SBarry Smith PetscErrorCode MatCreate_SeqAIJ_Inode(Mat B)
754c1414c8SBarry Smith {
764c1414c8SBarry Smith   Mat_SeqAIJ     *b=(Mat_SeqAIJ*)B->data;
774c1414c8SBarry Smith   PetscErrorCode ierr;
78ace3abfcSBarry Smith   PetscBool      no_inode,no_unroll;
794c1414c8SBarry Smith 
804c1414c8SBarry Smith   PetscFunctionBegin;
81da9f3051SSatish Balay   no_inode             = PETSC_FALSE;
82da9f3051SSatish Balay   no_unroll            = PETSC_FALSE;
834c1414c8SBarry Smith   b->inode.node_count  = 0;
844c1414c8SBarry Smith   b->inode.size        = 0;
854c1414c8SBarry Smith   b->inode.limit       = 5;
864c1414c8SBarry Smith   b->inode.max_limit   = 5;
8771f1c65dSBarry Smith   b->inode.ibdiagvalid = PETSC_FALSE;
8871f1c65dSBarry Smith   b->inode.ibdiag      = 0;
8971f1c65dSBarry Smith   b->inode.bdiag       = 0;
904c1414c8SBarry Smith 
917adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for SEQAIJ matrix","Mat");CHKERRQ(ierr);
92acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);CHKERRQ(ierr);
93da9f3051SSatish Balay     if (no_unroll) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);}
94acfcf0e5SJed Brown     ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes -slower-",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);CHKERRQ(ierr);
95da9f3051SSatish Balay     if (no_inode) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);}
96d74fe821SBarry Smith     ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",PETSC_NULL,b->inode.limit,&b->inode.limit,PETSC_NULL);CHKERRQ(ierr);
97d74fe821SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
98ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
99d74fe821SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
100d74fe821SBarry Smith 
1014c1414c8SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeAdjustForInodes_C",
1024108e4d5SBarry Smith                                      "MatInodeAdjustForInodes_SeqAIJ_Inode",
1034108e4d5SBarry Smith                                       MatInodeAdjustForInodes_SeqAIJ_Inode);CHKERRQ(ierr);
1044c1414c8SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatInodeGetInodeSizes_C",
1054108e4d5SBarry Smith                                      "MatInodeGetInodeSizes_SeqAIJ_Inode",
1064108e4d5SBarry Smith                                       MatInodeGetInodeSizes_SeqAIJ_Inode);CHKERRQ(ierr);
1074c1414c8SBarry Smith   PetscFunctionReturn(0);
1084c1414c8SBarry Smith }
1094c1414c8SBarry Smith 
1104c1414c8SBarry Smith #undef __FUNCT__
1114108e4d5SBarry Smith #define __FUNCT__ "MatSetOption_SeqAIJ_Inode"
112ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqAIJ_Inode(Mat A,MatOption op,PetscBool  flg)
1134c1414c8SBarry Smith {
1144c1414c8SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
115d74fe821SBarry Smith 
1164c1414c8SBarry Smith   PetscFunctionBegin;
1174c1414c8SBarry Smith   switch(op) {
1184c1414c8SBarry Smith     case MAT_USE_INODES:
1194e0d8c25SBarry Smith       a->inode.use         = flg;
1204c1414c8SBarry Smith       break;
1214c1414c8SBarry Smith     default:
1224c1414c8SBarry Smith       break;
1234c1414c8SBarry Smith   }
1244c1414c8SBarry Smith   PetscFunctionReturn(0);
1254c1414c8SBarry Smith }
1264c1414c8SBarry Smith 
1274c1414c8SBarry Smith #undef __FUNCT__
1284108e4d5SBarry Smith #define __FUNCT__ "MatDuplicate_SeqAIJ_Inode"
1294108e4d5SBarry Smith PetscErrorCode MatDuplicate_SeqAIJ_Inode(Mat A,MatDuplicateOption cpvalues,Mat *C)
1304c1414c8SBarry Smith {
1314c1414c8SBarry Smith   Mat            B=*C;
1324c1414c8SBarry Smith   Mat_SeqAIJ     *c=(Mat_SeqAIJ*)B->data,*a=(Mat_SeqAIJ*)A->data;
1334c1414c8SBarry Smith   PetscErrorCode ierr;
134d0f46423SBarry Smith   PetscInt       m=A->rmap->n;
1354c1414c8SBarry Smith 
1364c1414c8SBarry Smith   PetscFunctionBegin;
1374c1414c8SBarry Smith 
1384c1414c8SBarry Smith   c->inode.use          = a->inode.use;
1394c1414c8SBarry Smith   c->inode.limit        = a->inode.limit;
1404c1414c8SBarry Smith   c->inode.max_limit    = a->inode.max_limit;
1414c1414c8SBarry Smith   if (a->inode.size){
1424c1414c8SBarry Smith     ierr                = PetscMalloc((m+1)*sizeof(PetscInt),&c->inode.size);CHKERRQ(ierr);
1434c1414c8SBarry Smith     c->inode.node_count = a->inode.node_count;
1444c1414c8SBarry Smith     ierr                = PetscMemcpy(c->inode.size,a->inode.size,(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
1454c1414c8SBarry Smith   } else {
1464c1414c8SBarry Smith     c->inode.size       = 0;
1474c1414c8SBarry Smith     c->inode.node_count = 0;
1484c1414c8SBarry Smith   }
14971f1c65dSBarry Smith   c->inode.ibdiagvalid = PETSC_FALSE;
15071f1c65dSBarry Smith   c->inode.ibdiag      = 0;
15171f1c65dSBarry Smith   c->inode.bdiag       = 0;
1524c1414c8SBarry Smith   PetscFunctionReturn(0);
1534c1414c8SBarry Smith }
1544c1414c8SBarry Smith 
1554c1414c8SBarry Smith 
1564c1414c8SBarry Smith 
157