1be1d678aSKris Buschelman #define PETSCMAT_DLL 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 4a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/inline/spops.h" 9aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1049b5e25fSSatish Balay 1149b5e25fSSatish Balay #define CHUNKSIZE 10 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay /* 1449b5e25fSSatish Balay Checks for missing diagonals 1549b5e25fSSatish Balay */ 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 18dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A) 1949b5e25fSSatish Balay { 20045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 216849ba73SBarry Smith PetscErrorCode ierr; 2213f74950SBarry Smith PetscInt *diag,*jj = a->j,i; 2349b5e25fSSatish Balay 2449b5e25fSSatish Balay PetscFunctionBegin; 25045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2649b5e25fSSatish Balay diag = a->diag; 2749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2877431f27SBarry Smith if (jj[diag[i]] != i) SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Matrix is missing diagonal number %D",i); 2949b5e25fSSatish Balay } 3049b5e25fSSatish Balay PetscFunctionReturn(0); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay 334a2ae208SSatish Balay #undef __FUNCT__ 344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 35dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 3649b5e25fSSatish Balay { 37045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 386849ba73SBarry Smith PetscErrorCode ierr; 3909f38230SBarry Smith PetscInt i; 4049b5e25fSSatish Balay 4149b5e25fSSatish Balay PetscFunctionBegin; 4209f38230SBarry Smith if (!a->diag) { 4309f38230SBarry Smith ierr = PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 4409f38230SBarry Smith } 4509f38230SBarry Smith for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i]; 4649b5e25fSSatish Balay PetscFunctionReturn(0); 4749b5e25fSSatish Balay } 4849b5e25fSSatish Balay 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 518f7157efSSatish Balay static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 5249b5e25fSSatish Balay { 53a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 548f7157efSSatish Balay PetscInt i,j,n = a->mbs,nz = a->i[n],bs = A->rmap.bs; 558f7157efSSatish Balay PetscErrorCode ierr; 5649b5e25fSSatish Balay 5749b5e25fSSatish Balay PetscFunctionBegin; 58d3e5a4abSHong Zhang *nn = n; 59a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 608f7157efSSatish Balay if (!blockcompressed) { 618f7157efSSatish Balay /* malloc & create the natural set of indices */ 62f1d0d59dSSatish Balay ierr = PetscMalloc2((n+1)*bs,PetscInt,ia,nz*bs,PetscInt,ja);CHKERRQ(ierr); 638f7157efSSatish Balay for (i=0; i<n+1; i++) { 648f7157efSSatish Balay for (j=0; j<bs; j++) { 658f7157efSSatish Balay *ia[i*bs+j] = a->i[i]*bs+j+oshift; 668f7157efSSatish Balay } 678f7157efSSatish Balay } 688f7157efSSatish Balay for (i=0; i<nz; i++) { 698f7157efSSatish Balay for (j=0; j<bs; j++) { 708f7157efSSatish Balay *ja[i*bs+j] = a->j[i]*bs+j+oshift; 718f7157efSSatish Balay } 728f7157efSSatish Balay } 738f7157efSSatish Balay } else { /* blockcompressed */ 74a6ece127SHong Zhang if (oshift == 1) { 75a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 766c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 77a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 788f7157efSSatish Balay } 79a1373b80SHong Zhang *ia = a->i; *ja = a->j; 80a6ece127SHong Zhang } 818f7157efSSatish Balay 8249b5e25fSSatish Balay PetscFunctionReturn(0); 8349b5e25fSSatish Balay } 8449b5e25fSSatish Balay 854a2ae208SSatish Balay #undef __FUNCT__ 864a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 878f7157efSSatish Balay static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 8849b5e25fSSatish Balay { 89b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 908f7157efSSatish Balay PetscInt i,n = a->mbs,nz = a->i[n]; 918f7157efSSatish Balay PetscErrorCode ierr; 92a6ece127SHong Zhang 9349b5e25fSSatish Balay PetscFunctionBegin; 9449b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 95a6ece127SHong Zhang 968f7157efSSatish Balay if (!blockcompressed) { 978f7157efSSatish Balay ierr = PetscFree2(*ia,*ja);CHKERRQ(ierr); 988f7157efSSatish Balay } else if (oshift == 1) { /* blockcompressed */ 996c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 100a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 101a6ece127SHong Zhang } 1028f7157efSSatish Balay 103a6ece127SHong Zhang PetscFunctionReturn(0); 10449b5e25fSSatish Balay } 10549b5e25fSSatish Balay 1064a2ae208SSatish Balay #undef __FUNCT__ 1074a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 108dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 10949b5e25fSSatish Balay { 11049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 111dfbe8321SBarry Smith PetscErrorCode ierr; 11249b5e25fSSatish Balay 11349b5e25fSSatish Balay PetscFunctionBegin; 114a9f03627SSatish Balay #if defined(PETSC_USE_LOG) 115899cda47SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap.N,a->nz); 116a9f03627SSatish Balay #endif 117e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr); 1189bfd6278SHong Zhang if (a->row) {ierr = ISDestroy(a->row);CHKERRQ(ierr);} 1199bfd6278SHong Zhang if (a->col){ierr = ISDestroy(a->col);CHKERRQ(ierr);} 1209bfd6278SHong Zhang if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 12105b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 12205b42c5fSBarry Smith ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr); 12305b42c5fSBarry Smith ierr = PetscFree(a->solve_work);CHKERRQ(ierr); 124e2ee2a47SBarry Smith ierr = PetscFree(a->relax_work);CHKERRQ(ierr); 12505b42c5fSBarry Smith ierr = PetscFree(a->solves_work);CHKERRQ(ierr); 12605b42c5fSBarry Smith ierr = PetscFree(a->mult_work);CHKERRQ(ierr); 12705b42c5fSBarry Smith ierr = PetscFree(a->saved_values);CHKERRQ(ierr); 12805b42c5fSBarry Smith ierr = PetscFree(a->xtoy);CHKERRQ(ierr); 1291a3463dfSHong Zhang 1301a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 13149b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 132901853e0SKris Buschelman 133dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 134901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 135901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 136901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 137901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 138901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 139901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 14049b5e25fSSatish Balay PetscFunctionReturn(0); 14149b5e25fSSatish Balay } 14249b5e25fSSatish Balay 1434a2ae208SSatish Balay #undef __FUNCT__ 1444a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 1454e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscTruth flg) 14649b5e25fSSatish Balay { 147045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 14863ba0a88SBarry Smith PetscErrorCode ierr; 14949b5e25fSSatish Balay 15049b5e25fSSatish Balay PetscFunctionBegin; 1514d9d31abSKris Buschelman switch (op) { 1524d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1534e0d8c25SBarry Smith a->roworiented = flg; 1544d9d31abSKris Buschelman break; 1554d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1564e0d8c25SBarry Smith a->keepzeroedrows = flg; 1574d9d31abSKris Buschelman break; 158*512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 159*512a5fc5SBarry Smith a->nonew = (flg ? 0 : 1); 1604d9d31abSKris Buschelman break; 1614d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1624e0d8c25SBarry Smith a->nonew = (flg ? -1 : 0); 1634d9d31abSKris Buschelman break; 1644d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1654e0d8c25SBarry Smith a->nonew = (flg ? -2 : 0); 1664d9d31abSKris Buschelman break; 1674e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1684d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1694d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 170290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 1714d9d31abSKris Buschelman break; 1729a4540c5SBarry Smith case MAT_HERMITIAN: 1734e0d8c25SBarry Smith if (flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 17477e54ba9SKris Buschelman case MAT_SYMMETRIC: 17577e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1769a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1774e0d8c25SBarry Smith if (!flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 178290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);CHKERRQ(ierr); 179290bbb0aSBarry Smith break; 180941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 1814e0d8c25SBarry Smith a->ignore_ltriangular = flg; 182941593c8SHong Zhang break; 183941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 1844e0d8c25SBarry Smith a->ignore_ltriangular = flg; 18577e54ba9SKris Buschelman break; 186f5edf698SHong Zhang case MAT_GETROW_UPPERTRIANGULAR: 1874e0d8c25SBarry Smith a->getrow_utriangular = flg; 188f5edf698SHong Zhang break; 1894d9d31abSKris Buschelman default: 190ad86a440SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 19149b5e25fSSatish Balay } 19249b5e25fSSatish Balay PetscFunctionReturn(0); 19349b5e25fSSatish Balay } 19449b5e25fSSatish Balay 1954a2ae208SSatish Balay #undef __FUNCT__ 1964a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 19713f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 19849b5e25fSSatish Balay { 19949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2006849ba73SBarry Smith PetscErrorCode ierr; 20113f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 20249b5e25fSSatish Balay MatScalar *aa,*aa_i; 20387828ca2SBarry Smith PetscScalar *v_i; 20449b5e25fSSatish Balay 20549b5e25fSSatish Balay PetscFunctionBegin; 2064e0d8c25SBarry Smith if (A && !a->getrow_utriangular) SETERRQ(PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()"); 207f5edf698SHong Zhang /* Get the upper triangular part of the row */ 208899cda47SBarry Smith bs = A->rmap.bs; 20949b5e25fSSatish Balay ai = a->i; 21049b5e25fSSatish Balay aj = a->j; 21149b5e25fSSatish Balay aa = a->a; 21249b5e25fSSatish Balay bs2 = a->bs2; 21349b5e25fSSatish Balay 214899cda47SBarry Smith if (row < 0 || row >= A->rmap.N) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 21549b5e25fSSatish Balay 21649b5e25fSSatish Balay bn = row/bs; /* Block number */ 21749b5e25fSSatish Balay bp = row % bs; /* Block position */ 21849b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 21949b5e25fSSatish Balay *ncols = bs*M; 22049b5e25fSSatish Balay 22149b5e25fSSatish Balay if (v) { 22249b5e25fSSatish Balay *v = 0; 22349b5e25fSSatish Balay if (*ncols) { 22487828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 22549b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22649b5e25fSSatish Balay v_i = *v + i*bs; 22749b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 22849b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 22949b5e25fSSatish Balay } 23049b5e25fSSatish Balay } 23149b5e25fSSatish Balay } 23249b5e25fSSatish Balay 23349b5e25fSSatish Balay if (cols) { 23449b5e25fSSatish Balay *cols = 0; 23549b5e25fSSatish Balay if (*ncols) { 23613f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 23749b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23849b5e25fSSatish Balay cols_i = *cols + i*bs; 23949b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 24049b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24149b5e25fSSatish Balay } 24249b5e25fSSatish Balay } 24349b5e25fSSatish Balay } 24449b5e25fSSatish Balay 24549b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2465ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2475ddb2528SHong Zhang #ifdef column_search 24849b5e25fSSatish Balay v_i = *v + M*bs; 24949b5e25fSSatish Balay cols_i = *cols + M*bs; 25049b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25149b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25249b5e25fSSatish Balay for (j=0; j<M; j++){ 25349b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 25449b5e25fSSatish Balay if (itmp == bn){ 25549b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 25649b5e25fSSatish Balay for (k=0; k<bs; k++) { 25749b5e25fSSatish Balay *cols_i++ = i*bs+k; 25849b5e25fSSatish Balay *v_i++ = aa_i[k]; 25949b5e25fSSatish Balay } 26049b5e25fSSatish Balay *ncols += bs; 26149b5e25fSSatish Balay break; 26249b5e25fSSatish Balay } 26349b5e25fSSatish Balay } 26449b5e25fSSatish Balay } 2655ddb2528SHong Zhang #endif 26649b5e25fSSatish Balay PetscFunctionReturn(0); 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay 2694a2ae208SSatish Balay #undef __FUNCT__ 2704a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 27113f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27249b5e25fSSatish Balay { 273dfbe8321SBarry Smith PetscErrorCode ierr; 27449b5e25fSSatish Balay 27549b5e25fSSatish Balay PetscFunctionBegin; 27605b42c5fSBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 27705b42c5fSBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 27849b5e25fSSatish Balay PetscFunctionReturn(0); 27949b5e25fSSatish Balay } 28049b5e25fSSatish Balay 2814a2ae208SSatish Balay #undef __FUNCT__ 282f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ" 283f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 284f5edf698SHong Zhang { 285f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 286f5edf698SHong Zhang 287f5edf698SHong Zhang PetscFunctionBegin; 288f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 289f5edf698SHong Zhang PetscFunctionReturn(0); 290f5edf698SHong Zhang } 291f5edf698SHong Zhang #undef __FUNCT__ 292f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ" 293f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 294f5edf698SHong Zhang { 295f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 296f5edf698SHong Zhang 297f5edf698SHong Zhang PetscFunctionBegin; 298f5edf698SHong Zhang a->getrow_utriangular = PETSC_FALSE; 299f5edf698SHong Zhang PetscFunctionReturn(0); 300f5edf698SHong Zhang } 301f5edf698SHong Zhang 302f5edf698SHong Zhang #undef __FUNCT__ 3034a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 304dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 30549b5e25fSSatish Balay { 306dfbe8321SBarry Smith PetscErrorCode ierr; 30749b5e25fSSatish Balay PetscFunctionBegin; 308999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 3098115998fSBarry Smith PetscFunctionReturn(0); 31049b5e25fSSatish Balay } 31149b5e25fSSatish Balay 3124a2ae208SSatish Balay #undef __FUNCT__ 3134a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3146849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 31549b5e25fSSatish Balay { 31649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 317dfbe8321SBarry Smith PetscErrorCode ierr; 318899cda47SBarry Smith PetscInt i,j,bs = A->rmap.bs,k,l,bs2=a->bs2; 3192dcb1b2aSMatthew Knepley const char *name; 320f3ef73ceSBarry Smith PetscViewerFormat format; 32149b5e25fSSatish Balay 32249b5e25fSSatish Balay PetscFunctionBegin; 32380fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 324b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 325456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 32677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 327fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 32870d5e725SHong Zhang if (A->factor && bs>1){ 32970d5e725SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr); 33070d5e725SHong Zhang PetscFunctionReturn(0); 33170d5e725SHong Zhang } 332c9f458caSMatthew Knepley Mat aij; 333c9f458caSMatthew Knepley ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr); 334c9f458caSMatthew Knepley ierr = MatView(aij,viewer);CHKERRQ(ierr); 335c9f458caSMatthew Knepley ierr = MatDestroy(aij);CHKERRQ(ierr); 336fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 337b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 33849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33949b5e25fSSatish Balay for (j=0; j<bs; j++) { 34077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 34149b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 34249b5e25fSSatish Balay for (l=0; l<bs; l++) { 34349b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 34449b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 345a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 34649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34749b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 348a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 34949b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35049b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 351a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35249b5e25fSSatish Balay } 35349b5e25fSSatish Balay #else 35449b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 355a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35649b5e25fSSatish Balay } 35749b5e25fSSatish Balay #endif 35849b5e25fSSatish Balay } 35949b5e25fSSatish Balay } 360b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 36149b5e25fSSatish Balay } 36249b5e25fSSatish Balay } 363b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 364c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 365c1490034SHong Zhang PetscFunctionReturn(0); 36649b5e25fSSatish Balay } else { 3678608aa04SHong Zhang if (A->factor && bs>1){ 3688608aa04SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored. MatView_SeqSBAIJ_ASCII() may not display complete or logically correct entries!\n");CHKERRQ(ierr); 3698608aa04SHong Zhang } 370b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 37149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 37249b5e25fSSatish Balay for (j=0; j<bs; j++) { 37377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 37449b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 37549b5e25fSSatish Balay for (l=0; l<bs; l++) { 37649b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37749b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 378a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 37949b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38049b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 381a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 38249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38349b5e25fSSatish Balay } else { 384a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38549b5e25fSSatish Balay } 38649b5e25fSSatish Balay #else 387e9f7bc9eSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38849b5e25fSSatish Balay #endif 38949b5e25fSSatish Balay } 39049b5e25fSSatish Balay } 391b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 39249b5e25fSSatish Balay } 39349b5e25fSSatish Balay } 394b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 39549b5e25fSSatish Balay } 396b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 39749b5e25fSSatish Balay PetscFunctionReturn(0); 39849b5e25fSSatish Balay } 39949b5e25fSSatish Balay 4004a2ae208SSatish Balay #undef __FUNCT__ 4014a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 4026849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 40349b5e25fSSatish Balay { 40449b5e25fSSatish Balay Mat A = (Mat) Aa; 40549b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 4066849ba73SBarry Smith PetscErrorCode ierr; 407899cda47SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap.bs,bs2=a->bs2; 40813f74950SBarry Smith PetscMPIInt rank; 40949b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 41049b5e25fSSatish Balay MatScalar *aa; 41149b5e25fSSatish Balay MPI_Comm comm; 412b0a32e0cSBarry Smith PetscViewer viewer; 41349b5e25fSSatish Balay 41449b5e25fSSatish Balay PetscFunctionBegin; 41549b5e25fSSatish Balay /* 41649b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 41749b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 41849b5e25fSSatish Balay rest should return immediately. 41949b5e25fSSatish Balay */ 42049b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 42149b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 42249b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 42349b5e25fSSatish Balay 42449b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 42549b5e25fSSatish Balay 426b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 427b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 42849b5e25fSSatish Balay 42949b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 430b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 43149b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 43249b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 433899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 43449b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 43549b5e25fSSatish Balay aa = a->a + j*bs2; 43649b5e25fSSatish Balay for (k=0; k<bs; k++) { 43749b5e25fSSatish Balay for (l=0; l<bs; l++) { 43849b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 439b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44049b5e25fSSatish Balay } 44149b5e25fSSatish Balay } 44249b5e25fSSatish Balay } 44349b5e25fSSatish Balay } 444b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 44549b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 44649b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 447899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 44849b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 44949b5e25fSSatish Balay aa = a->a + j*bs2; 45049b5e25fSSatish Balay for (k=0; k<bs; k++) { 45149b5e25fSSatish Balay for (l=0; l<bs; l++) { 45249b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 453b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 45449b5e25fSSatish Balay } 45549b5e25fSSatish Balay } 45649b5e25fSSatish Balay } 45749b5e25fSSatish Balay } 45849b5e25fSSatish Balay 459b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 46049b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46149b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 462899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 46349b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46449b5e25fSSatish Balay aa = a->a + j*bs2; 46549b5e25fSSatish Balay for (k=0; k<bs; k++) { 46649b5e25fSSatish Balay for (l=0; l<bs; l++) { 46749b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 468b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46949b5e25fSSatish Balay } 47049b5e25fSSatish Balay } 47149b5e25fSSatish Balay } 47249b5e25fSSatish Balay } 47349b5e25fSSatish Balay PetscFunctionReturn(0); 47449b5e25fSSatish Balay } 47549b5e25fSSatish Balay 4764a2ae208SSatish Balay #undef __FUNCT__ 4774a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4786849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 47949b5e25fSSatish Balay { 480dfbe8321SBarry Smith PetscErrorCode ierr; 48149b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 482b0a32e0cSBarry Smith PetscDraw draw; 48349b5e25fSSatish Balay PetscTruth isnull; 48449b5e25fSSatish Balay 48549b5e25fSSatish Balay PetscFunctionBegin; 486b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 487b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 48849b5e25fSSatish Balay 48949b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 490899cda47SBarry Smith xr = A->rmap.N; yr = A->rmap.N; h = yr/10.0; w = xr/10.0; 49149b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 492b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 493b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 49449b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 49549b5e25fSSatish Balay PetscFunctionReturn(0); 49649b5e25fSSatish Balay } 49749b5e25fSSatish Balay 4984a2ae208SSatish Balay #undef __FUNCT__ 4994a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 500dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 50149b5e25fSSatish Balay { 502dfbe8321SBarry Smith PetscErrorCode ierr; 50332077d6dSBarry Smith PetscTruth iascii,isdraw; 50449b5e25fSSatish Balay 50549b5e25fSSatish Balay PetscFunctionBegin; 50632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 507fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 50832077d6dSBarry Smith if (iascii){ 50949b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 51049b5e25fSSatish Balay } else if (isdraw) { 51149b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 51249b5e25fSSatish Balay } else { 513a5e6ed63SBarry Smith Mat B; 514ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 515a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 516a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 51749b5e25fSSatish Balay } 51849b5e25fSSatish Balay PetscFunctionReturn(0); 51949b5e25fSSatish Balay } 52049b5e25fSSatish Balay 52149b5e25fSSatish Balay 5224a2ae208SSatish Balay #undef __FUNCT__ 5234a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 52413f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 52549b5e25fSSatish Balay { 526045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 52713f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 52813f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 529899cda47SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->rmap.bs,bs2=a->bs2; 53097e567efSBarry Smith MatScalar *ap,*aa = a->a; 53149b5e25fSSatish Balay 53249b5e25fSSatish Balay PetscFunctionBegin; 53349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 53449b5e25fSSatish Balay row = im[k]; brow = row/bs; 53597e567efSBarry Smith if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */ 536899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 53749b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 53849b5e25fSSatish Balay nrow = ailen[brow]; 53949b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 54097e567efSBarry Smith if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */ 541899cda47SBarry Smith if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1); 54249b5e25fSSatish Balay col = in[l] ; 54349b5e25fSSatish Balay bcol = col/bs; 54449b5e25fSSatish Balay cidx = col%bs; 54549b5e25fSSatish Balay ridx = row%bs; 54649b5e25fSSatish Balay high = nrow; 54749b5e25fSSatish Balay low = 0; /* assume unsorted */ 54849b5e25fSSatish Balay while (high-low > 5) { 54949b5e25fSSatish Balay t = (low+high)/2; 55049b5e25fSSatish Balay if (rp[t] > bcol) high = t; 55149b5e25fSSatish Balay else low = t; 55249b5e25fSSatish Balay } 55349b5e25fSSatish Balay for (i=low; i<high; i++) { 55449b5e25fSSatish Balay if (rp[i] > bcol) break; 55549b5e25fSSatish Balay if (rp[i] == bcol) { 55649b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 55749b5e25fSSatish Balay goto finished; 55849b5e25fSSatish Balay } 55949b5e25fSSatish Balay } 56097e567efSBarry Smith *v++ = 0.0; 56149b5e25fSSatish Balay finished:; 56249b5e25fSSatish Balay } 56349b5e25fSSatish Balay } 56449b5e25fSSatish Balay PetscFunctionReturn(0); 56549b5e25fSSatish Balay } 56649b5e25fSSatish Balay 56749b5e25fSSatish Balay 5684a2ae208SSatish Balay #undef __FUNCT__ 5694a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 57013f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 57149b5e25fSSatish Balay { 5720880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5736849ba73SBarry Smith PetscErrorCode ierr; 574e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 57513f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 576899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap.bs,stepval; 5770880e062SHong Zhang PetscTruth roworiented=a->roworiented; 578f15d580aSBarry Smith const MatScalar *value = v; 579f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5800880e062SHong Zhang 58149b5e25fSSatish Balay PetscFunctionBegin; 5820880e062SHong Zhang if (roworiented) { 5830880e062SHong Zhang stepval = (n-1)*bs; 5840880e062SHong Zhang } else { 5850880e062SHong Zhang stepval = (m-1)*bs; 5860880e062SHong Zhang } 5870880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5880880e062SHong Zhang row = im[k]; 5890880e062SHong Zhang if (row < 0) continue; 5902515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 59177431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5920880e062SHong Zhang #endif 5930880e062SHong Zhang rp = aj + ai[row]; 5940880e062SHong Zhang ap = aa + bs2*ai[row]; 5950880e062SHong Zhang rmax = imax[row]; 5960880e062SHong Zhang nrow = ailen[row]; 5970880e062SHong Zhang low = 0; 598818f2c47SBarry Smith high = nrow; 5990880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6000880e062SHong Zhang if (in[l] < 0) continue; 6010880e062SHong Zhang col = in[l]; 6022515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 60377431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 604b1823623SSatish Balay #endif 6050880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6060880e062SHong Zhang if (roworiented) { 6070880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6080880e062SHong Zhang } else { 6090880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6100880e062SHong Zhang } 6117cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 612e2ee6c50SBarry Smith lastcol = col; 6130880e062SHong Zhang while (high-low > 7) { 6140880e062SHong Zhang t = (low+high)/2; 6150880e062SHong Zhang if (rp[t] > col) high = t; 6160880e062SHong Zhang else low = t; 6170880e062SHong Zhang } 6180880e062SHong Zhang for (i=low; i<high; i++) { 6190880e062SHong Zhang if (rp[i] > col) break; 6200880e062SHong Zhang if (rp[i] == col) { 6210880e062SHong Zhang bap = ap + bs2*i; 6220880e062SHong Zhang if (roworiented) { 6230880e062SHong Zhang if (is == ADD_VALUES) { 6240880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6250880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6260880e062SHong Zhang bap[jj] += *value++; 6270880e062SHong Zhang } 6280880e062SHong Zhang } 6290880e062SHong Zhang } else { 6300880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6310880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6320880e062SHong Zhang bap[jj] = *value++; 6330880e062SHong Zhang } 6340880e062SHong Zhang } 6350880e062SHong Zhang } 6360880e062SHong Zhang } else { 6370880e062SHong Zhang if (is == ADD_VALUES) { 6380880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6390880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6400880e062SHong Zhang *bap++ += *value++; 6410880e062SHong Zhang } 6420880e062SHong Zhang } 6430880e062SHong Zhang } else { 6440880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6450880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6460880e062SHong Zhang *bap++ = *value++; 6470880e062SHong Zhang } 6480880e062SHong Zhang } 6490880e062SHong Zhang } 6500880e062SHong Zhang } 6510880e062SHong Zhang goto noinsert2; 6520880e062SHong Zhang } 6530880e062SHong Zhang } 6540880e062SHong Zhang if (nonew == 1) goto noinsert2; 655085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 656421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 657c03d1d03SSatish Balay N = nrow++ - 1; high++; 6580880e062SHong Zhang /* shift up all the later entries in this row */ 6590880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6600880e062SHong Zhang rp[ii+1] = rp[ii]; 6610880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6620880e062SHong Zhang } 6630880e062SHong Zhang if (N >= i) { 6640880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6650880e062SHong Zhang } 6660880e062SHong Zhang rp[i] = col; 6670880e062SHong Zhang bap = ap + bs2*i; 6680880e062SHong Zhang if (roworiented) { 6690880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6700880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6710880e062SHong Zhang bap[jj] = *value++; 6720880e062SHong Zhang } 6730880e062SHong Zhang } 6740880e062SHong Zhang } else { 6750880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6760880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6770880e062SHong Zhang *bap++ = *value++; 6780880e062SHong Zhang } 6790880e062SHong Zhang } 6800880e062SHong Zhang } 6810880e062SHong Zhang noinsert2:; 6820880e062SHong Zhang low = i; 6830880e062SHong Zhang } 6840880e062SHong Zhang ailen[row] = nrow; 6850880e062SHong Zhang } 6860880e062SHong Zhang PetscFunctionReturn(0); 68749b5e25fSSatish Balay } 68849b5e25fSSatish Balay 6894a2ae208SSatish Balay #undef __FUNCT__ 6904a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 691dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 69249b5e25fSSatish Balay { 69349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6946849ba73SBarry Smith PetscErrorCode ierr; 69513f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 696899cda47SBarry Smith PetscInt m = A->rmap.N,*ip,N,*ailen = a->ilen; 69713f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 69849b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 69949b5e25fSSatish Balay 70049b5e25fSSatish Balay PetscFunctionBegin; 70149b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 70249b5e25fSSatish Balay 70349b5e25fSSatish Balay if (m) rmax = ailen[0]; 70449b5e25fSSatish Balay for (i=1; i<mbs; i++) { 70549b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 70649b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 70749b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 70849b5e25fSSatish Balay if (fshift) { 70949b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 71049b5e25fSSatish Balay N = ailen[i]; 71149b5e25fSSatish Balay for (j=0; j<N; j++) { 71249b5e25fSSatish Balay ip[j-fshift] = ip[j]; 71349b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 71449b5e25fSSatish Balay } 71549b5e25fSSatish Balay } 71649b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 71749b5e25fSSatish Balay } 71849b5e25fSSatish Balay if (mbs) { 71949b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 72049b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 72149b5e25fSSatish Balay } 72249b5e25fSSatish Balay /* reset ilen and imax for each row */ 72349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 72449b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 72549b5e25fSSatish Balay } 7266c6c5352SBarry Smith a->nz = ai[mbs]; 72749b5e25fSSatish Balay 728b424e231SHong Zhang /* diagonals may have moved, reset it */ 729b424e231SHong Zhang if (a->diag) { 73013f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 73149b5e25fSSatish Balay } 732899cda47SBarry Smith ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap.N,A->rmap.bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr); 733ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 734ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 73549b5e25fSSatish Balay a->reallocs = 0; 73649b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 73749b5e25fSSatish Balay PetscFunctionReturn(0); 73849b5e25fSSatish Balay } 73949b5e25fSSatish Balay 74049b5e25fSSatish Balay /* 74149b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 74249b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 74349b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 74449b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 74549b5e25fSSatish Balay */ 7464a2ae208SSatish Balay #undef __FUNCT__ 7474a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 74813f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 74949b5e25fSSatish Balay { 75013f74950SBarry Smith PetscInt i,j,k,row; 75149b5e25fSSatish Balay PetscTruth flg; 75249b5e25fSSatish Balay 75349b5e25fSSatish Balay PetscFunctionBegin; 75449b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 75549b5e25fSSatish Balay row = idx[i]; 75649b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 75749b5e25fSSatish Balay sizes[j] = 1; 75849b5e25fSSatish Balay i++; 75949b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 76049b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 76149b5e25fSSatish Balay i++; 76249b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 76349b5e25fSSatish Balay flg = PETSC_TRUE; 76449b5e25fSSatish Balay for (k=1; k<bs; k++) { 76549b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 76649b5e25fSSatish Balay flg = PETSC_FALSE; 76749b5e25fSSatish Balay break; 76849b5e25fSSatish Balay } 76949b5e25fSSatish Balay } 770abc0a331SBarry Smith if (flg) { /* No break in the bs */ 77149b5e25fSSatish Balay sizes[j] = bs; 77249b5e25fSSatish Balay i+= bs; 77349b5e25fSSatish Balay } else { 77449b5e25fSSatish Balay sizes[j] = 1; 77549b5e25fSSatish Balay i++; 77649b5e25fSSatish Balay } 77749b5e25fSSatish Balay } 77849b5e25fSSatish Balay } 77949b5e25fSSatish Balay *bs_max = j; 78049b5e25fSSatish Balay PetscFunctionReturn(0); 78149b5e25fSSatish Balay } 78249b5e25fSSatish Balay 78349b5e25fSSatish Balay 78449b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 78549b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 78649b5e25fSSatish Balay */ 78749b5e25fSSatish Balay 7884a2ae208SSatish Balay #undef __FUNCT__ 7894a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 79013f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 79149b5e25fSSatish Balay { 79249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7936849ba73SBarry Smith PetscErrorCode ierr; 794e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 79513f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 796899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 79713f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 79849b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 79949b5e25fSSatish Balay 80049b5e25fSSatish Balay PetscFunctionBegin; 80149b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 80249b5e25fSSatish Balay row = im[k]; /* row number */ 80349b5e25fSSatish Balay brow = row/bs; /* block row number */ 80449b5e25fSSatish Balay if (row < 0) continue; 8052515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 806899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 80749b5e25fSSatish Balay #endif 80849b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 80949b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 81049b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 81149b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 81249b5e25fSSatish Balay low = 0; 81349b5e25fSSatish Balay 81449b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 81549b5e25fSSatish Balay if (in[l] < 0) continue; 8162515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 817899cda47SBarry Smith if (in[l] >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap.N-1); 81849b5e25fSSatish Balay #endif 81949b5e25fSSatish Balay col = in[l]; 82049b5e25fSSatish Balay bcol = col/bs; /* block col number */ 82149b5e25fSSatish Balay 822941593c8SHong Zhang if (brow > bcol) { 823941593c8SHong Zhang if (a->ignore_ltriangular){ 824941593c8SHong Zhang continue; /* ignore lower triangular values */ 825941593c8SHong Zhang } else { 8264e0d8c25SBarry Smith SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)"); 827941593c8SHong Zhang } 828941593c8SHong Zhang } 829f4989cb3SHong Zhang 83049b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8318549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 83249b5e25fSSatish Balay /* element value a(k,l) */ 83349b5e25fSSatish Balay if (roworiented) { 83449b5e25fSSatish Balay value = v[l + k*n]; 83549b5e25fSSatish Balay } else { 83649b5e25fSSatish Balay value = v[k + l*m]; 83749b5e25fSSatish Balay } 83849b5e25fSSatish Balay 83949b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 8407cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 841e2ee6c50SBarry Smith lastcol = col; 84249b5e25fSSatish Balay while (high-low > 7) { 84349b5e25fSSatish Balay t = (low+high)/2; 84449b5e25fSSatish Balay if (rp[t] > bcol) high = t; 84549b5e25fSSatish Balay else low = t; 84649b5e25fSSatish Balay } 84749b5e25fSSatish Balay for (i=low; i<high; i++) { 84849b5e25fSSatish Balay if (rp[i] > bcol) break; 84949b5e25fSSatish Balay if (rp[i] == bcol) { 85049b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 85149b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 85249b5e25fSSatish Balay else *bap = value; 8538549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8548549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8558549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8568549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8578549e402SHong Zhang else *bap = value; 8588549e402SHong Zhang } 85949b5e25fSSatish Balay goto noinsert1; 86049b5e25fSSatish Balay } 86149b5e25fSSatish Balay } 86249b5e25fSSatish Balay 86349b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 864085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 865421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 86649b5e25fSSatish Balay 867c03d1d03SSatish Balay N = nrow++ - 1; high++; 86849b5e25fSSatish Balay /* shift up all the later entries in this row */ 86949b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 87049b5e25fSSatish Balay rp[ii+1] = rp[ii]; 87149b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 87249b5e25fSSatish Balay } 87349b5e25fSSatish Balay if (N>=i) { 87449b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 87549b5e25fSSatish Balay } 87649b5e25fSSatish Balay rp[i] = bcol; 87749b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 87849b5e25fSSatish Balay noinsert1:; 87949b5e25fSSatish Balay low = i; 8808549e402SHong Zhang } 88149b5e25fSSatish Balay } /* end of loop over added columns */ 88249b5e25fSSatish Balay ailen[brow] = nrow; 88349b5e25fSSatish Balay } /* end of loop over added rows */ 88449b5e25fSSatish Balay PetscFunctionReturn(0); 88549b5e25fSSatish Balay } 88649b5e25fSSatish Balay 8874a2ae208SSatish Balay #undef __FUNCT__ 8884d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 889dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 89049b5e25fSSatish Balay { 8914ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 89249b5e25fSSatish Balay Mat outA; 893dfbe8321SBarry Smith PetscErrorCode ierr; 894c84f5b01SHong Zhang PetscTruth row_identity; 89549b5e25fSSatish Balay 89649b5e25fSSatish Balay PetscFunctionBegin; 897c84f5b01SHong Zhang if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 898c84f5b01SHong Zhang ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 899c84f5b01SHong Zhang if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 900c84f5b01SHong Zhang if (inA->rmap.bs != 1) SETERRQ1(PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap.bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */ 901c84f5b01SHong Zhang 90249b5e25fSSatish Balay outA = inA; 9031260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 90449b5e25fSSatish Balay 9051a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 90649b5e25fSSatish Balay /* 907c84f5b01SHong Zhang Blocksize < 8 have a special faster factorization/solver 908c84f5b01SHong Zhang for ICC(0) factorization with natural ordering 90949b5e25fSSatish Balay */ 910c84f5b01SHong Zhang switch (inA->rmap.bs){ /* Note: row_identity = PETSC_TRUE! */ 91149b5e25fSSatish Balay case 1: 912c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 9130fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 91407f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 915d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 916759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 917759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 918c84f5b01SHong Zhang ierr = PetscInfo(inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr); 919c84f5b01SHong Zhang break; 92049b5e25fSSatish Balay case 2: 921c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9221a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9231a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 924759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 925759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 926ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 92749b5e25fSSatish Balay break; 92849b5e25fSSatish Balay case 3: 929c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9301a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9311a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 932759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering; 933759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering; 934ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 93549b5e25fSSatish Balay break; 93649b5e25fSSatish Balay case 4: 937c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 9381a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 9391a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 940759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering; 941759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering; 942ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 94349b5e25fSSatish Balay break; 94449b5e25fSSatish Balay case 5: 945c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 9461a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 9471a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 948759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering; 949759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering; 950ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 95149b5e25fSSatish Balay break; 95249b5e25fSSatish Balay case 6: 953c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 9541a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 9551a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 956759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering; 957759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering; 958ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 95949b5e25fSSatish Balay break; 96049b5e25fSSatish Balay case 7: 961c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 9621a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 9631a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 964759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering; 965759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering; 966ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 96749b5e25fSSatish Balay break; 96849b5e25fSSatish Balay default: 969c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 970c84f5b01SHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; 971c84f5b01SHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 972759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 973759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 974c84f5b01SHong Zhang break; 975c84f5b01SHong Zhang } 97649b5e25fSSatish Balay 977c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 978c3122656SLisandro Dalcin if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } 979c84f5b01SHong Zhang a->row = row; 980c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 981c3122656SLisandro Dalcin if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } 982c84f5b01SHong Zhang a->col = row; 983c84f5b01SHong Zhang 984c84f5b01SHong Zhang /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 985c84f5b01SHong Zhang if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 98652e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 98749b5e25fSSatish Balay 98849b5e25fSSatish Balay if (!a->solve_work) { 989c84f5b01SHong Zhang ierr = PetscMalloc((inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 990c84f5b01SHong Zhang ierr = PetscLogObjectMemory(inA,(inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar));CHKERRQ(ierr); 99149b5e25fSSatish Balay } 99249b5e25fSSatish Balay 993af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 99449b5e25fSSatish Balay PetscFunctionReturn(0); 99549b5e25fSSatish Balay } 996950f1e5bSHong Zhang 99749b5e25fSSatish Balay EXTERN_C_BEGIN 9984a2ae208SSatish Balay #undef __FUNCT__ 9994a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1000be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 100149b5e25fSSatish Balay { 1002045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 100313f74950SBarry Smith PetscInt i,nz,n; 100449b5e25fSSatish Balay 100549b5e25fSSatish Balay PetscFunctionBegin; 10066c6c5352SBarry Smith nz = baij->maxnz; 1007899cda47SBarry Smith n = mat->cmap.n; 100849b5e25fSSatish Balay for (i=0; i<nz; i++) { 100949b5e25fSSatish Balay baij->j[i] = indices[i]; 101049b5e25fSSatish Balay } 10116c6c5352SBarry Smith baij->nz = nz; 101249b5e25fSSatish Balay for (i=0; i<n; i++) { 101349b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 101449b5e25fSSatish Balay } 101549b5e25fSSatish Balay PetscFunctionReturn(0); 101649b5e25fSSatish Balay } 101749b5e25fSSatish Balay EXTERN_C_END 101849b5e25fSSatish Balay 10194a2ae208SSatish Balay #undef __FUNCT__ 10204a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 102149b5e25fSSatish Balay /*@ 102219585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 102349b5e25fSSatish Balay in the matrix. 102449b5e25fSSatish Balay 102549b5e25fSSatish Balay Input Parameters: 102619585528SSatish Balay + mat - the SeqSBAIJ matrix 102749b5e25fSSatish Balay - indices - the column indices 102849b5e25fSSatish Balay 102949b5e25fSSatish Balay Level: advanced 103049b5e25fSSatish Balay 103149b5e25fSSatish Balay Notes: 103249b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 103349b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 103449b5e25fSSatish Balay of the MatSetValues() operation. 103549b5e25fSSatish Balay 103649b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1037d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 103849b5e25fSSatish Balay 1039ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 104049b5e25fSSatish Balay 1041ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 104249b5e25fSSatish Balay @*/ 1043be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 104449b5e25fSSatish Balay { 104513f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 104649b5e25fSSatish Balay 104749b5e25fSSatish Balay PetscFunctionBegin; 10484482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 10494482741eSBarry Smith PetscValidPointer(indices,2); 1050c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 105149b5e25fSSatish Balay if (f) { 105249b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 105349b5e25fSSatish Balay } else { 1054e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 105549b5e25fSSatish Balay } 105649b5e25fSSatish Balay PetscFunctionReturn(0); 105749b5e25fSSatish Balay } 105849b5e25fSSatish Balay 10594a2ae208SSatish Balay #undef __FUNCT__ 10603c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 10613c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 10623c896bc6SHong Zhang { 10633c896bc6SHong Zhang PetscErrorCode ierr; 10643c896bc6SHong Zhang 10653c896bc6SHong Zhang PetscFunctionBegin; 10663c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 10673c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 10683c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10693c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 10703c896bc6SHong Zhang 1071899cda47SBarry Smith if (a->i[A->rmap.N] != b->i[B->rmap.N]) { 10723c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 10733c896bc6SHong Zhang } 1074899cda47SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.N])*sizeof(PetscScalar));CHKERRQ(ierr); 10753c896bc6SHong Zhang } else { 1076f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 10773c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1078f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 10793c896bc6SHong Zhang } 10803c896bc6SHong Zhang PetscFunctionReturn(0); 10813c896bc6SHong Zhang } 10823c896bc6SHong Zhang 10833c896bc6SHong Zhang #undef __FUNCT__ 10844a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1085dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1086273d9f13SBarry Smith { 1087dfbe8321SBarry Smith PetscErrorCode ierr; 1088273d9f13SBarry Smith 1089273d9f13SBarry Smith PetscFunctionBegin; 10907edd0491SSatish Balay ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1091273d9f13SBarry Smith PetscFunctionReturn(0); 1092273d9f13SBarry Smith } 1093273d9f13SBarry Smith 1094a6ece127SHong Zhang #undef __FUNCT__ 1095a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1096dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1097a6ece127SHong Zhang { 1098a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1099a6ece127SHong Zhang PetscFunctionBegin; 1100a6ece127SHong Zhang *array = a->a; 1101a6ece127SHong Zhang PetscFunctionReturn(0); 1102a6ece127SHong Zhang } 1103a6ece127SHong Zhang 1104a6ece127SHong Zhang #undef __FUNCT__ 1105a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1106dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1107a6ece127SHong Zhang { 1108a6ece127SHong Zhang PetscFunctionBegin; 1109a6ece127SHong Zhang PetscFunctionReturn(0); 1110a6ece127SHong Zhang } 1111a6ece127SHong Zhang 111242ee4b1aSHong Zhang #include "petscblaslapack.h" 111342ee4b1aSHong Zhang #undef __FUNCT__ 111442ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1115f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 111642ee4b1aSHong Zhang { 111742ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1118dfbe8321SBarry Smith PetscErrorCode ierr; 1119899cda47SBarry Smith PetscInt i,bs=Y->rmap.bs,bs2,j; 11204ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 112142ee4b1aSHong Zhang 112242ee4b1aSHong Zhang PetscFunctionBegin; 112342ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1124f4df32b1SMatthew Knepley PetscScalar alpha = a; 1125f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1126c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1127c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1128c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1129c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1130c537a176SHong Zhang } 1131c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1132c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1133c4319e64SHong Zhang y->XtoY = X; 1134c537a176SHong Zhang } 1135c4319e64SHong Zhang bs2 = bs*bs; 11366c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1137c4319e64SHong Zhang j = 0; 1138c4319e64SHong Zhang while (j < bs2){ 1139f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1140c4319e64SHong Zhang j++; 1141c537a176SHong Zhang } 1142c4319e64SHong Zhang } 1143ae15b995SBarry Smith ierr = PetscInfo3(0,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %G\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));CHKERRQ(ierr); 114442ee4b1aSHong Zhang } else { 1145f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1146f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1147f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 114842ee4b1aSHong Zhang } 114942ee4b1aSHong Zhang PetscFunctionReturn(0); 115042ee4b1aSHong Zhang } 115142ee4b1aSHong Zhang 1152efcf0fc3SBarry Smith #undef __FUNCT__ 1153efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1154dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1155efcf0fc3SBarry Smith { 1156efcf0fc3SBarry Smith PetscFunctionBegin; 1157efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1158efcf0fc3SBarry Smith PetscFunctionReturn(0); 1159efcf0fc3SBarry Smith } 1160efcf0fc3SBarry Smith 1161efcf0fc3SBarry Smith #undef __FUNCT__ 1162efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1163dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1164efcf0fc3SBarry Smith { 1165efcf0fc3SBarry Smith PetscFunctionBegin; 1166efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1167efcf0fc3SBarry Smith PetscFunctionReturn(0); 1168efcf0fc3SBarry Smith } 1169efcf0fc3SBarry Smith 1170efcf0fc3SBarry Smith #undef __FUNCT__ 1171efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1172ab5e4463SMatthew Knepley PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1173efcf0fc3SBarry Smith { 1174efcf0fc3SBarry Smith PetscFunctionBegin; 1175efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1176efcf0fc3SBarry Smith PetscFunctionReturn(0); 1177efcf0fc3SBarry Smith } 1178efcf0fc3SBarry Smith 117999cafbc1SBarry Smith #undef __FUNCT__ 118099cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 118199cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 118299cafbc1SBarry Smith { 118399cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 118499cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 118599cafbc1SBarry Smith PetscScalar *aa = a->a; 118699cafbc1SBarry Smith 118799cafbc1SBarry Smith PetscFunctionBegin; 118899cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 118999cafbc1SBarry Smith PetscFunctionReturn(0); 119099cafbc1SBarry Smith } 119199cafbc1SBarry Smith 119299cafbc1SBarry Smith #undef __FUNCT__ 119399cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 119499cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 119599cafbc1SBarry Smith { 119699cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 119799cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 119899cafbc1SBarry Smith PetscScalar *aa = a->a; 119999cafbc1SBarry Smith 120099cafbc1SBarry Smith PetscFunctionBegin; 120199cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 120299cafbc1SBarry Smith PetscFunctionReturn(0); 120399cafbc1SBarry Smith } 120499cafbc1SBarry Smith 120549b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 120649b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 120749b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 120849b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 120949b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 121097304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1211431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1212e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 121349b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 121449b5e25fSSatish Balay 0, 121549b5e25fSSatish Balay 0, 121697304618SKris Buschelman /*10*/ 0, 121749b5e25fSSatish Balay 0, 12185f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1219d06b337dSHong Zhang MatRelax_SeqSBAIJ, 122049b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 122197304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 122249b5e25fSSatish Balay MatEqual_SeqSBAIJ, 122349b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 122449b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 122549b5e25fSSatish Balay MatNorm_SeqSBAIJ, 122697304618SKris Buschelman /*20*/ 0, 122749b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 122849b5e25fSSatish Balay 0, 122949b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 123049b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1231dcf5cc72SBarry Smith /*25*/ 0, 123249b5e25fSSatish Balay 0, 123349b5e25fSSatish Balay 0, 12345f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12355f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 123697304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1237c464158bSHong Zhang 0, 12384d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1239a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1240a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 124197304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 12429495ac64SHong Zhang MatForwardSolve_SeqSBAIJ_N, 12439495ac64SHong Zhang MatBackwardSolve_SeqSBAIJ_N, 124449b5e25fSSatish Balay 0, 1245c84f5b01SHong Zhang MatICCFactor_SeqSBAIJ, 124697304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 124749b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 124849b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 124949b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 12503c896bc6SHong Zhang MatCopy_SeqSBAIJ, 12518c07d4e3SBarry Smith /*45*/ 0, 125249b5e25fSSatish Balay MatScale_SeqSBAIJ, 125349b5e25fSSatish Balay 0, 125449b5e25fSSatish Balay 0, 125549b5e25fSSatish Balay 0, 1256521d7252SBarry Smith /*50*/ 0, 125749b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 125849b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 125949b5e25fSSatish Balay 0, 126049b5e25fSSatish Balay 0, 126197304618SKris Buschelman /*55*/ 0, 126249b5e25fSSatish Balay 0, 126349b5e25fSSatish Balay 0, 126449b5e25fSSatish Balay 0, 126549b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 126697304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 126749b5e25fSSatish Balay 0, 126849b5e25fSSatish Balay 0, 1269357abbc8SBarry Smith 0, 1270d959ec07SHong Zhang 0, 127197304618SKris Buschelman /*65*/ 0, 1272d959ec07SHong Zhang 0, 1273d959ec07SHong Zhang 0, 1274d959ec07SHong Zhang 0, 1275d959ec07SHong Zhang 0, 1276985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqSBAIJ, 12773e0d88b5SBarry Smith 0, 12783e0d88b5SBarry Smith 0, 12793e0d88b5SBarry Smith 0, 12803e0d88b5SBarry Smith 0, 128197304618SKris Buschelman /*75*/ 0, 12823e0d88b5SBarry Smith 0, 12833e0d88b5SBarry Smith 0, 12843e0d88b5SBarry Smith 0, 12853e0d88b5SBarry Smith 0, 128697304618SKris Buschelman /*80*/ 0, 12873e0d88b5SBarry Smith 0, 12883e0d88b5SBarry Smith 0, 12893e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 129097304618SKris Buschelman MatGetInertia_SeqSBAIJ, 12913e0d88b5SBarry Smith #else 129297304618SKris Buschelman 0, 12933e0d88b5SBarry Smith #endif 1294865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1295865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1296865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1297efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1298865e5f61SKris Buschelman 0, 1299865e5f61SKris Buschelman 0, 1300865e5f61SKris Buschelman /*90*/ 0, 1301865e5f61SKris Buschelman 0, 1302865e5f61SKris Buschelman 0, 1303865e5f61SKris Buschelman 0, 1304865e5f61SKris Buschelman 0, 1305865e5f61SKris Buschelman /*95*/ 0, 1306865e5f61SKris Buschelman 0, 1307865e5f61SKris Buschelman 0, 130899cafbc1SBarry Smith 0, 130999cafbc1SBarry Smith 0, 131099cafbc1SBarry Smith /*100*/0, 131199cafbc1SBarry Smith 0, 131299cafbc1SBarry Smith 0, 131399cafbc1SBarry Smith 0, 131499cafbc1SBarry Smith 0, 131599cafbc1SBarry Smith /*105*/0, 131699cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1317f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1318f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 1319f5edf698SHong Zhang MatRestoreRowUpperTriangular_SeqSBAIJ 132099cafbc1SBarry Smith }; 1321be1d678aSKris Buschelman 132249b5e25fSSatish Balay EXTERN_C_BEGIN 13234a2ae208SSatish Balay #undef __FUNCT__ 13244a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1325be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 132649b5e25fSSatish Balay { 13274afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1328899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 1329dfbe8321SBarry Smith PetscErrorCode ierr; 133049b5e25fSSatish Balay 133149b5e25fSSatish Balay PetscFunctionBegin; 133249b5e25fSSatish Balay if (aij->nonew != 1) { 1333*512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 133449b5e25fSSatish Balay } 133549b5e25fSSatish Balay 133649b5e25fSSatish Balay /* allocate space for values if not already there */ 133749b5e25fSSatish Balay if (!aij->saved_values) { 133887828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 133949b5e25fSSatish Balay } 134049b5e25fSSatish Balay 134149b5e25fSSatish Balay /* copy values over */ 134287828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 134349b5e25fSSatish Balay PetscFunctionReturn(0); 134449b5e25fSSatish Balay } 134549b5e25fSSatish Balay EXTERN_C_END 134649b5e25fSSatish Balay 134749b5e25fSSatish Balay EXTERN_C_BEGIN 13484a2ae208SSatish Balay #undef __FUNCT__ 13494a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1350be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 135149b5e25fSSatish Balay { 13524afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13536849ba73SBarry Smith PetscErrorCode ierr; 1354899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 135549b5e25fSSatish Balay 135649b5e25fSSatish Balay PetscFunctionBegin; 135749b5e25fSSatish Balay if (aij->nonew != 1) { 1358*512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 135949b5e25fSSatish Balay } 136049b5e25fSSatish Balay if (!aij->saved_values) { 1361e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 136249b5e25fSSatish Balay } 136349b5e25fSSatish Balay 136449b5e25fSSatish Balay /* copy values over */ 136587828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 136649b5e25fSSatish Balay PetscFunctionReturn(0); 136749b5e25fSSatish Balay } 136849b5e25fSSatish Balay EXTERN_C_END 136949b5e25fSSatish Balay 13708549e402SHong Zhang EXTERN_C_BEGIN 13714a2ae208SSatish Balay #undef __FUNCT__ 1372a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1373be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 137449b5e25fSSatish Balay { 1375c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13766849ba73SBarry Smith PetscErrorCode ierr; 1377ab93d7beSBarry Smith PetscInt i,mbs,bs2; 1378ab93d7beSBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg; 137949b5e25fSSatish Balay 138049b5e25fSSatish Balay PetscFunctionBegin; 1381273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1382e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1383899cda47SBarry Smith B->rmap.bs = B->cmap.bs = bs; 13846148ca0dSBarry Smith ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr); 13856148ca0dSBarry Smith ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr); 1386899cda47SBarry Smith 1387899cda47SBarry Smith mbs = B->rmap.N/bs; 138849b5e25fSSatish Balay bs2 = bs*bs; 138949b5e25fSSatish Balay 1390899cda47SBarry Smith if (mbs*bs != B->rmap.N) { 139129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 139249b5e25fSSatish Balay } 139349b5e25fSSatish Balay 1394ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1395ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1396ab93d7beSBarry Smith nz = 0; 1397ab93d7beSBarry Smith } 1398ab93d7beSBarry Smith 1399435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 140077431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 140149b5e25fSSatish Balay if (nnz) { 140249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 140377431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 140477431f27SBarry Smith if (nnz[i] > mbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D rowlength %D",i,nnz[i],mbs); 140549b5e25fSSatish Balay } 140649b5e25fSSatish Balay } 140749b5e25fSSatish Balay 1408e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 140949b5e25fSSatish Balay if (!flg) { 141049b5e25fSSatish Balay switch (bs) { 141149b5e25fSSatish Balay case 1: 14124ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 141349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1414d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1415e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 141649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 141749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1418431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1419431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 142049b5e25fSSatish Balay break; 142149b5e25fSSatish Balay case 2: 14224ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 142349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1424e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 142549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 142649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1427431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1428431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 142949b5e25fSSatish Balay break; 143049b5e25fSSatish Balay case 3: 14315f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 143249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1433e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 143449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 143549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1436431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1437431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 143849b5e25fSSatish Balay break; 143949b5e25fSSatish Balay case 4: 14405f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 144149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1442e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 144349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 144449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1445431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1446431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 144749b5e25fSSatish Balay break; 144849b5e25fSSatish Balay case 5: 14495f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 145049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1451e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 145249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 145349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1454431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1455431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 145649b5e25fSSatish Balay break; 145749b5e25fSSatish Balay case 6: 14585f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 145949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1460e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 146149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 146249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1463431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1464431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 146549b5e25fSSatish Balay break; 146649b5e25fSSatish Balay case 7: 1467de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 146849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1469e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1470de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 147149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1472431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1473431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 147449b5e25fSSatish Balay break; 147549b5e25fSSatish Balay } 147649b5e25fSSatish Balay } 147749b5e25fSSatish Balay 147849b5e25fSSatish Balay b->mbs = mbs; 14794afc71dfSHong Zhang b->nbs = mbs; 1480ab93d7beSBarry Smith if (!skipallocation) { 1481ab93d7beSBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1482ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1483ab93d7beSBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 148449b5e25fSSatish Balay if (!nnz) { 1485435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 148649b5e25fSSatish Balay else if (nz <= 0) nz = 1; 148749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14888cef66ccSHong Zhang b->imax[i] = nz; 148949b5e25fSSatish Balay } 1490153ea458SHong Zhang nz = nz*mbs; /* total nz */ 149149b5e25fSSatish Balay } else { 149249b5e25fSSatish Balay nz = 0; 14938cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 149449b5e25fSSatish Balay } 14956c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 149649b5e25fSSatish Balay 149749b5e25fSSatish Balay /* allocate the matrix space */ 1498899cda47SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.N+1,PetscInt,&b->i);CHKERRQ(ierr); 14996c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 150013f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 150149b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 150249b5e25fSSatish Balay 150349b5e25fSSatish Balay /* pointer to beginning of each row */ 150449b5e25fSSatish Balay b->i[0] = 0; 150549b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 150649b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 150749b5e25fSSatish Balay } 1508e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1509e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1510e811da20SHong Zhang } else { 1511e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1512e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1513ab93d7beSBarry Smith } 151449b5e25fSSatish Balay 1515899cda47SBarry Smith B->rmap.bs = bs; 151649b5e25fSSatish Balay b->bs2 = bs2; 15176c6c5352SBarry Smith b->nz = 0; 15186c6c5352SBarry Smith b->maxnz = nz*bs2; 1519153ea458SHong Zhang 152016cdd363SHong Zhang b->inew = 0; 152116cdd363SHong Zhang b->jnew = 0; 152216cdd363SHong Zhang b->anew = 0; 152316cdd363SHong Zhang b->a2anew = 0; 15241a3463dfSHong Zhang b->permute = PETSC_FALSE; 1525c464158bSHong Zhang PetscFunctionReturn(0); 1526c464158bSHong Zhang } 1527a23d5eceSKris Buschelman EXTERN_C_END 1528153ea458SHong Zhang 1529d769727bSBarry Smith EXTERN_C_BEGIN 1530f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1531f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1532d769727bSBarry Smith EXTERN_C_END 1533d769727bSBarry Smith 15340bad9183SKris Buschelman /*MC 1535fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15360bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15370bad9183SKris Buschelman 15380bad9183SKris Buschelman Options Database Keys: 15390bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15400bad9183SKris Buschelman 15410bad9183SKris Buschelman Level: beginner 15420bad9183SKris Buschelman 15430bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15440bad9183SKris Buschelman M*/ 15450bad9183SKris Buschelman 1546a23d5eceSKris Buschelman EXTERN_C_BEGIN 1547a23d5eceSKris Buschelman #undef __FUNCT__ 1548a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1549be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1550a23d5eceSKris Buschelman { 1551a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1552dfbe8321SBarry Smith PetscErrorCode ierr; 155313f74950SBarry Smith PetscMPIInt size; 1554941593c8SHong Zhang PetscTruth flg; 1555a23d5eceSKris Buschelman 1556a23d5eceSKris Buschelman PetscFunctionBegin; 1557a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1558a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1559a23d5eceSKris Buschelman 156038f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1561a23d5eceSKris Buschelman B->data = (void*)b; 1562a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1563a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1564a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1565a23d5eceSKris Buschelman B->factor = 0; 1566a23d5eceSKris Buschelman B->mapping = 0; 1567a23d5eceSKris Buschelman b->row = 0; 1568a23d5eceSKris Buschelman b->icol = 0; 1569a23d5eceSKris Buschelman b->reallocs = 0; 1570a23d5eceSKris Buschelman b->saved_values = 0; 1571a23d5eceSKris Buschelman 1572a23d5eceSKris Buschelman 1573a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1574a23d5eceSKris Buschelman b->nonew = 0; 1575a23d5eceSKris Buschelman b->diag = 0; 1576a23d5eceSKris Buschelman b->solve_work = 0; 1577a23d5eceSKris Buschelman b->mult_work = 0; 1578a23d5eceSKris Buschelman B->spptr = 0; 1579a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1580a23d5eceSKris Buschelman b->xtoy = 0; 1581a23d5eceSKris Buschelman b->XtoY = 0; 1582a23d5eceSKris Buschelman 1583a23d5eceSKris Buschelman b->inew = 0; 1584a23d5eceSKris Buschelman b->jnew = 0; 1585a23d5eceSKris Buschelman b->anew = 0; 1586a23d5eceSKris Buschelman b->a2anew = 0; 1587a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1588a23d5eceSKris Buschelman 1589941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 1590941593c8SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr); 1591941593c8SHong Zhang if (flg) b->ignore_ltriangular = PETSC_TRUE; 1592941593c8SHong Zhang 1593f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 1594f5edf698SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr); 1595f5edf698SHong Zhang if (flg) b->getrow_utriangular = PETSC_TRUE; 1596f5edf698SHong Zhang 1597a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1598a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1599a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1600a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1601a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1602a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1603a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1604a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1605a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16064e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16074e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16084e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1609a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1610a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1611a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1612a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1613a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1614a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 161523ce1328SBarry Smith 161623ce1328SBarry Smith B->symmetric = PETSC_TRUE; 161723ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 161823ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 161923ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 162017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 1621a23d5eceSKris Buschelman PetscFunctionReturn(0); 1622a23d5eceSKris Buschelman } 1623a23d5eceSKris Buschelman EXTERN_C_END 1624a23d5eceSKris Buschelman 1625a23d5eceSKris Buschelman #undef __FUNCT__ 1626a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1627a23d5eceSKris Buschelman /*@C 1628a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1629a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1630a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1631a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1632a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1633a23d5eceSKris Buschelman 1634a23d5eceSKris Buschelman Collective on Mat 1635a23d5eceSKris Buschelman 1636a23d5eceSKris Buschelman Input Parameters: 1637a23d5eceSKris Buschelman + A - the symmetric matrix 1638a23d5eceSKris Buschelman . bs - size of block 1639a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1640a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1641a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1642a23d5eceSKris Buschelman 1643a23d5eceSKris Buschelman Options Database Keys: 1644a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1645a23d5eceSKris Buschelman block calculations (much slower) 1646a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1647a23d5eceSKris Buschelman 1648a23d5eceSKris Buschelman Level: intermediate 1649a23d5eceSKris Buschelman 1650a23d5eceSKris Buschelman Notes: 1651a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1652a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1653a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1654a23d5eceSKris Buschelman matrices. 1655a23d5eceSKris Buschelman 165649a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 165749a6f317SBarry Smith 165849a6f317SBarry Smith 1659a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1660a23d5eceSKris Buschelman @*/ 1661be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 166213f74950SBarry Smith { 166313f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1664a23d5eceSKris Buschelman 1665a23d5eceSKris Buschelman PetscFunctionBegin; 1666a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1667a23d5eceSKris Buschelman if (f) { 1668a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1669a23d5eceSKris Buschelman } 1670a23d5eceSKris Buschelman PetscFunctionReturn(0); 1671a23d5eceSKris Buschelman } 167249b5e25fSSatish Balay 16734a2ae208SSatish Balay #undef __FUNCT__ 16744a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1675c464158bSHong Zhang /*@C 1676c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1677c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1678c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1679c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1680c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 168149b5e25fSSatish Balay 1682c464158bSHong Zhang Collective on MPI_Comm 1683c464158bSHong Zhang 1684c464158bSHong Zhang Input Parameters: 1685c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1686c464158bSHong Zhang . bs - size of block 1687c464158bSHong Zhang . m - number of rows, or number of columns 1688c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1689744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1690744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1691c464158bSHong Zhang 1692c464158bSHong Zhang Output Parameter: 1693c464158bSHong Zhang . A - the symmetric matrix 1694c464158bSHong Zhang 1695c464158bSHong Zhang Options Database Keys: 1696c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1697c464158bSHong Zhang block calculations (much slower) 1698c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1699c464158bSHong Zhang 1700c464158bSHong Zhang Level: intermediate 1701c464158bSHong Zhang 1702c464158bSHong Zhang Notes: 17036d6d819aSHong Zhang The number of rows and columns must be divisible by blocksize. 17046d6d819aSHong Zhang This matrix type does not support complex Hermitian operation. 1705c464158bSHong Zhang 1706c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1707c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1708c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1709c464158bSHong Zhang matrices. 1710c464158bSHong Zhang 171149a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 171249a6f317SBarry Smith 1713c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1714c464158bSHong Zhang @*/ 1715be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1716c464158bSHong Zhang { 1717dfbe8321SBarry Smith PetscErrorCode ierr; 1718c464158bSHong Zhang 1719c464158bSHong Zhang PetscFunctionBegin; 1720f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1721f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1722c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1723ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 172449b5e25fSSatish Balay PetscFunctionReturn(0); 172549b5e25fSSatish Balay } 172649b5e25fSSatish Balay 17274a2ae208SSatish Balay #undef __FUNCT__ 17284a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1729dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 173049b5e25fSSatish Balay { 173149b5e25fSSatish Balay Mat C; 173249b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17336849ba73SBarry Smith PetscErrorCode ierr; 1734b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 173549b5e25fSSatish Balay 173649b5e25fSSatish Balay PetscFunctionBegin; 173729bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 173849b5e25fSSatish Balay 173949b5e25fSSatish Balay *B = 0; 1740f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); 1741899cda47SBarry Smith ierr = MatSetSizes(C,A->rmap.N,A->cmap.n,A->rmap.N,A->cmap.n);CHKERRQ(ierr); 1742be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17431d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1744692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1745692f9cbeSHong Zhang 1746273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 174749b5e25fSSatish Balay C->factor = A->factor; 174849b5e25fSSatish Balay c->row = 0; 174949b5e25fSSatish Balay c->icol = 0; 175049b5e25fSSatish Balay c->saved_values = 0; 175149b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 175249b5e25fSSatish Balay C->assembled = PETSC_TRUE; 175349b5e25fSSatish Balay 1754899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->rmap,&C->rmap);CHKERRQ(ierr); 1755899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->cmap,&C->cmap);CHKERRQ(ierr); 175649b5e25fSSatish Balay c->bs2 = a->bs2; 175749b5e25fSSatish Balay c->mbs = a->mbs; 175849b5e25fSSatish Balay c->nbs = a->nbs; 175949b5e25fSSatish Balay 17608777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 176149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 176249b5e25fSSatish Balay c->imax[i] = a->imax[i]; 176349b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 176449b5e25fSSatish Balay } 176549b5e25fSSatish Balay 176649b5e25fSSatish Balay /* allocate the matrix space */ 1767b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 176849b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 176913f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1770b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 177149b5e25fSSatish Balay if (mbs > 0) { 177213f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 177349b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 177449b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 177549b5e25fSSatish Balay } else { 177649b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 177749b5e25fSSatish Balay } 177849b5e25fSSatish Balay } 177949b5e25fSSatish Balay 178049b5e25fSSatish Balay c->roworiented = a->roworiented; 178149b5e25fSSatish Balay c->nonew = a->nonew; 178249b5e25fSSatish Balay 178349b5e25fSSatish Balay if (a->diag) { 178413f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 178552e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 178649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 178749b5e25fSSatish Balay c->diag[i] = a->diag[i]; 178849b5e25fSSatish Balay } 178949b5e25fSSatish Balay } else c->diag = 0; 17906c6c5352SBarry Smith c->nz = a->nz; 17916c6c5352SBarry Smith c->maxnz = a->maxnz; 179249b5e25fSSatish Balay c->solve_work = 0; 179349b5e25fSSatish Balay c->mult_work = 0; 1794e6b907acSBarry Smith c->free_a = PETSC_TRUE; 1795e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 179649b5e25fSSatish Balay *B = C; 1797b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 179849b5e25fSSatish Balay PetscFunctionReturn(0); 179949b5e25fSSatish Balay } 180049b5e25fSSatish Balay 18014a2ae208SSatish Balay #undef __FUNCT__ 18024a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1803f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A) 180449b5e25fSSatish Balay { 180549b5e25fSSatish Balay Mat_SeqSBAIJ *a; 180649b5e25fSSatish Balay Mat B; 18076849ba73SBarry Smith PetscErrorCode ierr; 180813f74950SBarry Smith int fd; 180913f74950SBarry Smith PetscMPIInt size; 181013f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 181113f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 181213f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 181313f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 181487828ca2SBarry Smith PetscScalar *aa; 181549b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 181649b5e25fSSatish Balay 181749b5e25fSSatish Balay PetscFunctionBegin; 1818b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 181949b5e25fSSatish Balay bs2 = bs*bs; 182049b5e25fSSatish Balay 182149b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 182229bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1823b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 182449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1825552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 182649b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 182749b5e25fSSatish Balay 182849b5e25fSSatish Balay if (header[3] < 0) { 182929bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 183049b5e25fSSatish Balay } 183149b5e25fSSatish Balay 183229bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 183349b5e25fSSatish Balay 183449b5e25fSSatish Balay /* 183549b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 183649b5e25fSSatish Balay divisible by the blocksize 183749b5e25fSSatish Balay */ 183849b5e25fSSatish Balay mbs = M/bs; 183949b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 184049b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 184149b5e25fSSatish Balay else mbs++; 184249b5e25fSSatish Balay if (extra_rows) { 1843ae15b995SBarry Smith ierr = PetscInfo(0,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 184449b5e25fSSatish Balay } 184549b5e25fSSatish Balay 184649b5e25fSSatish Balay /* read in row lengths */ 184713f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 184849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 184949b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 185049b5e25fSSatish Balay 185149b5e25fSSatish Balay /* read in column indices */ 185213f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 185349b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 185449b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 185549b5e25fSSatish Balay 185649b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 185713f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 185813f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 185913f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 186013f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 186149b5e25fSSatish Balay masked = mask + mbs; 186249b5e25fSSatish Balay rowcount = 0; nzcount = 0; 186349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 186449b5e25fSSatish Balay nmask = 0; 186549b5e25fSSatish Balay for (j=0; j<bs; j++) { 186649b5e25fSSatish Balay kmax = rowlengths[rowcount]; 186749b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18682d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 186903630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 187049b5e25fSSatish Balay } 187149b5e25fSSatish Balay rowcount++; 187249b5e25fSSatish Balay } 1873574b2666SHong Zhang s_browlengths[i] += nmask; 1874574b2666SHong Zhang 187549b5e25fSSatish Balay /* zero out the mask elements we set */ 187649b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 187749b5e25fSSatish Balay } 187849b5e25fSSatish Balay 187949b5e25fSSatish Balay /* create our matrix */ 1880f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 1881f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 18829abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 1883ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 188449b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 188549b5e25fSSatish Balay 188649b5e25fSSatish Balay /* set matrix "i" values */ 188749b5e25fSSatish Balay a->i[0] = 0; 188849b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1889574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1890574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 189149b5e25fSSatish Balay } 18926c6c5352SBarry Smith a->nz = a->i[mbs]; 189349b5e25fSSatish Balay 189449b5e25fSSatish Balay /* read in nonzero values */ 189587828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 189649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 189749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 189849b5e25fSSatish Balay 189949b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 190049b5e25fSSatish Balay nzcount = 0; jcount = 0; 190149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 190249b5e25fSSatish Balay nzcountb = nzcount; 190349b5e25fSSatish Balay nmask = 0; 190449b5e25fSSatish Balay for (j=0; j<bs; j++) { 190549b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 190649b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19072d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 190803630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19092d703238SHong Zhang } 19102d703238SHong Zhang } 19112d703238SHong Zhang /* sort the masked values */ 19122d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19132d703238SHong Zhang 19142d703238SHong Zhang /* set "j" values into matrix */ 19152d703238SHong Zhang maskcount = 1; 19162d703238SHong Zhang for (j=0; j<nmask; j++) { 191749b5e25fSSatish Balay a->j[jcount++] = masked[j]; 191849b5e25fSSatish Balay mask[masked[j]] = maskcount++; 191949b5e25fSSatish Balay } 1920574b2666SHong Zhang 192149b5e25fSSatish Balay /* set "a" values into matrix */ 192249b5e25fSSatish Balay ishift = bs2*a->i[i]; 192349b5e25fSSatish Balay for (j=0; j<bs; j++) { 192449b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 192549b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1926574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1927574b2666SHong Zhang if (tmp >= i){ 192849b5e25fSSatish Balay block = mask[tmp] - 1; 192949b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 193049b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1931574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1932574b2666SHong Zhang } 1933574b2666SHong Zhang nzcountb++; 193449b5e25fSSatish Balay } 193549b5e25fSSatish Balay } 193649b5e25fSSatish Balay /* zero out the mask elements we set */ 193749b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 193849b5e25fSSatish Balay } 19396c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 194049b5e25fSSatish Balay 194149b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1942574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 194349b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 194449b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 194549b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 194649b5e25fSSatish Balay 19479abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19489abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 194949b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19509abb65ffSKris Buschelman *A = B; 195149b5e25fSSatish Balay PetscFunctionReturn(0); 195249b5e25fSSatish Balay } 1953574b2666SHong Zhang 1954d06b337dSHong Zhang #undef __FUNCT__ 1955d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 195613f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1957d06b337dSHong Zhang { 1958d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1959d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1960d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19616849ba73SBarry Smith PetscErrorCode ierr; 1962899cda47SBarry Smith PetscInt m=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j; 1963521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 1964d06b337dSHong Zhang 1965d06b337dSHong Zhang PetscFunctionBegin; 196651f519a2SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 196751f519a2SBarry Smith 1968b965ef7fSBarry Smith its = its*lits; 196977431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1970d06b337dSHong Zhang 1971d06b337dSHong Zhang if (bs > 1) 1972d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1973d06b337dSHong Zhang 19741ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1975d06b337dSHong Zhang if (xx != bb) { 19761ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1977d06b337dSHong Zhang } else { 1978d06b337dSHong Zhang b = x; 1979d06b337dSHong Zhang } 1980d06b337dSHong Zhang 1981e2ee2a47SBarry Smith if (!a->relax_work) { 1982e2ee2a47SBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->relax_work);CHKERRQ(ierr); 1983e2ee2a47SBarry Smith } 1984e2ee2a47SBarry Smith t = a->relax_work; 1985d06b337dSHong Zhang 1986d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19872798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1988290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 1989d06b337dSHong Zhang 1990d06b337dSHong Zhang for (i=0; i<m; i++){ 199144706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1992d06b337dSHong Zhang v = aa + ai[i] + 1; 1993d06b337dSHong Zhang vj = aj + ai[i] + 1; 1994d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1995e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 1996d06b337dSHong Zhang x[i] = omega*t[i]/d; 1997d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 1998d06b337dSHong Zhang } 1999d06b337dSHong Zhang } 2000d06b337dSHong Zhang 20012798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 200295d750ceSBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){ 200395d750ceSBarry Smith t = b; 2004d06b337dSHong Zhang } 200595d750ceSBarry Smith 2006d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2007d06b337dSHong Zhang d = *(aa + ai[i]); 2008d06b337dSHong Zhang v = aa + ai[i] + 1; 2009d06b337dSHong Zhang vj = aj + ai[i] + 1; 2010d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2011e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2012d06b337dSHong Zhang sum = t[i]; 2013d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2014d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2015d06b337dSHong Zhang } 201695d750ceSBarry Smith t = a->relax_work; 2017d06b337dSHong Zhang } 2018d06b337dSHong Zhang its--; 2019d06b337dSHong Zhang } 2020d06b337dSHong Zhang 2021d06b337dSHong Zhang while (its--) { 202244706875SHong Zhang /* 202344706875SHong Zhang forward sweep: 202444706875SHong Zhang for i=0,...,m-1: 202544706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 202644706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 202744706875SHong Zhang b = b - x[i]*U^T(i,:); 2028d06b337dSHong Zhang 202944706875SHong Zhang */ 20302798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2031290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2032d06b337dSHong Zhang 2033d06b337dSHong Zhang for (i=0; i<m; i++){ 203444706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2035d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2036d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2037d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2038d06b337dSHong Zhang sum = t[i]; 2039e6b907acSBarry Smith ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2040d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2041d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2042d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2043d06b337dSHong Zhang } 2044d06b337dSHong Zhang } 2045d06b337dSHong Zhang 20462798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 204744706875SHong Zhang /* 204844706875SHong Zhang backward sweep: 204944706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 205044706875SHong Zhang for i=m-1,...,0: 205144706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 205244706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 205344706875SHong Zhang */ 205495d750ceSBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2055290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2056d06b337dSHong Zhang 2057d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2058d06b337dSHong Zhang v = aa + ai[i] + 1; 2059d06b337dSHong Zhang vj = aj + ai[i] + 1; 2060d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2061efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2062e6b907acSBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 2063d06b337dSHong Zhang } 2064d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2065d06b337dSHong Zhang d = *(aa + ai[i]); 2066d06b337dSHong Zhang v = aa + ai[i] + 1; 2067d06b337dSHong Zhang vj = aj + ai[i] + 1; 2068d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2069e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2070d06b337dSHong Zhang sum = t[i]; 2071d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2072d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2073d06b337dSHong Zhang } 2074d06b337dSHong Zhang } 2075d06b337dSHong Zhang } 2076d06b337dSHong Zhang 20771ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2078d06b337dSHong Zhang if (bb != xx) { 20791ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2080d06b337dSHong Zhang } 2081d06b337dSHong Zhang PetscFunctionReturn(0); 2082d06b337dSHong Zhang } 2083d06b337dSHong Zhang 2084c75a6043SHong Zhang #undef __FUNCT__ 2085c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2086c75a6043SHong Zhang /*@ 2087c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2088c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2089c75a6043SHong Zhang 2090c75a6043SHong Zhang Collective on MPI_Comm 2091c75a6043SHong Zhang 2092c75a6043SHong Zhang Input Parameters: 2093c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2094c75a6043SHong Zhang . bs - size of block 2095c75a6043SHong Zhang . m - number of rows 2096c75a6043SHong Zhang . n - number of columns 2097c75a6043SHong Zhang . i - row indices 2098c75a6043SHong Zhang . j - column indices 2099c75a6043SHong Zhang - a - matrix values 2100c75a6043SHong Zhang 2101c75a6043SHong Zhang Output Parameter: 2102c75a6043SHong Zhang . mat - the matrix 2103c75a6043SHong Zhang 2104c75a6043SHong Zhang Level: intermediate 2105c75a6043SHong Zhang 2106c75a6043SHong Zhang Notes: 2107c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2108c75a6043SHong Zhang once the matrix is destroyed 2109c75a6043SHong Zhang 2110c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2111c75a6043SHong Zhang 2112c75a6043SHong Zhang The i and j indices are 0 based 2113c75a6043SHong Zhang 2114c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2115c75a6043SHong Zhang 2116c75a6043SHong Zhang @*/ 2117c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2118c75a6043SHong Zhang { 2119c75a6043SHong Zhang PetscErrorCode ierr; 2120c75a6043SHong Zhang PetscInt ii; 2121c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2122c75a6043SHong Zhang 2123c75a6043SHong Zhang PetscFunctionBegin; 2124c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2125c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2126c75a6043SHong Zhang 2127c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2128c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2129c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2130c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2131c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2132c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2133c75a6043SHong Zhang 2134c75a6043SHong Zhang sbaij->i = i; 2135c75a6043SHong Zhang sbaij->j = j; 2136c75a6043SHong Zhang sbaij->a = a; 2137c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2138c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2139e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2140e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2141c75a6043SHong Zhang 2142c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2143c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2144c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2145c75a6043SHong Zhang if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]); 2146c75a6043SHong Zhang #endif 2147c75a6043SHong Zhang } 2148c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2149c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2150c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2151c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2152c75a6043SHong Zhang } 2153c75a6043SHong Zhang #endif 2154c75a6043SHong Zhang 2155c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2156c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2157c75a6043SHong Zhang PetscFunctionReturn(0); 2158c75a6043SHong Zhang } 2159d06b337dSHong Zhang 2160d06b337dSHong Zhang 2161d06b337dSHong Zhang 216249b5e25fSSatish Balay 216349b5e25fSSatish Balay 2164