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) { 347c9f458caSMatthew Knepley Mat aij; 348c9f458caSMatthew Knepley ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr); 349c9f458caSMatthew Knepley ierr = MatView(aij,viewer);CHKERRQ(ierr); 350c9f458caSMatthew Knepley ierr = MatDestroy(aij);CHKERRQ(ierr); 351fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 352b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 35349b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 35449b5e25fSSatish Balay for (j=0; j<bs; j++) { 35577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 35649b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 35749b5e25fSSatish Balay for (l=0; l<bs; l++) { 35849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 35949b5e25fSSatish Balay 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 (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 363a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 36449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36549b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 366a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36749b5e25fSSatish Balay } 36849b5e25fSSatish Balay #else 36949b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 370a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 37149b5e25fSSatish Balay } 37249b5e25fSSatish Balay #endif 37349b5e25fSSatish Balay } 37449b5e25fSSatish Balay } 375b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 37649b5e25fSSatish Balay } 37749b5e25fSSatish Balay } 378b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 379c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 380c1490034SHong Zhang PetscFunctionReturn(0); 38149b5e25fSSatish Balay } else { 3828608aa04SHong Zhang if (A->factor && bs>1){ 3838608aa04SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored. MatView_SeqSBAIJ_ASCII() may not display complete or logically correct entries!\n");CHKERRQ(ierr); 3848608aa04SHong Zhang } 385b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 38649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 38749b5e25fSSatish Balay for (j=0; j<bs; j++) { 38877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 38949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 39049b5e25fSSatish Balay for (l=0; l<bs; l++) { 39149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 39249b5e25fSSatish Balay 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 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 396a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 39749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39849b5e25fSSatish Balay } else { 399a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40049b5e25fSSatish Balay } 40149b5e25fSSatish Balay #else 402e9f7bc9eSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 40349b5e25fSSatish Balay #endif 40449b5e25fSSatish Balay } 40549b5e25fSSatish Balay } 406b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 40749b5e25fSSatish Balay } 40849b5e25fSSatish Balay } 409b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 41049b5e25fSSatish Balay } 411b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 41249b5e25fSSatish Balay PetscFunctionReturn(0); 41349b5e25fSSatish Balay } 41449b5e25fSSatish Balay 4154a2ae208SSatish Balay #undef __FUNCT__ 4164a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 4176849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 41849b5e25fSSatish Balay { 41949b5e25fSSatish Balay Mat A = (Mat) Aa; 42049b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 4216849ba73SBarry Smith PetscErrorCode ierr; 422899cda47SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap.bs,bs2=a->bs2; 42313f74950SBarry Smith PetscMPIInt rank; 42449b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 42549b5e25fSSatish Balay MatScalar *aa; 42649b5e25fSSatish Balay MPI_Comm comm; 427b0a32e0cSBarry Smith PetscViewer viewer; 42849b5e25fSSatish Balay 42949b5e25fSSatish Balay PetscFunctionBegin; 43049b5e25fSSatish Balay /* 43149b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 43249b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 43349b5e25fSSatish Balay rest should return immediately. 43449b5e25fSSatish Balay */ 43549b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 43649b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 43749b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 43849b5e25fSSatish Balay 43949b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 44049b5e25fSSatish Balay 441b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 442b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 44349b5e25fSSatish Balay 44449b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 445b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 44649b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 44749b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 448899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 44949b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45049b5e25fSSatish Balay aa = a->a + j*bs2; 45149b5e25fSSatish Balay for (k=0; k<bs; k++) { 45249b5e25fSSatish Balay for (l=0; l<bs; l++) { 45349b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 454b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 45549b5e25fSSatish Balay } 45649b5e25fSSatish Balay } 45749b5e25fSSatish Balay } 45849b5e25fSSatish Balay } 459b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 46049b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46149b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 462899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 46349b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46449b5e25fSSatish Balay aa = a->a + j*bs2; 46549b5e25fSSatish Balay for (k=0; k<bs; k++) { 46649b5e25fSSatish Balay for (l=0; l<bs; l++) { 46749b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 468b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46949b5e25fSSatish Balay } 47049b5e25fSSatish Balay } 47149b5e25fSSatish Balay } 47249b5e25fSSatish Balay } 47349b5e25fSSatish Balay 474b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 47549b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 47649b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 477899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 47849b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47949b5e25fSSatish Balay aa = a->a + j*bs2; 48049b5e25fSSatish Balay for (k=0; k<bs; k++) { 48149b5e25fSSatish Balay for (l=0; l<bs; l++) { 48249b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 483b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 48449b5e25fSSatish Balay } 48549b5e25fSSatish Balay } 48649b5e25fSSatish Balay } 48749b5e25fSSatish Balay } 48849b5e25fSSatish Balay PetscFunctionReturn(0); 48949b5e25fSSatish Balay } 49049b5e25fSSatish Balay 4914a2ae208SSatish Balay #undef __FUNCT__ 4924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4936849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 49449b5e25fSSatish Balay { 495dfbe8321SBarry Smith PetscErrorCode ierr; 49649b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 497b0a32e0cSBarry Smith PetscDraw draw; 49849b5e25fSSatish Balay PetscTruth isnull; 49949b5e25fSSatish Balay 50049b5e25fSSatish Balay PetscFunctionBegin; 501b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 502b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 50349b5e25fSSatish Balay 50449b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 505899cda47SBarry Smith xr = A->rmap.N; yr = A->rmap.N; h = yr/10.0; w = xr/10.0; 50649b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 507b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 508b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 50949b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 51049b5e25fSSatish Balay PetscFunctionReturn(0); 51149b5e25fSSatish Balay } 51249b5e25fSSatish Balay 5134a2ae208SSatish Balay #undef __FUNCT__ 5144a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 515dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 51649b5e25fSSatish Balay { 517dfbe8321SBarry Smith PetscErrorCode ierr; 51832077d6dSBarry Smith PetscTruth iascii,isdraw; 51949b5e25fSSatish Balay 52049b5e25fSSatish Balay PetscFunctionBegin; 52132077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 522fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 52332077d6dSBarry Smith if (iascii){ 52449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 52549b5e25fSSatish Balay } else if (isdraw) { 52649b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 52749b5e25fSSatish Balay } else { 528a5e6ed63SBarry Smith Mat B; 529ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 530a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 531a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 53249b5e25fSSatish Balay } 53349b5e25fSSatish Balay PetscFunctionReturn(0); 53449b5e25fSSatish Balay } 53549b5e25fSSatish Balay 53649b5e25fSSatish Balay 5374a2ae208SSatish Balay #undef __FUNCT__ 5384a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 53913f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 54049b5e25fSSatish Balay { 541045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 54213f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 54313f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 544899cda47SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->rmap.bs,bs2=a->bs2; 54549b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 54649b5e25fSSatish Balay 54749b5e25fSSatish Balay PetscFunctionBegin; 54849b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 54949b5e25fSSatish Balay row = im[k]; brow = row/bs; 55077431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 551899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 55249b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 55349b5e25fSSatish Balay nrow = ailen[brow]; 55449b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 55577431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 556899cda47SBarry 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); 55749b5e25fSSatish Balay col = in[l] ; 55849b5e25fSSatish Balay bcol = col/bs; 55949b5e25fSSatish Balay cidx = col%bs; 56049b5e25fSSatish Balay ridx = row%bs; 56149b5e25fSSatish Balay high = nrow; 56249b5e25fSSatish Balay low = 0; /* assume unsorted */ 56349b5e25fSSatish Balay while (high-low > 5) { 56449b5e25fSSatish Balay t = (low+high)/2; 56549b5e25fSSatish Balay if (rp[t] > bcol) high = t; 56649b5e25fSSatish Balay else low = t; 56749b5e25fSSatish Balay } 56849b5e25fSSatish Balay for (i=low; i<high; i++) { 56949b5e25fSSatish Balay if (rp[i] > bcol) break; 57049b5e25fSSatish Balay if (rp[i] == bcol) { 57149b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 57249b5e25fSSatish Balay goto finished; 57349b5e25fSSatish Balay } 57449b5e25fSSatish Balay } 57549b5e25fSSatish Balay *v++ = zero; 57649b5e25fSSatish Balay finished:; 57749b5e25fSSatish Balay } 57849b5e25fSSatish Balay } 57949b5e25fSSatish Balay PetscFunctionReturn(0); 58049b5e25fSSatish Balay } 58149b5e25fSSatish Balay 58249b5e25fSSatish Balay 5834a2ae208SSatish Balay #undef __FUNCT__ 5844a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 58513f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 58649b5e25fSSatish Balay { 5870880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5886849ba73SBarry Smith PetscErrorCode ierr; 589e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 59013f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 591899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap.bs,stepval; 5920880e062SHong Zhang PetscTruth roworiented=a->roworiented; 593f15d580aSBarry Smith const MatScalar *value = v; 594f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5950880e062SHong Zhang 59649b5e25fSSatish Balay PetscFunctionBegin; 5970880e062SHong Zhang if (roworiented) { 5980880e062SHong Zhang stepval = (n-1)*bs; 5990880e062SHong Zhang } else { 6000880e062SHong Zhang stepval = (m-1)*bs; 6010880e062SHong Zhang } 6020880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 6030880e062SHong Zhang row = im[k]; 6040880e062SHong Zhang if (row < 0) continue; 6052515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 60677431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 6070880e062SHong Zhang #endif 6080880e062SHong Zhang rp = aj + ai[row]; 6090880e062SHong Zhang ap = aa + bs2*ai[row]; 6100880e062SHong Zhang rmax = imax[row]; 6110880e062SHong Zhang nrow = ailen[row]; 6120880e062SHong Zhang low = 0; 613818f2c47SBarry Smith high = nrow; 6140880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6150880e062SHong Zhang if (in[l] < 0) continue; 6160880e062SHong Zhang col = in[l]; 6172515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 61877431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 619b1823623SSatish Balay #endif 6200880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6210880e062SHong Zhang if (roworiented) { 6220880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6230880e062SHong Zhang } else { 6240880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6250880e062SHong Zhang } 6267cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 627e2ee6c50SBarry Smith lastcol = col; 6280880e062SHong Zhang while (high-low > 7) { 6290880e062SHong Zhang t = (low+high)/2; 6300880e062SHong Zhang if (rp[t] > col) high = t; 6310880e062SHong Zhang else low = t; 6320880e062SHong Zhang } 6330880e062SHong Zhang for (i=low; i<high; i++) { 6340880e062SHong Zhang if (rp[i] > col) break; 6350880e062SHong Zhang if (rp[i] == col) { 6360880e062SHong Zhang bap = ap + bs2*i; 6370880e062SHong Zhang if (roworiented) { 6380880e062SHong Zhang if (is == ADD_VALUES) { 6390880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6400880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6410880e062SHong Zhang bap[jj] += *value++; 6420880e062SHong Zhang } 6430880e062SHong Zhang } 6440880e062SHong Zhang } else { 6450880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6460880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6470880e062SHong Zhang bap[jj] = *value++; 6480880e062SHong Zhang } 6490880e062SHong Zhang } 6500880e062SHong Zhang } 6510880e062SHong Zhang } else { 6520880e062SHong Zhang if (is == ADD_VALUES) { 6530880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6540880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6550880e062SHong Zhang *bap++ += *value++; 6560880e062SHong Zhang } 6570880e062SHong Zhang } 6580880e062SHong Zhang } else { 6590880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6600880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6610880e062SHong Zhang *bap++ = *value++; 6620880e062SHong Zhang } 6630880e062SHong Zhang } 6640880e062SHong Zhang } 6650880e062SHong Zhang } 6660880e062SHong Zhang goto noinsert2; 6670880e062SHong Zhang } 6680880e062SHong Zhang } 6690880e062SHong Zhang if (nonew == 1) goto noinsert2; 670085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 671421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 672c03d1d03SSatish Balay N = nrow++ - 1; high++; 6730880e062SHong Zhang /* shift up all the later entries in this row */ 6740880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6750880e062SHong Zhang rp[ii+1] = rp[ii]; 6760880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6770880e062SHong Zhang } 6780880e062SHong Zhang if (N >= i) { 6790880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6800880e062SHong Zhang } 6810880e062SHong Zhang rp[i] = col; 6820880e062SHong Zhang bap = ap + bs2*i; 6830880e062SHong Zhang if (roworiented) { 6840880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6850880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6860880e062SHong Zhang bap[jj] = *value++; 6870880e062SHong Zhang } 6880880e062SHong Zhang } 6890880e062SHong Zhang } else { 6900880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6910880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6920880e062SHong Zhang *bap++ = *value++; 6930880e062SHong Zhang } 6940880e062SHong Zhang } 6950880e062SHong Zhang } 6960880e062SHong Zhang noinsert2:; 6970880e062SHong Zhang low = i; 6980880e062SHong Zhang } 6990880e062SHong Zhang ailen[row] = nrow; 7000880e062SHong Zhang } 7010880e062SHong Zhang PetscFunctionReturn(0); 70249b5e25fSSatish Balay } 70349b5e25fSSatish Balay 7044a2ae208SSatish Balay #undef __FUNCT__ 7054a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 706dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 70749b5e25fSSatish Balay { 70849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7096849ba73SBarry Smith PetscErrorCode ierr; 71013f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 711899cda47SBarry Smith PetscInt m = A->rmap.N,*ip,N,*ailen = a->ilen; 71213f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 71349b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 71449b5e25fSSatish Balay 71549b5e25fSSatish Balay PetscFunctionBegin; 71649b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 71749b5e25fSSatish Balay 71849b5e25fSSatish Balay if (m) rmax = ailen[0]; 71949b5e25fSSatish Balay for (i=1; i<mbs; i++) { 72049b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 72149b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 72249b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 72349b5e25fSSatish Balay if (fshift) { 72449b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 72549b5e25fSSatish Balay N = ailen[i]; 72649b5e25fSSatish Balay for (j=0; j<N; j++) { 72749b5e25fSSatish Balay ip[j-fshift] = ip[j]; 72849b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 72949b5e25fSSatish Balay } 73049b5e25fSSatish Balay } 73149b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 73249b5e25fSSatish Balay } 73349b5e25fSSatish Balay if (mbs) { 73449b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 73549b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 73649b5e25fSSatish Balay } 73749b5e25fSSatish Balay /* reset ilen and imax for each row */ 73849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 73949b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 74049b5e25fSSatish Balay } 7416c6c5352SBarry Smith a->nz = ai[mbs]; 74249b5e25fSSatish Balay 743b424e231SHong Zhang /* diagonals may have moved, reset it */ 744b424e231SHong Zhang if (a->diag) { 74513f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 74649b5e25fSSatish Balay } 747899cda47SBarry 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); 748ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 749ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 75049b5e25fSSatish Balay a->reallocs = 0; 75149b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 75249b5e25fSSatish Balay PetscFunctionReturn(0); 75349b5e25fSSatish Balay } 75449b5e25fSSatish Balay 75549b5e25fSSatish Balay /* 75649b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 75749b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 75849b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 75949b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 76049b5e25fSSatish Balay */ 7614a2ae208SSatish Balay #undef __FUNCT__ 7624a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 76313f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 76449b5e25fSSatish Balay { 76513f74950SBarry Smith PetscInt i,j,k,row; 76649b5e25fSSatish Balay PetscTruth flg; 76749b5e25fSSatish Balay 76849b5e25fSSatish Balay PetscFunctionBegin; 76949b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 77049b5e25fSSatish Balay row = idx[i]; 77149b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 77249b5e25fSSatish Balay sizes[j] = 1; 77349b5e25fSSatish Balay i++; 77449b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 77549b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 77649b5e25fSSatish Balay i++; 77749b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 77849b5e25fSSatish Balay flg = PETSC_TRUE; 77949b5e25fSSatish Balay for (k=1; k<bs; k++) { 78049b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 78149b5e25fSSatish Balay flg = PETSC_FALSE; 78249b5e25fSSatish Balay break; 78349b5e25fSSatish Balay } 78449b5e25fSSatish Balay } 785abc0a331SBarry Smith if (flg) { /* No break in the bs */ 78649b5e25fSSatish Balay sizes[j] = bs; 78749b5e25fSSatish Balay i+= bs; 78849b5e25fSSatish Balay } else { 78949b5e25fSSatish Balay sizes[j] = 1; 79049b5e25fSSatish Balay i++; 79149b5e25fSSatish Balay } 79249b5e25fSSatish Balay } 79349b5e25fSSatish Balay } 79449b5e25fSSatish Balay *bs_max = j; 79549b5e25fSSatish Balay PetscFunctionReturn(0); 79649b5e25fSSatish Balay } 79749b5e25fSSatish Balay 79849b5e25fSSatish Balay 79949b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 80049b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 80149b5e25fSSatish Balay */ 80249b5e25fSSatish Balay 8034a2ae208SSatish Balay #undef __FUNCT__ 8044a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 80513f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 80649b5e25fSSatish Balay { 80749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8086849ba73SBarry Smith PetscErrorCode ierr; 809e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 81013f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 811899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 81213f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 81349b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 81449b5e25fSSatish Balay 81549b5e25fSSatish Balay PetscFunctionBegin; 81649b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 81749b5e25fSSatish Balay row = im[k]; /* row number */ 81849b5e25fSSatish Balay brow = row/bs; /* block row number */ 81949b5e25fSSatish Balay if (row < 0) continue; 8202515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 821899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 82249b5e25fSSatish Balay #endif 82349b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 82449b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 82549b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 82649b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 82749b5e25fSSatish Balay low = 0; 82849b5e25fSSatish Balay 82949b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 83049b5e25fSSatish Balay if (in[l] < 0) continue; 8312515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 832899cda47SBarry 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); 83349b5e25fSSatish Balay #endif 83449b5e25fSSatish Balay col = in[l]; 83549b5e25fSSatish Balay bcol = col/bs; /* block col number */ 83649b5e25fSSatish Balay 837941593c8SHong Zhang if (brow > bcol) { 838941593c8SHong Zhang if (a->ignore_ltriangular){ 839941593c8SHong Zhang continue; /* ignore lower triangular values */ 840941593c8SHong Zhang } else { 841941593c8SHong 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)"); 842941593c8SHong Zhang } 843941593c8SHong Zhang } 844f4989cb3SHong Zhang 84549b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8468549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 84749b5e25fSSatish Balay /* element value a(k,l) */ 84849b5e25fSSatish Balay if (roworiented) { 84949b5e25fSSatish Balay value = v[l + k*n]; 85049b5e25fSSatish Balay } else { 85149b5e25fSSatish Balay value = v[k + l*m]; 85249b5e25fSSatish Balay } 85349b5e25fSSatish Balay 85449b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 8557cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 856e2ee6c50SBarry Smith lastcol = col; 85749b5e25fSSatish Balay while (high-low > 7) { 85849b5e25fSSatish Balay t = (low+high)/2; 85949b5e25fSSatish Balay if (rp[t] > bcol) high = t; 86049b5e25fSSatish Balay else low = t; 86149b5e25fSSatish Balay } 86249b5e25fSSatish Balay for (i=low; i<high; i++) { 86349b5e25fSSatish Balay if (rp[i] > bcol) break; 86449b5e25fSSatish Balay if (rp[i] == bcol) { 86549b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 86649b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 86749b5e25fSSatish Balay else *bap = value; 8688549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8698549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8708549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8718549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8728549e402SHong Zhang else *bap = value; 8738549e402SHong Zhang } 87449b5e25fSSatish Balay goto noinsert1; 87549b5e25fSSatish Balay } 87649b5e25fSSatish Balay } 87749b5e25fSSatish Balay 87849b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 879085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 880421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 88149b5e25fSSatish Balay 882c03d1d03SSatish Balay N = nrow++ - 1; high++; 88349b5e25fSSatish Balay /* shift up all the later entries in this row */ 88449b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 88549b5e25fSSatish Balay rp[ii+1] = rp[ii]; 88649b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 88749b5e25fSSatish Balay } 88849b5e25fSSatish Balay if (N>=i) { 88949b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 89049b5e25fSSatish Balay } 89149b5e25fSSatish Balay rp[i] = bcol; 89249b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 89349b5e25fSSatish Balay noinsert1:; 89449b5e25fSSatish Balay low = i; 8958549e402SHong Zhang } 89649b5e25fSSatish Balay } /* end of loop over added columns */ 89749b5e25fSSatish Balay ailen[brow] = nrow; 89849b5e25fSSatish Balay } /* end of loop over added rows */ 89949b5e25fSSatish Balay PetscFunctionReturn(0); 90049b5e25fSSatish Balay } 90149b5e25fSSatish Balay 9024a2ae208SSatish Balay #undef __FUNCT__ 9034d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 904dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 90549b5e25fSSatish Balay { 9064ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 90749b5e25fSSatish Balay Mat outA; 908dfbe8321SBarry Smith PetscErrorCode ierr; 909c84f5b01SHong Zhang PetscTruth row_identity; 91049b5e25fSSatish Balay 91149b5e25fSSatish Balay PetscFunctionBegin; 912c84f5b01SHong Zhang if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 913c84f5b01SHong Zhang ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 914c84f5b01SHong Zhang if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 915c84f5b01SHong 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()! */ 916c84f5b01SHong Zhang 91749b5e25fSSatish Balay outA = inA; 9181260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 91949b5e25fSSatish Balay 9201a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 92149b5e25fSSatish Balay /* 922c84f5b01SHong Zhang Blocksize < 8 have a special faster factorization/solver 923c84f5b01SHong Zhang for ICC(0) factorization with natural ordering 92449b5e25fSSatish Balay */ 925c84f5b01SHong Zhang switch (inA->rmap.bs){ /* Note: row_identity = PETSC_TRUE! */ 92649b5e25fSSatish Balay case 1: 927c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 9280fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 92907f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 930d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 931759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 932759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 933c84f5b01SHong Zhang ierr = PetscInfo(inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr); 934c84f5b01SHong Zhang break; 93549b5e25fSSatish Balay case 2: 936c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9371a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9381a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 939759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 940759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 941ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 94249b5e25fSSatish Balay break; 94349b5e25fSSatish Balay case 3: 944c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9451a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9461a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 947759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering; 948759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering; 949ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 95049b5e25fSSatish Balay break; 95149b5e25fSSatish Balay case 4: 952c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 9531a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 9541a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 955759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering; 956759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering; 957ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 95849b5e25fSSatish Balay break; 95949b5e25fSSatish Balay case 5: 960c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 9611a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 9621a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 963759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering; 964759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering; 965ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 96649b5e25fSSatish Balay break; 96749b5e25fSSatish Balay case 6: 968c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 9691a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 9701a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 971759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering; 972759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering; 973ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 97449b5e25fSSatish Balay break; 97549b5e25fSSatish Balay case 7: 976c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 9771a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 9781a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 979759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering; 980759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering; 981ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 98249b5e25fSSatish Balay break; 98349b5e25fSSatish Balay default: 984c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 985c84f5b01SHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; 986c84f5b01SHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 987759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 988759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 989c84f5b01SHong Zhang break; 990c84f5b01SHong Zhang } 99149b5e25fSSatish Balay 992c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 993c3122656SLisandro Dalcin if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } 994c84f5b01SHong Zhang a->row = row; 995c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 996c3122656SLisandro Dalcin if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } 997c84f5b01SHong Zhang a->col = row; 998c84f5b01SHong Zhang 999c84f5b01SHong Zhang /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 1000c84f5b01SHong Zhang if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 100152e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 100249b5e25fSSatish Balay 100349b5e25fSSatish Balay if (!a->solve_work) { 1004c84f5b01SHong Zhang ierr = PetscMalloc((inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 1005c84f5b01SHong Zhang ierr = PetscLogObjectMemory(inA,(inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar));CHKERRQ(ierr); 100649b5e25fSSatish Balay } 100749b5e25fSSatish Balay 1008af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 100949b5e25fSSatish Balay PetscFunctionReturn(0); 101049b5e25fSSatish Balay } 1011950f1e5bSHong Zhang 101249b5e25fSSatish Balay EXTERN_C_BEGIN 10134a2ae208SSatish Balay #undef __FUNCT__ 10144a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1015be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 101649b5e25fSSatish Balay { 1017045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 101813f74950SBarry Smith PetscInt i,nz,n; 101949b5e25fSSatish Balay 102049b5e25fSSatish Balay PetscFunctionBegin; 10216c6c5352SBarry Smith nz = baij->maxnz; 1022899cda47SBarry Smith n = mat->cmap.n; 102349b5e25fSSatish Balay for (i=0; i<nz; i++) { 102449b5e25fSSatish Balay baij->j[i] = indices[i]; 102549b5e25fSSatish Balay } 10266c6c5352SBarry Smith baij->nz = nz; 102749b5e25fSSatish Balay for (i=0; i<n; i++) { 102849b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 102949b5e25fSSatish Balay } 103049b5e25fSSatish Balay PetscFunctionReturn(0); 103149b5e25fSSatish Balay } 103249b5e25fSSatish Balay EXTERN_C_END 103349b5e25fSSatish Balay 10344a2ae208SSatish Balay #undef __FUNCT__ 10354a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 103649b5e25fSSatish Balay /*@ 103719585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 103849b5e25fSSatish Balay in the matrix. 103949b5e25fSSatish Balay 104049b5e25fSSatish Balay Input Parameters: 104119585528SSatish Balay + mat - the SeqSBAIJ matrix 104249b5e25fSSatish Balay - indices - the column indices 104349b5e25fSSatish Balay 104449b5e25fSSatish Balay Level: advanced 104549b5e25fSSatish Balay 104649b5e25fSSatish Balay Notes: 104749b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 104849b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 104949b5e25fSSatish Balay of the MatSetValues() operation. 105049b5e25fSSatish Balay 105149b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1052d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 105349b5e25fSSatish Balay 1054ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 105549b5e25fSSatish Balay 1056ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 105749b5e25fSSatish Balay @*/ 1058be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 105949b5e25fSSatish Balay { 106013f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 106149b5e25fSSatish Balay 106249b5e25fSSatish Balay PetscFunctionBegin; 10634482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 10644482741eSBarry Smith PetscValidPointer(indices,2); 1065c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 106649b5e25fSSatish Balay if (f) { 106749b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 106849b5e25fSSatish Balay } else { 1069e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 107049b5e25fSSatish Balay } 107149b5e25fSSatish Balay PetscFunctionReturn(0); 107249b5e25fSSatish Balay } 107349b5e25fSSatish Balay 10744a2ae208SSatish Balay #undef __FUNCT__ 10753c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 10763c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 10773c896bc6SHong Zhang { 10783c896bc6SHong Zhang PetscErrorCode ierr; 10793c896bc6SHong Zhang 10803c896bc6SHong Zhang PetscFunctionBegin; 10813c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 10823c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 10833c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10843c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 10853c896bc6SHong Zhang 1086899cda47SBarry Smith if (a->i[A->rmap.N] != b->i[B->rmap.N]) { 10873c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 10883c896bc6SHong Zhang } 1089899cda47SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.N])*sizeof(PetscScalar));CHKERRQ(ierr); 10903c896bc6SHong Zhang } else { 1091f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 10923c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1093f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 10943c896bc6SHong Zhang } 10953c896bc6SHong Zhang PetscFunctionReturn(0); 10963c896bc6SHong Zhang } 10973c896bc6SHong Zhang 10983c896bc6SHong Zhang #undef __FUNCT__ 10994a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1100dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1101273d9f13SBarry Smith { 1102dfbe8321SBarry Smith PetscErrorCode ierr; 1103273d9f13SBarry Smith 1104273d9f13SBarry Smith PetscFunctionBegin; 11057edd0491SSatish Balay ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1106273d9f13SBarry Smith PetscFunctionReturn(0); 1107273d9f13SBarry Smith } 1108273d9f13SBarry Smith 1109a6ece127SHong Zhang #undef __FUNCT__ 1110a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1111dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1112a6ece127SHong Zhang { 1113a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1114a6ece127SHong Zhang PetscFunctionBegin; 1115a6ece127SHong Zhang *array = a->a; 1116a6ece127SHong Zhang PetscFunctionReturn(0); 1117a6ece127SHong Zhang } 1118a6ece127SHong Zhang 1119a6ece127SHong Zhang #undef __FUNCT__ 1120a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1121dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1122a6ece127SHong Zhang { 1123a6ece127SHong Zhang PetscFunctionBegin; 1124a6ece127SHong Zhang PetscFunctionReturn(0); 1125a6ece127SHong Zhang } 1126a6ece127SHong Zhang 112742ee4b1aSHong Zhang #include "petscblaslapack.h" 112842ee4b1aSHong Zhang #undef __FUNCT__ 112942ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1130f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 113142ee4b1aSHong Zhang { 113242ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1133dfbe8321SBarry Smith PetscErrorCode ierr; 1134899cda47SBarry Smith PetscInt i,bs=Y->rmap.bs,bs2,j; 11354ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 113642ee4b1aSHong Zhang 113742ee4b1aSHong Zhang PetscFunctionBegin; 113842ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1139f4df32b1SMatthew Knepley PetscScalar alpha = a; 1140f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1141c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1142c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1143c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1144c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1145c537a176SHong Zhang } 1146c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1147c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1148c4319e64SHong Zhang y->XtoY = X; 1149c537a176SHong Zhang } 1150c4319e64SHong Zhang bs2 = bs*bs; 11516c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1152c4319e64SHong Zhang j = 0; 1153c4319e64SHong Zhang while (j < bs2){ 1154f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1155c4319e64SHong Zhang j++; 1156c537a176SHong Zhang } 1157c4319e64SHong Zhang } 1158ae15b995SBarry 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); 115942ee4b1aSHong Zhang } else { 1160f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1161f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1162f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 116342ee4b1aSHong Zhang } 116442ee4b1aSHong Zhang PetscFunctionReturn(0); 116542ee4b1aSHong Zhang } 116642ee4b1aSHong Zhang 1167efcf0fc3SBarry Smith #undef __FUNCT__ 1168efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1169dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1170efcf0fc3SBarry Smith { 1171efcf0fc3SBarry Smith PetscFunctionBegin; 1172efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1173efcf0fc3SBarry Smith PetscFunctionReturn(0); 1174efcf0fc3SBarry Smith } 1175efcf0fc3SBarry Smith 1176efcf0fc3SBarry Smith #undef __FUNCT__ 1177efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1178dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1179efcf0fc3SBarry Smith { 1180efcf0fc3SBarry Smith PetscFunctionBegin; 1181efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1182efcf0fc3SBarry Smith PetscFunctionReturn(0); 1183efcf0fc3SBarry Smith } 1184efcf0fc3SBarry Smith 1185efcf0fc3SBarry Smith #undef __FUNCT__ 1186efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1187ab5e4463SMatthew Knepley PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1188efcf0fc3SBarry Smith { 1189efcf0fc3SBarry Smith PetscFunctionBegin; 1190efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1191efcf0fc3SBarry Smith PetscFunctionReturn(0); 1192efcf0fc3SBarry Smith } 1193efcf0fc3SBarry Smith 119499cafbc1SBarry Smith #undef __FUNCT__ 119599cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 119699cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 119799cafbc1SBarry Smith { 119899cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 119999cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 120099cafbc1SBarry Smith PetscScalar *aa = a->a; 120199cafbc1SBarry Smith 120299cafbc1SBarry Smith PetscFunctionBegin; 120399cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 120499cafbc1SBarry Smith PetscFunctionReturn(0); 120599cafbc1SBarry Smith } 120699cafbc1SBarry Smith 120799cafbc1SBarry Smith #undef __FUNCT__ 120899cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 120999cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 121099cafbc1SBarry Smith { 121199cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 121299cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 121399cafbc1SBarry Smith PetscScalar *aa = a->a; 121499cafbc1SBarry Smith 121599cafbc1SBarry Smith PetscFunctionBegin; 121699cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 121799cafbc1SBarry Smith PetscFunctionReturn(0); 121899cafbc1SBarry Smith } 121999cafbc1SBarry Smith 122049b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 122149b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 122249b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 122349b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 122449b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 122597304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1226431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1227e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 122849b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 122949b5e25fSSatish Balay 0, 123049b5e25fSSatish Balay 0, 123197304618SKris Buschelman /*10*/ 0, 123249b5e25fSSatish Balay 0, 12335f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1234d06b337dSHong Zhang MatRelax_SeqSBAIJ, 123549b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 123697304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 123749b5e25fSSatish Balay MatEqual_SeqSBAIJ, 123849b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 123949b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 124049b5e25fSSatish Balay MatNorm_SeqSBAIJ, 124197304618SKris Buschelman /*20*/ 0, 124249b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 124349b5e25fSSatish Balay 0, 124449b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 124549b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1246dcf5cc72SBarry Smith /*25*/ 0, 124749b5e25fSSatish Balay 0, 124849b5e25fSSatish Balay 0, 12495f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12505f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 125197304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1252c464158bSHong Zhang 0, 12534d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1254a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1255a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 125697304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 12579495ac64SHong Zhang MatForwardSolve_SeqSBAIJ_N, 12589495ac64SHong Zhang MatBackwardSolve_SeqSBAIJ_N, 125949b5e25fSSatish Balay 0, 1260c84f5b01SHong Zhang MatICCFactor_SeqSBAIJ, 126197304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 126249b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 126349b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 126449b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 12653c896bc6SHong Zhang MatCopy_SeqSBAIJ, 12668c07d4e3SBarry Smith /*45*/ 0, 126749b5e25fSSatish Balay MatScale_SeqSBAIJ, 126849b5e25fSSatish Balay 0, 126949b5e25fSSatish Balay 0, 127049b5e25fSSatish Balay 0, 1271521d7252SBarry Smith /*50*/ 0, 127249b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 127349b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 127449b5e25fSSatish Balay 0, 127549b5e25fSSatish Balay 0, 127697304618SKris Buschelman /*55*/ 0, 127749b5e25fSSatish Balay 0, 127849b5e25fSSatish Balay 0, 127949b5e25fSSatish Balay 0, 128049b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 128197304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 128249b5e25fSSatish Balay 0, 128349b5e25fSSatish Balay 0, 1284357abbc8SBarry Smith 0, 1285d959ec07SHong Zhang 0, 128697304618SKris Buschelman /*65*/ 0, 1287d959ec07SHong Zhang 0, 1288d959ec07SHong Zhang 0, 1289d959ec07SHong Zhang 0, 1290d959ec07SHong Zhang 0, 1291985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqSBAIJ, 12923e0d88b5SBarry Smith 0, 12933e0d88b5SBarry Smith 0, 12943e0d88b5SBarry Smith 0, 12953e0d88b5SBarry Smith 0, 129697304618SKris Buschelman /*75*/ 0, 12973e0d88b5SBarry Smith 0, 12983e0d88b5SBarry Smith 0, 12993e0d88b5SBarry Smith 0, 13003e0d88b5SBarry Smith 0, 130197304618SKris Buschelman /*80*/ 0, 13023e0d88b5SBarry Smith 0, 13033e0d88b5SBarry Smith 0, 13043e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 130597304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13063e0d88b5SBarry Smith #else 130797304618SKris Buschelman 0, 13083e0d88b5SBarry Smith #endif 1309865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1310865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1311865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1312efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1313865e5f61SKris Buschelman 0, 1314865e5f61SKris Buschelman 0, 1315865e5f61SKris Buschelman /*90*/ 0, 1316865e5f61SKris Buschelman 0, 1317865e5f61SKris Buschelman 0, 1318865e5f61SKris Buschelman 0, 1319865e5f61SKris Buschelman 0, 1320865e5f61SKris Buschelman /*95*/ 0, 1321865e5f61SKris Buschelman 0, 1322865e5f61SKris Buschelman 0, 132399cafbc1SBarry Smith 0, 132499cafbc1SBarry Smith 0, 132599cafbc1SBarry Smith /*100*/0, 132699cafbc1SBarry Smith 0, 132799cafbc1SBarry Smith 0, 132899cafbc1SBarry Smith 0, 132999cafbc1SBarry Smith 0, 133099cafbc1SBarry Smith /*105*/0, 133199cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1332f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1333f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 1334f5edf698SHong Zhang MatRestoreRowUpperTriangular_SeqSBAIJ 133599cafbc1SBarry Smith }; 1336be1d678aSKris Buschelman 133749b5e25fSSatish Balay EXTERN_C_BEGIN 13384a2ae208SSatish Balay #undef __FUNCT__ 13394a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1340be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 134149b5e25fSSatish Balay { 13424afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1343899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 1344dfbe8321SBarry Smith PetscErrorCode ierr; 134549b5e25fSSatish Balay 134649b5e25fSSatish Balay PetscFunctionBegin; 134749b5e25fSSatish Balay if (aij->nonew != 1) { 1348e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 134949b5e25fSSatish Balay } 135049b5e25fSSatish Balay 135149b5e25fSSatish Balay /* allocate space for values if not already there */ 135249b5e25fSSatish Balay if (!aij->saved_values) { 135387828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 135449b5e25fSSatish Balay } 135549b5e25fSSatish Balay 135649b5e25fSSatish Balay /* copy values over */ 135787828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 135849b5e25fSSatish Balay PetscFunctionReturn(0); 135949b5e25fSSatish Balay } 136049b5e25fSSatish Balay EXTERN_C_END 136149b5e25fSSatish Balay 136249b5e25fSSatish Balay EXTERN_C_BEGIN 13634a2ae208SSatish Balay #undef __FUNCT__ 13644a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1365be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 136649b5e25fSSatish Balay { 13674afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13686849ba73SBarry Smith PetscErrorCode ierr; 1369899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 137049b5e25fSSatish Balay 137149b5e25fSSatish Balay PetscFunctionBegin; 137249b5e25fSSatish Balay if (aij->nonew != 1) { 1373e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 137449b5e25fSSatish Balay } 137549b5e25fSSatish Balay if (!aij->saved_values) { 1376e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 137749b5e25fSSatish Balay } 137849b5e25fSSatish Balay 137949b5e25fSSatish Balay /* copy values over */ 138087828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 138149b5e25fSSatish Balay PetscFunctionReturn(0); 138249b5e25fSSatish Balay } 138349b5e25fSSatish Balay EXTERN_C_END 138449b5e25fSSatish Balay 13858549e402SHong Zhang EXTERN_C_BEGIN 13864a2ae208SSatish Balay #undef __FUNCT__ 1387a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1388be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 138949b5e25fSSatish Balay { 1390c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13916849ba73SBarry Smith PetscErrorCode ierr; 1392ab93d7beSBarry Smith PetscInt i,mbs,bs2; 1393ab93d7beSBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg; 139449b5e25fSSatish Balay 139549b5e25fSSatish Balay PetscFunctionBegin; 1396273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1397e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1398899cda47SBarry Smith B->rmap.bs = B->cmap.bs = bs; 13996148ca0dSBarry Smith ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr); 14006148ca0dSBarry Smith ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr); 1401899cda47SBarry Smith 1402899cda47SBarry Smith mbs = B->rmap.N/bs; 140349b5e25fSSatish Balay bs2 = bs*bs; 140449b5e25fSSatish Balay 1405899cda47SBarry Smith if (mbs*bs != B->rmap.N) { 140629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 140749b5e25fSSatish Balay } 140849b5e25fSSatish Balay 1409ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1410ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1411ab93d7beSBarry Smith nz = 0; 1412ab93d7beSBarry Smith } 1413ab93d7beSBarry Smith 1414435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 141577431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 141649b5e25fSSatish Balay if (nnz) { 141749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 141877431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 141977431f27SBarry 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); 142049b5e25fSSatish Balay } 142149b5e25fSSatish Balay } 142249b5e25fSSatish Balay 1423e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 142449b5e25fSSatish Balay if (!flg) { 142549b5e25fSSatish Balay switch (bs) { 142649b5e25fSSatish Balay case 1: 14274ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 142849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1429d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1430e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 143149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 143249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1433431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1434431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 143549b5e25fSSatish Balay break; 143649b5e25fSSatish Balay case 2: 14374ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 143849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1439e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 144049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 144149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1442431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1443431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 144449b5e25fSSatish Balay break; 144549b5e25fSSatish Balay case 3: 14465f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 144749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1448e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 144949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 145049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1451431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1452431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 145349b5e25fSSatish Balay break; 145449b5e25fSSatish Balay case 4: 14555f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 145649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1457e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 145849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 145949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1460431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1461431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 146249b5e25fSSatish Balay break; 146349b5e25fSSatish Balay case 5: 14645f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 146549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1466e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 146749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 146849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1469431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1470431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 147149b5e25fSSatish Balay break; 147249b5e25fSSatish Balay case 6: 14735f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 147449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1475e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 147649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 147749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1478431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1479431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 148049b5e25fSSatish Balay break; 148149b5e25fSSatish Balay case 7: 1482de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 148349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1484e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1485de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 148649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1487431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1488431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 148949b5e25fSSatish Balay break; 149049b5e25fSSatish Balay } 149149b5e25fSSatish Balay } 149249b5e25fSSatish Balay 149349b5e25fSSatish Balay b->mbs = mbs; 14944afc71dfSHong Zhang b->nbs = mbs; 1495ab93d7beSBarry Smith if (!skipallocation) { 1496ab93d7beSBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1497ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1498ab93d7beSBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 149949b5e25fSSatish Balay if (!nnz) { 1500435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 150149b5e25fSSatish Balay else if (nz <= 0) nz = 1; 150249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15038cef66ccSHong Zhang b->imax[i] = nz; 150449b5e25fSSatish Balay } 1505153ea458SHong Zhang nz = nz*mbs; /* total nz */ 150649b5e25fSSatish Balay } else { 150749b5e25fSSatish Balay nz = 0; 15088cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 150949b5e25fSSatish Balay } 15106c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 151149b5e25fSSatish Balay 151249b5e25fSSatish Balay /* allocate the matrix space */ 1513899cda47SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.N+1,PetscInt,&b->i);CHKERRQ(ierr); 15146c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 151513f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 151649b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 151749b5e25fSSatish Balay 151849b5e25fSSatish Balay /* pointer to beginning of each row */ 151949b5e25fSSatish Balay b->i[0] = 0; 152049b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 152149b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 152249b5e25fSSatish Balay } 1523e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1524e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1525e811da20SHong Zhang } else { 1526e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1527e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1528ab93d7beSBarry Smith } 152949b5e25fSSatish Balay 1530899cda47SBarry Smith B->rmap.bs = bs; 153149b5e25fSSatish Balay b->bs2 = bs2; 15326c6c5352SBarry Smith b->nz = 0; 15336c6c5352SBarry Smith b->maxnz = nz*bs2; 1534153ea458SHong Zhang 153516cdd363SHong Zhang b->inew = 0; 153616cdd363SHong Zhang b->jnew = 0; 153716cdd363SHong Zhang b->anew = 0; 153816cdd363SHong Zhang b->a2anew = 0; 15391a3463dfSHong Zhang b->permute = PETSC_FALSE; 1540c464158bSHong Zhang PetscFunctionReturn(0); 1541c464158bSHong Zhang } 1542a23d5eceSKris Buschelman EXTERN_C_END 1543153ea458SHong Zhang 1544d769727bSBarry Smith EXTERN_C_BEGIN 1545f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1546f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1547d769727bSBarry Smith EXTERN_C_END 1548d769727bSBarry Smith 15490bad9183SKris Buschelman /*MC 1550fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15510bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15520bad9183SKris Buschelman 15530bad9183SKris Buschelman Options Database Keys: 15540bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15550bad9183SKris Buschelman 15560bad9183SKris Buschelman Level: beginner 15570bad9183SKris Buschelman 15580bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15590bad9183SKris Buschelman M*/ 15600bad9183SKris Buschelman 1561a23d5eceSKris Buschelman EXTERN_C_BEGIN 1562a23d5eceSKris Buschelman #undef __FUNCT__ 1563a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1564be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1565a23d5eceSKris Buschelman { 1566a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1567dfbe8321SBarry Smith PetscErrorCode ierr; 156813f74950SBarry Smith PetscMPIInt size; 1569941593c8SHong Zhang PetscTruth flg; 1570a23d5eceSKris Buschelman 1571a23d5eceSKris Buschelman PetscFunctionBegin; 1572a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1573a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1574a23d5eceSKris Buschelman 1575*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1576a23d5eceSKris Buschelman B->data = (void*)b; 1577a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1578a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1579a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1580a23d5eceSKris Buschelman B->factor = 0; 1581a23d5eceSKris Buschelman B->mapping = 0; 1582a23d5eceSKris Buschelman b->row = 0; 1583a23d5eceSKris Buschelman b->icol = 0; 1584a23d5eceSKris Buschelman b->reallocs = 0; 1585a23d5eceSKris Buschelman b->saved_values = 0; 1586a23d5eceSKris Buschelman 1587a23d5eceSKris Buschelman 1588a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1589a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1590a23d5eceSKris Buschelman b->nonew = 0; 1591a23d5eceSKris Buschelman b->diag = 0; 1592a23d5eceSKris Buschelman b->solve_work = 0; 1593a23d5eceSKris Buschelman b->mult_work = 0; 1594a23d5eceSKris Buschelman B->spptr = 0; 1595a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1596a23d5eceSKris Buschelman b->xtoy = 0; 1597a23d5eceSKris Buschelman b->XtoY = 0; 1598a23d5eceSKris Buschelman 1599a23d5eceSKris Buschelman b->inew = 0; 1600a23d5eceSKris Buschelman b->jnew = 0; 1601a23d5eceSKris Buschelman b->anew = 0; 1602a23d5eceSKris Buschelman b->a2anew = 0; 1603a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1604a23d5eceSKris Buschelman 1605941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 1606941593c8SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr); 1607941593c8SHong Zhang if (flg) b->ignore_ltriangular = PETSC_TRUE; 1608941593c8SHong Zhang 1609f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 1610f5edf698SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr); 1611f5edf698SHong Zhang if (flg) b->getrow_utriangular = PETSC_TRUE; 1612f5edf698SHong Zhang 1613a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1614a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1615a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1616a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1617a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1618a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1619a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1620a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1621a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16224e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16234e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16244e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1625a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1626a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1627a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1628a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1629a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1630a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 163123ce1328SBarry Smith 163223ce1328SBarry Smith B->symmetric = PETSC_TRUE; 163323ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 163423ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 163523ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 163617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 1637a23d5eceSKris Buschelman PetscFunctionReturn(0); 1638a23d5eceSKris Buschelman } 1639a23d5eceSKris Buschelman EXTERN_C_END 1640a23d5eceSKris Buschelman 1641a23d5eceSKris Buschelman #undef __FUNCT__ 1642a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1643a23d5eceSKris Buschelman /*@C 1644a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1645a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1646a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1647a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1648a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1649a23d5eceSKris Buschelman 1650a23d5eceSKris Buschelman Collective on Mat 1651a23d5eceSKris Buschelman 1652a23d5eceSKris Buschelman Input Parameters: 1653a23d5eceSKris Buschelman + A - the symmetric matrix 1654a23d5eceSKris Buschelman . bs - size of block 1655a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1656a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1657a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1658a23d5eceSKris Buschelman 1659a23d5eceSKris Buschelman Options Database Keys: 1660a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1661a23d5eceSKris Buschelman block calculations (much slower) 1662a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1663a23d5eceSKris Buschelman 1664a23d5eceSKris Buschelman Level: intermediate 1665a23d5eceSKris Buschelman 1666a23d5eceSKris Buschelman Notes: 1667a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1668a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1669a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1670a23d5eceSKris Buschelman matrices. 1671a23d5eceSKris Buschelman 167249a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 167349a6f317SBarry Smith 167449a6f317SBarry Smith 1675a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1676a23d5eceSKris Buschelman @*/ 1677be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 167813f74950SBarry Smith { 167913f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1680a23d5eceSKris Buschelman 1681a23d5eceSKris Buschelman PetscFunctionBegin; 1682a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1683a23d5eceSKris Buschelman if (f) { 1684a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1685a23d5eceSKris Buschelman } 1686a23d5eceSKris Buschelman PetscFunctionReturn(0); 1687a23d5eceSKris Buschelman } 168849b5e25fSSatish Balay 16894a2ae208SSatish Balay #undef __FUNCT__ 16904a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1691c464158bSHong Zhang /*@C 1692c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1693c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1694c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1695c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1696c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 169749b5e25fSSatish Balay 1698c464158bSHong Zhang Collective on MPI_Comm 1699c464158bSHong Zhang 1700c464158bSHong Zhang Input Parameters: 1701c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1702c464158bSHong Zhang . bs - size of block 1703c464158bSHong Zhang . m - number of rows, or number of columns 1704c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1705744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1706744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1707c464158bSHong Zhang 1708c464158bSHong Zhang Output Parameter: 1709c464158bSHong Zhang . A - the symmetric matrix 1710c464158bSHong Zhang 1711c464158bSHong Zhang Options Database Keys: 1712c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1713c464158bSHong Zhang block calculations (much slower) 1714c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1715c464158bSHong Zhang 1716c464158bSHong Zhang Level: intermediate 1717c464158bSHong Zhang 1718c464158bSHong Zhang Notes: 17196d6d819aSHong Zhang The number of rows and columns must be divisible by blocksize. 17206d6d819aSHong Zhang This matrix type does not support complex Hermitian operation. 1721c464158bSHong Zhang 1722c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1723c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1724c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1725c464158bSHong Zhang matrices. 1726c464158bSHong Zhang 172749a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 172849a6f317SBarry Smith 1729c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1730c464158bSHong Zhang @*/ 1731be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1732c464158bSHong Zhang { 1733dfbe8321SBarry Smith PetscErrorCode ierr; 1734c464158bSHong Zhang 1735c464158bSHong Zhang PetscFunctionBegin; 1736f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1737f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1738c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1739ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 174049b5e25fSSatish Balay PetscFunctionReturn(0); 174149b5e25fSSatish Balay } 174249b5e25fSSatish Balay 17434a2ae208SSatish Balay #undef __FUNCT__ 17444a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1745dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 174649b5e25fSSatish Balay { 174749b5e25fSSatish Balay Mat C; 174849b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17496849ba73SBarry Smith PetscErrorCode ierr; 1750b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 175149b5e25fSSatish Balay 175249b5e25fSSatish Balay PetscFunctionBegin; 175329bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 175449b5e25fSSatish Balay 175549b5e25fSSatish Balay *B = 0; 1756f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); 1757899cda47SBarry Smith ierr = MatSetSizes(C,A->rmap.N,A->cmap.n,A->rmap.N,A->cmap.n);CHKERRQ(ierr); 1758be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17591d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1760692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1761692f9cbeSHong Zhang 1762273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 176349b5e25fSSatish Balay C->factor = A->factor; 176449b5e25fSSatish Balay c->row = 0; 176549b5e25fSSatish Balay c->icol = 0; 176649b5e25fSSatish Balay c->saved_values = 0; 176749b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 176849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 176949b5e25fSSatish Balay 1770899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->rmap,&C->rmap);CHKERRQ(ierr); 1771899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->cmap,&C->cmap);CHKERRQ(ierr); 177249b5e25fSSatish Balay c->bs2 = a->bs2; 177349b5e25fSSatish Balay c->mbs = a->mbs; 177449b5e25fSSatish Balay c->nbs = a->nbs; 177549b5e25fSSatish Balay 17768777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 177749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 177849b5e25fSSatish Balay c->imax[i] = a->imax[i]; 177949b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 178049b5e25fSSatish Balay } 178149b5e25fSSatish Balay 178249b5e25fSSatish Balay /* allocate the matrix space */ 1783b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 178449b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 178513f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1786b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 178749b5e25fSSatish Balay if (mbs > 0) { 178813f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 178949b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 179049b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 179149b5e25fSSatish Balay } else { 179249b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 179349b5e25fSSatish Balay } 179449b5e25fSSatish Balay } 179549b5e25fSSatish Balay 179649b5e25fSSatish Balay c->sorted = a->sorted; 179749b5e25fSSatish Balay c->roworiented = a->roworiented; 179849b5e25fSSatish Balay c->nonew = a->nonew; 179949b5e25fSSatish Balay 180049b5e25fSSatish Balay if (a->diag) { 180113f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 180252e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 180349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 180449b5e25fSSatish Balay c->diag[i] = a->diag[i]; 180549b5e25fSSatish Balay } 180649b5e25fSSatish Balay } else c->diag = 0; 18076c6c5352SBarry Smith c->nz = a->nz; 18086c6c5352SBarry Smith c->maxnz = a->maxnz; 180949b5e25fSSatish Balay c->solve_work = 0; 181049b5e25fSSatish Balay c->mult_work = 0; 1811e6b907acSBarry Smith c->free_a = PETSC_TRUE; 1812e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 181349b5e25fSSatish Balay *B = C; 1814b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 181549b5e25fSSatish Balay PetscFunctionReturn(0); 181649b5e25fSSatish Balay } 181749b5e25fSSatish Balay 18184a2ae208SSatish Balay #undef __FUNCT__ 18194a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1820f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A) 182149b5e25fSSatish Balay { 182249b5e25fSSatish Balay Mat_SeqSBAIJ *a; 182349b5e25fSSatish Balay Mat B; 18246849ba73SBarry Smith PetscErrorCode ierr; 182513f74950SBarry Smith int fd; 182613f74950SBarry Smith PetscMPIInt size; 182713f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 182813f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 182913f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 183013f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 183187828ca2SBarry Smith PetscScalar *aa; 183249b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 183349b5e25fSSatish Balay 183449b5e25fSSatish Balay PetscFunctionBegin; 1835b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 183649b5e25fSSatish Balay bs2 = bs*bs; 183749b5e25fSSatish Balay 183849b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 183929bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1840b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 184149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1842552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 184349b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 184449b5e25fSSatish Balay 184549b5e25fSSatish Balay if (header[3] < 0) { 184629bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 184749b5e25fSSatish Balay } 184849b5e25fSSatish Balay 184929bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 185049b5e25fSSatish Balay 185149b5e25fSSatish Balay /* 185249b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 185349b5e25fSSatish Balay divisible by the blocksize 185449b5e25fSSatish Balay */ 185549b5e25fSSatish Balay mbs = M/bs; 185649b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 185749b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 185849b5e25fSSatish Balay else mbs++; 185949b5e25fSSatish Balay if (extra_rows) { 1860ae15b995SBarry Smith ierr = PetscInfo(0,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 186149b5e25fSSatish Balay } 186249b5e25fSSatish Balay 186349b5e25fSSatish Balay /* read in row lengths */ 186413f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 186549b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 186649b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 186749b5e25fSSatish Balay 186849b5e25fSSatish Balay /* read in column indices */ 186913f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 187049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 187149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 187249b5e25fSSatish Balay 187349b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 187413f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 187513f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 187613f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 187713f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 187849b5e25fSSatish Balay masked = mask + mbs; 187949b5e25fSSatish Balay rowcount = 0; nzcount = 0; 188049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 188149b5e25fSSatish Balay nmask = 0; 188249b5e25fSSatish Balay for (j=0; j<bs; j++) { 188349b5e25fSSatish Balay kmax = rowlengths[rowcount]; 188449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18852d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 188603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 188749b5e25fSSatish Balay } 188849b5e25fSSatish Balay rowcount++; 188949b5e25fSSatish Balay } 1890574b2666SHong Zhang s_browlengths[i] += nmask; 1891574b2666SHong Zhang 189249b5e25fSSatish Balay /* zero out the mask elements we set */ 189349b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 189449b5e25fSSatish Balay } 189549b5e25fSSatish Balay 189649b5e25fSSatish Balay /* create our matrix */ 1897f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 1898f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 18999abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 1900ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 190149b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 190249b5e25fSSatish Balay 190349b5e25fSSatish Balay /* set matrix "i" values */ 190449b5e25fSSatish Balay a->i[0] = 0; 190549b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1906574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1907574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 190849b5e25fSSatish Balay } 19096c6c5352SBarry Smith a->nz = a->i[mbs]; 191049b5e25fSSatish Balay 191149b5e25fSSatish Balay /* read in nonzero values */ 191287828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 191349b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 191449b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 191549b5e25fSSatish Balay 191649b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 191749b5e25fSSatish Balay nzcount = 0; jcount = 0; 191849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 191949b5e25fSSatish Balay nzcountb = nzcount; 192049b5e25fSSatish Balay nmask = 0; 192149b5e25fSSatish Balay for (j=0; j<bs; j++) { 192249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 192349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19242d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 192503630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19262d703238SHong Zhang } 19272d703238SHong Zhang } 19282d703238SHong Zhang /* sort the masked values */ 19292d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19302d703238SHong Zhang 19312d703238SHong Zhang /* set "j" values into matrix */ 19322d703238SHong Zhang maskcount = 1; 19332d703238SHong Zhang for (j=0; j<nmask; j++) { 193449b5e25fSSatish Balay a->j[jcount++] = masked[j]; 193549b5e25fSSatish Balay mask[masked[j]] = maskcount++; 193649b5e25fSSatish Balay } 1937574b2666SHong Zhang 193849b5e25fSSatish Balay /* set "a" values into matrix */ 193949b5e25fSSatish Balay ishift = bs2*a->i[i]; 194049b5e25fSSatish Balay for (j=0; j<bs; j++) { 194149b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 194249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1943574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1944574b2666SHong Zhang if (tmp >= i){ 194549b5e25fSSatish Balay block = mask[tmp] - 1; 194649b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 194749b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1948574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1949574b2666SHong Zhang } 1950574b2666SHong Zhang nzcountb++; 195149b5e25fSSatish Balay } 195249b5e25fSSatish Balay } 195349b5e25fSSatish Balay /* zero out the mask elements we set */ 195449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 195549b5e25fSSatish Balay } 19566c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 195749b5e25fSSatish Balay 195849b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1959574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 196049b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 196149b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 196249b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 196349b5e25fSSatish Balay 19649abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19659abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196649b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19679abb65ffSKris Buschelman *A = B; 196849b5e25fSSatish Balay PetscFunctionReturn(0); 196949b5e25fSSatish Balay } 1970574b2666SHong Zhang 1971d06b337dSHong Zhang #undef __FUNCT__ 1972d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 197313f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1974d06b337dSHong Zhang { 1975d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1976d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1977d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19786849ba73SBarry Smith PetscErrorCode ierr; 1979899cda47SBarry Smith PetscInt m=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j; 1980521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 1981d06b337dSHong Zhang 1982d06b337dSHong Zhang PetscFunctionBegin; 198351f519a2SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 198451f519a2SBarry Smith 1985b965ef7fSBarry Smith its = its*lits; 198677431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1987d06b337dSHong Zhang 1988d06b337dSHong Zhang if (bs > 1) 1989d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1990d06b337dSHong Zhang 19911ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1992d06b337dSHong Zhang if (xx != bb) { 19931ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1994d06b337dSHong Zhang } else { 1995d06b337dSHong Zhang b = x; 1996d06b337dSHong Zhang } 1997d06b337dSHong Zhang 1998e2ee2a47SBarry Smith if (!a->relax_work) { 1999e2ee2a47SBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->relax_work);CHKERRQ(ierr); 2000e2ee2a47SBarry Smith } 2001e2ee2a47SBarry Smith t = a->relax_work; 2002d06b337dSHong Zhang 2003d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 20042798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2005290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2006d06b337dSHong Zhang 2007d06b337dSHong Zhang for (i=0; i<m; i++){ 200844706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2009d06b337dSHong Zhang v = aa + ai[i] + 1; 2010d06b337dSHong Zhang vj = aj + ai[i] + 1; 2011d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2012e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2013d06b337dSHong Zhang x[i] = omega*t[i]/d; 2014d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 2015d06b337dSHong Zhang } 2016d06b337dSHong Zhang } 2017d06b337dSHong Zhang 20182798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 201995d750ceSBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){ 202095d750ceSBarry Smith t = b; 2021d06b337dSHong Zhang } 202295d750ceSBarry Smith 2023d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2024d06b337dSHong Zhang d = *(aa + ai[i]); 2025d06b337dSHong Zhang v = aa + ai[i] + 1; 2026d06b337dSHong Zhang vj = aj + ai[i] + 1; 2027d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2028e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2029d06b337dSHong Zhang sum = t[i]; 2030d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2031d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2032d06b337dSHong Zhang } 203395d750ceSBarry Smith t = a->relax_work; 2034d06b337dSHong Zhang } 2035d06b337dSHong Zhang its--; 2036d06b337dSHong Zhang } 2037d06b337dSHong Zhang 2038d06b337dSHong Zhang while (its--) { 203944706875SHong Zhang /* 204044706875SHong Zhang forward sweep: 204144706875SHong Zhang for i=0,...,m-1: 204244706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 204344706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 204444706875SHong Zhang b = b - x[i]*U^T(i,:); 2045d06b337dSHong Zhang 204644706875SHong Zhang */ 20472798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2048290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2049d06b337dSHong Zhang 2050d06b337dSHong Zhang for (i=0; i<m; i++){ 205144706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2052d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2053d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2054d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2055d06b337dSHong Zhang sum = t[i]; 2056e6b907acSBarry Smith ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2057d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2058d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2059d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2060d06b337dSHong Zhang } 2061d06b337dSHong Zhang } 2062d06b337dSHong Zhang 20632798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 206444706875SHong Zhang /* 206544706875SHong Zhang backward sweep: 206644706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 206744706875SHong Zhang for i=m-1,...,0: 206844706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 206944706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 207044706875SHong Zhang */ 207195d750ceSBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2072290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2073d06b337dSHong Zhang 2074d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2075d06b337dSHong Zhang v = aa + ai[i] + 1; 2076d06b337dSHong Zhang vj = aj + ai[i] + 1; 2077d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2078efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2079e6b907acSBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 2080d06b337dSHong Zhang } 2081d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2082d06b337dSHong Zhang d = *(aa + ai[i]); 2083d06b337dSHong Zhang v = aa + ai[i] + 1; 2084d06b337dSHong Zhang vj = aj + ai[i] + 1; 2085d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2086e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2087d06b337dSHong Zhang sum = t[i]; 2088d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2089d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2090d06b337dSHong Zhang } 2091d06b337dSHong Zhang } 2092d06b337dSHong Zhang } 2093d06b337dSHong Zhang 20941ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2095d06b337dSHong Zhang if (bb != xx) { 20961ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2097d06b337dSHong Zhang } 2098d06b337dSHong Zhang PetscFunctionReturn(0); 2099d06b337dSHong Zhang } 2100d06b337dSHong Zhang 2101c75a6043SHong Zhang #undef __FUNCT__ 2102c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2103c75a6043SHong Zhang /*@ 2104c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2105c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2106c75a6043SHong Zhang 2107c75a6043SHong Zhang Collective on MPI_Comm 2108c75a6043SHong Zhang 2109c75a6043SHong Zhang Input Parameters: 2110c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2111c75a6043SHong Zhang . bs - size of block 2112c75a6043SHong Zhang . m - number of rows 2113c75a6043SHong Zhang . n - number of columns 2114c75a6043SHong Zhang . i - row indices 2115c75a6043SHong Zhang . j - column indices 2116c75a6043SHong Zhang - a - matrix values 2117c75a6043SHong Zhang 2118c75a6043SHong Zhang Output Parameter: 2119c75a6043SHong Zhang . mat - the matrix 2120c75a6043SHong Zhang 2121c75a6043SHong Zhang Level: intermediate 2122c75a6043SHong Zhang 2123c75a6043SHong Zhang Notes: 2124c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2125c75a6043SHong Zhang once the matrix is destroyed 2126c75a6043SHong Zhang 2127c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2128c75a6043SHong Zhang 2129c75a6043SHong Zhang The i and j indices are 0 based 2130c75a6043SHong Zhang 2131c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2132c75a6043SHong Zhang 2133c75a6043SHong Zhang @*/ 2134c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2135c75a6043SHong Zhang { 2136c75a6043SHong Zhang PetscErrorCode ierr; 2137c75a6043SHong Zhang PetscInt ii; 2138c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2139c75a6043SHong Zhang 2140c75a6043SHong Zhang PetscFunctionBegin; 2141c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2142c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2143c75a6043SHong Zhang 2144c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2145c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2146c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2147c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2148c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2149c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2150c75a6043SHong Zhang 2151c75a6043SHong Zhang sbaij->i = i; 2152c75a6043SHong Zhang sbaij->j = j; 2153c75a6043SHong Zhang sbaij->a = a; 2154c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2155c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2156e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2157e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2158c75a6043SHong Zhang 2159c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2160c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2161c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2162c75a6043SHong 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]); 2163c75a6043SHong Zhang #endif 2164c75a6043SHong Zhang } 2165c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2166c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2167c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2168c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2169c75a6043SHong Zhang } 2170c75a6043SHong Zhang #endif 2171c75a6043SHong Zhang 2172c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2173c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2174c75a6043SHong Zhang PetscFunctionReturn(0); 2175c75a6043SHong Zhang } 2176d06b337dSHong Zhang 2177d06b337dSHong Zhang 2178d06b337dSHong Zhang 217949b5e25fSSatish Balay 218049b5e25fSSatish Balay 2181