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" 145dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op) 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: 1534d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1544d9d31abSKris Buschelman break; 1554d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1564d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1574d9d31abSKris Buschelman break; 1584d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1594d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1604d9d31abSKris Buschelman break; 1614d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1624d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1634d9d31abSKris Buschelman break; 1644d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1654d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1664d9d31abSKris Buschelman break; 1674d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1684d9d31abSKris Buschelman a->nonew = 1; 1694d9d31abSKris Buschelman break; 1704d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1714d9d31abSKris Buschelman a->nonew = -1; 1724d9d31abSKris Buschelman break; 1734d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1744d9d31abSKris Buschelman a->nonew = -2; 1754d9d31abSKris Buschelman break; 1764d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1774d9d31abSKris Buschelman a->nonew = 0; 1784d9d31abSKris Buschelman break; 1794d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1804d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1814d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1824d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1834d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 184290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 1854d9d31abSKris Buschelman break; 1864d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 18729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1889a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1899a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1909a4540c5SBarry Smith case MAT_HERMITIAN: 1919a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 19277e54ba9SKris Buschelman case MAT_SYMMETRIC: 19377e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1949a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1959a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1969a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 197290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);CHKERRQ(ierr); 198290bbb0aSBarry Smith break; 199941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 200941593c8SHong Zhang a->ignore_ltriangular = PETSC_TRUE; 201941593c8SHong Zhang break; 202941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 203941593c8SHong Zhang a->ignore_ltriangular = PETSC_FALSE; 20477e54ba9SKris Buschelman break; 205f5edf698SHong Zhang case MAT_GETROW_UPPERTRIANGULAR: 206f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 207f5edf698SHong Zhang break; 2084d9d31abSKris Buschelman default: 209ad86a440SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 21049b5e25fSSatish Balay } 21149b5e25fSSatish Balay PetscFunctionReturn(0); 21249b5e25fSSatish Balay } 21349b5e25fSSatish Balay 2144a2ae208SSatish Balay #undef __FUNCT__ 2154a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 21613f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 21749b5e25fSSatish Balay { 21849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2196849ba73SBarry Smith PetscErrorCode ierr; 22013f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 22149b5e25fSSatish Balay MatScalar *aa,*aa_i; 22287828ca2SBarry Smith PetscScalar *v_i; 22349b5e25fSSatish Balay 22449b5e25fSSatish Balay PetscFunctionBegin; 225f5edf698SHong Zhang 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) or MatGetRowUpperTriangular()"); 226f5edf698SHong Zhang /* Get the upper triangular part of the row */ 227899cda47SBarry Smith bs = A->rmap.bs; 22849b5e25fSSatish Balay ai = a->i; 22949b5e25fSSatish Balay aj = a->j; 23049b5e25fSSatish Balay aa = a->a; 23149b5e25fSSatish Balay bs2 = a->bs2; 23249b5e25fSSatish Balay 233899cda47SBarry Smith if (row < 0 || row >= A->rmap.N) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 23449b5e25fSSatish Balay 23549b5e25fSSatish Balay bn = row/bs; /* Block number */ 23649b5e25fSSatish Balay bp = row % bs; /* Block position */ 23749b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 23849b5e25fSSatish Balay *ncols = bs*M; 23949b5e25fSSatish Balay 24049b5e25fSSatish Balay if (v) { 24149b5e25fSSatish Balay *v = 0; 24249b5e25fSSatish Balay if (*ncols) { 24387828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 24449b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 24549b5e25fSSatish Balay v_i = *v + i*bs; 24649b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 24749b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 24849b5e25fSSatish Balay } 24949b5e25fSSatish Balay } 25049b5e25fSSatish Balay } 25149b5e25fSSatish Balay 25249b5e25fSSatish Balay if (cols) { 25349b5e25fSSatish Balay *cols = 0; 25449b5e25fSSatish Balay if (*ncols) { 25513f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 25649b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 25749b5e25fSSatish Balay cols_i = *cols + i*bs; 25849b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 25949b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 26049b5e25fSSatish Balay } 26149b5e25fSSatish Balay } 26249b5e25fSSatish Balay } 26349b5e25fSSatish Balay 26449b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2655ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2665ddb2528SHong Zhang #ifdef column_search 26749b5e25fSSatish Balay v_i = *v + M*bs; 26849b5e25fSSatish Balay cols_i = *cols + M*bs; 26949b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 27049b5e25fSSatish Balay M = ai[i+1] - ai[i]; 27149b5e25fSSatish Balay for (j=0; j<M; j++){ 27249b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 27349b5e25fSSatish Balay if (itmp == bn){ 27449b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 27549b5e25fSSatish Balay for (k=0; k<bs; k++) { 27649b5e25fSSatish Balay *cols_i++ = i*bs+k; 27749b5e25fSSatish Balay *v_i++ = aa_i[k]; 27849b5e25fSSatish Balay } 27949b5e25fSSatish Balay *ncols += bs; 28049b5e25fSSatish Balay break; 28149b5e25fSSatish Balay } 28249b5e25fSSatish Balay } 28349b5e25fSSatish Balay } 2845ddb2528SHong Zhang #endif 28549b5e25fSSatish Balay PetscFunctionReturn(0); 28649b5e25fSSatish Balay } 28749b5e25fSSatish Balay 2884a2ae208SSatish Balay #undef __FUNCT__ 2894a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 29013f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 29149b5e25fSSatish Balay { 292dfbe8321SBarry Smith PetscErrorCode ierr; 29349b5e25fSSatish Balay 29449b5e25fSSatish Balay PetscFunctionBegin; 29505b42c5fSBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 29605b42c5fSBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 29749b5e25fSSatish Balay PetscFunctionReturn(0); 29849b5e25fSSatish Balay } 29949b5e25fSSatish Balay 3004a2ae208SSatish Balay #undef __FUNCT__ 301f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ" 302f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 303f5edf698SHong Zhang { 304f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 305f5edf698SHong Zhang 306f5edf698SHong Zhang PetscFunctionBegin; 307f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 308f5edf698SHong Zhang PetscFunctionReturn(0); 309f5edf698SHong Zhang } 310f5edf698SHong Zhang #undef __FUNCT__ 311f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ" 312f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 313f5edf698SHong Zhang { 314f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 315f5edf698SHong Zhang 316f5edf698SHong Zhang PetscFunctionBegin; 317f5edf698SHong Zhang a->getrow_utriangular = PETSC_FALSE; 318f5edf698SHong Zhang PetscFunctionReturn(0); 319f5edf698SHong Zhang } 320f5edf698SHong Zhang 321f5edf698SHong Zhang #undef __FUNCT__ 3224a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 323dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 32449b5e25fSSatish Balay { 325dfbe8321SBarry Smith PetscErrorCode ierr; 32649b5e25fSSatish Balay PetscFunctionBegin; 327999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 3288115998fSBarry Smith PetscFunctionReturn(0); 32949b5e25fSSatish Balay } 33049b5e25fSSatish Balay 3314a2ae208SSatish Balay #undef __FUNCT__ 3324a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3336849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 33449b5e25fSSatish Balay { 33549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 336dfbe8321SBarry Smith PetscErrorCode ierr; 337899cda47SBarry Smith PetscInt i,j,bs = A->rmap.bs,k,l,bs2=a->bs2; 3382dcb1b2aSMatthew Knepley const char *name; 339f3ef73ceSBarry Smith PetscViewerFormat format; 34049b5e25fSSatish Balay 34149b5e25fSSatish Balay PetscFunctionBegin; 34280fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 343b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 344456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 34577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 346fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 34729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 348fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 349b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 35049b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 35149b5e25fSSatish Balay for (j=0; j<bs; j++) { 35277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 35349b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 35449b5e25fSSatish Balay for (l=0; l<bs; l++) { 35549b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 35649b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 357a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 35849b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35949b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 360a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 36149b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36249b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 363a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36449b5e25fSSatish Balay } 36549b5e25fSSatish Balay #else 36649b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 367a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 36849b5e25fSSatish Balay } 36949b5e25fSSatish Balay #endif 37049b5e25fSSatish Balay } 37149b5e25fSSatish Balay } 372b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 37349b5e25fSSatish Balay } 37449b5e25fSSatish Balay } 375b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 376c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 377c1490034SHong Zhang PetscFunctionReturn(0); 37849b5e25fSSatish Balay } else { 379*8608aa04SHong Zhang if (A->factor && bs>1){ 380*8608aa04SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored. MatView_SeqSBAIJ_ASCII() may not display complete or logically correct entries!\n");CHKERRQ(ierr); 381*8608aa04SHong Zhang } 382b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 38349b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 38449b5e25fSSatish Balay for (j=0; j<bs; j++) { 38577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 38649b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 38749b5e25fSSatish Balay for (l=0; l<bs; l++) { 38849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 38949b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 390a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 39149b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39249b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 393a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 39449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39549b5e25fSSatish Balay } else { 396a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39749b5e25fSSatish Balay } 39849b5e25fSSatish Balay #else 399e9f7bc9eSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 40049b5e25fSSatish Balay #endif 40149b5e25fSSatish Balay } 40249b5e25fSSatish Balay } 403b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 40449b5e25fSSatish Balay } 40549b5e25fSSatish Balay } 406b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 40749b5e25fSSatish Balay } 408b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 40949b5e25fSSatish Balay PetscFunctionReturn(0); 41049b5e25fSSatish Balay } 41149b5e25fSSatish Balay 4124a2ae208SSatish Balay #undef __FUNCT__ 4134a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 4146849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 41549b5e25fSSatish Balay { 41649b5e25fSSatish Balay Mat A = (Mat) Aa; 41749b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 4186849ba73SBarry Smith PetscErrorCode ierr; 419899cda47SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap.bs,bs2=a->bs2; 42013f74950SBarry Smith PetscMPIInt rank; 42149b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 42249b5e25fSSatish Balay MatScalar *aa; 42349b5e25fSSatish Balay MPI_Comm comm; 424b0a32e0cSBarry Smith PetscViewer viewer; 42549b5e25fSSatish Balay 42649b5e25fSSatish Balay PetscFunctionBegin; 42749b5e25fSSatish Balay /* 42849b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 42949b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 43049b5e25fSSatish Balay rest should return immediately. 43149b5e25fSSatish Balay */ 43249b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 43349b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 43449b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 43549b5e25fSSatish Balay 43649b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 43749b5e25fSSatish Balay 438b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 439b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 44049b5e25fSSatish Balay 44149b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 442b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 44349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 44449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 445899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 44649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 44749b5e25fSSatish Balay aa = a->a + j*bs2; 44849b5e25fSSatish Balay for (k=0; k<bs; k++) { 44949b5e25fSSatish Balay for (l=0; l<bs; l++) { 45049b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 451b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 45249b5e25fSSatish Balay } 45349b5e25fSSatish Balay } 45449b5e25fSSatish Balay } 45549b5e25fSSatish Balay } 456b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 45749b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45849b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 459899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 46049b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46149b5e25fSSatish Balay aa = a->a + j*bs2; 46249b5e25fSSatish Balay for (k=0; k<bs; k++) { 46349b5e25fSSatish Balay for (l=0; l<bs; l++) { 46449b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 465b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46649b5e25fSSatish Balay } 46749b5e25fSSatish Balay } 46849b5e25fSSatish Balay } 46949b5e25fSSatish Balay } 47049b5e25fSSatish Balay 471b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 47249b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 47349b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 474899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 47549b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47649b5e25fSSatish Balay aa = a->a + j*bs2; 47749b5e25fSSatish Balay for (k=0; k<bs; k++) { 47849b5e25fSSatish Balay for (l=0; l<bs; l++) { 47949b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 480b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 48149b5e25fSSatish Balay } 48249b5e25fSSatish Balay } 48349b5e25fSSatish Balay } 48449b5e25fSSatish Balay } 48549b5e25fSSatish Balay PetscFunctionReturn(0); 48649b5e25fSSatish Balay } 48749b5e25fSSatish Balay 4884a2ae208SSatish Balay #undef __FUNCT__ 4894a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4906849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 49149b5e25fSSatish Balay { 492dfbe8321SBarry Smith PetscErrorCode ierr; 49349b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 494b0a32e0cSBarry Smith PetscDraw draw; 49549b5e25fSSatish Balay PetscTruth isnull; 49649b5e25fSSatish Balay 49749b5e25fSSatish Balay PetscFunctionBegin; 498b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 499b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 50049b5e25fSSatish Balay 50149b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 502899cda47SBarry Smith xr = A->rmap.N; yr = A->rmap.N; h = yr/10.0; w = xr/10.0; 50349b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 504b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 505b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 50649b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 50749b5e25fSSatish Balay PetscFunctionReturn(0); 50849b5e25fSSatish Balay } 50949b5e25fSSatish Balay 5104a2ae208SSatish Balay #undef __FUNCT__ 5114a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 512dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 51349b5e25fSSatish Balay { 514dfbe8321SBarry Smith PetscErrorCode ierr; 51532077d6dSBarry Smith PetscTruth iascii,isdraw; 51649b5e25fSSatish Balay 51749b5e25fSSatish Balay PetscFunctionBegin; 51832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 519fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 52032077d6dSBarry Smith if (iascii){ 52149b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 52249b5e25fSSatish Balay } else if (isdraw) { 52349b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 52449b5e25fSSatish Balay } else { 525a5e6ed63SBarry Smith Mat B; 526ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 527a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 528a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 52949b5e25fSSatish Balay } 53049b5e25fSSatish Balay PetscFunctionReturn(0); 53149b5e25fSSatish Balay } 53249b5e25fSSatish Balay 53349b5e25fSSatish Balay 5344a2ae208SSatish Balay #undef __FUNCT__ 5354a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 53613f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 53749b5e25fSSatish Balay { 538045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 53913f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 54013f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 541899cda47SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->rmap.bs,bs2=a->bs2; 54249b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 54349b5e25fSSatish Balay 54449b5e25fSSatish Balay PetscFunctionBegin; 54549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 54649b5e25fSSatish Balay row = im[k]; brow = row/bs; 54777431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 548899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 54949b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 55049b5e25fSSatish Balay nrow = ailen[brow]; 55149b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 55277431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 553899cda47SBarry 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); 55449b5e25fSSatish Balay col = in[l] ; 55549b5e25fSSatish Balay bcol = col/bs; 55649b5e25fSSatish Balay cidx = col%bs; 55749b5e25fSSatish Balay ridx = row%bs; 55849b5e25fSSatish Balay high = nrow; 55949b5e25fSSatish Balay low = 0; /* assume unsorted */ 56049b5e25fSSatish Balay while (high-low > 5) { 56149b5e25fSSatish Balay t = (low+high)/2; 56249b5e25fSSatish Balay if (rp[t] > bcol) high = t; 56349b5e25fSSatish Balay else low = t; 56449b5e25fSSatish Balay } 56549b5e25fSSatish Balay for (i=low; i<high; i++) { 56649b5e25fSSatish Balay if (rp[i] > bcol) break; 56749b5e25fSSatish Balay if (rp[i] == bcol) { 56849b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 56949b5e25fSSatish Balay goto finished; 57049b5e25fSSatish Balay } 57149b5e25fSSatish Balay } 57249b5e25fSSatish Balay *v++ = zero; 57349b5e25fSSatish Balay finished:; 57449b5e25fSSatish Balay } 57549b5e25fSSatish Balay } 57649b5e25fSSatish Balay PetscFunctionReturn(0); 57749b5e25fSSatish Balay } 57849b5e25fSSatish Balay 57949b5e25fSSatish Balay 5804a2ae208SSatish Balay #undef __FUNCT__ 5814a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 58213f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 58349b5e25fSSatish Balay { 5840880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5856849ba73SBarry Smith PetscErrorCode ierr; 586e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 58713f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 588899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap.bs,stepval; 5890880e062SHong Zhang PetscTruth roworiented=a->roworiented; 590f15d580aSBarry Smith const MatScalar *value = v; 591f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5920880e062SHong Zhang 59349b5e25fSSatish Balay PetscFunctionBegin; 5940880e062SHong Zhang if (roworiented) { 5950880e062SHong Zhang stepval = (n-1)*bs; 5960880e062SHong Zhang } else { 5970880e062SHong Zhang stepval = (m-1)*bs; 5980880e062SHong Zhang } 5990880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 6000880e062SHong Zhang row = im[k]; 6010880e062SHong Zhang if (row < 0) continue; 6022515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 60377431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 6040880e062SHong Zhang #endif 6050880e062SHong Zhang rp = aj + ai[row]; 6060880e062SHong Zhang ap = aa + bs2*ai[row]; 6070880e062SHong Zhang rmax = imax[row]; 6080880e062SHong Zhang nrow = ailen[row]; 6090880e062SHong Zhang low = 0; 610818f2c47SBarry Smith high = nrow; 6110880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6120880e062SHong Zhang if (in[l] < 0) continue; 6130880e062SHong Zhang col = in[l]; 6142515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 61577431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 616b1823623SSatish Balay #endif 6170880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6180880e062SHong Zhang if (roworiented) { 6190880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6200880e062SHong Zhang } else { 6210880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6220880e062SHong Zhang } 6237cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 624e2ee6c50SBarry Smith lastcol = col; 6250880e062SHong Zhang while (high-low > 7) { 6260880e062SHong Zhang t = (low+high)/2; 6270880e062SHong Zhang if (rp[t] > col) high = t; 6280880e062SHong Zhang else low = t; 6290880e062SHong Zhang } 6300880e062SHong Zhang for (i=low; i<high; i++) { 6310880e062SHong Zhang if (rp[i] > col) break; 6320880e062SHong Zhang if (rp[i] == col) { 6330880e062SHong Zhang bap = ap + bs2*i; 6340880e062SHong Zhang if (roworiented) { 6350880e062SHong Zhang if (is == ADD_VALUES) { 6360880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6370880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6380880e062SHong Zhang bap[jj] += *value++; 6390880e062SHong Zhang } 6400880e062SHong Zhang } 6410880e062SHong Zhang } else { 6420880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6430880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6440880e062SHong Zhang bap[jj] = *value++; 6450880e062SHong Zhang } 6460880e062SHong Zhang } 6470880e062SHong Zhang } 6480880e062SHong Zhang } else { 6490880e062SHong Zhang if (is == ADD_VALUES) { 6500880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6510880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6520880e062SHong Zhang *bap++ += *value++; 6530880e062SHong Zhang } 6540880e062SHong Zhang } 6550880e062SHong Zhang } else { 6560880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6570880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6580880e062SHong Zhang *bap++ = *value++; 6590880e062SHong Zhang } 6600880e062SHong Zhang } 6610880e062SHong Zhang } 6620880e062SHong Zhang } 6630880e062SHong Zhang goto noinsert2; 6640880e062SHong Zhang } 6650880e062SHong Zhang } 6660880e062SHong Zhang if (nonew == 1) goto noinsert2; 667085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 668421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 669c03d1d03SSatish Balay N = nrow++ - 1; high++; 6700880e062SHong Zhang /* shift up all the later entries in this row */ 6710880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6720880e062SHong Zhang rp[ii+1] = rp[ii]; 6730880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6740880e062SHong Zhang } 6750880e062SHong Zhang if (N >= i) { 6760880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6770880e062SHong Zhang } 6780880e062SHong Zhang rp[i] = col; 6790880e062SHong Zhang bap = ap + bs2*i; 6800880e062SHong Zhang if (roworiented) { 6810880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6820880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6830880e062SHong Zhang bap[jj] = *value++; 6840880e062SHong Zhang } 6850880e062SHong Zhang } 6860880e062SHong Zhang } else { 6870880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6880880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6890880e062SHong Zhang *bap++ = *value++; 6900880e062SHong Zhang } 6910880e062SHong Zhang } 6920880e062SHong Zhang } 6930880e062SHong Zhang noinsert2:; 6940880e062SHong Zhang low = i; 6950880e062SHong Zhang } 6960880e062SHong Zhang ailen[row] = nrow; 6970880e062SHong Zhang } 6980880e062SHong Zhang PetscFunctionReturn(0); 69949b5e25fSSatish Balay } 70049b5e25fSSatish Balay 7014a2ae208SSatish Balay #undef __FUNCT__ 7024a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 703dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 70449b5e25fSSatish Balay { 70549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7066849ba73SBarry Smith PetscErrorCode ierr; 70713f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 708899cda47SBarry Smith PetscInt m = A->rmap.N,*ip,N,*ailen = a->ilen; 70913f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 71049b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 71149b5e25fSSatish Balay 71249b5e25fSSatish Balay PetscFunctionBegin; 71349b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 71449b5e25fSSatish Balay 71549b5e25fSSatish Balay if (m) rmax = ailen[0]; 71649b5e25fSSatish Balay for (i=1; i<mbs; i++) { 71749b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 71849b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 71949b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 72049b5e25fSSatish Balay if (fshift) { 72149b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 72249b5e25fSSatish Balay N = ailen[i]; 72349b5e25fSSatish Balay for (j=0; j<N; j++) { 72449b5e25fSSatish Balay ip[j-fshift] = ip[j]; 72549b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 72649b5e25fSSatish Balay } 72749b5e25fSSatish Balay } 72849b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 72949b5e25fSSatish Balay } 73049b5e25fSSatish Balay if (mbs) { 73149b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 73249b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 73349b5e25fSSatish Balay } 73449b5e25fSSatish Balay /* reset ilen and imax for each row */ 73549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 73649b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 73749b5e25fSSatish Balay } 7386c6c5352SBarry Smith a->nz = ai[mbs]; 73949b5e25fSSatish Balay 740b424e231SHong Zhang /* diagonals may have moved, reset it */ 741b424e231SHong Zhang if (a->diag) { 74213f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 74349b5e25fSSatish Balay } 744899cda47SBarry 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); 745ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 746ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 74749b5e25fSSatish Balay a->reallocs = 0; 74849b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 74949b5e25fSSatish Balay PetscFunctionReturn(0); 75049b5e25fSSatish Balay } 75149b5e25fSSatish Balay 75249b5e25fSSatish Balay /* 75349b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 75449b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 75549b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 75649b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 75749b5e25fSSatish Balay */ 7584a2ae208SSatish Balay #undef __FUNCT__ 7594a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 76013f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 76149b5e25fSSatish Balay { 76213f74950SBarry Smith PetscInt i,j,k,row; 76349b5e25fSSatish Balay PetscTruth flg; 76449b5e25fSSatish Balay 76549b5e25fSSatish Balay PetscFunctionBegin; 76649b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 76749b5e25fSSatish Balay row = idx[i]; 76849b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 76949b5e25fSSatish Balay sizes[j] = 1; 77049b5e25fSSatish Balay i++; 77149b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 77249b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 77349b5e25fSSatish Balay i++; 77449b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 77549b5e25fSSatish Balay flg = PETSC_TRUE; 77649b5e25fSSatish Balay for (k=1; k<bs; k++) { 77749b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 77849b5e25fSSatish Balay flg = PETSC_FALSE; 77949b5e25fSSatish Balay break; 78049b5e25fSSatish Balay } 78149b5e25fSSatish Balay } 782abc0a331SBarry Smith if (flg) { /* No break in the bs */ 78349b5e25fSSatish Balay sizes[j] = bs; 78449b5e25fSSatish Balay i+= bs; 78549b5e25fSSatish Balay } else { 78649b5e25fSSatish Balay sizes[j] = 1; 78749b5e25fSSatish Balay i++; 78849b5e25fSSatish Balay } 78949b5e25fSSatish Balay } 79049b5e25fSSatish Balay } 79149b5e25fSSatish Balay *bs_max = j; 79249b5e25fSSatish Balay PetscFunctionReturn(0); 79349b5e25fSSatish Balay } 79449b5e25fSSatish Balay 79549b5e25fSSatish Balay 79649b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 79749b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 79849b5e25fSSatish Balay */ 79949b5e25fSSatish Balay 8004a2ae208SSatish Balay #undef __FUNCT__ 8014a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 80213f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 80349b5e25fSSatish Balay { 80449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8056849ba73SBarry Smith PetscErrorCode ierr; 806e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 80713f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 808899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 80913f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 81049b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 81149b5e25fSSatish Balay 81249b5e25fSSatish Balay PetscFunctionBegin; 81349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 81449b5e25fSSatish Balay row = im[k]; /* row number */ 81549b5e25fSSatish Balay brow = row/bs; /* block row number */ 81649b5e25fSSatish Balay if (row < 0) continue; 8172515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 818899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 81949b5e25fSSatish Balay #endif 82049b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 82149b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 82249b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 82349b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 82449b5e25fSSatish Balay low = 0; 82549b5e25fSSatish Balay 82649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 82749b5e25fSSatish Balay if (in[l] < 0) continue; 8282515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 829899cda47SBarry 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); 83049b5e25fSSatish Balay #endif 83149b5e25fSSatish Balay col = in[l]; 83249b5e25fSSatish Balay bcol = col/bs; /* block col number */ 83349b5e25fSSatish Balay 834941593c8SHong Zhang if (brow > bcol) { 835941593c8SHong Zhang if (a->ignore_ltriangular){ 836941593c8SHong Zhang continue; /* ignore lower triangular values */ 837941593c8SHong Zhang } else { 838941593c8SHong Zhang 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)"); 839941593c8SHong Zhang } 840941593c8SHong Zhang } 841f4989cb3SHong Zhang 84249b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8438549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 84449b5e25fSSatish Balay /* element value a(k,l) */ 84549b5e25fSSatish Balay if (roworiented) { 84649b5e25fSSatish Balay value = v[l + k*n]; 84749b5e25fSSatish Balay } else { 84849b5e25fSSatish Balay value = v[k + l*m]; 84949b5e25fSSatish Balay } 85049b5e25fSSatish Balay 85149b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 8527cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 853e2ee6c50SBarry Smith lastcol = col; 85449b5e25fSSatish Balay while (high-low > 7) { 85549b5e25fSSatish Balay t = (low+high)/2; 85649b5e25fSSatish Balay if (rp[t] > bcol) high = t; 85749b5e25fSSatish Balay else low = t; 85849b5e25fSSatish Balay } 85949b5e25fSSatish Balay for (i=low; i<high; i++) { 86049b5e25fSSatish Balay if (rp[i] > bcol) break; 86149b5e25fSSatish Balay if (rp[i] == bcol) { 86249b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 86349b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 86449b5e25fSSatish Balay else *bap = value; 8658549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8668549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8678549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8688549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8698549e402SHong Zhang else *bap = value; 8708549e402SHong Zhang } 87149b5e25fSSatish Balay goto noinsert1; 87249b5e25fSSatish Balay } 87349b5e25fSSatish Balay } 87449b5e25fSSatish Balay 87549b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 876085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 877421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 87849b5e25fSSatish Balay 879c03d1d03SSatish Balay N = nrow++ - 1; high++; 88049b5e25fSSatish Balay /* shift up all the later entries in this row */ 88149b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 88249b5e25fSSatish Balay rp[ii+1] = rp[ii]; 88349b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 88449b5e25fSSatish Balay } 88549b5e25fSSatish Balay if (N>=i) { 88649b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 88749b5e25fSSatish Balay } 88849b5e25fSSatish Balay rp[i] = bcol; 88949b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 89049b5e25fSSatish Balay noinsert1:; 89149b5e25fSSatish Balay low = i; 8928549e402SHong Zhang } 89349b5e25fSSatish Balay } /* end of loop over added columns */ 89449b5e25fSSatish Balay ailen[brow] = nrow; 89549b5e25fSSatish Balay } /* end of loop over added rows */ 89649b5e25fSSatish Balay PetscFunctionReturn(0); 89749b5e25fSSatish Balay } 89849b5e25fSSatish Balay 8994a2ae208SSatish Balay #undef __FUNCT__ 9004d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 901dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 90249b5e25fSSatish Balay { 9034ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 90449b5e25fSSatish Balay Mat outA; 905dfbe8321SBarry Smith PetscErrorCode ierr; 906c84f5b01SHong Zhang PetscTruth row_identity; 90749b5e25fSSatish Balay 90849b5e25fSSatish Balay PetscFunctionBegin; 909c84f5b01SHong Zhang if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 910c84f5b01SHong Zhang ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 911c84f5b01SHong Zhang if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 912c84f5b01SHong 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()! */ 913c84f5b01SHong Zhang 91449b5e25fSSatish Balay outA = inA; 9151260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 91649b5e25fSSatish Balay 9171a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 91849b5e25fSSatish Balay /* 919c84f5b01SHong Zhang Blocksize < 8 have a special faster factorization/solver 920c84f5b01SHong Zhang for ICC(0) factorization with natural ordering 92149b5e25fSSatish Balay */ 922c84f5b01SHong Zhang switch (inA->rmap.bs){ /* Note: row_identity = PETSC_TRUE! */ 92349b5e25fSSatish Balay case 1: 924c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 9250fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 92607f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 927d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 928759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 929759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 930c84f5b01SHong Zhang ierr = PetscInfo(inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr); 931c84f5b01SHong Zhang break; 93249b5e25fSSatish Balay case 2: 933c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9341a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9351a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 936759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 937759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 938ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 93949b5e25fSSatish Balay break; 94049b5e25fSSatish Balay case 3: 941c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9421a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9431a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 944759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering; 945759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering; 946ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 94749b5e25fSSatish Balay break; 94849b5e25fSSatish Balay case 4: 949c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 9501a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 9511a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 952759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering; 953759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering; 954ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 95549b5e25fSSatish Balay break; 95649b5e25fSSatish Balay case 5: 957c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 9581a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 9591a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 960759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering; 961759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering; 962ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 96349b5e25fSSatish Balay break; 96449b5e25fSSatish Balay case 6: 965c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 9661a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 9671a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 968759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering; 969759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering; 970ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 97149b5e25fSSatish Balay break; 97249b5e25fSSatish Balay case 7: 973c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 9741a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 9751a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 976759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering; 977759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering; 978ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 97949b5e25fSSatish Balay break; 98049b5e25fSSatish Balay default: 981c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 982c84f5b01SHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; 983c84f5b01SHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 984759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 985759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 986c84f5b01SHong Zhang break; 987c84f5b01SHong Zhang } 98849b5e25fSSatish Balay 989c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 990c3122656SLisandro Dalcin if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } 991c84f5b01SHong Zhang a->row = row; 992c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 993c3122656SLisandro Dalcin if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } 994c84f5b01SHong Zhang a->col = row; 995c84f5b01SHong Zhang 996c84f5b01SHong Zhang /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 997c84f5b01SHong Zhang if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 99852e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 99949b5e25fSSatish Balay 100049b5e25fSSatish Balay if (!a->solve_work) { 1001c84f5b01SHong Zhang ierr = PetscMalloc((inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 1002c84f5b01SHong Zhang ierr = PetscLogObjectMemory(inA,(inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar));CHKERRQ(ierr); 100349b5e25fSSatish Balay } 100449b5e25fSSatish Balay 1005af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 100649b5e25fSSatish Balay PetscFunctionReturn(0); 100749b5e25fSSatish Balay } 1008950f1e5bSHong Zhang 100949b5e25fSSatish Balay EXTERN_C_BEGIN 10104a2ae208SSatish Balay #undef __FUNCT__ 10114a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1012be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 101349b5e25fSSatish Balay { 1014045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 101513f74950SBarry Smith PetscInt i,nz,n; 101649b5e25fSSatish Balay 101749b5e25fSSatish Balay PetscFunctionBegin; 10186c6c5352SBarry Smith nz = baij->maxnz; 1019899cda47SBarry Smith n = mat->cmap.n; 102049b5e25fSSatish Balay for (i=0; i<nz; i++) { 102149b5e25fSSatish Balay baij->j[i] = indices[i]; 102249b5e25fSSatish Balay } 10236c6c5352SBarry Smith baij->nz = nz; 102449b5e25fSSatish Balay for (i=0; i<n; i++) { 102549b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 102649b5e25fSSatish Balay } 102749b5e25fSSatish Balay PetscFunctionReturn(0); 102849b5e25fSSatish Balay } 102949b5e25fSSatish Balay EXTERN_C_END 103049b5e25fSSatish Balay 10314a2ae208SSatish Balay #undef __FUNCT__ 10324a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 103349b5e25fSSatish Balay /*@ 103419585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 103549b5e25fSSatish Balay in the matrix. 103649b5e25fSSatish Balay 103749b5e25fSSatish Balay Input Parameters: 103819585528SSatish Balay + mat - the SeqSBAIJ matrix 103949b5e25fSSatish Balay - indices - the column indices 104049b5e25fSSatish Balay 104149b5e25fSSatish Balay Level: advanced 104249b5e25fSSatish Balay 104349b5e25fSSatish Balay Notes: 104449b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 104549b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 104649b5e25fSSatish Balay of the MatSetValues() operation. 104749b5e25fSSatish Balay 104849b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1049d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 105049b5e25fSSatish Balay 1051ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 105249b5e25fSSatish Balay 1053ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 105449b5e25fSSatish Balay @*/ 1055be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 105649b5e25fSSatish Balay { 105713f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 105849b5e25fSSatish Balay 105949b5e25fSSatish Balay PetscFunctionBegin; 10604482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 10614482741eSBarry Smith PetscValidPointer(indices,2); 1062c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 106349b5e25fSSatish Balay if (f) { 106449b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 106549b5e25fSSatish Balay } else { 1066e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 106749b5e25fSSatish Balay } 106849b5e25fSSatish Balay PetscFunctionReturn(0); 106949b5e25fSSatish Balay } 107049b5e25fSSatish Balay 10714a2ae208SSatish Balay #undef __FUNCT__ 10723c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 10733c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 10743c896bc6SHong Zhang { 10753c896bc6SHong Zhang PetscErrorCode ierr; 10763c896bc6SHong Zhang 10773c896bc6SHong Zhang PetscFunctionBegin; 10783c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 10793c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 10803c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10813c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 10823c896bc6SHong Zhang 1083899cda47SBarry Smith if (a->i[A->rmap.N] != b->i[B->rmap.N]) { 10843c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 10853c896bc6SHong Zhang } 1086899cda47SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.N])*sizeof(PetscScalar));CHKERRQ(ierr); 10873c896bc6SHong Zhang } else { 1088f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 10893c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1090f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 10913c896bc6SHong Zhang } 10923c896bc6SHong Zhang PetscFunctionReturn(0); 10933c896bc6SHong Zhang } 10943c896bc6SHong Zhang 10953c896bc6SHong Zhang #undef __FUNCT__ 10964a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1097dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1098273d9f13SBarry Smith { 1099dfbe8321SBarry Smith PetscErrorCode ierr; 1100273d9f13SBarry Smith 1101273d9f13SBarry Smith PetscFunctionBegin; 11027edd0491SSatish Balay ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1103273d9f13SBarry Smith PetscFunctionReturn(0); 1104273d9f13SBarry Smith } 1105273d9f13SBarry Smith 1106a6ece127SHong Zhang #undef __FUNCT__ 1107a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1108dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1109a6ece127SHong Zhang { 1110a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1111a6ece127SHong Zhang PetscFunctionBegin; 1112a6ece127SHong Zhang *array = a->a; 1113a6ece127SHong Zhang PetscFunctionReturn(0); 1114a6ece127SHong Zhang } 1115a6ece127SHong Zhang 1116a6ece127SHong Zhang #undef __FUNCT__ 1117a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1118dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1119a6ece127SHong Zhang { 1120a6ece127SHong Zhang PetscFunctionBegin; 1121a6ece127SHong Zhang PetscFunctionReturn(0); 1122a6ece127SHong Zhang } 1123a6ece127SHong Zhang 112442ee4b1aSHong Zhang #include "petscblaslapack.h" 112542ee4b1aSHong Zhang #undef __FUNCT__ 112642ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1127f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 112842ee4b1aSHong Zhang { 112942ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1130dfbe8321SBarry Smith PetscErrorCode ierr; 1131899cda47SBarry Smith PetscInt i,bs=Y->rmap.bs,bs2,j; 11324ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 113342ee4b1aSHong Zhang 113442ee4b1aSHong Zhang PetscFunctionBegin; 113542ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1136f4df32b1SMatthew Knepley PetscScalar alpha = a; 1137f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1138c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1139c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1140c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1141c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1142c537a176SHong Zhang } 1143c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1144c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1145c4319e64SHong Zhang y->XtoY = X; 1146c537a176SHong Zhang } 1147c4319e64SHong Zhang bs2 = bs*bs; 11486c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1149c4319e64SHong Zhang j = 0; 1150c4319e64SHong Zhang while (j < bs2){ 1151f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1152c4319e64SHong Zhang j++; 1153c537a176SHong Zhang } 1154c4319e64SHong Zhang } 1155ae15b995SBarry 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); 115642ee4b1aSHong Zhang } else { 1157f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1158f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1159f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 116042ee4b1aSHong Zhang } 116142ee4b1aSHong Zhang PetscFunctionReturn(0); 116242ee4b1aSHong Zhang } 116342ee4b1aSHong Zhang 1164efcf0fc3SBarry Smith #undef __FUNCT__ 1165efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1166dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1167efcf0fc3SBarry Smith { 1168efcf0fc3SBarry Smith PetscFunctionBegin; 1169efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1170efcf0fc3SBarry Smith PetscFunctionReturn(0); 1171efcf0fc3SBarry Smith } 1172efcf0fc3SBarry Smith 1173efcf0fc3SBarry Smith #undef __FUNCT__ 1174efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1175dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1176efcf0fc3SBarry Smith { 1177efcf0fc3SBarry Smith PetscFunctionBegin; 1178efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1179efcf0fc3SBarry Smith PetscFunctionReturn(0); 1180efcf0fc3SBarry Smith } 1181efcf0fc3SBarry Smith 1182efcf0fc3SBarry Smith #undef __FUNCT__ 1183efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1184dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1185efcf0fc3SBarry Smith { 1186efcf0fc3SBarry Smith PetscFunctionBegin; 1187efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1188efcf0fc3SBarry Smith PetscFunctionReturn(0); 1189efcf0fc3SBarry Smith } 1190efcf0fc3SBarry Smith 119199cafbc1SBarry Smith #undef __FUNCT__ 119299cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 119399cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 119499cafbc1SBarry Smith { 119599cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 119699cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 119799cafbc1SBarry Smith PetscScalar *aa = a->a; 119899cafbc1SBarry Smith 119999cafbc1SBarry Smith PetscFunctionBegin; 120099cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 120199cafbc1SBarry Smith PetscFunctionReturn(0); 120299cafbc1SBarry Smith } 120399cafbc1SBarry Smith 120499cafbc1SBarry Smith #undef __FUNCT__ 120599cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 120699cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 120799cafbc1SBarry Smith { 120899cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 120999cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 121099cafbc1SBarry Smith PetscScalar *aa = a->a; 121199cafbc1SBarry Smith 121299cafbc1SBarry Smith PetscFunctionBegin; 121399cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 121499cafbc1SBarry Smith PetscFunctionReturn(0); 121599cafbc1SBarry Smith } 121699cafbc1SBarry Smith 121749b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 121849b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 121949b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 122049b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 122149b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 122297304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1223431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1224e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 122549b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 122649b5e25fSSatish Balay 0, 122749b5e25fSSatish Balay 0, 122897304618SKris Buschelman /*10*/ 0, 122949b5e25fSSatish Balay 0, 12305f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1231d06b337dSHong Zhang MatRelax_SeqSBAIJ, 123249b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 123397304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 123449b5e25fSSatish Balay MatEqual_SeqSBAIJ, 123549b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 123649b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 123749b5e25fSSatish Balay MatNorm_SeqSBAIJ, 123897304618SKris Buschelman /*20*/ 0, 123949b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 124049b5e25fSSatish Balay 0, 124149b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 124249b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1243dcf5cc72SBarry Smith /*25*/ 0, 124449b5e25fSSatish Balay 0, 124549b5e25fSSatish Balay 0, 12465f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12475f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 124897304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1249c464158bSHong Zhang 0, 12504d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1251a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1252a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 125397304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 12549495ac64SHong Zhang MatForwardSolve_SeqSBAIJ_N, 12559495ac64SHong Zhang MatBackwardSolve_SeqSBAIJ_N, 125649b5e25fSSatish Balay 0, 1257c84f5b01SHong Zhang MatICCFactor_SeqSBAIJ, 125897304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 125949b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 126049b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 126149b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 12623c896bc6SHong Zhang MatCopy_SeqSBAIJ, 12638c07d4e3SBarry Smith /*45*/ 0, 126449b5e25fSSatish Balay MatScale_SeqSBAIJ, 126549b5e25fSSatish Balay 0, 126649b5e25fSSatish Balay 0, 126749b5e25fSSatish Balay 0, 1268521d7252SBarry Smith /*50*/ 0, 126949b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 127049b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 127149b5e25fSSatish Balay 0, 127249b5e25fSSatish Balay 0, 127397304618SKris Buschelman /*55*/ 0, 127449b5e25fSSatish Balay 0, 127549b5e25fSSatish Balay 0, 127649b5e25fSSatish Balay 0, 127749b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 127897304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 127949b5e25fSSatish Balay 0, 128049b5e25fSSatish Balay 0, 1281357abbc8SBarry Smith 0, 1282d959ec07SHong Zhang 0, 128397304618SKris Buschelman /*65*/ 0, 1284d959ec07SHong Zhang 0, 1285d959ec07SHong Zhang 0, 1286d959ec07SHong Zhang 0, 1287d959ec07SHong Zhang 0, 1288985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqSBAIJ, 12893e0d88b5SBarry Smith 0, 12903e0d88b5SBarry Smith 0, 12913e0d88b5SBarry Smith 0, 12923e0d88b5SBarry Smith 0, 129397304618SKris Buschelman /*75*/ 0, 12943e0d88b5SBarry Smith 0, 12953e0d88b5SBarry Smith 0, 12963e0d88b5SBarry Smith 0, 12973e0d88b5SBarry Smith 0, 129897304618SKris Buschelman /*80*/ 0, 12993e0d88b5SBarry Smith 0, 13003e0d88b5SBarry Smith 0, 13013e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 130297304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13033e0d88b5SBarry Smith #else 130497304618SKris Buschelman 0, 13053e0d88b5SBarry Smith #endif 1306865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1307865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1308865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1309efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1310865e5f61SKris Buschelman 0, 1311865e5f61SKris Buschelman 0, 1312865e5f61SKris Buschelman /*90*/ 0, 1313865e5f61SKris Buschelman 0, 1314865e5f61SKris Buschelman 0, 1315865e5f61SKris Buschelman 0, 1316865e5f61SKris Buschelman 0, 1317865e5f61SKris Buschelman /*95*/ 0, 1318865e5f61SKris Buschelman 0, 1319865e5f61SKris Buschelman 0, 132099cafbc1SBarry Smith 0, 132199cafbc1SBarry Smith 0, 132299cafbc1SBarry Smith /*100*/0, 132399cafbc1SBarry Smith 0, 132499cafbc1SBarry Smith 0, 132599cafbc1SBarry Smith 0, 132699cafbc1SBarry Smith 0, 132799cafbc1SBarry Smith /*105*/0, 132899cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1329f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1330f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 1331f5edf698SHong Zhang MatRestoreRowUpperTriangular_SeqSBAIJ 133299cafbc1SBarry Smith }; 1333be1d678aSKris Buschelman 133449b5e25fSSatish Balay EXTERN_C_BEGIN 13354a2ae208SSatish Balay #undef __FUNCT__ 13364a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1337be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 133849b5e25fSSatish Balay { 13394afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1340899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 1341dfbe8321SBarry Smith PetscErrorCode ierr; 134249b5e25fSSatish Balay 134349b5e25fSSatish Balay PetscFunctionBegin; 134449b5e25fSSatish Balay if (aij->nonew != 1) { 1345e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 134649b5e25fSSatish Balay } 134749b5e25fSSatish Balay 134849b5e25fSSatish Balay /* allocate space for values if not already there */ 134949b5e25fSSatish Balay if (!aij->saved_values) { 135087828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 135149b5e25fSSatish Balay } 135249b5e25fSSatish Balay 135349b5e25fSSatish Balay /* copy values over */ 135487828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 135549b5e25fSSatish Balay PetscFunctionReturn(0); 135649b5e25fSSatish Balay } 135749b5e25fSSatish Balay EXTERN_C_END 135849b5e25fSSatish Balay 135949b5e25fSSatish Balay EXTERN_C_BEGIN 13604a2ae208SSatish Balay #undef __FUNCT__ 13614a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1362be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 136349b5e25fSSatish Balay { 13644afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13656849ba73SBarry Smith PetscErrorCode ierr; 1366899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 136749b5e25fSSatish Balay 136849b5e25fSSatish Balay PetscFunctionBegin; 136949b5e25fSSatish Balay if (aij->nonew != 1) { 1370e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 137149b5e25fSSatish Balay } 137249b5e25fSSatish Balay if (!aij->saved_values) { 1373e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 137449b5e25fSSatish Balay } 137549b5e25fSSatish Balay 137649b5e25fSSatish Balay /* copy values over */ 137787828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 137849b5e25fSSatish Balay PetscFunctionReturn(0); 137949b5e25fSSatish Balay } 138049b5e25fSSatish Balay EXTERN_C_END 138149b5e25fSSatish Balay 13828549e402SHong Zhang EXTERN_C_BEGIN 13834a2ae208SSatish Balay #undef __FUNCT__ 1384a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1385be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 138649b5e25fSSatish Balay { 1387c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13886849ba73SBarry Smith PetscErrorCode ierr; 1389ab93d7beSBarry Smith PetscInt i,mbs,bs2; 1390ab93d7beSBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg; 139149b5e25fSSatish Balay 139249b5e25fSSatish Balay PetscFunctionBegin; 1393273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1394e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1395899cda47SBarry Smith B->rmap.bs = B->cmap.bs = bs; 13966148ca0dSBarry Smith ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr); 13976148ca0dSBarry Smith ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr); 1398899cda47SBarry Smith 1399899cda47SBarry Smith mbs = B->rmap.N/bs; 140049b5e25fSSatish Balay bs2 = bs*bs; 140149b5e25fSSatish Balay 1402899cda47SBarry Smith if (mbs*bs != B->rmap.N) { 140329bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 140449b5e25fSSatish Balay } 140549b5e25fSSatish Balay 1406ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1407ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1408ab93d7beSBarry Smith nz = 0; 1409ab93d7beSBarry Smith } 1410ab93d7beSBarry Smith 1411435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 141277431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 141349b5e25fSSatish Balay if (nnz) { 141449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 141577431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 141677431f27SBarry 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); 141749b5e25fSSatish Balay } 141849b5e25fSSatish Balay } 141949b5e25fSSatish Balay 1420e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 142149b5e25fSSatish Balay if (!flg) { 142249b5e25fSSatish Balay switch (bs) { 142349b5e25fSSatish Balay case 1: 14244ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 142549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1426d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1427e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 142849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 142949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1430431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1431431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 143249b5e25fSSatish Balay break; 143349b5e25fSSatish Balay case 2: 14344ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 143549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1436e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 143749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 143849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1439431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1440431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 144149b5e25fSSatish Balay break; 144249b5e25fSSatish Balay case 3: 14435f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 144449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1445e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 144649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 144749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1448431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1449431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 145049b5e25fSSatish Balay break; 145149b5e25fSSatish Balay case 4: 14525f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 145349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1454e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 145549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 145649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1457431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1458431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 145949b5e25fSSatish Balay break; 146049b5e25fSSatish Balay case 5: 14615f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 146249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1463e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 146449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 146549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1466431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1467431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 146849b5e25fSSatish Balay break; 146949b5e25fSSatish Balay case 6: 14705f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 147149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1472e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 147349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 147449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1475431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1476431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 147749b5e25fSSatish Balay break; 147849b5e25fSSatish Balay case 7: 1479de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 148049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1481e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1482de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 148349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1484431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1485431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 148649b5e25fSSatish Balay break; 148749b5e25fSSatish Balay } 148849b5e25fSSatish Balay } 148949b5e25fSSatish Balay 149049b5e25fSSatish Balay b->mbs = mbs; 14914afc71dfSHong Zhang b->nbs = mbs; 1492ab93d7beSBarry Smith if (!skipallocation) { 1493ab93d7beSBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1494ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1495ab93d7beSBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 149649b5e25fSSatish Balay if (!nnz) { 1497435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 149849b5e25fSSatish Balay else if (nz <= 0) nz = 1; 149949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15008cef66ccSHong Zhang b->imax[i] = nz; 150149b5e25fSSatish Balay } 1502153ea458SHong Zhang nz = nz*mbs; /* total nz */ 150349b5e25fSSatish Balay } else { 150449b5e25fSSatish Balay nz = 0; 15058cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 150649b5e25fSSatish Balay } 15076c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 150849b5e25fSSatish Balay 150949b5e25fSSatish Balay /* allocate the matrix space */ 1510899cda47SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.N+1,PetscInt,&b->i);CHKERRQ(ierr); 15116c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 151213f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 151349b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 151449b5e25fSSatish Balay 151549b5e25fSSatish Balay /* pointer to beginning of each row */ 151649b5e25fSSatish Balay b->i[0] = 0; 151749b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 151849b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 151949b5e25fSSatish Balay } 1520e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1521e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1522e811da20SHong Zhang } else { 1523e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1524e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1525ab93d7beSBarry Smith } 152649b5e25fSSatish Balay 1527899cda47SBarry Smith B->rmap.bs = bs; 152849b5e25fSSatish Balay b->bs2 = bs2; 15296c6c5352SBarry Smith b->nz = 0; 15306c6c5352SBarry Smith b->maxnz = nz*bs2; 1531153ea458SHong Zhang 153216cdd363SHong Zhang b->inew = 0; 153316cdd363SHong Zhang b->jnew = 0; 153416cdd363SHong Zhang b->anew = 0; 153516cdd363SHong Zhang b->a2anew = 0; 15361a3463dfSHong Zhang b->permute = PETSC_FALSE; 1537c464158bSHong Zhang PetscFunctionReturn(0); 1538c464158bSHong Zhang } 1539a23d5eceSKris Buschelman EXTERN_C_END 1540153ea458SHong Zhang 1541d769727bSBarry Smith EXTERN_C_BEGIN 1542f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1543f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1544d769727bSBarry Smith EXTERN_C_END 1545d769727bSBarry Smith 15460bad9183SKris Buschelman /*MC 1547fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15480bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15490bad9183SKris Buschelman 15500bad9183SKris Buschelman Options Database Keys: 15510bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15520bad9183SKris Buschelman 15530bad9183SKris Buschelman Level: beginner 15540bad9183SKris Buschelman 15550bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15560bad9183SKris Buschelman M*/ 15570bad9183SKris Buschelman 1558a23d5eceSKris Buschelman EXTERN_C_BEGIN 1559a23d5eceSKris Buschelman #undef __FUNCT__ 1560a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1561be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1562a23d5eceSKris Buschelman { 1563a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1564dfbe8321SBarry Smith PetscErrorCode ierr; 156513f74950SBarry Smith PetscMPIInt size; 1566941593c8SHong Zhang PetscTruth flg; 1567a23d5eceSKris Buschelman 1568a23d5eceSKris Buschelman PetscFunctionBegin; 1569a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1570a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1571a23d5eceSKris Buschelman 1572a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1573a23d5eceSKris Buschelman B->data = (void*)b; 1574a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1575a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1576a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1577a23d5eceSKris Buschelman B->factor = 0; 1578a23d5eceSKris Buschelman B->mapping = 0; 1579a23d5eceSKris Buschelman b->row = 0; 1580a23d5eceSKris Buschelman b->icol = 0; 1581a23d5eceSKris Buschelman b->reallocs = 0; 1582a23d5eceSKris Buschelman b->saved_values = 0; 1583a23d5eceSKris Buschelman 1584a23d5eceSKris Buschelman 1585a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1586a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1587a23d5eceSKris Buschelman b->nonew = 0; 1588a23d5eceSKris Buschelman b->diag = 0; 1589a23d5eceSKris Buschelman b->solve_work = 0; 1590a23d5eceSKris Buschelman b->mult_work = 0; 1591a23d5eceSKris Buschelman B->spptr = 0; 1592a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1593a23d5eceSKris Buschelman b->xtoy = 0; 1594a23d5eceSKris Buschelman b->XtoY = 0; 1595a23d5eceSKris Buschelman 1596a23d5eceSKris Buschelman b->inew = 0; 1597a23d5eceSKris Buschelman b->jnew = 0; 1598a23d5eceSKris Buschelman b->anew = 0; 1599a23d5eceSKris Buschelman b->a2anew = 0; 1600a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1601a23d5eceSKris Buschelman 1602941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 1603941593c8SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr); 1604941593c8SHong Zhang if (flg) b->ignore_ltriangular = PETSC_TRUE; 1605941593c8SHong Zhang 1606f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 1607f5edf698SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr); 1608f5edf698SHong Zhang if (flg) b->getrow_utriangular = PETSC_TRUE; 1609f5edf698SHong Zhang 1610a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1611a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1612a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1613a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1614a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1615a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1616a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1617a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1618a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16194e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16204e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16214e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1622a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1623a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1624a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1625a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1626a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1627a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 162823ce1328SBarry Smith 162923ce1328SBarry Smith B->symmetric = PETSC_TRUE; 163023ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 163123ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 163223ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 163317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 1634a23d5eceSKris Buschelman PetscFunctionReturn(0); 1635a23d5eceSKris Buschelman } 1636a23d5eceSKris Buschelman EXTERN_C_END 1637a23d5eceSKris Buschelman 1638a23d5eceSKris Buschelman #undef __FUNCT__ 1639a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1640a23d5eceSKris Buschelman /*@C 1641a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1642a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1643a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1644a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1645a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1646a23d5eceSKris Buschelman 1647a23d5eceSKris Buschelman Collective on Mat 1648a23d5eceSKris Buschelman 1649a23d5eceSKris Buschelman Input Parameters: 1650a23d5eceSKris Buschelman + A - the symmetric matrix 1651a23d5eceSKris Buschelman . bs - size of block 1652a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1653a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1654a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1655a23d5eceSKris Buschelman 1656a23d5eceSKris Buschelman Options Database Keys: 1657a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1658a23d5eceSKris Buschelman block calculations (much slower) 1659a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1660a23d5eceSKris Buschelman 1661a23d5eceSKris Buschelman Level: intermediate 1662a23d5eceSKris Buschelman 1663a23d5eceSKris Buschelman Notes: 1664a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1665a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1666a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1667a23d5eceSKris Buschelman matrices. 1668a23d5eceSKris Buschelman 166949a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 167049a6f317SBarry Smith 167149a6f317SBarry Smith 1672a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1673a23d5eceSKris Buschelman @*/ 1674be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 167513f74950SBarry Smith { 167613f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1677a23d5eceSKris Buschelman 1678a23d5eceSKris Buschelman PetscFunctionBegin; 1679a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1680a23d5eceSKris Buschelman if (f) { 1681a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1682a23d5eceSKris Buschelman } 1683a23d5eceSKris Buschelman PetscFunctionReturn(0); 1684a23d5eceSKris Buschelman } 168549b5e25fSSatish Balay 16864a2ae208SSatish Balay #undef __FUNCT__ 16874a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1688c464158bSHong Zhang /*@C 1689c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1690c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1691c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1692c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1693c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 169449b5e25fSSatish Balay 1695c464158bSHong Zhang Collective on MPI_Comm 1696c464158bSHong Zhang 1697c464158bSHong Zhang Input Parameters: 1698c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1699c464158bSHong Zhang . bs - size of block 1700c464158bSHong Zhang . m - number of rows, or number of columns 1701c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1702744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1703744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1704c464158bSHong Zhang 1705c464158bSHong Zhang Output Parameter: 1706c464158bSHong Zhang . A - the symmetric matrix 1707c464158bSHong Zhang 1708c464158bSHong Zhang Options Database Keys: 1709c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1710c464158bSHong Zhang block calculations (much slower) 1711c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1712c464158bSHong Zhang 1713c464158bSHong Zhang Level: intermediate 1714c464158bSHong Zhang 1715c464158bSHong Zhang Notes: 1716c464158bSHong Zhang 1717c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1718c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1719c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1720c464158bSHong Zhang matrices. 1721c464158bSHong Zhang 172249a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 172349a6f317SBarry Smith 1724c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1725c464158bSHong Zhang @*/ 1726be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1727c464158bSHong Zhang { 1728dfbe8321SBarry Smith PetscErrorCode ierr; 1729c464158bSHong Zhang 1730c464158bSHong Zhang PetscFunctionBegin; 1731f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1732f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1733c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1734ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 173549b5e25fSSatish Balay PetscFunctionReturn(0); 173649b5e25fSSatish Balay } 173749b5e25fSSatish Balay 17384a2ae208SSatish Balay #undef __FUNCT__ 17394a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1740dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 174149b5e25fSSatish Balay { 174249b5e25fSSatish Balay Mat C; 174349b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17446849ba73SBarry Smith PetscErrorCode ierr; 1745b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 174649b5e25fSSatish Balay 174749b5e25fSSatish Balay PetscFunctionBegin; 174829bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 174949b5e25fSSatish Balay 175049b5e25fSSatish Balay *B = 0; 1751f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); 1752899cda47SBarry Smith ierr = MatSetSizes(C,A->rmap.N,A->cmap.n,A->rmap.N,A->cmap.n);CHKERRQ(ierr); 1753be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17541d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1755692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1756692f9cbeSHong Zhang 1757273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 175849b5e25fSSatish Balay C->factor = A->factor; 175949b5e25fSSatish Balay c->row = 0; 176049b5e25fSSatish Balay c->icol = 0; 176149b5e25fSSatish Balay c->saved_values = 0; 176249b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 176349b5e25fSSatish Balay C->assembled = PETSC_TRUE; 176449b5e25fSSatish Balay 1765899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->rmap,&C->rmap);CHKERRQ(ierr); 1766899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->cmap,&C->cmap);CHKERRQ(ierr); 176749b5e25fSSatish Balay c->bs2 = a->bs2; 176849b5e25fSSatish Balay c->mbs = a->mbs; 176949b5e25fSSatish Balay c->nbs = a->nbs; 177049b5e25fSSatish Balay 17718777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 177249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 177349b5e25fSSatish Balay c->imax[i] = a->imax[i]; 177449b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 177549b5e25fSSatish Balay } 177649b5e25fSSatish Balay 177749b5e25fSSatish Balay /* allocate the matrix space */ 1778b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 177949b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 178013f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1781b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 178249b5e25fSSatish Balay if (mbs > 0) { 178313f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 178449b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 178549b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 178649b5e25fSSatish Balay } else { 178749b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 178849b5e25fSSatish Balay } 178949b5e25fSSatish Balay } 179049b5e25fSSatish Balay 179149b5e25fSSatish Balay c->sorted = a->sorted; 179249b5e25fSSatish Balay c->roworiented = a->roworiented; 179349b5e25fSSatish Balay c->nonew = a->nonew; 179449b5e25fSSatish Balay 179549b5e25fSSatish Balay if (a->diag) { 179613f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 179752e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 179849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 179949b5e25fSSatish Balay c->diag[i] = a->diag[i]; 180049b5e25fSSatish Balay } 180149b5e25fSSatish Balay } else c->diag = 0; 18026c6c5352SBarry Smith c->nz = a->nz; 18036c6c5352SBarry Smith c->maxnz = a->maxnz; 180449b5e25fSSatish Balay c->solve_work = 0; 180549b5e25fSSatish Balay c->mult_work = 0; 1806e6b907acSBarry Smith c->free_a = PETSC_TRUE; 1807e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 180849b5e25fSSatish Balay *B = C; 1809b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 181049b5e25fSSatish Balay PetscFunctionReturn(0); 181149b5e25fSSatish Balay } 181249b5e25fSSatish Balay 18134a2ae208SSatish Balay #undef __FUNCT__ 18144a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1815f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A) 181649b5e25fSSatish Balay { 181749b5e25fSSatish Balay Mat_SeqSBAIJ *a; 181849b5e25fSSatish Balay Mat B; 18196849ba73SBarry Smith PetscErrorCode ierr; 182013f74950SBarry Smith int fd; 182113f74950SBarry Smith PetscMPIInt size; 182213f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 182313f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 182413f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 182513f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 182687828ca2SBarry Smith PetscScalar *aa; 182749b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 182849b5e25fSSatish Balay 182949b5e25fSSatish Balay PetscFunctionBegin; 1830b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 183149b5e25fSSatish Balay bs2 = bs*bs; 183249b5e25fSSatish Balay 183349b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 183429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1835b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 183649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1837552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 183849b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 183949b5e25fSSatish Balay 184049b5e25fSSatish Balay if (header[3] < 0) { 184129bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 184249b5e25fSSatish Balay } 184349b5e25fSSatish Balay 184429bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 184549b5e25fSSatish Balay 184649b5e25fSSatish Balay /* 184749b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 184849b5e25fSSatish Balay divisible by the blocksize 184949b5e25fSSatish Balay */ 185049b5e25fSSatish Balay mbs = M/bs; 185149b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 185249b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 185349b5e25fSSatish Balay else mbs++; 185449b5e25fSSatish Balay if (extra_rows) { 1855ae15b995SBarry Smith ierr = PetscInfo(0,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 185649b5e25fSSatish Balay } 185749b5e25fSSatish Balay 185849b5e25fSSatish Balay /* read in row lengths */ 185913f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 186049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 186149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 186249b5e25fSSatish Balay 186349b5e25fSSatish Balay /* read in column indices */ 186413f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 186549b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 186649b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 186749b5e25fSSatish Balay 186849b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 186913f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 187013f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 187113f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 187213f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 187349b5e25fSSatish Balay masked = mask + mbs; 187449b5e25fSSatish Balay rowcount = 0; nzcount = 0; 187549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 187649b5e25fSSatish Balay nmask = 0; 187749b5e25fSSatish Balay for (j=0; j<bs; j++) { 187849b5e25fSSatish Balay kmax = rowlengths[rowcount]; 187949b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18802d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 188103630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 188249b5e25fSSatish Balay } 188349b5e25fSSatish Balay rowcount++; 188449b5e25fSSatish Balay } 1885574b2666SHong Zhang s_browlengths[i] += nmask; 1886574b2666SHong Zhang 188749b5e25fSSatish Balay /* zero out the mask elements we set */ 188849b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 188949b5e25fSSatish Balay } 189049b5e25fSSatish Balay 189149b5e25fSSatish Balay /* create our matrix */ 1892f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 1893f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 18949abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 1895ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 189649b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 189749b5e25fSSatish Balay 189849b5e25fSSatish Balay /* set matrix "i" values */ 189949b5e25fSSatish Balay a->i[0] = 0; 190049b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1901574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1902574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 190349b5e25fSSatish Balay } 19046c6c5352SBarry Smith a->nz = a->i[mbs]; 190549b5e25fSSatish Balay 190649b5e25fSSatish Balay /* read in nonzero values */ 190787828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 190849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 190949b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 191049b5e25fSSatish Balay 191149b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 191249b5e25fSSatish Balay nzcount = 0; jcount = 0; 191349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 191449b5e25fSSatish Balay nzcountb = nzcount; 191549b5e25fSSatish Balay nmask = 0; 191649b5e25fSSatish Balay for (j=0; j<bs; j++) { 191749b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 191849b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19192d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 192003630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19212d703238SHong Zhang } 19222d703238SHong Zhang } 19232d703238SHong Zhang /* sort the masked values */ 19242d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19252d703238SHong Zhang 19262d703238SHong Zhang /* set "j" values into matrix */ 19272d703238SHong Zhang maskcount = 1; 19282d703238SHong Zhang for (j=0; j<nmask; j++) { 192949b5e25fSSatish Balay a->j[jcount++] = masked[j]; 193049b5e25fSSatish Balay mask[masked[j]] = maskcount++; 193149b5e25fSSatish Balay } 1932574b2666SHong Zhang 193349b5e25fSSatish Balay /* set "a" values into matrix */ 193449b5e25fSSatish Balay ishift = bs2*a->i[i]; 193549b5e25fSSatish Balay for (j=0; j<bs; j++) { 193649b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 193749b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1938574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1939574b2666SHong Zhang if (tmp >= i){ 194049b5e25fSSatish Balay block = mask[tmp] - 1; 194149b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 194249b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1943574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1944574b2666SHong Zhang } 1945574b2666SHong Zhang nzcountb++; 194649b5e25fSSatish Balay } 194749b5e25fSSatish Balay } 194849b5e25fSSatish Balay /* zero out the mask elements we set */ 194949b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 195049b5e25fSSatish Balay } 19516c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 195249b5e25fSSatish Balay 195349b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1954574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 195549b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 195649b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 195749b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 195849b5e25fSSatish Balay 19599abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19609abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196149b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19629abb65ffSKris Buschelman *A = B; 196349b5e25fSSatish Balay PetscFunctionReturn(0); 196449b5e25fSSatish Balay } 1965574b2666SHong Zhang 1966d06b337dSHong Zhang #undef __FUNCT__ 1967d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 196813f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1969d06b337dSHong Zhang { 1970d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1971d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1972d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19736849ba73SBarry Smith PetscErrorCode ierr; 1974899cda47SBarry Smith PetscInt m=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j; 1975521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 1976d06b337dSHong Zhang 1977d06b337dSHong Zhang PetscFunctionBegin; 197851f519a2SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 197951f519a2SBarry Smith 1980b965ef7fSBarry Smith its = its*lits; 198177431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1982d06b337dSHong Zhang 1983d06b337dSHong Zhang if (bs > 1) 1984d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1985d06b337dSHong Zhang 19861ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1987d06b337dSHong Zhang if (xx != bb) { 19881ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1989d06b337dSHong Zhang } else { 1990d06b337dSHong Zhang b = x; 1991d06b337dSHong Zhang } 1992d06b337dSHong Zhang 1993e2ee2a47SBarry Smith if (!a->relax_work) { 1994e2ee2a47SBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->relax_work);CHKERRQ(ierr); 1995e2ee2a47SBarry Smith } 1996e2ee2a47SBarry Smith t = a->relax_work; 1997d06b337dSHong Zhang 1998d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19992798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2000290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2001d06b337dSHong Zhang 2002d06b337dSHong Zhang for (i=0; i<m; i++){ 200344706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2004d06b337dSHong Zhang v = aa + ai[i] + 1; 2005d06b337dSHong Zhang vj = aj + ai[i] + 1; 2006d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2007e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2008d06b337dSHong Zhang x[i] = omega*t[i]/d; 2009d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 2010d06b337dSHong Zhang } 2011d06b337dSHong Zhang } 2012d06b337dSHong Zhang 20132798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 201495d750ceSBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){ 201595d750ceSBarry Smith t = b; 2016d06b337dSHong Zhang } 201795d750ceSBarry Smith 2018d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2019d06b337dSHong Zhang d = *(aa + ai[i]); 2020d06b337dSHong Zhang v = aa + ai[i] + 1; 2021d06b337dSHong Zhang vj = aj + ai[i] + 1; 2022d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2023e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2024d06b337dSHong Zhang sum = t[i]; 2025d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2026d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2027d06b337dSHong Zhang } 202895d750ceSBarry Smith t = a->relax_work; 2029d06b337dSHong Zhang } 2030d06b337dSHong Zhang its--; 2031d06b337dSHong Zhang } 2032d06b337dSHong Zhang 2033d06b337dSHong Zhang while (its--) { 203444706875SHong Zhang /* 203544706875SHong Zhang forward sweep: 203644706875SHong Zhang for i=0,...,m-1: 203744706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 203844706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 203944706875SHong Zhang b = b - x[i]*U^T(i,:); 2040d06b337dSHong Zhang 204144706875SHong Zhang */ 20422798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2043290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2044d06b337dSHong Zhang 2045d06b337dSHong Zhang for (i=0; i<m; i++){ 204644706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2047d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2048d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2049d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2050d06b337dSHong Zhang sum = t[i]; 2051e6b907acSBarry Smith ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2052d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2053d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2054d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2055d06b337dSHong Zhang } 2056d06b337dSHong Zhang } 2057d06b337dSHong Zhang 20582798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 205944706875SHong Zhang /* 206044706875SHong Zhang backward sweep: 206144706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 206244706875SHong Zhang for i=m-1,...,0: 206344706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 206444706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 206544706875SHong Zhang */ 206695d750ceSBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2067290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2068d06b337dSHong Zhang 2069d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2070d06b337dSHong Zhang v = aa + ai[i] + 1; 2071d06b337dSHong Zhang vj = aj + ai[i] + 1; 2072d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2073efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2074e6b907acSBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 2075d06b337dSHong Zhang } 2076d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2077d06b337dSHong Zhang d = *(aa + ai[i]); 2078d06b337dSHong Zhang v = aa + ai[i] + 1; 2079d06b337dSHong Zhang vj = aj + ai[i] + 1; 2080d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2081e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2082d06b337dSHong Zhang sum = t[i]; 2083d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2084d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2085d06b337dSHong Zhang } 2086d06b337dSHong Zhang } 2087d06b337dSHong Zhang } 2088d06b337dSHong Zhang 20891ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2090d06b337dSHong Zhang if (bb != xx) { 20911ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2092d06b337dSHong Zhang } 2093d06b337dSHong Zhang PetscFunctionReturn(0); 2094d06b337dSHong Zhang } 2095d06b337dSHong Zhang 2096c75a6043SHong Zhang #undef __FUNCT__ 2097c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2098c75a6043SHong Zhang /*@ 2099c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2100c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2101c75a6043SHong Zhang 2102c75a6043SHong Zhang Collective on MPI_Comm 2103c75a6043SHong Zhang 2104c75a6043SHong Zhang Input Parameters: 2105c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2106c75a6043SHong Zhang . bs - size of block 2107c75a6043SHong Zhang . m - number of rows 2108c75a6043SHong Zhang . n - number of columns 2109c75a6043SHong Zhang . i - row indices 2110c75a6043SHong Zhang . j - column indices 2111c75a6043SHong Zhang - a - matrix values 2112c75a6043SHong Zhang 2113c75a6043SHong Zhang Output Parameter: 2114c75a6043SHong Zhang . mat - the matrix 2115c75a6043SHong Zhang 2116c75a6043SHong Zhang Level: intermediate 2117c75a6043SHong Zhang 2118c75a6043SHong Zhang Notes: 2119c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2120c75a6043SHong Zhang once the matrix is destroyed 2121c75a6043SHong Zhang 2122c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2123c75a6043SHong Zhang 2124c75a6043SHong Zhang The i and j indices are 0 based 2125c75a6043SHong Zhang 2126c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2127c75a6043SHong Zhang 2128c75a6043SHong Zhang @*/ 2129c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2130c75a6043SHong Zhang { 2131c75a6043SHong Zhang PetscErrorCode ierr; 2132c75a6043SHong Zhang PetscInt ii; 2133c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2134c75a6043SHong Zhang 2135c75a6043SHong Zhang PetscFunctionBegin; 2136c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2137c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2138c75a6043SHong Zhang 2139c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2140c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2141c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2142c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2143c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2144c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2145c75a6043SHong Zhang 2146c75a6043SHong Zhang sbaij->i = i; 2147c75a6043SHong Zhang sbaij->j = j; 2148c75a6043SHong Zhang sbaij->a = a; 2149c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2150c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2151e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2152e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2153c75a6043SHong Zhang 2154c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2155c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2156c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2157c75a6043SHong 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]); 2158c75a6043SHong Zhang #endif 2159c75a6043SHong Zhang } 2160c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2161c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2162c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2163c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2164c75a6043SHong Zhang } 2165c75a6043SHong Zhang #endif 2166c75a6043SHong Zhang 2167c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2168c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2169c75a6043SHong Zhang PetscFunctionReturn(0); 2170c75a6043SHong Zhang } 2171d06b337dSHong Zhang 2172d06b337dSHong Zhang 2173d06b337dSHong Zhang 217449b5e25fSSatish Balay 217549b5e25fSSatish Balay 2176