1be1d678aSKris Buschelman #define PETSCMAT_DLL 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 4a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/inline/spops.h" 9aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1049b5e25fSSatish Balay 1149b5e25fSSatish Balay #define CHUNKSIZE 10 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay /* 1449b5e25fSSatish Balay Checks for missing diagonals 1549b5e25fSSatish Balay */ 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 18dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A) 1949b5e25fSSatish Balay { 20045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 216849ba73SBarry Smith PetscErrorCode ierr; 2213f74950SBarry Smith PetscInt *diag,*jj = a->j,i; 2349b5e25fSSatish Balay 2449b5e25fSSatish Balay PetscFunctionBegin; 25045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2649b5e25fSSatish Balay diag = a->diag; 2749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2877431f27SBarry Smith if (jj[diag[i]] != i) SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Matrix is missing diagonal number %D",i); 2949b5e25fSSatish Balay } 3049b5e25fSSatish Balay PetscFunctionReturn(0); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay 334a2ae208SSatish Balay #undef __FUNCT__ 344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 35dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 3649b5e25fSSatish Balay { 37045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 386849ba73SBarry Smith PetscErrorCode ierr; 3909f38230SBarry Smith PetscInt i; 4049b5e25fSSatish Balay 4149b5e25fSSatish Balay PetscFunctionBegin; 4209f38230SBarry Smith if (!a->diag) { 4309f38230SBarry Smith ierr = PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 4409f38230SBarry Smith } 4509f38230SBarry Smith for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i]; 4649b5e25fSSatish Balay PetscFunctionReturn(0); 4749b5e25fSSatish Balay } 4849b5e25fSSatish Balay 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 5113f74950SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 5249b5e25fSSatish Balay { 53a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5413f74950SBarry Smith PetscInt n = a->mbs,i; 5549b5e25fSSatish Balay 5649b5e25fSSatish Balay PetscFunctionBegin; 57d3e5a4abSHong Zhang *nn = n; 58a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 59a1373b80SHong Zhang 60a6ece127SHong Zhang if (oshift == 1) { 61a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 6213f74950SBarry Smith PetscInt nz = a->i[n]; 636c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 64a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 65a1373b80SHong Zhang *ia = a->i; *ja = a->j; 66a1373b80SHong Zhang } else { 67a1373b80SHong Zhang *ia = a->i; *ja = a->j; 68a6ece127SHong Zhang } 6949b5e25fSSatish Balay PetscFunctionReturn(0); 7049b5e25fSSatish Balay } 7149b5e25fSSatish Balay 724a2ae208SSatish Balay #undef __FUNCT__ 734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 7413f74950SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 7549b5e25fSSatish Balay { 76b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7713f74950SBarry Smith PetscInt i,n = a->mbs; 78a6ece127SHong Zhang 7949b5e25fSSatish Balay PetscFunctionBegin; 8049b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 81a6ece127SHong Zhang 82a6ece127SHong Zhang if (oshift == 1) { 8313f74950SBarry Smith PetscInt nz = a->i[n]-1; 846c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 85a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 86a6ece127SHong Zhang } 87a6ece127SHong Zhang PetscFunctionReturn(0); 8849b5e25fSSatish Balay } 8949b5e25fSSatish Balay 904a2ae208SSatish Balay #undef __FUNCT__ 914a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 92dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 9349b5e25fSSatish Balay { 9449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 95dfbe8321SBarry Smith PetscErrorCode ierr; 9649b5e25fSSatish Balay 9749b5e25fSSatish Balay PetscFunctionBegin; 98a9f03627SSatish Balay #if defined(PETSC_USE_LOG) 99899cda47SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap.N,a->nz); 100a9f03627SSatish Balay #endif 101e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr); 10249b5e25fSSatish Balay if (a->row) { 10349b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 10449b5e25fSSatish Balay } 10505b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 10605b42c5fSBarry Smith ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr); 10749b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 10805b42c5fSBarry Smith ierr = PetscFree(a->solve_work);CHKERRQ(ierr); 109e2ee2a47SBarry Smith ierr = PetscFree(a->relax_work);CHKERRQ(ierr); 11005b42c5fSBarry Smith ierr = PetscFree(a->solves_work);CHKERRQ(ierr); 11105b42c5fSBarry Smith ierr = PetscFree(a->mult_work);CHKERRQ(ierr); 11205b42c5fSBarry Smith ierr = PetscFree(a->saved_values);CHKERRQ(ierr); 11305b42c5fSBarry Smith ierr = PetscFree(a->xtoy);CHKERRQ(ierr); 1141a3463dfSHong Zhang 1151a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 11649b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 117901853e0SKris Buschelman 118901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 119901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 120901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 121901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 122901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 123901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 12449b5e25fSSatish Balay PetscFunctionReturn(0); 12549b5e25fSSatish Balay } 12649b5e25fSSatish Balay 1274a2ae208SSatish Balay #undef __FUNCT__ 1284a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 129dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op) 13049b5e25fSSatish Balay { 131045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 13263ba0a88SBarry Smith PetscErrorCode ierr; 13349b5e25fSSatish Balay 13449b5e25fSSatish Balay PetscFunctionBegin; 1354d9d31abSKris Buschelman switch (op) { 1364d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1374d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1384d9d31abSKris Buschelman break; 1394d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1404d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1414d9d31abSKris Buschelman break; 1424d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1434d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1444d9d31abSKris Buschelman break; 1454d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1464d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1474d9d31abSKris Buschelman break; 1484d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1494d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1504d9d31abSKris Buschelman break; 1514d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1524d9d31abSKris Buschelman a->nonew = 1; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1554d9d31abSKris Buschelman a->nonew = -1; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1584d9d31abSKris Buschelman a->nonew = -2; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1614d9d31abSKris Buschelman a->nonew = 0; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1644d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1654d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1664d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1674d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 168ad86a440SBarry Smith ierr = PetscInfo1(A,"Option %d ignored\n",op);CHKERRQ(ierr); 1694d9d31abSKris Buschelman break; 1704d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 17129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1729a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1739a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1749a4540c5SBarry Smith case MAT_HERMITIAN: 1759a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 17677e54ba9SKris Buschelman case MAT_SYMMETRIC: 17777e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1789a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1799a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1809a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 181941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 182941593c8SHong Zhang a->ignore_ltriangular = PETSC_TRUE; 183941593c8SHong Zhang break; 184941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 185941593c8SHong Zhang a->ignore_ltriangular = PETSC_FALSE; 18677e54ba9SKris Buschelman break; 187f5edf698SHong Zhang case MAT_GETROW_UPPERTRIANGULAR: 188f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 189f5edf698SHong Zhang break; 1904d9d31abSKris Buschelman default: 191ad86a440SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 19249b5e25fSSatish Balay } 19349b5e25fSSatish Balay PetscFunctionReturn(0); 19449b5e25fSSatish Balay } 19549b5e25fSSatish Balay 1964a2ae208SSatish Balay #undef __FUNCT__ 1974a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 19813f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 19949b5e25fSSatish Balay { 20049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2016849ba73SBarry Smith PetscErrorCode ierr; 20213f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 20349b5e25fSSatish Balay MatScalar *aa,*aa_i; 20487828ca2SBarry Smith PetscScalar *v_i; 20549b5e25fSSatish Balay 20649b5e25fSSatish Balay PetscFunctionBegin; 207f5edf698SHong Zhang if (A && !a->getrow_utriangular) SETERRQ(PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR) or MatGetRowUpperTriangular()"); 208f5edf698SHong Zhang /* Get the upper triangular part of the row */ 209899cda47SBarry Smith bs = A->rmap.bs; 21049b5e25fSSatish Balay ai = a->i; 21149b5e25fSSatish Balay aj = a->j; 21249b5e25fSSatish Balay aa = a->a; 21349b5e25fSSatish Balay bs2 = a->bs2; 21449b5e25fSSatish Balay 215899cda47SBarry Smith if (row < 0 || row >= A->rmap.N) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 21649b5e25fSSatish Balay 21749b5e25fSSatish Balay bn = row/bs; /* Block number */ 21849b5e25fSSatish Balay bp = row % bs; /* Block position */ 21949b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 22049b5e25fSSatish Balay *ncols = bs*M; 22149b5e25fSSatish Balay 22249b5e25fSSatish Balay if (v) { 22349b5e25fSSatish Balay *v = 0; 22449b5e25fSSatish Balay if (*ncols) { 22587828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 22649b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22749b5e25fSSatish Balay v_i = *v + i*bs; 22849b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 22949b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 23049b5e25fSSatish Balay } 23149b5e25fSSatish Balay } 23249b5e25fSSatish Balay } 23349b5e25fSSatish Balay 23449b5e25fSSatish Balay if (cols) { 23549b5e25fSSatish Balay *cols = 0; 23649b5e25fSSatish Balay if (*ncols) { 23713f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 23849b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23949b5e25fSSatish Balay cols_i = *cols + i*bs; 24049b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 24149b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24249b5e25fSSatish Balay } 24349b5e25fSSatish Balay } 24449b5e25fSSatish Balay } 24549b5e25fSSatish Balay 24649b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2475ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2485ddb2528SHong Zhang #ifdef column_search 24949b5e25fSSatish Balay v_i = *v + M*bs; 25049b5e25fSSatish Balay cols_i = *cols + M*bs; 25149b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25249b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25349b5e25fSSatish Balay for (j=0; j<M; j++){ 25449b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 25549b5e25fSSatish Balay if (itmp == bn){ 25649b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 25749b5e25fSSatish Balay for (k=0; k<bs; k++) { 25849b5e25fSSatish Balay *cols_i++ = i*bs+k; 25949b5e25fSSatish Balay *v_i++ = aa_i[k]; 26049b5e25fSSatish Balay } 26149b5e25fSSatish Balay *ncols += bs; 26249b5e25fSSatish Balay break; 26349b5e25fSSatish Balay } 26449b5e25fSSatish Balay } 26549b5e25fSSatish Balay } 2665ddb2528SHong Zhang #endif 26749b5e25fSSatish Balay PetscFunctionReturn(0); 26849b5e25fSSatish Balay } 26949b5e25fSSatish Balay 2704a2ae208SSatish Balay #undef __FUNCT__ 2714a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 27213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27349b5e25fSSatish Balay { 274dfbe8321SBarry Smith PetscErrorCode ierr; 27549b5e25fSSatish Balay 27649b5e25fSSatish Balay PetscFunctionBegin; 27705b42c5fSBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 27805b42c5fSBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 27949b5e25fSSatish Balay PetscFunctionReturn(0); 28049b5e25fSSatish Balay } 28149b5e25fSSatish Balay 2824a2ae208SSatish Balay #undef __FUNCT__ 283f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ" 284f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 285f5edf698SHong Zhang { 286f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 287f5edf698SHong Zhang 288f5edf698SHong Zhang PetscFunctionBegin; 289f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 290f5edf698SHong Zhang PetscFunctionReturn(0); 291f5edf698SHong Zhang } 292f5edf698SHong Zhang #undef __FUNCT__ 293f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ" 294f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 295f5edf698SHong Zhang { 296f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 297f5edf698SHong Zhang 298f5edf698SHong Zhang PetscFunctionBegin; 299f5edf698SHong Zhang a->getrow_utriangular = PETSC_FALSE; 300f5edf698SHong Zhang PetscFunctionReturn(0); 301f5edf698SHong Zhang } 302f5edf698SHong Zhang 303f5edf698SHong Zhang #undef __FUNCT__ 3044a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 305dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 30649b5e25fSSatish Balay { 307dfbe8321SBarry Smith PetscErrorCode ierr; 30849b5e25fSSatish Balay PetscFunctionBegin; 309999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 3108115998fSBarry Smith PetscFunctionReturn(0); 31149b5e25fSSatish Balay } 31249b5e25fSSatish Balay 3134a2ae208SSatish Balay #undef __FUNCT__ 3144a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3156849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 31649b5e25fSSatish Balay { 31749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 318dfbe8321SBarry Smith PetscErrorCode ierr; 319899cda47SBarry Smith PetscInt i,j,bs = A->rmap.bs,k,l,bs2=a->bs2; 3202dcb1b2aSMatthew Knepley const char *name; 321f3ef73ceSBarry Smith PetscViewerFormat format; 32249b5e25fSSatish Balay 32349b5e25fSSatish Balay PetscFunctionBegin; 32480fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 325b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 326456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 32777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 328fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 32929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 330fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 331b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 33249b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33349b5e25fSSatish Balay for (j=0; j<bs; j++) { 33477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 33549b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33649b5e25fSSatish Balay for (l=0; l<bs; l++) { 33749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 33849b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 339a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 34049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34149b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 342a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 34349b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34449b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 345a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34649b5e25fSSatish Balay } 34749b5e25fSSatish Balay #else 34849b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 349a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35049b5e25fSSatish Balay } 35149b5e25fSSatish Balay #endif 35249b5e25fSSatish Balay } 35349b5e25fSSatish Balay } 354b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 35549b5e25fSSatish Balay } 35649b5e25fSSatish Balay } 357b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 358c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 359c1490034SHong Zhang PetscFunctionReturn(0); 36049b5e25fSSatish Balay } else { 361b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36249b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 36349b5e25fSSatish Balay for (j=0; j<bs; j++) { 36477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 36549b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 36649b5e25fSSatish Balay for (l=0; l<bs; l++) { 36749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 36849b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 369a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 37049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37149b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 372a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 37349b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37449b5e25fSSatish Balay } else { 375a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37649b5e25fSSatish Balay } 37749b5e25fSSatish Balay #else 378*e9f7bc9eSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 37949b5e25fSSatish Balay #endif 38049b5e25fSSatish Balay } 38149b5e25fSSatish Balay } 382b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 38349b5e25fSSatish Balay } 38449b5e25fSSatish Balay } 385b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 38649b5e25fSSatish Balay } 387b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 38849b5e25fSSatish Balay PetscFunctionReturn(0); 38949b5e25fSSatish Balay } 39049b5e25fSSatish Balay 3914a2ae208SSatish Balay #undef __FUNCT__ 3924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 3936849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 39449b5e25fSSatish Balay { 39549b5e25fSSatish Balay Mat A = (Mat) Aa; 39649b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 3976849ba73SBarry Smith PetscErrorCode ierr; 398899cda47SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap.bs,bs2=a->bs2; 39913f74950SBarry Smith PetscMPIInt rank; 40049b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 40149b5e25fSSatish Balay MatScalar *aa; 40249b5e25fSSatish Balay MPI_Comm comm; 403b0a32e0cSBarry Smith PetscViewer viewer; 40449b5e25fSSatish Balay 40549b5e25fSSatish Balay PetscFunctionBegin; 40649b5e25fSSatish Balay /* 40749b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 40849b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 40949b5e25fSSatish Balay rest should return immediately. 41049b5e25fSSatish Balay */ 41149b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 41249b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 41349b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 41449b5e25fSSatish Balay 41549b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 41649b5e25fSSatish Balay 417b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 418b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 41949b5e25fSSatish Balay 42049b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 421b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 42249b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 42349b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 424899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 42549b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42649b5e25fSSatish Balay aa = a->a + j*bs2; 42749b5e25fSSatish Balay for (k=0; k<bs; k++) { 42849b5e25fSSatish Balay for (l=0; l<bs; l++) { 42949b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 430b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 43149b5e25fSSatish Balay } 43249b5e25fSSatish Balay } 43349b5e25fSSatish Balay } 43449b5e25fSSatish Balay } 435b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 43649b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 43749b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 438899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 43949b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 44049b5e25fSSatish Balay aa = a->a + j*bs2; 44149b5e25fSSatish Balay for (k=0; k<bs; k++) { 44249b5e25fSSatish Balay for (l=0; l<bs; l++) { 44349b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 444b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44549b5e25fSSatish Balay } 44649b5e25fSSatish Balay } 44749b5e25fSSatish Balay } 44849b5e25fSSatish Balay } 44949b5e25fSSatish Balay 450b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 45149b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45249b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 453899cda47SBarry Smith y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0; 45449b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45549b5e25fSSatish Balay aa = a->a + j*bs2; 45649b5e25fSSatish Balay for (k=0; k<bs; k++) { 45749b5e25fSSatish Balay for (l=0; l<bs; l++) { 45849b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 459b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46049b5e25fSSatish Balay } 46149b5e25fSSatish Balay } 46249b5e25fSSatish Balay } 46349b5e25fSSatish Balay } 46449b5e25fSSatish Balay PetscFunctionReturn(0); 46549b5e25fSSatish Balay } 46649b5e25fSSatish Balay 4674a2ae208SSatish Balay #undef __FUNCT__ 4684a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4696849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 47049b5e25fSSatish Balay { 471dfbe8321SBarry Smith PetscErrorCode ierr; 47249b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 473b0a32e0cSBarry Smith PetscDraw draw; 47449b5e25fSSatish Balay PetscTruth isnull; 47549b5e25fSSatish Balay 47649b5e25fSSatish Balay PetscFunctionBegin; 47749b5e25fSSatish Balay 478b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 479b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 48049b5e25fSSatish Balay 48149b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 482899cda47SBarry Smith xr = A->rmap.N; yr = A->rmap.N; h = yr/10.0; w = xr/10.0; 48349b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 484b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 485b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 48649b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 48749b5e25fSSatish Balay PetscFunctionReturn(0); 48849b5e25fSSatish Balay } 48949b5e25fSSatish Balay 4904a2ae208SSatish Balay #undef __FUNCT__ 4914a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 492dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 49349b5e25fSSatish Balay { 494dfbe8321SBarry Smith PetscErrorCode ierr; 49532077d6dSBarry Smith PetscTruth iascii,isdraw; 49649b5e25fSSatish Balay 49749b5e25fSSatish Balay PetscFunctionBegin; 49832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 499fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 50032077d6dSBarry Smith if (iascii){ 50149b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 50249b5e25fSSatish Balay } else if (isdraw) { 50349b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 50449b5e25fSSatish Balay } else { 505a5e6ed63SBarry Smith Mat B; 506ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 507a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 508a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 50949b5e25fSSatish Balay } 51049b5e25fSSatish Balay PetscFunctionReturn(0); 51149b5e25fSSatish Balay } 51249b5e25fSSatish Balay 51349b5e25fSSatish Balay 5144a2ae208SSatish Balay #undef __FUNCT__ 5154a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 51613f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 51749b5e25fSSatish Balay { 518045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 51913f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 52013f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 521899cda47SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->rmap.bs,bs2=a->bs2; 52249b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 52349b5e25fSSatish Balay 52449b5e25fSSatish Balay PetscFunctionBegin; 52549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 52649b5e25fSSatish Balay row = im[k]; brow = row/bs; 52777431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 528899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 52949b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 53049b5e25fSSatish Balay nrow = ailen[brow]; 53149b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 53277431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 533899cda47SBarry 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); 53449b5e25fSSatish Balay col = in[l] ; 53549b5e25fSSatish Balay bcol = col/bs; 53649b5e25fSSatish Balay cidx = col%bs; 53749b5e25fSSatish Balay ridx = row%bs; 53849b5e25fSSatish Balay high = nrow; 53949b5e25fSSatish Balay low = 0; /* assume unsorted */ 54049b5e25fSSatish Balay while (high-low > 5) { 54149b5e25fSSatish Balay t = (low+high)/2; 54249b5e25fSSatish Balay if (rp[t] > bcol) high = t; 54349b5e25fSSatish Balay else low = t; 54449b5e25fSSatish Balay } 54549b5e25fSSatish Balay for (i=low; i<high; i++) { 54649b5e25fSSatish Balay if (rp[i] > bcol) break; 54749b5e25fSSatish Balay if (rp[i] == bcol) { 54849b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 54949b5e25fSSatish Balay goto finished; 55049b5e25fSSatish Balay } 55149b5e25fSSatish Balay } 55249b5e25fSSatish Balay *v++ = zero; 55349b5e25fSSatish Balay finished:; 55449b5e25fSSatish Balay } 55549b5e25fSSatish Balay } 55649b5e25fSSatish Balay PetscFunctionReturn(0); 55749b5e25fSSatish Balay } 55849b5e25fSSatish Balay 55949b5e25fSSatish Balay 5604a2ae208SSatish Balay #undef __FUNCT__ 5614a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 56213f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 56349b5e25fSSatish Balay { 5640880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5656849ba73SBarry Smith PetscErrorCode ierr; 566e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 56713f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 568899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap.bs,stepval; 5690880e062SHong Zhang PetscTruth roworiented=a->roworiented; 570f15d580aSBarry Smith const MatScalar *value = v; 571f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5720880e062SHong Zhang 57349b5e25fSSatish Balay PetscFunctionBegin; 5740880e062SHong Zhang if (roworiented) { 5750880e062SHong Zhang stepval = (n-1)*bs; 5760880e062SHong Zhang } else { 5770880e062SHong Zhang stepval = (m-1)*bs; 5780880e062SHong Zhang } 5790880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5800880e062SHong Zhang row = im[k]; 5810880e062SHong Zhang if (row < 0) continue; 5822515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 58377431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5840880e062SHong Zhang #endif 5850880e062SHong Zhang rp = aj + ai[row]; 5860880e062SHong Zhang ap = aa + bs2*ai[row]; 5870880e062SHong Zhang rmax = imax[row]; 5880880e062SHong Zhang nrow = ailen[row]; 5890880e062SHong Zhang low = 0; 590818f2c47SBarry Smith high = nrow; 5910880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5920880e062SHong Zhang if (in[l] < 0) continue; 5930880e062SHong Zhang col = in[l]; 5942515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 59577431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 596b1823623SSatish Balay #endif 5970880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 5980880e062SHong Zhang if (roworiented) { 5990880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6000880e062SHong Zhang } else { 6010880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6020880e062SHong Zhang } 6037cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 604e2ee6c50SBarry Smith lastcol = col; 6050880e062SHong Zhang while (high-low > 7) { 6060880e062SHong Zhang t = (low+high)/2; 6070880e062SHong Zhang if (rp[t] > col) high = t; 6080880e062SHong Zhang else low = t; 6090880e062SHong Zhang } 6100880e062SHong Zhang for (i=low; i<high; i++) { 6110880e062SHong Zhang if (rp[i] > col) break; 6120880e062SHong Zhang if (rp[i] == col) { 6130880e062SHong Zhang bap = ap + bs2*i; 6140880e062SHong Zhang if (roworiented) { 6150880e062SHong Zhang if (is == ADD_VALUES) { 6160880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6170880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6180880e062SHong Zhang bap[jj] += *value++; 6190880e062SHong Zhang } 6200880e062SHong Zhang } 6210880e062SHong Zhang } else { 6220880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6230880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6240880e062SHong Zhang bap[jj] = *value++; 6250880e062SHong Zhang } 6260880e062SHong Zhang } 6270880e062SHong Zhang } 6280880e062SHong Zhang } else { 6290880e062SHong Zhang if (is == ADD_VALUES) { 6300880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6310880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6320880e062SHong Zhang *bap++ += *value++; 6330880e062SHong Zhang } 6340880e062SHong Zhang } 6350880e062SHong Zhang } else { 6360880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6370880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6380880e062SHong Zhang *bap++ = *value++; 6390880e062SHong Zhang } 6400880e062SHong Zhang } 6410880e062SHong Zhang } 6420880e062SHong Zhang } 6430880e062SHong Zhang goto noinsert2; 6440880e062SHong Zhang } 6450880e062SHong Zhang } 6460880e062SHong Zhang if (nonew == 1) goto noinsert2; 647085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 648e6b907acSBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew); 649c03d1d03SSatish Balay N = nrow++ - 1; high++; 6500880e062SHong Zhang /* shift up all the later entries in this row */ 6510880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6520880e062SHong Zhang rp[ii+1] = rp[ii]; 6530880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6540880e062SHong Zhang } 6550880e062SHong Zhang if (N >= i) { 6560880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6570880e062SHong Zhang } 6580880e062SHong Zhang rp[i] = col; 6590880e062SHong Zhang bap = ap + bs2*i; 6600880e062SHong Zhang if (roworiented) { 6610880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6620880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6630880e062SHong Zhang bap[jj] = *value++; 6640880e062SHong Zhang } 6650880e062SHong Zhang } 6660880e062SHong Zhang } else { 6670880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6680880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6690880e062SHong Zhang *bap++ = *value++; 6700880e062SHong Zhang } 6710880e062SHong Zhang } 6720880e062SHong Zhang } 6730880e062SHong Zhang noinsert2:; 6740880e062SHong Zhang low = i; 6750880e062SHong Zhang } 6760880e062SHong Zhang ailen[row] = nrow; 6770880e062SHong Zhang } 6780880e062SHong Zhang PetscFunctionReturn(0); 67949b5e25fSSatish Balay } 68049b5e25fSSatish Balay 6814a2ae208SSatish Balay #undef __FUNCT__ 6824a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 683dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 68449b5e25fSSatish Balay { 68549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6866849ba73SBarry Smith PetscErrorCode ierr; 68713f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 688899cda47SBarry Smith PetscInt m = A->rmap.N,*ip,N,*ailen = a->ilen; 68913f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 69049b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 69149b5e25fSSatish Balay 69249b5e25fSSatish Balay PetscFunctionBegin; 69349b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 69449b5e25fSSatish Balay 69549b5e25fSSatish Balay if (m) rmax = ailen[0]; 69649b5e25fSSatish Balay for (i=1; i<mbs; i++) { 69749b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 69849b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 69949b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 70049b5e25fSSatish Balay if (fshift) { 70149b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 70249b5e25fSSatish Balay N = ailen[i]; 70349b5e25fSSatish Balay for (j=0; j<N; j++) { 70449b5e25fSSatish Balay ip[j-fshift] = ip[j]; 70549b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 70649b5e25fSSatish Balay } 70749b5e25fSSatish Balay } 70849b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 70949b5e25fSSatish Balay } 71049b5e25fSSatish Balay if (mbs) { 71149b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 71249b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 71349b5e25fSSatish Balay } 71449b5e25fSSatish Balay /* reset ilen and imax for each row */ 71549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 71649b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 71749b5e25fSSatish Balay } 7186c6c5352SBarry Smith a->nz = ai[mbs]; 71949b5e25fSSatish Balay 720b424e231SHong Zhang /* diagonals may have moved, reset it */ 721b424e231SHong Zhang if (a->diag) { 72213f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 72349b5e25fSSatish Balay } 724899cda47SBarry 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); 725ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 726ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 72749b5e25fSSatish Balay a->reallocs = 0; 72849b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 72949b5e25fSSatish Balay PetscFunctionReturn(0); 73049b5e25fSSatish Balay } 73149b5e25fSSatish Balay 73249b5e25fSSatish Balay /* 73349b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 73449b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 73549b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 73649b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 73749b5e25fSSatish Balay */ 7384a2ae208SSatish Balay #undef __FUNCT__ 7394a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 74013f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 74149b5e25fSSatish Balay { 74213f74950SBarry Smith PetscInt i,j,k,row; 74349b5e25fSSatish Balay PetscTruth flg; 74449b5e25fSSatish Balay 74549b5e25fSSatish Balay PetscFunctionBegin; 74649b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 74749b5e25fSSatish Balay row = idx[i]; 74849b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 74949b5e25fSSatish Balay sizes[j] = 1; 75049b5e25fSSatish Balay i++; 75149b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 75249b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 75349b5e25fSSatish Balay i++; 75449b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 75549b5e25fSSatish Balay flg = PETSC_TRUE; 75649b5e25fSSatish Balay for (k=1; k<bs; k++) { 75749b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 75849b5e25fSSatish Balay flg = PETSC_FALSE; 75949b5e25fSSatish Balay break; 76049b5e25fSSatish Balay } 76149b5e25fSSatish Balay } 762abc0a331SBarry Smith if (flg) { /* No break in the bs */ 76349b5e25fSSatish Balay sizes[j] = bs; 76449b5e25fSSatish Balay i+= bs; 76549b5e25fSSatish Balay } else { 76649b5e25fSSatish Balay sizes[j] = 1; 76749b5e25fSSatish Balay i++; 76849b5e25fSSatish Balay } 76949b5e25fSSatish Balay } 77049b5e25fSSatish Balay } 77149b5e25fSSatish Balay *bs_max = j; 77249b5e25fSSatish Balay PetscFunctionReturn(0); 77349b5e25fSSatish Balay } 77449b5e25fSSatish Balay 77549b5e25fSSatish Balay 77649b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 77749b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 77849b5e25fSSatish Balay */ 77949b5e25fSSatish Balay 7804a2ae208SSatish Balay #undef __FUNCT__ 7814a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 78213f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 78349b5e25fSSatish Balay { 78449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7856849ba73SBarry Smith PetscErrorCode ierr; 786e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 78713f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 788899cda47SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol; 78913f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 79049b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 79149b5e25fSSatish Balay 79249b5e25fSSatish Balay PetscFunctionBegin; 79349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 79449b5e25fSSatish Balay row = im[k]; /* row number */ 79549b5e25fSSatish Balay brow = row/bs; /* block row number */ 79649b5e25fSSatish Balay if (row < 0) continue; 7972515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 798899cda47SBarry Smith if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1); 79949b5e25fSSatish Balay #endif 80049b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 80149b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 80249b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 80349b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 80449b5e25fSSatish Balay low = 0; 80549b5e25fSSatish Balay 80649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 80749b5e25fSSatish Balay if (in[l] < 0) continue; 8082515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 809899cda47SBarry 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); 81049b5e25fSSatish Balay #endif 81149b5e25fSSatish Balay col = in[l]; 81249b5e25fSSatish Balay bcol = col/bs; /* block col number */ 81349b5e25fSSatish Balay 814941593c8SHong Zhang if (brow > bcol) { 815941593c8SHong Zhang if (a->ignore_ltriangular){ 816941593c8SHong Zhang continue; /* ignore lower triangular values */ 817941593c8SHong Zhang } else { 818941593c8SHong Zhang SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR)"); 819941593c8SHong Zhang } 820941593c8SHong Zhang } 821f4989cb3SHong Zhang 82249b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8238549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 82449b5e25fSSatish Balay /* element value a(k,l) */ 82549b5e25fSSatish Balay if (roworiented) { 82649b5e25fSSatish Balay value = v[l + k*n]; 82749b5e25fSSatish Balay } else { 82849b5e25fSSatish Balay value = v[k + l*m]; 82949b5e25fSSatish Balay } 83049b5e25fSSatish Balay 83149b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 8327cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 833e2ee6c50SBarry Smith lastcol = col; 83449b5e25fSSatish Balay while (high-low > 7) { 83549b5e25fSSatish Balay t = (low+high)/2; 83649b5e25fSSatish Balay if (rp[t] > bcol) high = t; 83749b5e25fSSatish Balay else low = t; 83849b5e25fSSatish Balay } 83949b5e25fSSatish Balay for (i=low; i<high; i++) { 84049b5e25fSSatish Balay if (rp[i] > bcol) break; 84149b5e25fSSatish Balay if (rp[i] == bcol) { 84249b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 84349b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 84449b5e25fSSatish Balay else *bap = value; 8458549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8468549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8478549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8488549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8498549e402SHong Zhang else *bap = value; 8508549e402SHong Zhang } 85149b5e25fSSatish Balay goto noinsert1; 85249b5e25fSSatish Balay } 85349b5e25fSSatish Balay } 85449b5e25fSSatish Balay 85549b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 856085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 857e6b907acSBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew); 85849b5e25fSSatish Balay 859c03d1d03SSatish Balay N = nrow++ - 1; high++; 86049b5e25fSSatish Balay /* shift up all the later entries in this row */ 86149b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 86249b5e25fSSatish Balay rp[ii+1] = rp[ii]; 86349b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 86449b5e25fSSatish Balay } 86549b5e25fSSatish Balay if (N>=i) { 86649b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 86749b5e25fSSatish Balay } 86849b5e25fSSatish Balay rp[i] = bcol; 86949b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 87049b5e25fSSatish Balay noinsert1:; 87149b5e25fSSatish Balay low = i; 8728549e402SHong Zhang } 87349b5e25fSSatish Balay } /* end of loop over added columns */ 87449b5e25fSSatish Balay ailen[brow] = nrow; 87549b5e25fSSatish Balay } /* end of loop over added rows */ 87649b5e25fSSatish Balay PetscFunctionReturn(0); 87749b5e25fSSatish Balay } 87849b5e25fSSatish Balay 8796849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 8806849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 88113f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt); 88213f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 88313f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]); 884f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar); 8856849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 8866849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 8876849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec); 8886849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 8896849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 8906849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat); 8916849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec); 89249b5e25fSSatish Balay 8936849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 8946849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 8956849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 8966849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 8976849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 8986849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 8996849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 9006849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 90149b5e25fSSatish Balay 9026849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 903d59c15a7SBarry Smith 9046849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 9056849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 9066849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 9076849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 9086849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 9096849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 9106849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 9116849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 91207f98182SSatish Balay 913af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,MatFactorInfo*,Mat*); 914af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,MatFactorInfo*,Mat*); 915af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,MatFactorInfo*,Mat*); 916af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,MatFactorInfo*,Mat*); 917af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,MatFactorInfo*,Mat*); 918af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,MatFactorInfo*,Mat*); 919af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,MatFactorInfo*,Mat*); 920af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,MatFactorInfo*,Mat*); 92113f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*); 92249b5e25fSSatish Balay 9236849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 9246849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 9256849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 9266849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 9276849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 9286849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 9296849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 9306849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 93149b5e25fSSatish Balay 9326849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 9336849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 9346849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 9356849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 9366849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 9376849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 9386849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 9396849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 94049b5e25fSSatish Balay 9414d101231SSatish Balay #ifdef HAVE_MatICCFactor 9426849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work! */ 9434a2ae208SSatish Balay #undef __FUNCT__ 9444d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 945dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 94649b5e25fSSatish Balay { 9474ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 94849b5e25fSSatish Balay Mat outA; 949dfbe8321SBarry Smith PetscErrorCode ierr; 95049b5e25fSSatish Balay PetscTruth row_identity,col_identity; 95149b5e25fSSatish Balay 95249b5e25fSSatish Balay PetscFunctionBegin; 95349b5e25fSSatish Balay outA = inA; 9541260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 95549b5e25fSSatish Balay 9561a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 95749b5e25fSSatish Balay /* 95849b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 95949b5e25fSSatish Balay for ILU(0) factorization with natural ordering 96049b5e25fSSatish Balay */ 961899cda47SBarry Smith switch (a->rmap.bs) { 96249b5e25fSSatish Balay case 1: 9630fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 96407f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 965d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 966ae15b995SBarry Smith ierr = PetscInfo((inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr); 96749b5e25fSSatish Balay case 2: 9681a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9691a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9701a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 971ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 97249b5e25fSSatish Balay break; 97349b5e25fSSatish Balay case 3: 9741a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9751a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9761a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 977ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 97849b5e25fSSatish Balay break; 97949b5e25fSSatish Balay case 4: 9801a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 9811a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 9821a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 983ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 98449b5e25fSSatish Balay break; 98549b5e25fSSatish Balay case 5: 9861a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 9871a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 9881a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 989ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 99049b5e25fSSatish Balay break; 99149b5e25fSSatish Balay case 6: 9921a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 9931a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 9941a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 995ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 99649b5e25fSSatish Balay break; 99749b5e25fSSatish Balay case 7: 9981a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 9991a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10001a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 1001ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 100249b5e25fSSatish Balay break; 100349b5e25fSSatish Balay default: 100449b5e25fSSatish Balay a->row = row; 10051a3463dfSHong Zhang a->icol = col; 100649b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 100749b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 100849b5e25fSSatish Balay 100949b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 101049b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 101152e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 101249b5e25fSSatish Balay 101349b5e25fSSatish Balay if (!a->solve_work) { 1014899cda47SBarry Smith ierr = PetscMalloc((A->rmap.N+a->rmap.bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 1015899cda47SBarry Smith ierr = PetscLogObjectMemory(inA,(A->rmap.N+a->rmap.bs)*sizeof(PetscScalar));CHKERRQ(ierr); 101649b5e25fSSatish Balay } 101749b5e25fSSatish Balay } 101849b5e25fSSatish Balay 1019af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 102049b5e25fSSatish Balay PetscFunctionReturn(0); 102149b5e25fSSatish Balay } 1022950f1e5bSHong Zhang #endif 1023950f1e5bSHong Zhang 10244a2ae208SSatish Balay #undef __FUNCT__ 10254a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 1026dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A) 102749b5e25fSSatish Balay { 102849b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 102949b5e25fSSatish Balay MPI_Comm comm = A->comm; 1030dfbe8321SBarry Smith PetscErrorCode ierr; 103149b5e25fSSatish Balay 103249b5e25fSSatish Balay PetscFunctionBegin; 103349b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 103449b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 103549b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 1036f5edf698SHong Zhang ierr = (*PetscHelpPrintf)(comm," -mat_ignore_lower_triangular: Ignore lower triangular values set by user\n");CHKERRQ(ierr); 1037f5edf698SHong Zhang ierr = (*PetscHelpPrintf)(comm," -mat_getrow_uppertriangular: Enable MatGetRow() for retrieving the upper triangular part of the row\n");CHKERRQ(ierr); 103849b5e25fSSatish Balay PetscFunctionReturn(0); 103949b5e25fSSatish Balay } 104049b5e25fSSatish Balay 104149b5e25fSSatish Balay EXTERN_C_BEGIN 10424a2ae208SSatish Balay #undef __FUNCT__ 10434a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1044be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 104549b5e25fSSatish Balay { 1046045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 104713f74950SBarry Smith PetscInt i,nz,n; 104849b5e25fSSatish Balay 104949b5e25fSSatish Balay PetscFunctionBegin; 10506c6c5352SBarry Smith nz = baij->maxnz; 1051899cda47SBarry Smith n = mat->cmap.n; 105249b5e25fSSatish Balay for (i=0; i<nz; i++) { 105349b5e25fSSatish Balay baij->j[i] = indices[i]; 105449b5e25fSSatish Balay } 10556c6c5352SBarry Smith baij->nz = nz; 105649b5e25fSSatish Balay for (i=0; i<n; i++) { 105749b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 105849b5e25fSSatish Balay } 105949b5e25fSSatish Balay PetscFunctionReturn(0); 106049b5e25fSSatish Balay } 106149b5e25fSSatish Balay EXTERN_C_END 106249b5e25fSSatish Balay 10634a2ae208SSatish Balay #undef __FUNCT__ 10644a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 106549b5e25fSSatish Balay /*@ 106619585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 106749b5e25fSSatish Balay in the matrix. 106849b5e25fSSatish Balay 106949b5e25fSSatish Balay Input Parameters: 107019585528SSatish Balay + mat - the SeqSBAIJ matrix 107149b5e25fSSatish Balay - indices - the column indices 107249b5e25fSSatish Balay 107349b5e25fSSatish Balay Level: advanced 107449b5e25fSSatish Balay 107549b5e25fSSatish Balay Notes: 107649b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 107749b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 107849b5e25fSSatish Balay of the MatSetValues() operation. 107949b5e25fSSatish Balay 108049b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1081d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 108249b5e25fSSatish Balay 1083ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 108449b5e25fSSatish Balay 1085ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 108649b5e25fSSatish Balay @*/ 1087be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 108849b5e25fSSatish Balay { 108913f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 109049b5e25fSSatish Balay 109149b5e25fSSatish Balay PetscFunctionBegin; 10924482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 10934482741eSBarry Smith PetscValidPointer(indices,2); 1094c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 109549b5e25fSSatish Balay if (f) { 109649b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 109749b5e25fSSatish Balay } else { 1098e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 109949b5e25fSSatish Balay } 110049b5e25fSSatish Balay PetscFunctionReturn(0); 110149b5e25fSSatish Balay } 110249b5e25fSSatish Balay 11034a2ae208SSatish Balay #undef __FUNCT__ 11043c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 11053c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 11063c896bc6SHong Zhang { 11073c896bc6SHong Zhang PetscErrorCode ierr; 11083c896bc6SHong Zhang 11093c896bc6SHong Zhang PetscFunctionBegin; 11103c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 11113c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 11123c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 11133c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 11143c896bc6SHong Zhang 1115899cda47SBarry Smith if (a->i[A->rmap.N] != b->i[B->rmap.N]) { 11163c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 11173c896bc6SHong Zhang } 1118899cda47SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.N])*sizeof(PetscScalar));CHKERRQ(ierr); 11193c896bc6SHong Zhang } else { 1120f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 11213c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1122f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 11233c896bc6SHong Zhang } 11243c896bc6SHong Zhang PetscFunctionReturn(0); 11253c896bc6SHong Zhang } 11263c896bc6SHong Zhang 11273c896bc6SHong Zhang #undef __FUNCT__ 11284a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1129dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1130273d9f13SBarry Smith { 1131dfbe8321SBarry Smith PetscErrorCode ierr; 1132273d9f13SBarry Smith 1133273d9f13SBarry Smith PetscFunctionBegin; 11347edd0491SSatish Balay ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1135273d9f13SBarry Smith PetscFunctionReturn(0); 1136273d9f13SBarry Smith } 1137273d9f13SBarry Smith 1138a6ece127SHong Zhang #undef __FUNCT__ 1139a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1140dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1141a6ece127SHong Zhang { 1142a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1143a6ece127SHong Zhang PetscFunctionBegin; 1144a6ece127SHong Zhang *array = a->a; 1145a6ece127SHong Zhang PetscFunctionReturn(0); 1146a6ece127SHong Zhang } 1147a6ece127SHong Zhang 1148a6ece127SHong Zhang #undef __FUNCT__ 1149a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1150dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1151a6ece127SHong Zhang { 1152a6ece127SHong Zhang PetscFunctionBegin; 1153a6ece127SHong Zhang PetscFunctionReturn(0); 1154a6ece127SHong Zhang } 1155a6ece127SHong Zhang 115642ee4b1aSHong Zhang #include "petscblaslapack.h" 115742ee4b1aSHong Zhang #undef __FUNCT__ 115842ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1159f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 116042ee4b1aSHong Zhang { 116142ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1162dfbe8321SBarry Smith PetscErrorCode ierr; 1163899cda47SBarry Smith PetscInt i,bs=Y->rmap.bs,bs2,j; 11644ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 116542ee4b1aSHong Zhang 116642ee4b1aSHong Zhang PetscFunctionBegin; 116742ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1168f4df32b1SMatthew Knepley PetscScalar alpha = a; 1169f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1170c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1171c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1172c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1173c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1174c537a176SHong Zhang } 1175c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1176c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1177c4319e64SHong Zhang y->XtoY = X; 1178c537a176SHong Zhang } 1179c4319e64SHong Zhang bs2 = bs*bs; 11806c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1181c4319e64SHong Zhang j = 0; 1182c4319e64SHong Zhang while (j < bs2){ 1183f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1184c4319e64SHong Zhang j++; 1185c537a176SHong Zhang } 1186c4319e64SHong Zhang } 1187ae15b995SBarry Smith ierr = PetscInfo3(0,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %G\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));CHKERRQ(ierr); 118842ee4b1aSHong Zhang } else { 1189f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1190f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1191f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 119242ee4b1aSHong Zhang } 119342ee4b1aSHong Zhang PetscFunctionReturn(0); 119442ee4b1aSHong Zhang } 119542ee4b1aSHong Zhang 1196efcf0fc3SBarry Smith #undef __FUNCT__ 1197efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1198dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1199efcf0fc3SBarry Smith { 1200efcf0fc3SBarry Smith PetscFunctionBegin; 1201efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1202efcf0fc3SBarry Smith PetscFunctionReturn(0); 1203efcf0fc3SBarry Smith } 1204efcf0fc3SBarry Smith 1205efcf0fc3SBarry Smith #undef __FUNCT__ 1206efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1207dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1208efcf0fc3SBarry Smith { 1209efcf0fc3SBarry Smith PetscFunctionBegin; 1210efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1211efcf0fc3SBarry Smith PetscFunctionReturn(0); 1212efcf0fc3SBarry Smith } 1213efcf0fc3SBarry Smith 1214efcf0fc3SBarry Smith #undef __FUNCT__ 1215efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1216dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1217efcf0fc3SBarry Smith { 1218efcf0fc3SBarry Smith PetscFunctionBegin; 1219efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1220efcf0fc3SBarry Smith PetscFunctionReturn(0); 1221efcf0fc3SBarry Smith } 1222efcf0fc3SBarry Smith 122399cafbc1SBarry Smith #undef __FUNCT__ 122499cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 122599cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 122699cafbc1SBarry Smith { 122799cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 122899cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 122999cafbc1SBarry Smith PetscScalar *aa = a->a; 123099cafbc1SBarry Smith 123199cafbc1SBarry Smith PetscFunctionBegin; 123299cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 123399cafbc1SBarry Smith PetscFunctionReturn(0); 123499cafbc1SBarry Smith } 123599cafbc1SBarry Smith 123699cafbc1SBarry Smith #undef __FUNCT__ 123799cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 123899cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 123999cafbc1SBarry Smith { 124099cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 124199cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 124299cafbc1SBarry Smith PetscScalar *aa = a->a; 124399cafbc1SBarry Smith 124499cafbc1SBarry Smith PetscFunctionBegin; 124599cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 124699cafbc1SBarry Smith PetscFunctionReturn(0); 124799cafbc1SBarry Smith } 124899cafbc1SBarry Smith 124949b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 125049b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 125149b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 125249b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 125349b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 125497304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1255431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1256e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 125749b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 125849b5e25fSSatish Balay 0, 125949b5e25fSSatish Balay 0, 126097304618SKris Buschelman /*10*/ 0, 126149b5e25fSSatish Balay 0, 12625f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1263d06b337dSHong Zhang MatRelax_SeqSBAIJ, 126449b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 126597304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 126649b5e25fSSatish Balay MatEqual_SeqSBAIJ, 126749b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 126849b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 126949b5e25fSSatish Balay MatNorm_SeqSBAIJ, 127097304618SKris Buschelman /*20*/ 0, 127149b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 127249b5e25fSSatish Balay 0, 127349b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 127449b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1275dcf5cc72SBarry Smith /*25*/ 0, 127649b5e25fSSatish Balay 0, 127749b5e25fSSatish Balay 0, 12785f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12795f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 128097304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1281c464158bSHong Zhang 0, 12824d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1283a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1284a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 128597304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 128649b5e25fSSatish Balay 0, 128749b5e25fSSatish Balay 0, 128849b5e25fSSatish Balay 0, 1289950f1e5bSHong Zhang 0, 129097304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 129149b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 129249b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 129349b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 12943c896bc6SHong Zhang MatCopy_SeqSBAIJ, 129597304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 129649b5e25fSSatish Balay MatScale_SeqSBAIJ, 129749b5e25fSSatish Balay 0, 129849b5e25fSSatish Balay 0, 129949b5e25fSSatish Balay 0, 1300521d7252SBarry Smith /*50*/ 0, 130149b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 130249b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 130349b5e25fSSatish Balay 0, 130449b5e25fSSatish Balay 0, 130597304618SKris Buschelman /*55*/ 0, 130649b5e25fSSatish Balay 0, 130749b5e25fSSatish Balay 0, 130849b5e25fSSatish Balay 0, 130949b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 131097304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 131149b5e25fSSatish Balay 0, 131249b5e25fSSatish Balay 0, 1313357abbc8SBarry Smith 0, 1314d959ec07SHong Zhang 0, 131597304618SKris Buschelman /*65*/ 0, 1316d959ec07SHong Zhang 0, 1317d959ec07SHong Zhang 0, 1318d959ec07SHong Zhang 0, 1319d959ec07SHong Zhang 0, 132097304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13213e0d88b5SBarry Smith 0, 13223e0d88b5SBarry Smith 0, 13233e0d88b5SBarry Smith 0, 13243e0d88b5SBarry Smith 0, 132597304618SKris Buschelman /*75*/ 0, 13263e0d88b5SBarry Smith 0, 13273e0d88b5SBarry Smith 0, 13283e0d88b5SBarry Smith 0, 13293e0d88b5SBarry Smith 0, 133097304618SKris Buschelman /*80*/ 0, 13313e0d88b5SBarry Smith 0, 13323e0d88b5SBarry Smith 0, 13333e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 133497304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13353e0d88b5SBarry Smith #else 133697304618SKris Buschelman 0, 13373e0d88b5SBarry Smith #endif 1338865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1339865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1340865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1341efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1342865e5f61SKris Buschelman 0, 1343865e5f61SKris Buschelman 0, 1344865e5f61SKris Buschelman /*90*/ 0, 1345865e5f61SKris Buschelman 0, 1346865e5f61SKris Buschelman 0, 1347865e5f61SKris Buschelman 0, 1348865e5f61SKris Buschelman 0, 1349865e5f61SKris Buschelman /*95*/ 0, 1350865e5f61SKris Buschelman 0, 1351865e5f61SKris Buschelman 0, 135299cafbc1SBarry Smith 0, 135399cafbc1SBarry Smith 0, 135499cafbc1SBarry Smith /*100*/0, 135599cafbc1SBarry Smith 0, 135699cafbc1SBarry Smith 0, 135799cafbc1SBarry Smith 0, 135899cafbc1SBarry Smith 0, 135999cafbc1SBarry Smith /*105*/0, 136099cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1361f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1362f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 1363f5edf698SHong Zhang MatRestoreRowUpperTriangular_SeqSBAIJ 136499cafbc1SBarry Smith }; 1365be1d678aSKris Buschelman 136649b5e25fSSatish Balay EXTERN_C_BEGIN 13674a2ae208SSatish Balay #undef __FUNCT__ 13684a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1369be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 137049b5e25fSSatish Balay { 13714afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1372899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 1373dfbe8321SBarry Smith PetscErrorCode ierr; 137449b5e25fSSatish Balay 137549b5e25fSSatish Balay PetscFunctionBegin; 137649b5e25fSSatish Balay if (aij->nonew != 1) { 1377e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 137849b5e25fSSatish Balay } 137949b5e25fSSatish Balay 138049b5e25fSSatish Balay /* allocate space for values if not already there */ 138149b5e25fSSatish Balay if (!aij->saved_values) { 138287828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 138349b5e25fSSatish Balay } 138449b5e25fSSatish Balay 138549b5e25fSSatish Balay /* copy values over */ 138687828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 138749b5e25fSSatish Balay PetscFunctionReturn(0); 138849b5e25fSSatish Balay } 138949b5e25fSSatish Balay EXTERN_C_END 139049b5e25fSSatish Balay 139149b5e25fSSatish Balay EXTERN_C_BEGIN 13924a2ae208SSatish Balay #undef __FUNCT__ 13934a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1394be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 139549b5e25fSSatish Balay { 13964afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13976849ba73SBarry Smith PetscErrorCode ierr; 1398899cda47SBarry Smith PetscInt nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2; 139949b5e25fSSatish Balay 140049b5e25fSSatish Balay PetscFunctionBegin; 140149b5e25fSSatish Balay if (aij->nonew != 1) { 1402e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 140349b5e25fSSatish Balay } 140449b5e25fSSatish Balay if (!aij->saved_values) { 1405e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 140649b5e25fSSatish Balay } 140749b5e25fSSatish Balay 140849b5e25fSSatish Balay /* copy values over */ 140987828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 141049b5e25fSSatish Balay PetscFunctionReturn(0); 141149b5e25fSSatish Balay } 141249b5e25fSSatish Balay EXTERN_C_END 141349b5e25fSSatish Balay 14148549e402SHong Zhang EXTERN_C_BEGIN 14154a2ae208SSatish Balay #undef __FUNCT__ 1416a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1417be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 141849b5e25fSSatish Balay { 1419c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14206849ba73SBarry Smith PetscErrorCode ierr; 1421ab93d7beSBarry Smith PetscInt i,mbs,bs2; 1422ab93d7beSBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg; 142349b5e25fSSatish Balay 142449b5e25fSSatish Balay PetscFunctionBegin; 1425273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1426e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1427899cda47SBarry Smith B->rmap.bs = B->cmap.bs = bs; 1428899cda47SBarry Smith ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr); 1429899cda47SBarry Smith ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr); 1430899cda47SBarry Smith 1431899cda47SBarry Smith mbs = B->rmap.N/bs; 143249b5e25fSSatish Balay bs2 = bs*bs; 143349b5e25fSSatish Balay 1434899cda47SBarry Smith if (mbs*bs != B->rmap.N) { 143529bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 143649b5e25fSSatish Balay } 143749b5e25fSSatish Balay 1438ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1439ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1440ab93d7beSBarry Smith nz = 0; 1441ab93d7beSBarry Smith } 1442ab93d7beSBarry Smith 1443435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 144477431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 144549b5e25fSSatish Balay if (nnz) { 144649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 144777431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 144877431f27SBarry 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); 144949b5e25fSSatish Balay } 145049b5e25fSSatish Balay } 145149b5e25fSSatish Balay 1452e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 145349b5e25fSSatish Balay if (!flg) { 145449b5e25fSSatish Balay switch (bs) { 145549b5e25fSSatish Balay case 1: 14564ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 145749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1458d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1459e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 146049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 146149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1462431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1463431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 146449b5e25fSSatish Balay break; 146549b5e25fSSatish Balay case 2: 14664ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 146749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1468e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 146949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 147049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1471431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1472431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 147349b5e25fSSatish Balay break; 147449b5e25fSSatish Balay case 3: 14755f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 147649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1477e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 147849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 147949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1480431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1481431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 148249b5e25fSSatish Balay break; 148349b5e25fSSatish Balay case 4: 14845f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 148549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1486e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 148749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 148849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1489431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1490431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 149149b5e25fSSatish Balay break; 149249b5e25fSSatish Balay case 5: 14935f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 149449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1495e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 149649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 149749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1498431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1499431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 150049b5e25fSSatish Balay break; 150149b5e25fSSatish Balay case 6: 15025f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 150349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1504e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 150549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 150649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1507431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1508431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 150949b5e25fSSatish Balay break; 151049b5e25fSSatish Balay case 7: 1511de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 151249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1513e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1514de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 151549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1516431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1517431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 151849b5e25fSSatish Balay break; 151949b5e25fSSatish Balay } 152049b5e25fSSatish Balay } 152149b5e25fSSatish Balay 152249b5e25fSSatish Balay b->mbs = mbs; 15234afc71dfSHong Zhang b->nbs = mbs; 1524ab93d7beSBarry Smith if (!skipallocation) { 1525ab93d7beSBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1526ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1527ab93d7beSBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 152849b5e25fSSatish Balay if (!nnz) { 1529435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 153049b5e25fSSatish Balay else if (nz <= 0) nz = 1; 153149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15328cef66ccSHong Zhang b->imax[i] = nz; 153349b5e25fSSatish Balay } 1534153ea458SHong Zhang nz = nz*mbs; /* total nz */ 153549b5e25fSSatish Balay } else { 153649b5e25fSSatish Balay nz = 0; 15378cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 153849b5e25fSSatish Balay } 15396c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 154049b5e25fSSatish Balay 154149b5e25fSSatish Balay /* allocate the matrix space */ 1542899cda47SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.N+1,PetscInt,&b->i);CHKERRQ(ierr); 15436c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 154413f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 154549b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 154649b5e25fSSatish Balay 154749b5e25fSSatish Balay /* pointer to beginning of each row */ 154849b5e25fSSatish Balay b->i[0] = 0; 154949b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 155049b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 155149b5e25fSSatish Balay } 1552e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1553e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1554e811da20SHong Zhang } else { 1555e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1556e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1557ab93d7beSBarry Smith } 155849b5e25fSSatish Balay 1559899cda47SBarry Smith B->rmap.bs = bs; 156049b5e25fSSatish Balay b->bs2 = bs2; 15616c6c5352SBarry Smith b->nz = 0; 15626c6c5352SBarry Smith b->maxnz = nz*bs2; 1563153ea458SHong Zhang 156416cdd363SHong Zhang b->inew = 0; 156516cdd363SHong Zhang b->jnew = 0; 156616cdd363SHong Zhang b->anew = 0; 156716cdd363SHong Zhang b->a2anew = 0; 15681a3463dfSHong Zhang b->permute = PETSC_FALSE; 1569c464158bSHong Zhang PetscFunctionReturn(0); 1570c464158bSHong Zhang } 1571a23d5eceSKris Buschelman EXTERN_C_END 1572153ea458SHong Zhang 1573d769727bSBarry Smith EXTERN_C_BEGIN 1574f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1575f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1576d769727bSBarry Smith EXTERN_C_END 1577d769727bSBarry Smith 15780bad9183SKris Buschelman /*MC 1579fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15800bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15810bad9183SKris Buschelman 15820bad9183SKris Buschelman Options Database Keys: 15830bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15840bad9183SKris Buschelman 15850bad9183SKris Buschelman Level: beginner 15860bad9183SKris Buschelman 15870bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15880bad9183SKris Buschelman M*/ 15890bad9183SKris Buschelman 1590a23d5eceSKris Buschelman EXTERN_C_BEGIN 1591a23d5eceSKris Buschelman #undef __FUNCT__ 1592a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1593be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1594a23d5eceSKris Buschelman { 1595a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1596dfbe8321SBarry Smith PetscErrorCode ierr; 159713f74950SBarry Smith PetscMPIInt size; 1598941593c8SHong Zhang PetscTruth flg; 1599a23d5eceSKris Buschelman 1600a23d5eceSKris Buschelman PetscFunctionBegin; 1601a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1602a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1603a23d5eceSKris Buschelman 1604a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1605a23d5eceSKris Buschelman B->data = (void*)b; 1606a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1607a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1608a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1609a23d5eceSKris Buschelman B->factor = 0; 1610a23d5eceSKris Buschelman B->mapping = 0; 1611a23d5eceSKris Buschelman b->row = 0; 1612a23d5eceSKris Buschelman b->icol = 0; 1613a23d5eceSKris Buschelman b->reallocs = 0; 1614a23d5eceSKris Buschelman b->saved_values = 0; 1615a23d5eceSKris Buschelman 1616a23d5eceSKris Buschelman 1617a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1618a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1619a23d5eceSKris Buschelman b->nonew = 0; 1620a23d5eceSKris Buschelman b->diag = 0; 1621a23d5eceSKris Buschelman b->solve_work = 0; 1622a23d5eceSKris Buschelman b->mult_work = 0; 1623a23d5eceSKris Buschelman B->spptr = 0; 1624a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1625a23d5eceSKris Buschelman b->xtoy = 0; 1626a23d5eceSKris Buschelman b->XtoY = 0; 1627a23d5eceSKris Buschelman 1628a23d5eceSKris Buschelman b->inew = 0; 1629a23d5eceSKris Buschelman b->jnew = 0; 1630a23d5eceSKris Buschelman b->anew = 0; 1631a23d5eceSKris Buschelman b->a2anew = 0; 1632a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1633a23d5eceSKris Buschelman 1634941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 1635941593c8SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr); 1636941593c8SHong Zhang if (flg) b->ignore_ltriangular = PETSC_TRUE; 1637941593c8SHong Zhang 1638f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 1639f5edf698SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr); 1640f5edf698SHong Zhang if (flg) b->getrow_utriangular = PETSC_TRUE; 1641f5edf698SHong Zhang 1642a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1643a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1644a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1645a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1646a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1647a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1648a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1649a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1650a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16514e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16524e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16534e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1654a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1655a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1656a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1657a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1658a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1659a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 166023ce1328SBarry Smith 166123ce1328SBarry Smith B->symmetric = PETSC_TRUE; 166223ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 166323ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 166423ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 1665a23d5eceSKris Buschelman PetscFunctionReturn(0); 1666a23d5eceSKris Buschelman } 1667a23d5eceSKris Buschelman EXTERN_C_END 1668a23d5eceSKris Buschelman 1669a23d5eceSKris Buschelman #undef __FUNCT__ 1670a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1671a23d5eceSKris Buschelman /*@C 1672a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1673a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1674a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1675a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1676a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1677a23d5eceSKris Buschelman 1678a23d5eceSKris Buschelman Collective on Mat 1679a23d5eceSKris Buschelman 1680a23d5eceSKris Buschelman Input Parameters: 1681a23d5eceSKris Buschelman + A - the symmetric matrix 1682a23d5eceSKris Buschelman . bs - size of block 1683a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1684a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1685a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1686a23d5eceSKris Buschelman 1687a23d5eceSKris Buschelman Options Database Keys: 1688a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1689a23d5eceSKris Buschelman block calculations (much slower) 1690a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1691a23d5eceSKris Buschelman 1692a23d5eceSKris Buschelman Level: intermediate 1693a23d5eceSKris Buschelman 1694a23d5eceSKris Buschelman Notes: 1695a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1696a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1697a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1698a23d5eceSKris Buschelman matrices. 1699a23d5eceSKris Buschelman 170049a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 170149a6f317SBarry Smith 170249a6f317SBarry Smith 1703a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1704a23d5eceSKris Buschelman @*/ 1705be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 170613f74950SBarry Smith { 170713f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1708a23d5eceSKris Buschelman 1709a23d5eceSKris Buschelman PetscFunctionBegin; 1710a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1711a23d5eceSKris Buschelman if (f) { 1712a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1713a23d5eceSKris Buschelman } 1714a23d5eceSKris Buschelman PetscFunctionReturn(0); 1715a23d5eceSKris Buschelman } 171649b5e25fSSatish Balay 17174a2ae208SSatish Balay #undef __FUNCT__ 17184a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1719c464158bSHong Zhang /*@C 1720c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1721c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1722c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1723c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1724c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 172549b5e25fSSatish Balay 1726c464158bSHong Zhang Collective on MPI_Comm 1727c464158bSHong Zhang 1728c464158bSHong Zhang Input Parameters: 1729c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1730c464158bSHong Zhang . bs - size of block 1731c464158bSHong Zhang . m - number of rows, or number of columns 1732c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1733744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1734744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1735c464158bSHong Zhang 1736c464158bSHong Zhang Output Parameter: 1737c464158bSHong Zhang . A - the symmetric matrix 1738c464158bSHong Zhang 1739c464158bSHong Zhang Options Database Keys: 1740c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1741c464158bSHong Zhang block calculations (much slower) 1742c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1743c464158bSHong Zhang 1744c464158bSHong Zhang Level: intermediate 1745c464158bSHong Zhang 1746c464158bSHong Zhang Notes: 1747c464158bSHong Zhang 1748c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1749c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1750c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1751c464158bSHong Zhang matrices. 1752c464158bSHong Zhang 175349a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 175449a6f317SBarry Smith 1755c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1756c464158bSHong Zhang @*/ 1757be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1758c464158bSHong Zhang { 1759dfbe8321SBarry Smith PetscErrorCode ierr; 1760c464158bSHong Zhang 1761c464158bSHong Zhang PetscFunctionBegin; 1762f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1763f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1764c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1765ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 176649b5e25fSSatish Balay PetscFunctionReturn(0); 176749b5e25fSSatish Balay } 176849b5e25fSSatish Balay 17694a2ae208SSatish Balay #undef __FUNCT__ 17704a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1771dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 177249b5e25fSSatish Balay { 177349b5e25fSSatish Balay Mat C; 177449b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17756849ba73SBarry Smith PetscErrorCode ierr; 1776b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 177749b5e25fSSatish Balay 177849b5e25fSSatish Balay PetscFunctionBegin; 177929bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 178049b5e25fSSatish Balay 178149b5e25fSSatish Balay *B = 0; 1782f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); 1783899cda47SBarry Smith ierr = MatSetSizes(C,A->rmap.N,A->cmap.n,A->rmap.N,A->cmap.n);CHKERRQ(ierr); 1784be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17851d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1786692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1787692f9cbeSHong Zhang 1788273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 178949b5e25fSSatish Balay C->factor = A->factor; 179049b5e25fSSatish Balay c->row = 0; 179149b5e25fSSatish Balay c->icol = 0; 179249b5e25fSSatish Balay c->saved_values = 0; 179349b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 179449b5e25fSSatish Balay C->assembled = PETSC_TRUE; 179549b5e25fSSatish Balay 1796899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->rmap,&C->rmap);CHKERRQ(ierr); 1797899cda47SBarry Smith ierr = PetscMapCopy(A->comm,&A->cmap,&C->cmap);CHKERRQ(ierr); 179849b5e25fSSatish Balay c->bs2 = a->bs2; 179949b5e25fSSatish Balay c->mbs = a->mbs; 180049b5e25fSSatish Balay c->nbs = a->nbs; 180149b5e25fSSatish Balay 18028777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 180349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 180449b5e25fSSatish Balay c->imax[i] = a->imax[i]; 180549b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 180649b5e25fSSatish Balay } 180749b5e25fSSatish Balay 180849b5e25fSSatish Balay /* allocate the matrix space */ 1809b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 181049b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 181113f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1812b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 181349b5e25fSSatish Balay if (mbs > 0) { 181413f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 181549b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 181649b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 181749b5e25fSSatish Balay } else { 181849b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 181949b5e25fSSatish Balay } 182049b5e25fSSatish Balay } 182149b5e25fSSatish Balay 182249b5e25fSSatish Balay c->sorted = a->sorted; 182349b5e25fSSatish Balay c->roworiented = a->roworiented; 182449b5e25fSSatish Balay c->nonew = a->nonew; 182549b5e25fSSatish Balay 182649b5e25fSSatish Balay if (a->diag) { 182713f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 182852e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 182949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 183049b5e25fSSatish Balay c->diag[i] = a->diag[i]; 183149b5e25fSSatish Balay } 183249b5e25fSSatish Balay } else c->diag = 0; 18336c6c5352SBarry Smith c->nz = a->nz; 18346c6c5352SBarry Smith c->maxnz = a->maxnz; 183549b5e25fSSatish Balay c->solve_work = 0; 183649b5e25fSSatish Balay c->mult_work = 0; 1837e6b907acSBarry Smith c->free_a = PETSC_TRUE; 1838e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 183949b5e25fSSatish Balay *B = C; 1840b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 184149b5e25fSSatish Balay PetscFunctionReturn(0); 184249b5e25fSSatish Balay } 184349b5e25fSSatish Balay 18444a2ae208SSatish Balay #undef __FUNCT__ 18454a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1846f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A) 184749b5e25fSSatish Balay { 184849b5e25fSSatish Balay Mat_SeqSBAIJ *a; 184949b5e25fSSatish Balay Mat B; 18506849ba73SBarry Smith PetscErrorCode ierr; 185113f74950SBarry Smith int fd; 185213f74950SBarry Smith PetscMPIInt size; 185313f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 185413f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 185513f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 185613f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 185787828ca2SBarry Smith PetscScalar *aa; 185849b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 185949b5e25fSSatish Balay 186049b5e25fSSatish Balay PetscFunctionBegin; 1861b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 186249b5e25fSSatish Balay bs2 = bs*bs; 186349b5e25fSSatish Balay 186449b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 186529bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1866b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 186749b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1868552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 186949b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 187049b5e25fSSatish Balay 187149b5e25fSSatish Balay if (header[3] < 0) { 187229bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 187349b5e25fSSatish Balay } 187449b5e25fSSatish Balay 187529bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 187649b5e25fSSatish Balay 187749b5e25fSSatish Balay /* 187849b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 187949b5e25fSSatish Balay divisible by the blocksize 188049b5e25fSSatish Balay */ 188149b5e25fSSatish Balay mbs = M/bs; 188249b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 188349b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 188449b5e25fSSatish Balay else mbs++; 188549b5e25fSSatish Balay if (extra_rows) { 1886ae15b995SBarry Smith ierr = PetscInfo(0,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 188749b5e25fSSatish Balay } 188849b5e25fSSatish Balay 188949b5e25fSSatish Balay /* read in row lengths */ 189013f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 189149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 189249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 189349b5e25fSSatish Balay 189449b5e25fSSatish Balay /* read in column indices */ 189513f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 189649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 189749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 189849b5e25fSSatish Balay 189949b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 190013f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 190113f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 190213f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 190313f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 190449b5e25fSSatish Balay masked = mask + mbs; 190549b5e25fSSatish Balay rowcount = 0; nzcount = 0; 190649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 190749b5e25fSSatish Balay nmask = 0; 190849b5e25fSSatish Balay for (j=0; j<bs; j++) { 190949b5e25fSSatish Balay kmax = rowlengths[rowcount]; 191049b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19112d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 191203630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 191349b5e25fSSatish Balay } 191449b5e25fSSatish Balay rowcount++; 191549b5e25fSSatish Balay } 1916574b2666SHong Zhang s_browlengths[i] += nmask; 1917574b2666SHong Zhang 191849b5e25fSSatish Balay /* zero out the mask elements we set */ 191949b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 192049b5e25fSSatish Balay } 192149b5e25fSSatish Balay 192249b5e25fSSatish Balay /* create our matrix */ 1923f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 1924f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 19259abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 1926ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 192749b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 192849b5e25fSSatish Balay 192949b5e25fSSatish Balay /* set matrix "i" values */ 193049b5e25fSSatish Balay a->i[0] = 0; 193149b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1932574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1933574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 193449b5e25fSSatish Balay } 19356c6c5352SBarry Smith a->nz = a->i[mbs]; 193649b5e25fSSatish Balay 193749b5e25fSSatish Balay /* read in nonzero values */ 193887828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 193949b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 194049b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 194149b5e25fSSatish Balay 194249b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 194349b5e25fSSatish Balay nzcount = 0; jcount = 0; 194449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 194549b5e25fSSatish Balay nzcountb = nzcount; 194649b5e25fSSatish Balay nmask = 0; 194749b5e25fSSatish Balay for (j=0; j<bs; j++) { 194849b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 194949b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19502d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 195103630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19522d703238SHong Zhang } 19532d703238SHong Zhang } 19542d703238SHong Zhang /* sort the masked values */ 19552d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19562d703238SHong Zhang 19572d703238SHong Zhang /* set "j" values into matrix */ 19582d703238SHong Zhang maskcount = 1; 19592d703238SHong Zhang for (j=0; j<nmask; j++) { 196049b5e25fSSatish Balay a->j[jcount++] = masked[j]; 196149b5e25fSSatish Balay mask[masked[j]] = maskcount++; 196249b5e25fSSatish Balay } 1963574b2666SHong Zhang 196449b5e25fSSatish Balay /* set "a" values into matrix */ 196549b5e25fSSatish Balay ishift = bs2*a->i[i]; 196649b5e25fSSatish Balay for (j=0; j<bs; j++) { 196749b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 196849b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1969574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1970574b2666SHong Zhang if (tmp >= i){ 197149b5e25fSSatish Balay block = mask[tmp] - 1; 197249b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 197349b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1974574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1975574b2666SHong Zhang } 1976574b2666SHong Zhang nzcountb++; 197749b5e25fSSatish Balay } 197849b5e25fSSatish Balay } 197949b5e25fSSatish Balay /* zero out the mask elements we set */ 198049b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 198149b5e25fSSatish Balay } 19826c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 198349b5e25fSSatish Balay 198449b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1985574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 198649b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 198749b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 198849b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 198949b5e25fSSatish Balay 19909abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19919abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 199249b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19939abb65ffSKris Buschelman *A = B; 199449b5e25fSSatish Balay PetscFunctionReturn(0); 199549b5e25fSSatish Balay } 1996574b2666SHong Zhang 1997d06b337dSHong Zhang #undef __FUNCT__ 1998d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 199913f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 2000d06b337dSHong Zhang { 2001d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2002d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 2003d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 20046849ba73SBarry Smith PetscErrorCode ierr; 2005899cda47SBarry Smith PetscInt m=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j; 2006521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 2007d06b337dSHong Zhang 2008d06b337dSHong Zhang PetscFunctionBegin; 2009b965ef7fSBarry Smith its = its*lits; 201077431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 2011d06b337dSHong Zhang 2012d06b337dSHong Zhang if (bs > 1) 2013d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 2014d06b337dSHong Zhang 20151ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2016d06b337dSHong Zhang if (xx != bb) { 20171ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 2018d06b337dSHong Zhang } else { 2019d06b337dSHong Zhang b = x; 2020d06b337dSHong Zhang } 2021d06b337dSHong Zhang 2022e2ee2a47SBarry Smith if (!a->relax_work) { 2023e2ee2a47SBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->relax_work);CHKERRQ(ierr); 2024e2ee2a47SBarry Smith } 2025e2ee2a47SBarry Smith t = a->relax_work; 2026d06b337dSHong Zhang 2027d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 20282798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2029d06b337dSHong Zhang for (i=0; i<m; i++) 2030d06b337dSHong Zhang t[i] = b[i]; 2031d06b337dSHong Zhang 2032d06b337dSHong Zhang for (i=0; i<m; i++){ 203344706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2034d06b337dSHong Zhang v = aa + ai[i] + 1; 2035d06b337dSHong Zhang vj = aj + ai[i] + 1; 2036d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2037e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2038d06b337dSHong Zhang x[i] = omega*t[i]/d; 2039d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 2040d06b337dSHong Zhang } 2041d06b337dSHong Zhang } 2042d06b337dSHong Zhang 20432798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 204495d750ceSBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){ 204595d750ceSBarry Smith t = b; 2046d06b337dSHong Zhang } 204795d750ceSBarry Smith 2048d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2049d06b337dSHong Zhang d = *(aa + ai[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 sum = t[i]; 2055d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2056d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2057d06b337dSHong Zhang } 205895d750ceSBarry Smith t = a->relax_work; 2059d06b337dSHong Zhang } 2060d06b337dSHong Zhang its--; 2061d06b337dSHong Zhang } 2062d06b337dSHong Zhang 2063d06b337dSHong Zhang while (its--) { 206444706875SHong Zhang /* 206544706875SHong Zhang forward sweep: 206644706875SHong Zhang for i=0,...,m-1: 206744706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 206844706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 206944706875SHong Zhang b = b - x[i]*U^T(i,:); 2070d06b337dSHong Zhang 207144706875SHong Zhang */ 20722798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2073d06b337dSHong Zhang for (i=0; i<m; i++) 2074d06b337dSHong Zhang t[i] = b[i]; 2075d06b337dSHong Zhang 2076d06b337dSHong Zhang for (i=0; i<m; i++){ 207744706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2078d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2079d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2080d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2081d06b337dSHong Zhang sum = t[i]; 2082e6b907acSBarry Smith ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2083d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2084d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2085d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2086d06b337dSHong Zhang } 2087d06b337dSHong Zhang } 2088d06b337dSHong Zhang 20892798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 209044706875SHong Zhang /* 209144706875SHong Zhang backward sweep: 209244706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 209344706875SHong Zhang for i=m-1,...,0: 209444706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 209544706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 209644706875SHong Zhang */ 209795d750ceSBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2098d06b337dSHong Zhang for (i=0; i<m; i++) 2099d06b337dSHong Zhang t[i] = b[i]; 2100d06b337dSHong Zhang 2101d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2102d06b337dSHong Zhang v = aa + ai[i] + 1; 2103d06b337dSHong Zhang vj = aj + ai[i] + 1; 2104d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2105efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2106e6b907acSBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 2107d06b337dSHong Zhang } 2108d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2109d06b337dSHong Zhang d = *(aa + ai[i]); 2110d06b337dSHong Zhang v = aa + ai[i] + 1; 2111d06b337dSHong Zhang vj = aj + ai[i] + 1; 2112d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2113e6b907acSBarry Smith ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2114d06b337dSHong Zhang sum = t[i]; 2115d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2116d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2117d06b337dSHong Zhang } 2118d06b337dSHong Zhang } 2119d06b337dSHong Zhang } 2120d06b337dSHong Zhang 21211ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2122d06b337dSHong Zhang if (bb != xx) { 21231ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2124d06b337dSHong Zhang } 2125d06b337dSHong Zhang PetscFunctionReturn(0); 2126d06b337dSHong Zhang } 2127d06b337dSHong Zhang 2128c75a6043SHong Zhang #undef __FUNCT__ 2129c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2130c75a6043SHong Zhang /*@ 2131c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2132c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2133c75a6043SHong Zhang 2134c75a6043SHong Zhang Collective on MPI_Comm 2135c75a6043SHong Zhang 2136c75a6043SHong Zhang Input Parameters: 2137c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2138c75a6043SHong Zhang . bs - size of block 2139c75a6043SHong Zhang . m - number of rows 2140c75a6043SHong Zhang . n - number of columns 2141c75a6043SHong Zhang . i - row indices 2142c75a6043SHong Zhang . j - column indices 2143c75a6043SHong Zhang - a - matrix values 2144c75a6043SHong Zhang 2145c75a6043SHong Zhang Output Parameter: 2146c75a6043SHong Zhang . mat - the matrix 2147c75a6043SHong Zhang 2148c75a6043SHong Zhang Level: intermediate 2149c75a6043SHong Zhang 2150c75a6043SHong Zhang Notes: 2151c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2152c75a6043SHong Zhang once the matrix is destroyed 2153c75a6043SHong Zhang 2154c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2155c75a6043SHong Zhang 2156c75a6043SHong Zhang The i and j indices are 0 based 2157c75a6043SHong Zhang 2158c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2159c75a6043SHong Zhang 2160c75a6043SHong Zhang @*/ 2161c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2162c75a6043SHong Zhang { 2163c75a6043SHong Zhang PetscErrorCode ierr; 2164c75a6043SHong Zhang PetscInt ii; 2165c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2166c75a6043SHong Zhang 2167c75a6043SHong Zhang PetscFunctionBegin; 2168c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2169c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2170c75a6043SHong Zhang 2171c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2172c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2173c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2174c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2175c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2176c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2177c75a6043SHong Zhang 2178c75a6043SHong Zhang sbaij->i = i; 2179c75a6043SHong Zhang sbaij->j = j; 2180c75a6043SHong Zhang sbaij->a = a; 2181c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2182c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2183e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2184e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2185c75a6043SHong Zhang 2186c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2187c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2188c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2189c75a6043SHong 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]); 2190c75a6043SHong Zhang #endif 2191c75a6043SHong Zhang } 2192c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2193c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2194c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2195c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2196c75a6043SHong Zhang } 2197c75a6043SHong Zhang #endif 2198c75a6043SHong Zhang 2199c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2200c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2201c75a6043SHong Zhang PetscFunctionReturn(0); 2202c75a6043SHong Zhang } 2203d06b337dSHong Zhang 2204d06b337dSHong Zhang 2205d06b337dSHong Zhang 220649b5e25fSSatish Balay 220749b5e25fSSatish Balay 2208