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" 182af78befSBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscTruth *missing,PetscInt *dd) 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; 272af78befSBarry Smith *missing = PETSC_FALSE; 2849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 292af78befSBarry Smith if (jj[diag[i]] != i) { 302af78befSBarry Smith *missing = PETSC_TRUE; 312af78befSBarry Smith if (dd) *dd = i; 322af78befSBarry Smith break; 332af78befSBarry Smith } 3449b5e25fSSatish Balay } 3549b5e25fSSatish Balay PetscFunctionReturn(0); 3649b5e25fSSatish Balay } 3749b5e25fSSatish Balay 384a2ae208SSatish Balay #undef __FUNCT__ 394a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 40dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 4149b5e25fSSatish Balay { 42045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 436849ba73SBarry Smith PetscErrorCode ierr; 4409f38230SBarry Smith PetscInt i; 4549b5e25fSSatish Balay 4649b5e25fSSatish Balay PetscFunctionBegin; 4709f38230SBarry Smith if (!a->diag) { 4809f38230SBarry Smith ierr = PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 4909f38230SBarry Smith } 5009f38230SBarry Smith for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i]; 5149b5e25fSSatish Balay PetscFunctionReturn(0); 5249b5e25fSSatish Balay } 5349b5e25fSSatish Balay 544a2ae208SSatish Balay #undef __FUNCT__ 554a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 568f7157efSSatish Balay static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 5749b5e25fSSatish Balay { 58a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 598f7157efSSatish Balay PetscInt i,j,n = a->mbs,nz = a->i[n],bs = A->rmap.bs; 608f7157efSSatish Balay PetscErrorCode ierr; 6149b5e25fSSatish Balay 6249b5e25fSSatish Balay PetscFunctionBegin; 63d3e5a4abSHong Zhang *nn = n; 64a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 658f7157efSSatish Balay if (!blockcompressed) { 668f7157efSSatish Balay /* malloc & create the natural set of indices */ 67f1d0d59dSSatish Balay ierr = PetscMalloc2((n+1)*bs,PetscInt,ia,nz*bs,PetscInt,ja);CHKERRQ(ierr); 688f7157efSSatish Balay for (i=0; i<n+1; i++) { 698f7157efSSatish Balay for (j=0; j<bs; j++) { 708f7157efSSatish Balay *ia[i*bs+j] = a->i[i]*bs+j+oshift; 718f7157efSSatish Balay } 728f7157efSSatish Balay } 738f7157efSSatish Balay for (i=0; i<nz; i++) { 748f7157efSSatish Balay for (j=0; j<bs; j++) { 758f7157efSSatish Balay *ja[i*bs+j] = a->j[i]*bs+j+oshift; 768f7157efSSatish Balay } 778f7157efSSatish Balay } 788f7157efSSatish Balay } else { /* blockcompressed */ 79a6ece127SHong Zhang if (oshift == 1) { 80a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 816c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 82a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 838f7157efSSatish Balay } 84a1373b80SHong Zhang *ia = a->i; *ja = a->j; 85a6ece127SHong Zhang } 868f7157efSSatish Balay 8749b5e25fSSatish Balay PetscFunctionReturn(0); 8849b5e25fSSatish Balay } 8949b5e25fSSatish Balay 904a2ae208SSatish Balay #undef __FUNCT__ 914a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 928f7157efSSatish Balay static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 9349b5e25fSSatish Balay { 94b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 958f7157efSSatish Balay PetscInt i,n = a->mbs,nz = a->i[n]; 968f7157efSSatish Balay PetscErrorCode ierr; 97a6ece127SHong Zhang 9849b5e25fSSatish Balay PetscFunctionBegin; 9949b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 100a6ece127SHong Zhang 1018f7157efSSatish Balay if (!blockcompressed) { 1028f7157efSSatish Balay ierr = PetscFree2(*ia,*ja);CHKERRQ(ierr); 1038f7157efSSatish Balay } else if (oshift == 1) { /* blockcompressed */ 1046c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 105a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 106a6ece127SHong Zhang } 1078f7157efSSatish Balay 108a6ece127SHong Zhang PetscFunctionReturn(0); 10949b5e25fSSatish Balay } 11049b5e25fSSatish Balay 1114a2ae208SSatish Balay #undef __FUNCT__ 1124a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 113dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 11449b5e25fSSatish Balay { 11549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 116dfbe8321SBarry Smith PetscErrorCode ierr; 11749b5e25fSSatish Balay 11849b5e25fSSatish Balay PetscFunctionBegin; 119a9f03627SSatish Balay #if defined(PETSC_USE_LOG) 120899cda47SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap.N,a->nz); 121a9f03627SSatish Balay #endif 122e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr); 1239bfd6278SHong Zhang if (a->row) {ierr = ISDestroy(a->row);CHKERRQ(ierr);} 1249bfd6278SHong Zhang if (a->col){ierr = ISDestroy(a->col);CHKERRQ(ierr);} 1259bfd6278SHong Zhang if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 12605b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 12705b42c5fSBarry Smith ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr); 12805b42c5fSBarry Smith ierr = PetscFree(a->solve_work);CHKERRQ(ierr); 129e2ee2a47SBarry Smith ierr = PetscFree(a->relax_work);CHKERRQ(ierr); 13005b42c5fSBarry Smith ierr = PetscFree(a->solves_work);CHKERRQ(ierr); 13105b42c5fSBarry Smith ierr = PetscFree(a->mult_work);CHKERRQ(ierr); 13205b42c5fSBarry Smith ierr = PetscFree(a->saved_values);CHKERRQ(ierr); 13305b42c5fSBarry Smith ierr = PetscFree(a->xtoy);CHKERRQ(ierr); 1341a3463dfSHong Zhang 1351a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 13649b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 137901853e0SKris Buschelman 138dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 139901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 140901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 141901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 142901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 143901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 144901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 14549b5e25fSSatish Balay PetscFunctionReturn(0); 14649b5e25fSSatish Balay } 14749b5e25fSSatish Balay 1484a2ae208SSatish Balay #undef __FUNCT__ 1494a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 1504e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscTruth flg) 15149b5e25fSSatish Balay { 152045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 15363ba0a88SBarry Smith PetscErrorCode ierr; 15449b5e25fSSatish Balay 15549b5e25fSSatish Balay PetscFunctionBegin; 1564d9d31abSKris Buschelman switch (op) { 1574d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1584e0d8c25SBarry Smith a->roworiented = flg; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1614e0d8c25SBarry Smith a->keepzeroedrows = flg; 1624d9d31abSKris Buschelman break; 163512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 164512a5fc5SBarry Smith a->nonew = (flg ? 0 : 1); 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1674e0d8c25SBarry Smith a->nonew = (flg ? -1 : 0); 1684d9d31abSKris Buschelman break; 1694d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1704e0d8c25SBarry Smith a->nonew = (flg ? -2 : 0); 1714d9d31abSKris Buschelman break; 1724e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1734d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1744d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 175290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 1764d9d31abSKris Buschelman break; 1779a4540c5SBarry Smith case MAT_HERMITIAN: 1784e0d8c25SBarry Smith if (flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 17977e54ba9SKris Buschelman case MAT_SYMMETRIC: 18077e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1819a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1824e0d8c25SBarry Smith if (!flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 183290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);CHKERRQ(ierr); 184290bbb0aSBarry Smith break; 185941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 1864e0d8c25SBarry Smith a->ignore_ltriangular = flg; 187941593c8SHong Zhang break; 188941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 1894e0d8c25SBarry Smith a->ignore_ltriangular = flg; 19077e54ba9SKris Buschelman break; 191f5edf698SHong Zhang case MAT_GETROW_UPPERTRIANGULAR: 1924e0d8c25SBarry Smith a->getrow_utriangular = flg; 193f5edf698SHong Zhang break; 1944d9d31abSKris Buschelman default: 195ad86a440SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 19649b5e25fSSatish Balay } 19749b5e25fSSatish Balay PetscFunctionReturn(0); 19849b5e25fSSatish Balay } 19949b5e25fSSatish Balay 2004a2ae208SSatish Balay #undef __FUNCT__ 2014a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 20213f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 20349b5e25fSSatish Balay { 20449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2056849ba73SBarry Smith PetscErrorCode ierr; 20613f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 20749b5e25fSSatish Balay MatScalar *aa,*aa_i; 20887828ca2SBarry Smith PetscScalar *v_i; 20949b5e25fSSatish Balay 21049b5e25fSSatish Balay PetscFunctionBegin; 2114e0d8c25SBarry Smith if (A && !a->getrow_utriangular) SETERRQ(PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()"); 212f5edf698SHong Zhang /* Get the upper triangular part of the row */ 213899cda47SBarry Smith bs = A->rmap.bs; 21449b5e25fSSatish Balay ai = a->i; 21549b5e25fSSatish Balay aj = a->j; 21649b5e25fSSatish Balay aa = a->a; 21749b5e25fSSatish Balay bs2 = a->bs2; 21849b5e25fSSatish Balay 219899cda47SBarry Smith if (row < 0 || row >= A->rmap.N) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 22049b5e25fSSatish Balay 22149b5e25fSSatish Balay bn = row/bs; /* Block number */ 22249b5e25fSSatish Balay bp = row % bs; /* Block position */ 22349b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 22449b5e25fSSatish Balay *ncols = bs*M; 22549b5e25fSSatish Balay 22649b5e25fSSatish Balay if (v) { 22749b5e25fSSatish Balay *v = 0; 22849b5e25fSSatish Balay if (*ncols) { 22987828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 23049b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23149b5e25fSSatish Balay v_i = *v + i*bs; 23249b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 23349b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 23449b5e25fSSatish Balay } 23549b5e25fSSatish Balay } 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay 23849b5e25fSSatish Balay if (cols) { 23949b5e25fSSatish Balay *cols = 0; 24049b5e25fSSatish Balay if (*ncols) { 24113f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 24249b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 24349b5e25fSSatish Balay cols_i = *cols + i*bs; 24449b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 24549b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24649b5e25fSSatish Balay } 24749b5e25fSSatish Balay } 24849b5e25fSSatish Balay } 24949b5e25fSSatish Balay 25049b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2515ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2525ddb2528SHong Zhang #ifdef column_search 25349b5e25fSSatish Balay v_i = *v + M*bs; 25449b5e25fSSatish Balay cols_i = *cols + M*bs; 25549b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25649b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25749b5e25fSSatish Balay for (j=0; j<M; j++){ 25849b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 25949b5e25fSSatish Balay if (itmp == bn){ 26049b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 26149b5e25fSSatish Balay for (k=0; k<bs; k++) { 26249b5e25fSSatish Balay *cols_i++ = i*bs+k; 26349b5e25fSSatish Balay *v_i++ = aa_i[k]; 26449b5e25fSSatish Balay } 26549b5e25fSSatish Balay *ncols += bs; 26649b5e25fSSatish Balay break; 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay } 26949b5e25fSSatish Balay } 2705ddb2528SHong Zhang #endif 27149b5e25fSSatish Balay PetscFunctionReturn(0); 27249b5e25fSSatish Balay } 27349b5e25fSSatish Balay 2744a2ae208SSatish Balay #undef __FUNCT__ 2754a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 27613f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27749b5e25fSSatish Balay { 278dfbe8321SBarry Smith PetscErrorCode ierr; 27949b5e25fSSatish Balay 28049b5e25fSSatish Balay PetscFunctionBegin; 28105b42c5fSBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 28205b42c5fSBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 28349b5e25fSSatish Balay PetscFunctionReturn(0); 28449b5e25fSSatish Balay } 28549b5e25fSSatish Balay 2864a2ae208SSatish Balay #undef __FUNCT__ 287f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ" 288f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 289f5edf698SHong Zhang { 290f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 291f5edf698SHong Zhang 292f5edf698SHong Zhang PetscFunctionBegin; 293f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 294f5edf698SHong Zhang PetscFunctionReturn(0); 295f5edf698SHong Zhang } 296f5edf698SHong Zhang #undef __FUNCT__ 297f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ" 298f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 299f5edf698SHong Zhang { 300f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 301f5edf698SHong Zhang 302f5edf698SHong Zhang PetscFunctionBegin; 303f5edf698SHong Zhang a->getrow_utriangular = PETSC_FALSE; 304f5edf698SHong Zhang PetscFunctionReturn(0); 305f5edf698SHong Zhang } 306f5edf698SHong Zhang 307f5edf698SHong Zhang #undef __FUNCT__ 3084a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 309fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B) 31049b5e25fSSatish Balay { 311dfbe8321SBarry Smith PetscErrorCode ierr; 31249b5e25fSSatish Balay PetscFunctionBegin; 313815cbec1SBarry Smith if (reuse == MAT_INITIAL_MATRIX || *B != A) { 314999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 315fc4dec0aSBarry Smith } 3168115998fSBarry Smith PetscFunctionReturn(0); 31749b5e25fSSatish Balay } 31849b5e25fSSatish Balay 3194a2ae208SSatish Balay #undef __FUNCT__ 3204a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3216849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 32249b5e25fSSatish Balay { 32349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 324dfbe8321SBarry Smith PetscErrorCode ierr; 325899cda47SBarry Smith PetscInt i,j,bs = A->rmap.bs,k,l,bs2=a->bs2; 3262dcb1b2aSMatthew Knepley const char *name; 327f3ef73ceSBarry Smith PetscViewerFormat format; 32849b5e25fSSatish Balay 32949b5e25fSSatish Balay PetscFunctionBegin; 33080fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 331b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 332456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 33377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 334fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 335d2507d54SMatthew Knepley Mat aij; 336d2507d54SMatthew Knepley 33770d5e725SHong Zhang if (A->factor && bs>1){ 33870d5e725SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr); 33970d5e725SHong Zhang PetscFunctionReturn(0); 34070d5e725SHong Zhang } 341c9f458caSMatthew Knepley ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr); 342c9f458caSMatthew Knepley ierr = MatView(aij,viewer);CHKERRQ(ierr); 343c9f458caSMatthew Knepley ierr = MatDestroy(aij);CHKERRQ(ierr); 344fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 345b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 34649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 34749b5e25fSSatish Balay for (j=0; j<bs; j++) { 34877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 34949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 35049b5e25fSSatish Balay for (l=0; l<bs; l++) { 35149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 35249b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 353a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 35449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35549b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 356a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 35749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35849b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 359a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36049b5e25fSSatish Balay } 36149b5e25fSSatish Balay #else 36249b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 363a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 36449b5e25fSSatish Balay } 36549b5e25fSSatish Balay #endif 36649b5e25fSSatish Balay } 36749b5e25fSSatish Balay } 368b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 36949b5e25fSSatish Balay } 37049b5e25fSSatish Balay } 371b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 372c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 373c1490034SHong Zhang PetscFunctionReturn(0); 37449b5e25fSSatish Balay } else { 3758608aa04SHong Zhang if (A->factor && bs>1){ 3768608aa04SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored. MatView_SeqSBAIJ_ASCII() may not display complete or logically correct entries!\n");CHKERRQ(ierr); 3778608aa04SHong Zhang } 378b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 37949b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 38049b5e25fSSatish Balay for (j=0; j<bs; j++) { 38177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 38249b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 38349b5e25fSSatish Balay for (l=0; l<bs; l++) { 38449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 38549b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 386a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 38749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38849b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 389a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 39049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39149b5e25fSSatish Balay } else { 392a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39349b5e25fSSatish Balay } 39449b5e25fSSatish Balay #else 395e9f7bc9eSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 39649b5e25fSSatish Balay #endif 39749b5e25fSSatish Balay } 39849b5e25fSSatish Balay } 399b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 40049b5e25fSSatish Balay } 40149b5e25fSSatish Balay } 402b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 40349b5e25fSSatish Balay } 404b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 40549b5e25fSSatish Balay PetscFunctionReturn(0); 40649b5e25fSSatish Balay } 40749b5e25fSSatish Balay 4084a2ae208SSatish Balay #undef __FUNCT__ 4094a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 4106849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 41149b5e25fSSatish Balay { 41249b5e25fSSatish Balay Mat A = (Mat) Aa; 41349b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 4146849ba73SBarry Smith PetscErrorCode ierr; 415899cda47SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap.bs,bs2=a->bs2; 41613f74950SBarry Smith PetscMPIInt rank; 41749b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 41849b5e25fSSatish Balay MatScalar *aa; 41949b5e25fSSatish Balay MPI_Comm comm; 420b0a32e0cSBarry Smith PetscViewer viewer; 42149b5e25fSSatish Balay 42249b5e25fSSatish Balay PetscFunctionBegin; 42349b5e25fSSatish Balay /* 42449b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 42549b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 42649b5e25fSSatish Balay rest should return immediately. 42749b5e25fSSatish Balay */ 42849b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 42949b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 43049b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 43149b5e25fSSatish Balay 43249b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 43349b5e25fSSatish Balay 434b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 435b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 43649b5e25fSSatish Balay 43749b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 438b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 43949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 44049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 441899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 44249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 44349b5e25fSSatish Balay aa = a->a + j*bs2; 44449b5e25fSSatish Balay for (k=0; k<bs; k++) { 44549b5e25fSSatish Balay for (l=0; l<bs; l++) { 44649b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 447b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44849b5e25fSSatish Balay } 44949b5e25fSSatish Balay } 45049b5e25fSSatish Balay } 45149b5e25fSSatish Balay } 452b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 45349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 455899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 45649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45749b5e25fSSatish Balay aa = a->a + j*bs2; 45849b5e25fSSatish Balay for (k=0; k<bs; k++) { 45949b5e25fSSatish Balay for (l=0; l<bs; l++) { 46049b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 461b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46249b5e25fSSatish Balay } 46349b5e25fSSatish Balay } 46449b5e25fSSatish Balay } 46549b5e25fSSatish Balay } 46649b5e25fSSatish Balay 467b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 46849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 470899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 47149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47249b5e25fSSatish Balay aa = a->a + j*bs2; 47349b5e25fSSatish Balay for (k=0; k<bs; k++) { 47449b5e25fSSatish Balay for (l=0; l<bs; l++) { 47549b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 476b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 47749b5e25fSSatish Balay } 47849b5e25fSSatish Balay } 47949b5e25fSSatish Balay } 48049b5e25fSSatish Balay } 48149b5e25fSSatish Balay PetscFunctionReturn(0); 48249b5e25fSSatish Balay } 48349b5e25fSSatish Balay 4844a2ae208SSatish Balay #undef __FUNCT__ 4854a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4866849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 48749b5e25fSSatish Balay { 488dfbe8321SBarry Smith PetscErrorCode ierr; 48949b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 490b0a32e0cSBarry Smith PetscDraw draw; 49149b5e25fSSatish Balay PetscTruth isnull; 49249b5e25fSSatish Balay 49349b5e25fSSatish Balay PetscFunctionBegin; 494b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 495b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 49649b5e25fSSatish Balay 49749b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 498899cda47SBarry Smith xr = A->rmap.N; yr = A->rmap.N; h = yr/10.0; w = xr/10.0; 49949b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 500b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 501b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 50249b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 50349b5e25fSSatish Balay PetscFunctionReturn(0); 50449b5e25fSSatish Balay } 50549b5e25fSSatish Balay 5064a2ae208SSatish Balay #undef __FUNCT__ 5074a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 508dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 50949b5e25fSSatish Balay { 510dfbe8321SBarry Smith PetscErrorCode ierr; 51132077d6dSBarry Smith PetscTruth iascii,isdraw; 51249b5e25fSSatish Balay 51349b5e25fSSatish Balay PetscFunctionBegin; 51432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 515fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 51632077d6dSBarry Smith if (iascii){ 51749b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 51849b5e25fSSatish Balay } else if (isdraw) { 51949b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 52049b5e25fSSatish Balay } else { 521a5e6ed63SBarry Smith Mat B; 522ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 523a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 524a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 52549b5e25fSSatish Balay } 52649b5e25fSSatish Balay PetscFunctionReturn(0); 52749b5e25fSSatish Balay } 52849b5e25fSSatish Balay 52949b5e25fSSatish Balay 5304a2ae208SSatish Balay #undef __FUNCT__ 5314a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 53213f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 53349b5e25fSSatish Balay { 534045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 53513f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 53613f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 537899cda47SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->rmap.bs,bs2=a->bs2; 53897e567efSBarry Smith MatScalar *ap,*aa = a->a; 53949b5e25fSSatish Balay 54049b5e25fSSatish Balay PetscFunctionBegin; 54149b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 54249b5e25fSSatish Balay row = im[k]; brow = row/bs; 54397e567efSBarry Smith if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */ 544899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 54549b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 54649b5e25fSSatish Balay nrow = ailen[brow]; 54749b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 54897e567efSBarry Smith if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */ 549899cda47SBarry 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); 55049b5e25fSSatish Balay col = in[l] ; 55149b5e25fSSatish Balay bcol = col/bs; 55249b5e25fSSatish Balay cidx = col%bs; 55349b5e25fSSatish Balay ridx = row%bs; 55449b5e25fSSatish Balay high = nrow; 55549b5e25fSSatish Balay low = 0; /* assume unsorted */ 55649b5e25fSSatish Balay while (high-low > 5) { 55749b5e25fSSatish Balay t = (low+high)/2; 55849b5e25fSSatish Balay if (rp[t] > bcol) high = t; 55949b5e25fSSatish Balay else low = t; 56049b5e25fSSatish Balay } 56149b5e25fSSatish Balay for (i=low; i<high; i++) { 56249b5e25fSSatish Balay if (rp[i] > bcol) break; 56349b5e25fSSatish Balay if (rp[i] == bcol) { 56449b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 56549b5e25fSSatish Balay goto finished; 56649b5e25fSSatish Balay } 56749b5e25fSSatish Balay } 56897e567efSBarry Smith *v++ = 0.0; 56949b5e25fSSatish Balay finished:; 57049b5e25fSSatish Balay } 57149b5e25fSSatish Balay } 57249b5e25fSSatish Balay PetscFunctionReturn(0); 57349b5e25fSSatish Balay } 57449b5e25fSSatish Balay 57549b5e25fSSatish Balay 5764a2ae208SSatish Balay #undef __FUNCT__ 5774a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 57813f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 57949b5e25fSSatish Balay { 5800880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5816849ba73SBarry Smith PetscErrorCode ierr; 582e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 58313f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 584899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap.bs,stepval; 5850880e062SHong Zhang PetscTruth roworiented=a->roworiented; 586dd6ea824SBarry Smith const PetscScalar *value = v; 587f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5880880e062SHong Zhang 58949b5e25fSSatish Balay PetscFunctionBegin; 5900880e062SHong Zhang if (roworiented) { 5910880e062SHong Zhang stepval = (n-1)*bs; 5920880e062SHong Zhang } else { 5930880e062SHong Zhang stepval = (m-1)*bs; 5940880e062SHong Zhang } 5950880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5960880e062SHong Zhang row = im[k]; 5970880e062SHong Zhang if (row < 0) continue; 5982515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 59977431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 6000880e062SHong Zhang #endif 6010880e062SHong Zhang rp = aj + ai[row]; 6020880e062SHong Zhang ap = aa + bs2*ai[row]; 6030880e062SHong Zhang rmax = imax[row]; 6040880e062SHong Zhang nrow = ailen[row]; 6050880e062SHong Zhang low = 0; 606818f2c47SBarry Smith high = nrow; 6070880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6080880e062SHong Zhang if (in[l] < 0) continue; 6090880e062SHong Zhang col = in[l]; 6102515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 61177431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 612b1823623SSatish Balay #endif 6130880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6140880e062SHong Zhang if (roworiented) { 6150880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6160880e062SHong Zhang } else { 6170880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6180880e062SHong Zhang } 6197cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 620e2ee6c50SBarry Smith lastcol = col; 6210880e062SHong Zhang while (high-low > 7) { 6220880e062SHong Zhang t = (low+high)/2; 6230880e062SHong Zhang if (rp[t] > col) high = t; 6240880e062SHong Zhang else low = t; 6250880e062SHong Zhang } 6260880e062SHong Zhang for (i=low; i<high; i++) { 6270880e062SHong Zhang if (rp[i] > col) break; 6280880e062SHong Zhang if (rp[i] == col) { 6290880e062SHong Zhang bap = ap + bs2*i; 6300880e062SHong Zhang if (roworiented) { 6310880e062SHong Zhang if (is == ADD_VALUES) { 6320880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6330880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6340880e062SHong Zhang bap[jj] += *value++; 6350880e062SHong Zhang } 6360880e062SHong Zhang } 6370880e062SHong Zhang } else { 6380880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6390880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6400880e062SHong Zhang bap[jj] = *value++; 6410880e062SHong Zhang } 6420880e062SHong Zhang } 6430880e062SHong Zhang } 6440880e062SHong Zhang } else { 6450880e062SHong Zhang if (is == ADD_VALUES) { 6460880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6470880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6480880e062SHong Zhang *bap++ += *value++; 6490880e062SHong Zhang } 6500880e062SHong Zhang } 6510880e062SHong Zhang } else { 6520880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6530880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6540880e062SHong Zhang *bap++ = *value++; 6550880e062SHong Zhang } 6560880e062SHong Zhang } 6570880e062SHong Zhang } 6580880e062SHong Zhang } 6590880e062SHong Zhang goto noinsert2; 6600880e062SHong Zhang } 6610880e062SHong Zhang } 6620880e062SHong Zhang if (nonew == 1) goto noinsert2; 663085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 664421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 665c03d1d03SSatish Balay N = nrow++ - 1; high++; 6660880e062SHong Zhang /* shift up all the later entries in this row */ 6670880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6680880e062SHong Zhang rp[ii+1] = rp[ii]; 6690880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6700880e062SHong Zhang } 6710880e062SHong Zhang if (N >= i) { 6720880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6730880e062SHong Zhang } 6740880e062SHong Zhang rp[i] = col; 6750880e062SHong Zhang bap = ap + bs2*i; 6760880e062SHong Zhang if (roworiented) { 6770880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6780880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6790880e062SHong Zhang bap[jj] = *value++; 6800880e062SHong Zhang } 6810880e062SHong Zhang } 6820880e062SHong Zhang } else { 6830880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6840880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6850880e062SHong Zhang *bap++ = *value++; 6860880e062SHong Zhang } 6870880e062SHong Zhang } 6880880e062SHong Zhang } 6890880e062SHong Zhang noinsert2:; 6900880e062SHong Zhang low = i; 6910880e062SHong Zhang } 6920880e062SHong Zhang ailen[row] = nrow; 6930880e062SHong Zhang } 6940880e062SHong Zhang PetscFunctionReturn(0); 69549b5e25fSSatish Balay } 69649b5e25fSSatish Balay 6974a2ae208SSatish Balay #undef __FUNCT__ 6984a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 699dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 70049b5e25fSSatish Balay { 70149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7026849ba73SBarry Smith PetscErrorCode ierr; 70313f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 704899cda47SBarry Smith PetscInt m = A->rmap.N,*ip,N,*ailen = a->ilen; 70513f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 70649b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 70749b5e25fSSatish Balay 70849b5e25fSSatish Balay PetscFunctionBegin; 70949b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 71049b5e25fSSatish Balay 71149b5e25fSSatish Balay if (m) rmax = ailen[0]; 71249b5e25fSSatish Balay for (i=1; i<mbs; i++) { 71349b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 71449b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 71549b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 71649b5e25fSSatish Balay if (fshift) { 71749b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 71849b5e25fSSatish Balay N = ailen[i]; 71949b5e25fSSatish Balay for (j=0; j<N; j++) { 72049b5e25fSSatish Balay ip[j-fshift] = ip[j]; 72149b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 72249b5e25fSSatish Balay } 72349b5e25fSSatish Balay } 72449b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 72549b5e25fSSatish Balay } 72649b5e25fSSatish Balay if (mbs) { 72749b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 72849b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 72949b5e25fSSatish Balay } 73049b5e25fSSatish Balay /* reset ilen and imax for each row */ 73149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 73249b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 73349b5e25fSSatish Balay } 7346c6c5352SBarry Smith a->nz = ai[mbs]; 73549b5e25fSSatish Balay 736b424e231SHong Zhang /* diagonals may have moved, reset it */ 737b424e231SHong Zhang if (a->diag) { 73813f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 73949b5e25fSSatish Balay } 740899cda47SBarry 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); 741ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 742ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 74349b5e25fSSatish Balay a->reallocs = 0; 74449b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 74549b5e25fSSatish Balay PetscFunctionReturn(0); 74649b5e25fSSatish Balay } 74749b5e25fSSatish Balay 74849b5e25fSSatish Balay /* 74949b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 75049b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 75149b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 75249b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 75349b5e25fSSatish Balay */ 7544a2ae208SSatish Balay #undef __FUNCT__ 7554a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 75613f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 75749b5e25fSSatish Balay { 75813f74950SBarry Smith PetscInt i,j,k,row; 75949b5e25fSSatish Balay PetscTruth flg; 76049b5e25fSSatish Balay 76149b5e25fSSatish Balay PetscFunctionBegin; 76249b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 76349b5e25fSSatish Balay row = idx[i]; 76449b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 76549b5e25fSSatish Balay sizes[j] = 1; 76649b5e25fSSatish Balay i++; 76749b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 76849b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 76949b5e25fSSatish Balay i++; 77049b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 77149b5e25fSSatish Balay flg = PETSC_TRUE; 77249b5e25fSSatish Balay for (k=1; k<bs; k++) { 77349b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 77449b5e25fSSatish Balay flg = PETSC_FALSE; 77549b5e25fSSatish Balay break; 77649b5e25fSSatish Balay } 77749b5e25fSSatish Balay } 778abc0a331SBarry Smith if (flg) { /* No break in the bs */ 77949b5e25fSSatish Balay sizes[j] = bs; 78049b5e25fSSatish Balay i+= bs; 78149b5e25fSSatish Balay } else { 78249b5e25fSSatish Balay sizes[j] = 1; 78349b5e25fSSatish Balay i++; 78449b5e25fSSatish Balay } 78549b5e25fSSatish Balay } 78649b5e25fSSatish Balay } 78749b5e25fSSatish Balay *bs_max = j; 78849b5e25fSSatish Balay PetscFunctionReturn(0); 78949b5e25fSSatish Balay } 79049b5e25fSSatish Balay 79149b5e25fSSatish Balay 79249b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 79349b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 79449b5e25fSSatish Balay */ 79549b5e25fSSatish Balay 7964a2ae208SSatish Balay #undef __FUNCT__ 7974a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 79813f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 79949b5e25fSSatish Balay { 80049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8016849ba73SBarry Smith PetscErrorCode ierr; 802e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 80313f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 804899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 80513f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 80649b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 80749b5e25fSSatish Balay 80849b5e25fSSatish Balay PetscFunctionBegin; 80949b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 81049b5e25fSSatish Balay row = im[k]; /* row number */ 81149b5e25fSSatish Balay brow = row/bs; /* block row number */ 81249b5e25fSSatish Balay if (row < 0) continue; 8132515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 814899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 81549b5e25fSSatish Balay #endif 81649b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 81749b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 81849b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 81949b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 82049b5e25fSSatish Balay low = 0; 82149b5e25fSSatish Balay 82249b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 82349b5e25fSSatish Balay if (in[l] < 0) continue; 8242515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 825899cda47SBarry 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); 82649b5e25fSSatish Balay #endif 82749b5e25fSSatish Balay col = in[l]; 82849b5e25fSSatish Balay bcol = col/bs; /* block col number */ 82949b5e25fSSatish Balay 830941593c8SHong Zhang if (brow > bcol) { 831941593c8SHong Zhang if (a->ignore_ltriangular){ 832941593c8SHong Zhang continue; /* ignore lower triangular values */ 833941593c8SHong Zhang } else { 8344e0d8c25SBarry Smith SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)"); 835941593c8SHong Zhang } 836941593c8SHong Zhang } 837f4989cb3SHong Zhang 83849b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8398549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 84049b5e25fSSatish Balay /* element value a(k,l) */ 84149b5e25fSSatish Balay if (roworiented) { 84249b5e25fSSatish Balay value = v[l + k*n]; 84349b5e25fSSatish Balay } else { 84449b5e25fSSatish Balay value = v[k + l*m]; 84549b5e25fSSatish Balay } 84649b5e25fSSatish Balay 84749b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 8487cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 849e2ee6c50SBarry Smith lastcol = col; 85049b5e25fSSatish Balay while (high-low > 7) { 85149b5e25fSSatish Balay t = (low+high)/2; 85249b5e25fSSatish Balay if (rp[t] > bcol) high = t; 85349b5e25fSSatish Balay else low = t; 85449b5e25fSSatish Balay } 85549b5e25fSSatish Balay for (i=low; i<high; i++) { 85649b5e25fSSatish Balay if (rp[i] > bcol) break; 85749b5e25fSSatish Balay if (rp[i] == bcol) { 85849b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 85949b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 86049b5e25fSSatish Balay else *bap = value; 8618549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8628549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8638549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8648549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8658549e402SHong Zhang else *bap = value; 8668549e402SHong Zhang } 86749b5e25fSSatish Balay goto noinsert1; 86849b5e25fSSatish Balay } 86949b5e25fSSatish Balay } 87049b5e25fSSatish Balay 87149b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 872085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 873421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 87449b5e25fSSatish Balay 875c03d1d03SSatish Balay N = nrow++ - 1; high++; 87649b5e25fSSatish Balay /* shift up all the later entries in this row */ 87749b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 87849b5e25fSSatish Balay rp[ii+1] = rp[ii]; 87949b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 88049b5e25fSSatish Balay } 88149b5e25fSSatish Balay if (N>=i) { 88249b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 88349b5e25fSSatish Balay } 88449b5e25fSSatish Balay rp[i] = bcol; 88549b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 88649b5e25fSSatish Balay noinsert1:; 88749b5e25fSSatish Balay low = i; 8888549e402SHong Zhang } 88949b5e25fSSatish Balay } /* end of loop over added columns */ 89049b5e25fSSatish Balay ailen[brow] = nrow; 89149b5e25fSSatish Balay } /* end of loop over added rows */ 89249b5e25fSSatish Balay PetscFunctionReturn(0); 89349b5e25fSSatish Balay } 89449b5e25fSSatish Balay 8954a2ae208SSatish Balay #undef __FUNCT__ 8964d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 897dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 89849b5e25fSSatish Balay { 8994ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 90049b5e25fSSatish Balay Mat outA; 901dfbe8321SBarry Smith PetscErrorCode ierr; 902c84f5b01SHong Zhang PetscTruth row_identity; 90349b5e25fSSatish Balay 90449b5e25fSSatish Balay PetscFunctionBegin; 905c84f5b01SHong Zhang if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 906c84f5b01SHong Zhang ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 907c84f5b01SHong Zhang if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 908c84f5b01SHong 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()! */ 909c84f5b01SHong Zhang 91049b5e25fSSatish Balay outA = inA; 911*5c9eb25fSBarry Smith inA->factor = MAT_FACTOR_CHOLESKY; 91249b5e25fSSatish Balay 9131a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 91449b5e25fSSatish Balay /* 915c84f5b01SHong Zhang Blocksize < 8 have a special faster factorization/solver 916c84f5b01SHong Zhang for ICC(0) factorization with natural ordering 91749b5e25fSSatish Balay */ 918c84f5b01SHong Zhang switch (inA->rmap.bs){ /* Note: row_identity = PETSC_TRUE! */ 91949b5e25fSSatish Balay case 1: 920c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 9210fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 92207f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 923d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 924759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 925759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 926c84f5b01SHong Zhang ierr = PetscInfo(inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr); 927c84f5b01SHong Zhang break; 92849b5e25fSSatish Balay case 2: 929c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9301a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9311a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 932759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 933759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 934ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 93549b5e25fSSatish Balay break; 93649b5e25fSSatish Balay case 3: 937c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9381a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9391a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 940759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering; 941759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering; 942ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 94349b5e25fSSatish Balay break; 94449b5e25fSSatish Balay case 4: 945c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 9461a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 9471a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 948759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering; 949759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering; 950ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 95149b5e25fSSatish Balay break; 95249b5e25fSSatish Balay case 5: 953c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 9541a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 9551a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 956759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering; 957759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering; 958ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 95949b5e25fSSatish Balay break; 96049b5e25fSSatish Balay case 6: 961c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 9621a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 9631a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 964759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering; 965759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering; 966ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 96749b5e25fSSatish Balay break; 96849b5e25fSSatish Balay case 7: 969c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 9701a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 9711a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 972759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering; 973759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering; 974ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 97549b5e25fSSatish Balay break; 97649b5e25fSSatish Balay default: 977c84f5b01SHong Zhang inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 978c84f5b01SHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering; 979c84f5b01SHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 980759322afSHong Zhang inA->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 981759322afSHong Zhang inA->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 982c84f5b01SHong Zhang break; 983c84f5b01SHong Zhang } 98449b5e25fSSatish Balay 985c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 986c3122656SLisandro Dalcin if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } 987c84f5b01SHong Zhang a->row = row; 988c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 989c3122656SLisandro Dalcin if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } 990c84f5b01SHong Zhang a->col = row; 991c84f5b01SHong Zhang 992c84f5b01SHong Zhang /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 993c84f5b01SHong Zhang if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 99452e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 99549b5e25fSSatish Balay 99649b5e25fSSatish Balay if (!a->solve_work) { 997c84f5b01SHong Zhang ierr = PetscMalloc((inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 998c84f5b01SHong Zhang ierr = PetscLogObjectMemory(inA,(inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar));CHKERRQ(ierr); 99949b5e25fSSatish Balay } 100049b5e25fSSatish Balay 1001af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 100249b5e25fSSatish Balay PetscFunctionReturn(0); 100349b5e25fSSatish Balay } 1004950f1e5bSHong Zhang 100549b5e25fSSatish Balay EXTERN_C_BEGIN 10064a2ae208SSatish Balay #undef __FUNCT__ 10074a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1008be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 100949b5e25fSSatish Balay { 1010045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 101113f74950SBarry Smith PetscInt i,nz,n; 101249b5e25fSSatish Balay 101349b5e25fSSatish Balay PetscFunctionBegin; 10146c6c5352SBarry Smith nz = baij->maxnz; 1015899cda47SBarry Smith n = mat->cmap.n; 101649b5e25fSSatish Balay for (i=0; i<nz; i++) { 101749b5e25fSSatish Balay baij->j[i] = indices[i]; 101849b5e25fSSatish Balay } 10196c6c5352SBarry Smith baij->nz = nz; 102049b5e25fSSatish Balay for (i=0; i<n; i++) { 102149b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 102249b5e25fSSatish Balay } 102349b5e25fSSatish Balay PetscFunctionReturn(0); 102449b5e25fSSatish Balay } 102549b5e25fSSatish Balay EXTERN_C_END 102649b5e25fSSatish Balay 10274a2ae208SSatish Balay #undef __FUNCT__ 10284a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 102949b5e25fSSatish Balay /*@ 103019585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 103149b5e25fSSatish Balay in the matrix. 103249b5e25fSSatish Balay 103349b5e25fSSatish Balay Input Parameters: 103419585528SSatish Balay + mat - the SeqSBAIJ matrix 103549b5e25fSSatish Balay - indices - the column indices 103649b5e25fSSatish Balay 103749b5e25fSSatish Balay Level: advanced 103849b5e25fSSatish Balay 103949b5e25fSSatish Balay Notes: 104049b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 104149b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 104249b5e25fSSatish Balay of the MatSetValues() operation. 104349b5e25fSSatish Balay 104449b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1045d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 104649b5e25fSSatish Balay 1047ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 104849b5e25fSSatish Balay 1049ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 105049b5e25fSSatish Balay @*/ 1051be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 105249b5e25fSSatish Balay { 105313f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 105449b5e25fSSatish Balay 105549b5e25fSSatish Balay PetscFunctionBegin; 10564482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 10574482741eSBarry Smith PetscValidPointer(indices,2); 1058c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 105949b5e25fSSatish Balay if (f) { 106049b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 106149b5e25fSSatish Balay } else { 1062e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 106349b5e25fSSatish Balay } 106449b5e25fSSatish Balay PetscFunctionReturn(0); 106549b5e25fSSatish Balay } 106649b5e25fSSatish Balay 10674a2ae208SSatish Balay #undef __FUNCT__ 10683c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 10693c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 10703c896bc6SHong Zhang { 10713c896bc6SHong Zhang PetscErrorCode ierr; 10723c896bc6SHong Zhang 10733c896bc6SHong Zhang PetscFunctionBegin; 10743c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 10753c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 10763c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10773c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 10783c896bc6SHong Zhang 1079899cda47SBarry Smith if (a->i[A->rmap.N] != b->i[B->rmap.N]) { 10803c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 10813c896bc6SHong Zhang } 1082899cda47SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.N])*sizeof(PetscScalar));CHKERRQ(ierr); 10833c896bc6SHong Zhang } else { 1084f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 10853c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1086f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 10873c896bc6SHong Zhang } 10883c896bc6SHong Zhang PetscFunctionReturn(0); 10893c896bc6SHong Zhang } 10903c896bc6SHong Zhang 10913c896bc6SHong Zhang #undef __FUNCT__ 10924a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1093dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1094273d9f13SBarry Smith { 1095dfbe8321SBarry Smith PetscErrorCode ierr; 1096273d9f13SBarry Smith 1097273d9f13SBarry Smith PetscFunctionBegin; 10987edd0491SSatish Balay ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1099273d9f13SBarry Smith PetscFunctionReturn(0); 1100273d9f13SBarry Smith } 1101273d9f13SBarry Smith 1102a6ece127SHong Zhang #undef __FUNCT__ 1103a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1104dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1105a6ece127SHong Zhang { 1106a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1107a6ece127SHong Zhang PetscFunctionBegin; 1108a6ece127SHong Zhang *array = a->a; 1109a6ece127SHong Zhang PetscFunctionReturn(0); 1110a6ece127SHong Zhang } 1111a6ece127SHong Zhang 1112a6ece127SHong Zhang #undef __FUNCT__ 1113a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1114dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1115a6ece127SHong Zhang { 1116a6ece127SHong Zhang PetscFunctionBegin; 1117a6ece127SHong Zhang PetscFunctionReturn(0); 1118a6ece127SHong Zhang } 1119a6ece127SHong Zhang 112042ee4b1aSHong Zhang #include "petscblaslapack.h" 112142ee4b1aSHong Zhang #undef __FUNCT__ 112242ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1123f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 112442ee4b1aSHong Zhang { 112542ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1126dfbe8321SBarry Smith PetscErrorCode ierr; 1127899cda47SBarry Smith PetscInt i,bs=Y->rmap.bs,bs2,j; 11280805154bSBarry Smith PetscBLASInt one = 1,bnz = PetscBLASIntCast(x->nz); 112942ee4b1aSHong Zhang 113042ee4b1aSHong Zhang PetscFunctionBegin; 113142ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1132f4df32b1SMatthew Knepley PetscScalar alpha = a; 1133f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1134c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1135c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1136c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1137c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1138c537a176SHong Zhang } 1139c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1140c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1141c4319e64SHong Zhang y->XtoY = X; 1142c537a176SHong Zhang } 1143c4319e64SHong Zhang bs2 = bs*bs; 11446c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1145c4319e64SHong Zhang j = 0; 1146c4319e64SHong Zhang while (j < bs2){ 1147f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1148c4319e64SHong Zhang j++; 1149c537a176SHong Zhang } 1150c4319e64SHong Zhang } 11511e2582c4SBarry Smith ierr = PetscInfo3(Y,"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); 115242ee4b1aSHong Zhang } else { 1153f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1154f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1155f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 115642ee4b1aSHong Zhang } 115742ee4b1aSHong Zhang PetscFunctionReturn(0); 115842ee4b1aSHong Zhang } 115942ee4b1aSHong Zhang 1160efcf0fc3SBarry Smith #undef __FUNCT__ 1161efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1162dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1163efcf0fc3SBarry Smith { 1164efcf0fc3SBarry Smith PetscFunctionBegin; 1165efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1166efcf0fc3SBarry Smith PetscFunctionReturn(0); 1167efcf0fc3SBarry Smith } 1168efcf0fc3SBarry Smith 1169efcf0fc3SBarry Smith #undef __FUNCT__ 1170efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1171dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1172efcf0fc3SBarry Smith { 1173efcf0fc3SBarry Smith PetscFunctionBegin; 1174efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1175efcf0fc3SBarry Smith PetscFunctionReturn(0); 1176efcf0fc3SBarry Smith } 1177efcf0fc3SBarry Smith 1178efcf0fc3SBarry Smith #undef __FUNCT__ 1179efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1180ab5e4463SMatthew Knepley PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1181efcf0fc3SBarry Smith { 1182efcf0fc3SBarry Smith PetscFunctionBegin; 1183efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1184efcf0fc3SBarry Smith PetscFunctionReturn(0); 1185efcf0fc3SBarry Smith } 1186efcf0fc3SBarry Smith 118799cafbc1SBarry Smith #undef __FUNCT__ 118899cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 118999cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 119099cafbc1SBarry Smith { 119199cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 119299cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 1193dd6ea824SBarry Smith MatScalar *aa = a->a; 119499cafbc1SBarry Smith 119599cafbc1SBarry Smith PetscFunctionBegin; 119699cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 119799cafbc1SBarry Smith PetscFunctionReturn(0); 119899cafbc1SBarry Smith } 119999cafbc1SBarry Smith 120099cafbc1SBarry Smith #undef __FUNCT__ 120199cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 120299cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 120399cafbc1SBarry Smith { 120499cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 120599cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 1206dd6ea824SBarry Smith MatScalar *aa = a->a; 120799cafbc1SBarry Smith 120899cafbc1SBarry Smith PetscFunctionBegin; 120999cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 121099cafbc1SBarry Smith PetscFunctionReturn(0); 121199cafbc1SBarry Smith } 121299cafbc1SBarry Smith 121349b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 121449b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 121549b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 121649b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 121749b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 121897304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1219431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1220e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 122149b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 122249b5e25fSSatish Balay 0, 122349b5e25fSSatish Balay 0, 122497304618SKris Buschelman /*10*/ 0, 122549b5e25fSSatish Balay 0, 12265f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1227d06b337dSHong Zhang MatRelax_SeqSBAIJ, 122849b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 122997304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 123049b5e25fSSatish Balay MatEqual_SeqSBAIJ, 123149b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 123249b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 123349b5e25fSSatish Balay MatNorm_SeqSBAIJ, 123497304618SKris Buschelman /*20*/ 0, 123549b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 123649b5e25fSSatish Balay 0, 123749b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 123849b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1239dcf5cc72SBarry Smith /*25*/ 0, 124049b5e25fSSatish Balay 0, 124149b5e25fSSatish Balay 0, 12425f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12435f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 124497304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1245c464158bSHong Zhang 0, 12464d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1247a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1248a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 124997304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 12509495ac64SHong Zhang MatForwardSolve_SeqSBAIJ_N, 12519495ac64SHong Zhang MatBackwardSolve_SeqSBAIJ_N, 125249b5e25fSSatish Balay 0, 1253c84f5b01SHong Zhang MatICCFactor_SeqSBAIJ, 125497304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 125549b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 125649b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 125749b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 12583c896bc6SHong Zhang MatCopy_SeqSBAIJ, 12598c07d4e3SBarry Smith /*45*/ 0, 126049b5e25fSSatish Balay MatScale_SeqSBAIJ, 126149b5e25fSSatish Balay 0, 126249b5e25fSSatish Balay 0, 126349b5e25fSSatish Balay 0, 1264521d7252SBarry Smith /*50*/ 0, 126549b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 126649b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 126749b5e25fSSatish Balay 0, 126849b5e25fSSatish Balay 0, 126997304618SKris Buschelman /*55*/ 0, 127049b5e25fSSatish Balay 0, 127149b5e25fSSatish Balay 0, 127249b5e25fSSatish Balay 0, 127349b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 127497304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 127549b5e25fSSatish Balay 0, 127649b5e25fSSatish Balay 0, 1277357abbc8SBarry Smith 0, 1278d959ec07SHong Zhang 0, 127997304618SKris Buschelman /*65*/ 0, 1280d959ec07SHong Zhang 0, 1281d959ec07SHong Zhang 0, 1282d959ec07SHong Zhang 0, 1283d959ec07SHong Zhang 0, 1284985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqSBAIJ, 12853e0d88b5SBarry Smith 0, 12863e0d88b5SBarry Smith 0, 12873e0d88b5SBarry Smith 0, 12883e0d88b5SBarry Smith 0, 128997304618SKris Buschelman /*75*/ 0, 12903e0d88b5SBarry Smith 0, 12913e0d88b5SBarry Smith 0, 12923e0d88b5SBarry Smith 0, 12933e0d88b5SBarry Smith 0, 129497304618SKris Buschelman /*80*/ 0, 12953e0d88b5SBarry Smith 0, 12963e0d88b5SBarry Smith 0, 12973e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 129897304618SKris Buschelman MatGetInertia_SeqSBAIJ, 12993e0d88b5SBarry Smith #else 130097304618SKris Buschelman 0, 13013e0d88b5SBarry Smith #endif 1302865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1303865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1304865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1305efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1306865e5f61SKris Buschelman 0, 1307865e5f61SKris Buschelman 0, 1308865e5f61SKris Buschelman /*90*/ 0, 1309865e5f61SKris Buschelman 0, 1310865e5f61SKris Buschelman 0, 1311865e5f61SKris Buschelman 0, 1312865e5f61SKris Buschelman 0, 1313865e5f61SKris Buschelman /*95*/ 0, 1314865e5f61SKris Buschelman 0, 1315865e5f61SKris Buschelman 0, 131699cafbc1SBarry Smith 0, 131799cafbc1SBarry Smith 0, 131899cafbc1SBarry Smith /*100*/0, 131999cafbc1SBarry Smith 0, 132099cafbc1SBarry Smith 0, 132199cafbc1SBarry Smith 0, 132299cafbc1SBarry Smith 0, 132399cafbc1SBarry Smith /*105*/0, 132499cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1325f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1326f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 13272af78befSBarry Smith MatRestoreRowUpperTriangular_SeqSBAIJ, 13282af78befSBarry Smith /*110*/0, 13292af78befSBarry Smith 0, 13302af78befSBarry Smith 0, 13312af78befSBarry Smith 0, 13322af78befSBarry Smith MatMissingDiagonal_SeqSBAIJ 133399cafbc1SBarry Smith }; 1334be1d678aSKris Buschelman 133549b5e25fSSatish Balay EXTERN_C_BEGIN 13364a2ae208SSatish Balay #undef __FUNCT__ 13374a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1338be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 133949b5e25fSSatish Balay { 13404afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1341899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 1342dfbe8321SBarry Smith PetscErrorCode ierr; 134349b5e25fSSatish Balay 134449b5e25fSSatish Balay PetscFunctionBegin; 134549b5e25fSSatish Balay if (aij->nonew != 1) { 1346512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 134749b5e25fSSatish Balay } 134849b5e25fSSatish Balay 134949b5e25fSSatish Balay /* allocate space for values if not already there */ 135049b5e25fSSatish Balay if (!aij->saved_values) { 135187828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 135249b5e25fSSatish Balay } 135349b5e25fSSatish Balay 135449b5e25fSSatish Balay /* copy values over */ 135587828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 135649b5e25fSSatish Balay PetscFunctionReturn(0); 135749b5e25fSSatish Balay } 135849b5e25fSSatish Balay EXTERN_C_END 135949b5e25fSSatish Balay 136049b5e25fSSatish Balay EXTERN_C_BEGIN 13614a2ae208SSatish Balay #undef __FUNCT__ 13624a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1363be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 136449b5e25fSSatish Balay { 13654afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13666849ba73SBarry Smith PetscErrorCode ierr; 1367899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 136849b5e25fSSatish Balay 136949b5e25fSSatish Balay PetscFunctionBegin; 137049b5e25fSSatish Balay if (aij->nonew != 1) { 1371512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 137249b5e25fSSatish Balay } 137349b5e25fSSatish Balay if (!aij->saved_values) { 1374e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 137549b5e25fSSatish Balay } 137649b5e25fSSatish Balay 137749b5e25fSSatish Balay /* copy values over */ 137887828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 137949b5e25fSSatish Balay PetscFunctionReturn(0); 138049b5e25fSSatish Balay } 138149b5e25fSSatish Balay EXTERN_C_END 138249b5e25fSSatish Balay 13838549e402SHong Zhang EXTERN_C_BEGIN 13844a2ae208SSatish Balay #undef __FUNCT__ 1385a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1386be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 138749b5e25fSSatish Balay { 1388c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13896849ba73SBarry Smith PetscErrorCode ierr; 1390ab93d7beSBarry Smith PetscInt i,mbs,bs2; 1391ab93d7beSBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg; 139249b5e25fSSatish Balay 139349b5e25fSSatish Balay PetscFunctionBegin; 1394273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 13957adad957SLisandro Dalcin ierr = PetscOptionsGetInt(((PetscObject)B)->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1396899cda47SBarry Smith B->rmap.bs = B->cmap.bs = bs; 13976148ca0dSBarry Smith ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr); 13986148ca0dSBarry Smith ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr); 1399899cda47SBarry Smith 1400899cda47SBarry Smith mbs = B->rmap.N/bs; 140149b5e25fSSatish Balay bs2 = bs*bs; 140249b5e25fSSatish Balay 1403899cda47SBarry Smith if (mbs*bs != B->rmap.N) { 140429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 140549b5e25fSSatish Balay } 140649b5e25fSSatish Balay 1407ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1408ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1409ab93d7beSBarry Smith nz = 0; 1410ab93d7beSBarry Smith } 1411ab93d7beSBarry Smith 1412435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 141377431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 141449b5e25fSSatish Balay if (nnz) { 141549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 141677431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 141777431f27SBarry 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); 141849b5e25fSSatish Balay } 141949b5e25fSSatish Balay } 142049b5e25fSSatish Balay 14217adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)B)->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 142249b5e25fSSatish Balay if (!flg) { 142349b5e25fSSatish Balay switch (bs) { 142449b5e25fSSatish Balay case 1: 14254ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 142649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1427d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1428e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 142949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 143049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1431431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1432431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 143349b5e25fSSatish Balay break; 143449b5e25fSSatish Balay case 2: 14354ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 143649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1437e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 143849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 143949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1440431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1441431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 144249b5e25fSSatish Balay break; 144349b5e25fSSatish Balay case 3: 14445f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 144549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1446e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 144749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 144849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1449431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1450431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 145149b5e25fSSatish Balay break; 145249b5e25fSSatish Balay case 4: 14535f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 145449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1455e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 145649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 145749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1458431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1459431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 146049b5e25fSSatish Balay break; 146149b5e25fSSatish Balay case 5: 14625f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 146349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1464e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 146549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 146649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1467431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1468431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 146949b5e25fSSatish Balay break; 147049b5e25fSSatish Balay case 6: 14715f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 147249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1473e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 147449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 147549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1476431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1477431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 147849b5e25fSSatish Balay break; 147949b5e25fSSatish Balay case 7: 1480de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 148149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1482e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1483de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 148449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1485431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1486431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 148749b5e25fSSatish Balay break; 148849b5e25fSSatish Balay } 148949b5e25fSSatish Balay } 149049b5e25fSSatish Balay 149149b5e25fSSatish Balay b->mbs = mbs; 14924afc71dfSHong Zhang b->nbs = mbs; 1493ab93d7beSBarry Smith if (!skipallocation) { 14942ee49352SLisandro Dalcin if (!b->imax) { 1495ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 14962ee49352SLisandro Dalcin ierr = PetscLogObjectMemory(B,2*mbs*sizeof(PetscInt));CHKERRQ(ierr); 14972ee49352SLisandro Dalcin } 149849b5e25fSSatish Balay if (!nnz) { 1499435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 150049b5e25fSSatish Balay else if (nz <= 0) nz = 1; 150149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15028cef66ccSHong Zhang b->imax[i] = nz; 150349b5e25fSSatish Balay } 1504153ea458SHong Zhang nz = nz*mbs; /* total nz */ 150549b5e25fSSatish Balay } else { 150649b5e25fSSatish Balay nz = 0; 15078cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 150849b5e25fSSatish Balay } 15092ee49352SLisandro Dalcin /* b->ilen will count nonzeros in each block row so far. */ 15102ee49352SLisandro Dalcin for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 15116c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 151249b5e25fSSatish Balay 151349b5e25fSSatish Balay /* allocate the matrix space */ 15142ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);CHKERRQ(ierr); 1515899cda47SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.N+1,PetscInt,&b->i);CHKERRQ(ierr); 15162ee49352SLisandro Dalcin ierr = PetscLogObjectMemory(B,(B->rmap.N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 15176c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 151813f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 151949b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 152049b5e25fSSatish Balay 152149b5e25fSSatish Balay /* pointer to beginning of each row */ 152249b5e25fSSatish Balay b->i[0] = 0; 152349b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 152449b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 152549b5e25fSSatish Balay } 1526e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1527e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1528e811da20SHong Zhang } else { 1529e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1530e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1531ab93d7beSBarry Smith } 153249b5e25fSSatish Balay 1533899cda47SBarry Smith B->rmap.bs = bs; 153449b5e25fSSatish Balay b->bs2 = bs2; 15356c6c5352SBarry Smith b->nz = 0; 15366c6c5352SBarry Smith b->maxnz = nz*bs2; 1537153ea458SHong Zhang 153816cdd363SHong Zhang b->inew = 0; 153916cdd363SHong Zhang b->jnew = 0; 154016cdd363SHong Zhang b->anew = 0; 154116cdd363SHong Zhang b->a2anew = 0; 15421a3463dfSHong Zhang b->permute = PETSC_FALSE; 1543c464158bSHong Zhang PetscFunctionReturn(0); 1544c464158bSHong Zhang } 1545a23d5eceSKris Buschelman EXTERN_C_END 1546153ea458SHong Zhang 1547d769727bSBarry Smith EXTERN_C_BEGIN 1548f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1549f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1550d769727bSBarry Smith EXTERN_C_END 1551d769727bSBarry Smith 1552*5c9eb25fSBarry Smith #undef __FUNCT__ 1553*5c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_seqsbaij_petsc" 1554*5c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 1555*5c9eb25fSBarry Smith { 1556*5c9eb25fSBarry Smith PetscInt n = A->rmap.n; 1557*5c9eb25fSBarry Smith PetscErrorCode ierr; 1558*5c9eb25fSBarry Smith 1559*5c9eb25fSBarry Smith PetscFunctionBegin; 1560*5c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 1561*5c9eb25fSBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 1562*5c9eb25fSBarry Smith if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 1563*5c9eb25fSBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 1564*5c9eb25fSBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1565*5c9eb25fSBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 1566*5c9eb25fSBarry Smith PetscFunctionReturn(0); 1567*5c9eb25fSBarry Smith } 1568*5c9eb25fSBarry Smith 1569*5c9eb25fSBarry Smith EXTERN_C_BEGIN 1570*5c9eb25fSBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat,MatFactorType,Mat*); 1571*5c9eb25fSBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_spooles(Mat,MatFactorType,Mat*); 1572*5c9eb25fSBarry Smith EXTERN_C_END 1573*5c9eb25fSBarry Smith 15740bad9183SKris Buschelman /*MC 1575fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15760bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15770bad9183SKris Buschelman 15780bad9183SKris Buschelman Options Database Keys: 15790bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15800bad9183SKris Buschelman 15810bad9183SKris Buschelman Level: beginner 15820bad9183SKris Buschelman 15830bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15840bad9183SKris Buschelman M*/ 15850bad9183SKris Buschelman 1586a23d5eceSKris Buschelman EXTERN_C_BEGIN 1587a23d5eceSKris Buschelman #undef __FUNCT__ 1588a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1589be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1590a23d5eceSKris Buschelman { 1591a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1592dfbe8321SBarry Smith PetscErrorCode ierr; 159313f74950SBarry Smith PetscMPIInt size; 1594941593c8SHong Zhang PetscTruth flg; 1595a23d5eceSKris Buschelman 1596a23d5eceSKris Buschelman PetscFunctionBegin; 15977adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 1598a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1599a23d5eceSKris Buschelman 160038f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1601a23d5eceSKris Buschelman B->data = (void*)b; 1602a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1603a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1604a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1605a23d5eceSKris Buschelman B->mapping = 0; 1606a23d5eceSKris Buschelman b->row = 0; 1607a23d5eceSKris Buschelman b->icol = 0; 1608a23d5eceSKris Buschelman b->reallocs = 0; 1609a23d5eceSKris Buschelman b->saved_values = 0; 1610a23d5eceSKris Buschelman 1611a23d5eceSKris Buschelman 1612a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1613a23d5eceSKris Buschelman b->nonew = 0; 1614a23d5eceSKris Buschelman b->diag = 0; 1615a23d5eceSKris Buschelman b->solve_work = 0; 1616a23d5eceSKris Buschelman b->mult_work = 0; 1617a23d5eceSKris Buschelman B->spptr = 0; 1618a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1619a23d5eceSKris Buschelman b->xtoy = 0; 1620a23d5eceSKris Buschelman b->XtoY = 0; 1621a23d5eceSKris Buschelman 1622a23d5eceSKris Buschelman b->inew = 0; 1623a23d5eceSKris Buschelman b->jnew = 0; 1624a23d5eceSKris Buschelman b->anew = 0; 1625a23d5eceSKris Buschelman b->a2anew = 0; 1626a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1627a23d5eceSKris Buschelman 1628941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 1629941593c8SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr); 1630941593c8SHong Zhang if (flg) b->ignore_ltriangular = PETSC_TRUE; 1631941593c8SHong Zhang 1632f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 1633f5edf698SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr); 1634f5edf698SHong Zhang if (flg) b->getrow_utriangular = PETSC_TRUE; 1635f5edf698SHong Zhang 1636*5c9eb25fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqsbaij_spooles_C", 1637*5c9eb25fSBarry Smith "MatGetFactor_seqsbaij_spooles", 1638*5c9eb25fSBarry Smith MatGetFactor_seqsbaij_spooles);CHKERRQ(ierr); 1639*5c9eb25fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqsbaij_mumps_C", 1640*5c9eb25fSBarry Smith "MatGetFactor_seqsbaij_mumps", 1641*5c9eb25fSBarry Smith MatGetFactor_seqsbaij_mumps);CHKERRQ(ierr); 1642*5c9eb25fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqsbaij_petsc_C", 1643*5c9eb25fSBarry Smith "MatGetFactor_seqsbaij_petsc", 1644*5c9eb25fSBarry Smith MatGetFactor_seqsbaij_petsc);CHKERRQ(ierr); 1645a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1646a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1647a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1648a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1649a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1650a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1651a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1652a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1653a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16544e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16554e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16564e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1657a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1658a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1659a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1660a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1661a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1662a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 166323ce1328SBarry Smith 166423ce1328SBarry Smith B->symmetric = PETSC_TRUE; 166523ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 166623ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 166723ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 166817667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 1669a23d5eceSKris Buschelman PetscFunctionReturn(0); 1670a23d5eceSKris Buschelman } 1671a23d5eceSKris Buschelman EXTERN_C_END 1672a23d5eceSKris Buschelman 1673a23d5eceSKris Buschelman #undef __FUNCT__ 1674a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1675a23d5eceSKris Buschelman /*@C 1676a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1677a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1678a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1679a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1680a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1681a23d5eceSKris Buschelman 1682a23d5eceSKris Buschelman Collective on Mat 1683a23d5eceSKris Buschelman 1684a23d5eceSKris Buschelman Input Parameters: 1685a23d5eceSKris Buschelman + A - the symmetric matrix 1686a23d5eceSKris Buschelman . bs - size of block 1687a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1688a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1689a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1690a23d5eceSKris Buschelman 1691a23d5eceSKris Buschelman Options Database Keys: 1692a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1693a23d5eceSKris Buschelman block calculations (much slower) 1694a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1695a23d5eceSKris Buschelman 1696a23d5eceSKris Buschelman Level: intermediate 1697a23d5eceSKris Buschelman 1698a23d5eceSKris Buschelman Notes: 1699a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1700a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1701a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1702a23d5eceSKris Buschelman matrices. 1703a23d5eceSKris Buschelman 1704aa95bbe8SBarry Smith You can call MatGetInfo() to get information on how effective the preallocation was; 1705aa95bbe8SBarry Smith for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 1706aa95bbe8SBarry Smith You can also run with the option -info and look for messages with the string 1707aa95bbe8SBarry Smith malloc in them to see if additional memory allocation was needed. 1708aa95bbe8SBarry Smith 170949a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 171049a6f317SBarry Smith 171149a6f317SBarry Smith 1712a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1713a23d5eceSKris Buschelman @*/ 1714be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 171513f74950SBarry Smith { 171613f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1717a23d5eceSKris Buschelman 1718a23d5eceSKris Buschelman PetscFunctionBegin; 1719a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1720a23d5eceSKris Buschelman if (f) { 1721a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1722a23d5eceSKris Buschelman } 1723a23d5eceSKris Buschelman PetscFunctionReturn(0); 1724a23d5eceSKris Buschelman } 172549b5e25fSSatish Balay 17264a2ae208SSatish Balay #undef __FUNCT__ 17274a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1728c464158bSHong Zhang /*@C 1729c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1730c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1731c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1732c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1733c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 173449b5e25fSSatish Balay 1735c464158bSHong Zhang Collective on MPI_Comm 1736c464158bSHong Zhang 1737c464158bSHong Zhang Input Parameters: 1738c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1739c464158bSHong Zhang . bs - size of block 1740c464158bSHong Zhang . m - number of rows, or number of columns 1741c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1742744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1743744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1744c464158bSHong Zhang 1745c464158bSHong Zhang Output Parameter: 1746c464158bSHong Zhang . A - the symmetric matrix 1747c464158bSHong Zhang 1748c464158bSHong Zhang Options Database Keys: 1749c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1750c464158bSHong Zhang block calculations (much slower) 1751c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1752c464158bSHong Zhang 1753c464158bSHong Zhang Level: intermediate 1754c464158bSHong Zhang 1755175b88e8SBarry Smith It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 1756175b88e8SBarry Smith MatXXXXSetPreallocation() paradgm instead of this routine directly. This is definitely 1757175b88e8SBarry Smith true if you plan to use the external direct solvers such as SuperLU, MUMPS or Spooles. 1758175b88e8SBarry Smith [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 1759175b88e8SBarry Smith 1760c464158bSHong Zhang Notes: 17616d6d819aSHong Zhang The number of rows and columns must be divisible by blocksize. 17626d6d819aSHong Zhang This matrix type does not support complex Hermitian operation. 1763c464158bSHong Zhang 1764c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1765c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1766c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1767c464158bSHong Zhang matrices. 1768c464158bSHong Zhang 176949a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 177049a6f317SBarry Smith 1771c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1772c464158bSHong Zhang @*/ 1773be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1774c464158bSHong Zhang { 1775dfbe8321SBarry Smith PetscErrorCode ierr; 1776c464158bSHong Zhang 1777c464158bSHong Zhang PetscFunctionBegin; 1778f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1779f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1780c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1781ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 178249b5e25fSSatish Balay PetscFunctionReturn(0); 178349b5e25fSSatish Balay } 178449b5e25fSSatish Balay 17854a2ae208SSatish Balay #undef __FUNCT__ 17864a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1787dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 178849b5e25fSSatish Balay { 178949b5e25fSSatish Balay Mat C; 179049b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17916849ba73SBarry Smith PetscErrorCode ierr; 1792b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 179349b5e25fSSatish Balay 179449b5e25fSSatish Balay PetscFunctionBegin; 179529bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 179649b5e25fSSatish Balay 179749b5e25fSSatish Balay *B = 0; 17987adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 1799899cda47SBarry Smith ierr = MatSetSizes(C,A->rmap.N,A->cmap.n,A->rmap.N,A->cmap.n);CHKERRQ(ierr); 18007adad957SLisandro Dalcin ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 18011d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1802692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1803692f9cbeSHong Zhang 1804273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 180549b5e25fSSatish Balay C->factor = A->factor; 180649b5e25fSSatish Balay c->row = 0; 180749b5e25fSSatish Balay c->icol = 0; 180849b5e25fSSatish Balay c->saved_values = 0; 180949b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 181049b5e25fSSatish Balay C->assembled = PETSC_TRUE; 181149b5e25fSSatish Balay 18127adad957SLisandro Dalcin ierr = PetscMapCopy(((PetscObject)A)->comm,&A->rmap,&C->rmap);CHKERRQ(ierr); 18137adad957SLisandro Dalcin ierr = PetscMapCopy(((PetscObject)A)->comm,&A->cmap,&C->cmap);CHKERRQ(ierr); 181449b5e25fSSatish Balay c->bs2 = a->bs2; 181549b5e25fSSatish Balay c->mbs = a->mbs; 181649b5e25fSSatish Balay c->nbs = a->nbs; 181749b5e25fSSatish Balay 18188777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 181949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 182049b5e25fSSatish Balay c->imax[i] = a->imax[i]; 182149b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 182249b5e25fSSatish Balay } 182349b5e25fSSatish Balay 182449b5e25fSSatish Balay /* allocate the matrix space */ 1825b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 182649b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 182713f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1828b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 182949b5e25fSSatish Balay if (mbs > 0) { 183013f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 183149b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 183249b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 183349b5e25fSSatish Balay } else { 183449b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 183549b5e25fSSatish Balay } 183649b5e25fSSatish Balay } 183749b5e25fSSatish Balay 183849b5e25fSSatish Balay c->roworiented = a->roworiented; 183949b5e25fSSatish Balay c->nonew = a->nonew; 184049b5e25fSSatish Balay 184149b5e25fSSatish Balay if (a->diag) { 184213f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 184352e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 184449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 184549b5e25fSSatish Balay c->diag[i] = a->diag[i]; 184649b5e25fSSatish Balay } 184749b5e25fSSatish Balay } else c->diag = 0; 18486c6c5352SBarry Smith c->nz = a->nz; 18496c6c5352SBarry Smith c->maxnz = a->maxnz; 185049b5e25fSSatish Balay c->solve_work = 0; 185149b5e25fSSatish Balay c->mult_work = 0; 1852e6b907acSBarry Smith c->free_a = PETSC_TRUE; 1853e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 185449b5e25fSSatish Balay *B = C; 18557adad957SLisandro Dalcin ierr = PetscFListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);CHKERRQ(ierr); 185649b5e25fSSatish Balay PetscFunctionReturn(0); 185749b5e25fSSatish Balay } 185849b5e25fSSatish Balay 18594a2ae208SSatish Balay #undef __FUNCT__ 18604a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1861f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A) 186249b5e25fSSatish Balay { 186349b5e25fSSatish Balay Mat_SeqSBAIJ *a; 186449b5e25fSSatish Balay Mat B; 18656849ba73SBarry Smith PetscErrorCode ierr; 186613f74950SBarry Smith int fd; 186713f74950SBarry Smith PetscMPIInt size; 186813f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 186913f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 187013f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 187113f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 187287828ca2SBarry Smith PetscScalar *aa; 187349b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 187449b5e25fSSatish Balay 187549b5e25fSSatish Balay PetscFunctionBegin; 1876b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 187749b5e25fSSatish Balay bs2 = bs*bs; 187849b5e25fSSatish Balay 187949b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 188029bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1881b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 188249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1883552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 188449b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 188549b5e25fSSatish Balay 188649b5e25fSSatish Balay if (header[3] < 0) { 188729bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 188849b5e25fSSatish Balay } 188949b5e25fSSatish Balay 189029bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 189149b5e25fSSatish Balay 189249b5e25fSSatish Balay /* 189349b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 189449b5e25fSSatish Balay divisible by the blocksize 189549b5e25fSSatish Balay */ 189649b5e25fSSatish Balay mbs = M/bs; 189749b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 189849b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 189949b5e25fSSatish Balay else mbs++; 190049b5e25fSSatish Balay if (extra_rows) { 19011e2582c4SBarry Smith ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 190249b5e25fSSatish Balay } 190349b5e25fSSatish Balay 190449b5e25fSSatish Balay /* read in row lengths */ 190513f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 190649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 190749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 190849b5e25fSSatish Balay 190949b5e25fSSatish Balay /* read in column indices */ 191013f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 191149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 191249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 191349b5e25fSSatish Balay 191449b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 191513f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 191613f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 191713f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 191813f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 191949b5e25fSSatish Balay masked = mask + mbs; 192049b5e25fSSatish Balay rowcount = 0; nzcount = 0; 192149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 192249b5e25fSSatish Balay nmask = 0; 192349b5e25fSSatish Balay for (j=0; j<bs; j++) { 192449b5e25fSSatish Balay kmax = rowlengths[rowcount]; 192549b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19262d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 192703630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 192849b5e25fSSatish Balay } 192949b5e25fSSatish Balay rowcount++; 193049b5e25fSSatish Balay } 1931574b2666SHong Zhang s_browlengths[i] += nmask; 1932574b2666SHong Zhang 193349b5e25fSSatish Balay /* zero out the mask elements we set */ 193449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 193549b5e25fSSatish Balay } 193649b5e25fSSatish Balay 193749b5e25fSSatish Balay /* create our matrix */ 1938f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 1939f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 19409abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 1941ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 194249b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 194349b5e25fSSatish Balay 194449b5e25fSSatish Balay /* set matrix "i" values */ 194549b5e25fSSatish Balay a->i[0] = 0; 194649b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1947574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1948574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 194949b5e25fSSatish Balay } 19506c6c5352SBarry Smith a->nz = a->i[mbs]; 195149b5e25fSSatish Balay 195249b5e25fSSatish Balay /* read in nonzero values */ 195387828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 195449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 195549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 195649b5e25fSSatish Balay 195749b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 195849b5e25fSSatish Balay nzcount = 0; jcount = 0; 195949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 196049b5e25fSSatish Balay nzcountb = nzcount; 196149b5e25fSSatish Balay nmask = 0; 196249b5e25fSSatish Balay for (j=0; j<bs; j++) { 196349b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 196449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19652d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 196603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19672d703238SHong Zhang } 19682d703238SHong Zhang } 19692d703238SHong Zhang /* sort the masked values */ 19702d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19712d703238SHong Zhang 19722d703238SHong Zhang /* set "j" values into matrix */ 19732d703238SHong Zhang maskcount = 1; 19742d703238SHong Zhang for (j=0; j<nmask; j++) { 197549b5e25fSSatish Balay a->j[jcount++] = masked[j]; 197649b5e25fSSatish Balay mask[masked[j]] = maskcount++; 197749b5e25fSSatish Balay } 1978574b2666SHong Zhang 197949b5e25fSSatish Balay /* set "a" values into matrix */ 198049b5e25fSSatish Balay ishift = bs2*a->i[i]; 198149b5e25fSSatish Balay for (j=0; j<bs; j++) { 198249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 198349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1984574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1985574b2666SHong Zhang if (tmp >= i){ 198649b5e25fSSatish Balay block = mask[tmp] - 1; 198749b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 198849b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1989574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1990574b2666SHong Zhang } 1991574b2666SHong Zhang nzcountb++; 199249b5e25fSSatish Balay } 199349b5e25fSSatish Balay } 199449b5e25fSSatish Balay /* zero out the mask elements we set */ 199549b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 199649b5e25fSSatish Balay } 19976c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 199849b5e25fSSatish Balay 199949b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2000574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 200149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 200249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 200349b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 200449b5e25fSSatish Balay 20059abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20069abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 200749b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 20089abb65ffSKris Buschelman *A = B; 200949b5e25fSSatish Balay PetscFunctionReturn(0); 201049b5e25fSSatish Balay } 2011574b2666SHong Zhang 2012d06b337dSHong Zhang #undef __FUNCT__ 2013d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 201413f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 2015d06b337dSHong Zhang { 2016d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2017d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 2018d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 20196849ba73SBarry Smith PetscErrorCode ierr; 2020899cda47SBarry Smith PetscInt m=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j; 2021521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 2022d06b337dSHong Zhang 2023d06b337dSHong Zhang PetscFunctionBegin; 202451f519a2SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat"); 202551f519a2SBarry Smith 2026b965ef7fSBarry Smith its = its*lits; 202777431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 2028d06b337dSHong Zhang 2029d06b337dSHong Zhang if (bs > 1) 2030d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 2031d06b337dSHong Zhang 20321ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2033d06b337dSHong Zhang if (xx != bb) { 20341ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 2035d06b337dSHong Zhang } else { 2036d06b337dSHong Zhang b = x; 2037d06b337dSHong Zhang } 2038d06b337dSHong Zhang 2039e2ee2a47SBarry Smith if (!a->relax_work) { 2040e2ee2a47SBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->relax_work);CHKERRQ(ierr); 2041e2ee2a47SBarry Smith } 2042e2ee2a47SBarry Smith t = a->relax_work; 2043d06b337dSHong Zhang 2044d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 20452798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2046290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2047d06b337dSHong Zhang 2048d06b337dSHong Zhang for (i=0; i<m; i++){ 204944706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2050d06b337dSHong Zhang v = aa + ai[i] + 1; 2051d06b337dSHong Zhang vj = aj + ai[i] + 1; 2052d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2053e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2054d06b337dSHong Zhang x[i] = omega*t[i]/d; 2055d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 2056d06b337dSHong Zhang } 2057d06b337dSHong Zhang } 2058d06b337dSHong Zhang 20592798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 206095d750ceSBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){ 206195d750ceSBarry Smith t = b; 2062d06b337dSHong Zhang } 206395d750ceSBarry Smith 2064d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2065d06b337dSHong Zhang d = *(aa + ai[i]); 2066d06b337dSHong Zhang v = aa + ai[i] + 1; 2067d06b337dSHong Zhang vj = aj + ai[i] + 1; 2068d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2069e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2070d06b337dSHong Zhang sum = t[i]; 2071d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2072d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2073d06b337dSHong Zhang } 207495d750ceSBarry Smith t = a->relax_work; 2075d06b337dSHong Zhang } 2076d06b337dSHong Zhang its--; 2077d06b337dSHong Zhang } 2078d06b337dSHong Zhang 2079d06b337dSHong Zhang while (its--) { 208044706875SHong Zhang /* 208144706875SHong Zhang forward sweep: 208244706875SHong Zhang for i=0,...,m-1: 208344706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 208444706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 208544706875SHong Zhang b = b - x[i]*U^T(i,:); 2086d06b337dSHong Zhang 208744706875SHong Zhang */ 20882798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2089290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2090d06b337dSHong Zhang 2091d06b337dSHong Zhang for (i=0; i<m; i++){ 209244706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2093d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2094d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2095d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2096d06b337dSHong Zhang sum = t[i]; 2097e6b907acSBarry Smith ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2098d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2099d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2100d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2101d06b337dSHong Zhang } 2102d06b337dSHong Zhang } 2103d06b337dSHong Zhang 21042798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 210544706875SHong Zhang /* 210644706875SHong Zhang backward sweep: 210744706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 210844706875SHong Zhang for i=m-1,...,0: 210944706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 211044706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 211144706875SHong Zhang */ 211295d750ceSBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2113290bbb0aSBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 2114d06b337dSHong Zhang 2115d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2116d06b337dSHong Zhang v = aa + ai[i] + 1; 2117d06b337dSHong Zhang vj = aj + ai[i] + 1; 2118d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2119efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2120e6b907acSBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 2121d06b337dSHong Zhang } 2122d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2123d06b337dSHong Zhang d = *(aa + ai[i]); 2124d06b337dSHong Zhang v = aa + ai[i] + 1; 2125d06b337dSHong Zhang vj = aj + ai[i] + 1; 2126d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2127e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2128d06b337dSHong Zhang sum = t[i]; 2129d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2130d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2131d06b337dSHong Zhang } 2132d06b337dSHong Zhang } 2133d06b337dSHong Zhang } 2134d06b337dSHong Zhang 21351ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2136d06b337dSHong Zhang if (bb != xx) { 21371ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2138d06b337dSHong Zhang } 2139d06b337dSHong Zhang PetscFunctionReturn(0); 2140d06b337dSHong Zhang } 2141d06b337dSHong Zhang 2142c75a6043SHong Zhang #undef __FUNCT__ 2143c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2144c75a6043SHong Zhang /*@ 2145c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2146c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2147c75a6043SHong Zhang 2148c75a6043SHong Zhang Collective on MPI_Comm 2149c75a6043SHong Zhang 2150c75a6043SHong Zhang Input Parameters: 2151c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2152c75a6043SHong Zhang . bs - size of block 2153c75a6043SHong Zhang . m - number of rows 2154c75a6043SHong Zhang . n - number of columns 2155c75a6043SHong Zhang . i - row indices 2156c75a6043SHong Zhang . j - column indices 2157c75a6043SHong Zhang - a - matrix values 2158c75a6043SHong Zhang 2159c75a6043SHong Zhang Output Parameter: 2160c75a6043SHong Zhang . mat - the matrix 2161c75a6043SHong Zhang 2162c75a6043SHong Zhang Level: intermediate 2163c75a6043SHong Zhang 2164c75a6043SHong Zhang Notes: 2165c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2166c75a6043SHong Zhang once the matrix is destroyed 2167c75a6043SHong Zhang 2168c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2169c75a6043SHong Zhang 2170c75a6043SHong Zhang The i and j indices are 0 based 2171c75a6043SHong Zhang 2172c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2173c75a6043SHong Zhang 2174c75a6043SHong Zhang @*/ 2175c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2176c75a6043SHong Zhang { 2177c75a6043SHong Zhang PetscErrorCode ierr; 2178c75a6043SHong Zhang PetscInt ii; 2179c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2180c75a6043SHong Zhang 2181c75a6043SHong Zhang PetscFunctionBegin; 2182c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2183c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2184c75a6043SHong Zhang 2185c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2186c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2187c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2188c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2189c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2190c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2191c75a6043SHong Zhang 2192c75a6043SHong Zhang sbaij->i = i; 2193c75a6043SHong Zhang sbaij->j = j; 2194c75a6043SHong Zhang sbaij->a = a; 2195c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2196c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2197e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2198e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2199c75a6043SHong Zhang 2200c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2201c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2202c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2203c75a6043SHong 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]); 2204c75a6043SHong Zhang #endif 2205c75a6043SHong Zhang } 2206c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2207c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2208c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2209c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2210c75a6043SHong Zhang } 2211c75a6043SHong Zhang #endif 2212c75a6043SHong Zhang 2213c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2214c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2215c75a6043SHong Zhang PetscFunctionReturn(0); 2216c75a6043SHong Zhang } 2217d06b337dSHong Zhang 2218d06b337dSHong Zhang 2219d06b337dSHong Zhang 222049b5e25fSSatish Balay 222149b5e25fSSatish Balay 2222