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; 3913f74950SBarry Smith PetscInt i,mbs = a->mbs; 4049b5e25fSSatish Balay 4149b5e25fSSatish Balay PetscFunctionBegin; 4249b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4349b5e25fSSatish Balay 4413f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 4552e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 46b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4749b5e25fSSatish Balay PetscFunctionReturn(0); 4849b5e25fSSatish Balay } 4949b5e25fSSatish Balay 504a2ae208SSatish Balay #undef __FUNCT__ 514a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 5213f74950SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 5349b5e25fSSatish Balay { 54a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5513f74950SBarry Smith PetscInt n = a->mbs,i; 5649b5e25fSSatish Balay 5749b5e25fSSatish Balay PetscFunctionBegin; 58d3e5a4abSHong Zhang *nn = n; 59a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 60a1373b80SHong Zhang 61a6ece127SHong Zhang if (oshift == 1) { 62a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 6313f74950SBarry Smith PetscInt nz = a->i[n]; 646c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 65a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 66a1373b80SHong Zhang *ia = a->i; *ja = a->j; 67a1373b80SHong Zhang } else { 68a1373b80SHong Zhang *ia = a->i; *ja = a->j; 69a6ece127SHong Zhang } 7049b5e25fSSatish Balay PetscFunctionReturn(0); 7149b5e25fSSatish Balay } 7249b5e25fSSatish Balay 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 7513f74950SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 7649b5e25fSSatish Balay { 77b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7813f74950SBarry Smith PetscInt i,n = a->mbs; 79a6ece127SHong Zhang 8049b5e25fSSatish Balay PetscFunctionBegin; 8149b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 82a6ece127SHong Zhang 83a6ece127SHong Zhang if (oshift == 1) { 8413f74950SBarry Smith PetscInt nz = a->i[n]-1; 856c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 86a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 87a6ece127SHong Zhang } 88a6ece127SHong Zhang PetscFunctionReturn(0); 8949b5e25fSSatish Balay } 9049b5e25fSSatish Balay 914a2ae208SSatish Balay #undef __FUNCT__ 924a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 93dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 9449b5e25fSSatish Balay { 9549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 96dfbe8321SBarry Smith PetscErrorCode ierr; 9749b5e25fSSatish Balay 9849b5e25fSSatish Balay PetscFunctionBegin; 99a9f03627SSatish Balay #if defined(PETSC_USE_LOG) 10077431f27SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->m,a->nz); 101a9f03627SSatish Balay #endif 102085a36d4SBarry Smith ierr = MatSeqXAIJFreeAIJ(a->singlemalloc,&a->a,&a->j,&a->i);CHKERRQ(ierr); 10349b5e25fSSatish Balay if (a->row) { 10449b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 10549b5e25fSSatish Balay } 10649b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 107ab93d7beSBarry Smith if (a->imax) {ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);} 10849b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 109d59c15a7SBarry Smith if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 110d59c15a7SBarry Smith if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);} 111d59c15a7SBarry Smith if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 11249b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 113c4319e64SHong Zhang if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 1141a3463dfSHong Zhang 1151a3463dfSHong Zhang if (a->inew){ 1161a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1171a3463dfSHong Zhang a->inew = 0; 1181a3463dfSHong Zhang } 11949b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 120901853e0SKris Buschelman 121901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 122901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 123901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 124901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 125901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 126901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 12749b5e25fSSatish Balay PetscFunctionReturn(0); 12849b5e25fSSatish Balay } 12949b5e25fSSatish Balay 1304a2ae208SSatish Balay #undef __FUNCT__ 1314a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 132dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op) 13349b5e25fSSatish Balay { 134045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 13563ba0a88SBarry Smith PetscErrorCode ierr; 13649b5e25fSSatish Balay 13749b5e25fSSatish Balay PetscFunctionBegin; 1384d9d31abSKris Buschelman switch (op) { 1394d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1404d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1414d9d31abSKris Buschelman break; 1424d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1434d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1444d9d31abSKris Buschelman break; 1454d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1464d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1474d9d31abSKris Buschelman break; 1484d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1494d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1504d9d31abSKris Buschelman break; 1514d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1524d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1554d9d31abSKris Buschelman a->nonew = 1; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1584d9d31abSKris Buschelman a->nonew = -1; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1614d9d31abSKris Buschelman a->nonew = -2; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1644d9d31abSKris Buschelman a->nonew = 0; 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1674d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1684d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1694d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1704d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 171*ae15b995SBarry Smith ierr = PetscInfo(A,"Option ignored\n");CHKERRQ(ierr); 1724d9d31abSKris Buschelman break; 1734d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 17429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1759a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1769a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1779a4540c5SBarry Smith case MAT_HERMITIAN: 1789a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 17977e54ba9SKris Buschelman case MAT_SYMMETRIC: 18077e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1819a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1829a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1839a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 184941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 185941593c8SHong Zhang a->ignore_ltriangular = PETSC_TRUE; 186941593c8SHong Zhang break; 187941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 188941593c8SHong Zhang a->ignore_ltriangular = PETSC_FALSE; 18977e54ba9SKris Buschelman break; 190f5edf698SHong Zhang case MAT_GETROW_UPPERTRIANGULAR: 191f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 192f5edf698SHong Zhang break; 1934d9d31abSKris Buschelman default: 19429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 19549b5e25fSSatish Balay } 19649b5e25fSSatish Balay PetscFunctionReturn(0); 19749b5e25fSSatish Balay } 19849b5e25fSSatish Balay 1994a2ae208SSatish Balay #undef __FUNCT__ 2004a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 20113f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 20249b5e25fSSatish Balay { 20349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2046849ba73SBarry Smith PetscErrorCode ierr; 20513f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 20649b5e25fSSatish Balay MatScalar *aa,*aa_i; 20787828ca2SBarry Smith PetscScalar *v_i; 20849b5e25fSSatish Balay 20949b5e25fSSatish Balay PetscFunctionBegin; 210f5edf698SHong 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()"); 211f5edf698SHong Zhang /* Get the upper triangular part of the row */ 212521d7252SBarry Smith bs = A->bs; 21349b5e25fSSatish Balay ai = a->i; 21449b5e25fSSatish Balay aj = a->j; 21549b5e25fSSatish Balay aa = a->a; 21649b5e25fSSatish Balay bs2 = a->bs2; 21749b5e25fSSatish Balay 21877431f27SBarry Smith if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 21949b5e25fSSatish Balay 22049b5e25fSSatish Balay bn = row/bs; /* Block number */ 22149b5e25fSSatish Balay bp = row % bs; /* Block position */ 22249b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 22349b5e25fSSatish Balay *ncols = bs*M; 22449b5e25fSSatish Balay 22549b5e25fSSatish Balay if (v) { 22649b5e25fSSatish Balay *v = 0; 22749b5e25fSSatish Balay if (*ncols) { 22887828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 22949b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23049b5e25fSSatish Balay v_i = *v + i*bs; 23149b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 23249b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 23349b5e25fSSatish Balay } 23449b5e25fSSatish Balay } 23549b5e25fSSatish Balay } 23649b5e25fSSatish Balay 23749b5e25fSSatish Balay if (cols) { 23849b5e25fSSatish Balay *cols = 0; 23949b5e25fSSatish Balay if (*ncols) { 24013f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 24149b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 24249b5e25fSSatish Balay cols_i = *cols + i*bs; 24349b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 24449b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24549b5e25fSSatish Balay } 24649b5e25fSSatish Balay } 24749b5e25fSSatish Balay } 24849b5e25fSSatish Balay 24949b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2505ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2515ddb2528SHong Zhang #ifdef column_search 25249b5e25fSSatish Balay v_i = *v + M*bs; 25349b5e25fSSatish Balay cols_i = *cols + M*bs; 25449b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25549b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25649b5e25fSSatish Balay for (j=0; j<M; j++){ 25749b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 25849b5e25fSSatish Balay if (itmp == bn){ 25949b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 26049b5e25fSSatish Balay for (k=0; k<bs; k++) { 26149b5e25fSSatish Balay *cols_i++ = i*bs+k; 26249b5e25fSSatish Balay *v_i++ = aa_i[k]; 26349b5e25fSSatish Balay } 26449b5e25fSSatish Balay *ncols += bs; 26549b5e25fSSatish Balay break; 26649b5e25fSSatish Balay } 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay } 2695ddb2528SHong Zhang #endif 27049b5e25fSSatish Balay PetscFunctionReturn(0); 27149b5e25fSSatish Balay } 27249b5e25fSSatish Balay 2734a2ae208SSatish Balay #undef __FUNCT__ 2744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 27513f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27649b5e25fSSatish Balay { 277dfbe8321SBarry Smith PetscErrorCode ierr; 27849b5e25fSSatish Balay 27949b5e25fSSatish Balay PetscFunctionBegin; 28049b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 28149b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 28249b5e25fSSatish Balay PetscFunctionReturn(0); 28349b5e25fSSatish Balay } 28449b5e25fSSatish Balay 2854a2ae208SSatish Balay #undef __FUNCT__ 286f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ" 287f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 288f5edf698SHong Zhang { 289f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 290f5edf698SHong Zhang 291f5edf698SHong Zhang PetscFunctionBegin; 292f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 293f5edf698SHong Zhang PetscFunctionReturn(0); 294f5edf698SHong Zhang } 295f5edf698SHong Zhang #undef __FUNCT__ 296f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ" 297f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 298f5edf698SHong Zhang { 299f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 300f5edf698SHong Zhang 301f5edf698SHong Zhang PetscFunctionBegin; 302f5edf698SHong Zhang a->getrow_utriangular = PETSC_FALSE; 303f5edf698SHong Zhang PetscFunctionReturn(0); 304f5edf698SHong Zhang } 305f5edf698SHong Zhang 306f5edf698SHong Zhang #undef __FUNCT__ 3074a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 308dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 30949b5e25fSSatish Balay { 310dfbe8321SBarry Smith PetscErrorCode ierr; 31149b5e25fSSatish Balay PetscFunctionBegin; 312999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 3138115998fSBarry Smith PetscFunctionReturn(0); 31449b5e25fSSatish Balay } 31549b5e25fSSatish Balay 3164a2ae208SSatish Balay #undef __FUNCT__ 3174a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3186849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 31949b5e25fSSatish Balay { 32049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 321dfbe8321SBarry Smith PetscErrorCode ierr; 322521d7252SBarry Smith PetscInt i,j,bs = A->bs,k,l,bs2=a->bs2; 3232dcb1b2aSMatthew Knepley const char *name; 324f3ef73ceSBarry Smith PetscViewerFormat format; 32549b5e25fSSatish Balay 32649b5e25fSSatish Balay PetscFunctionBegin; 32780fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 328b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 329456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 33077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 331fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 33229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 333fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 334b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 33549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33649b5e25fSSatish Balay for (j=0; j<bs; j++) { 33777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 33849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33949b5e25fSSatish Balay for (l=0; l<bs; l++) { 34049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 34149b5e25fSSatish Balay 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 (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 345a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 34649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34749b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 348a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34949b5e25fSSatish Balay } 35049b5e25fSSatish Balay #else 35149b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 352a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35349b5e25fSSatish Balay } 35449b5e25fSSatish Balay #endif 35549b5e25fSSatish Balay } 35649b5e25fSSatish Balay } 357b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 35849b5e25fSSatish Balay } 35949b5e25fSSatish Balay } 360b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 361c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 362c1490034SHong Zhang PetscFunctionReturn(0); 36349b5e25fSSatish Balay } else { 364b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 36649b5e25fSSatish Balay for (j=0; j<bs; j++) { 36777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 36849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 36949b5e25fSSatish Balay for (l=0; l<bs; l++) { 37049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37149b5e25fSSatish Balay 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 if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 375a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 37649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37749b5e25fSSatish Balay } else { 378a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37949b5e25fSSatish Balay } 38049b5e25fSSatish Balay #else 381a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %G ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38249b5e25fSSatish Balay #endif 38349b5e25fSSatish Balay } 38449b5e25fSSatish Balay } 385b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 38649b5e25fSSatish Balay } 38749b5e25fSSatish Balay } 388b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 38949b5e25fSSatish Balay } 390b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 39149b5e25fSSatish Balay PetscFunctionReturn(0); 39249b5e25fSSatish Balay } 39349b5e25fSSatish Balay 3944a2ae208SSatish Balay #undef __FUNCT__ 3954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 3966849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 39749b5e25fSSatish Balay { 39849b5e25fSSatish Balay Mat A = (Mat) Aa; 39949b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 4006849ba73SBarry Smith PetscErrorCode ierr; 40113f74950SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2; 40213f74950SBarry Smith PetscMPIInt rank; 40349b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 40449b5e25fSSatish Balay MatScalar *aa; 40549b5e25fSSatish Balay MPI_Comm comm; 406b0a32e0cSBarry Smith PetscViewer viewer; 40749b5e25fSSatish Balay 40849b5e25fSSatish Balay PetscFunctionBegin; 40949b5e25fSSatish Balay /* 41049b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 41149b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 41249b5e25fSSatish Balay rest should return immediately. 41349b5e25fSSatish Balay */ 41449b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 41549b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 41649b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 41749b5e25fSSatish Balay 41849b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 41949b5e25fSSatish Balay 420b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 421b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 42249b5e25fSSatish Balay 42349b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 424b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 42549b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 42649b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 427c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 42849b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42949b5e25fSSatish Balay aa = a->a + j*bs2; 43049b5e25fSSatish Balay for (k=0; k<bs; k++) { 43149b5e25fSSatish Balay for (l=0; l<bs; l++) { 43249b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 433b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 43449b5e25fSSatish Balay } 43549b5e25fSSatish Balay } 43649b5e25fSSatish Balay } 43749b5e25fSSatish Balay } 438b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 43949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 44049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 441c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 44249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 44349b5e25fSSatish Balay aa = a->a + j*bs2; 44449b5e25fSSatish Balay for (k=0; k<bs; k++) { 44549b5e25fSSatish Balay for (l=0; l<bs; l++) { 44649b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 447b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44849b5e25fSSatish Balay } 44949b5e25fSSatish Balay } 45049b5e25fSSatish Balay } 45149b5e25fSSatish Balay } 45249b5e25fSSatish Balay 453b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 45449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 456c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45849b5e25fSSatish Balay aa = a->a + j*bs2; 45949b5e25fSSatish Balay for (k=0; k<bs; k++) { 46049b5e25fSSatish Balay for (l=0; l<bs; l++) { 46149b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 462b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46349b5e25fSSatish Balay } 46449b5e25fSSatish Balay } 46549b5e25fSSatish Balay } 46649b5e25fSSatish Balay } 46749b5e25fSSatish Balay PetscFunctionReturn(0); 46849b5e25fSSatish Balay } 46949b5e25fSSatish Balay 4704a2ae208SSatish Balay #undef __FUNCT__ 4714a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4726849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 47349b5e25fSSatish Balay { 474dfbe8321SBarry Smith PetscErrorCode ierr; 47549b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 476b0a32e0cSBarry Smith PetscDraw draw; 47749b5e25fSSatish Balay PetscTruth isnull; 47849b5e25fSSatish Balay 47949b5e25fSSatish Balay PetscFunctionBegin; 48049b5e25fSSatish Balay 481b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 482b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 48349b5e25fSSatish Balay 48449b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 485c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 48649b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 487b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 488b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 48949b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 49049b5e25fSSatish Balay PetscFunctionReturn(0); 49149b5e25fSSatish Balay } 49249b5e25fSSatish Balay 4934a2ae208SSatish Balay #undef __FUNCT__ 4944a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 495dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 49649b5e25fSSatish Balay { 497dfbe8321SBarry Smith PetscErrorCode ierr; 49832077d6dSBarry Smith PetscTruth iascii,isdraw; 49949b5e25fSSatish Balay 50049b5e25fSSatish Balay PetscFunctionBegin; 50132077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 502fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 50332077d6dSBarry Smith if (iascii){ 50449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 50549b5e25fSSatish Balay } else if (isdraw) { 50649b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 50749b5e25fSSatish Balay } else { 508a5e6ed63SBarry Smith Mat B; 509ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 510a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 511a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 51249b5e25fSSatish Balay } 51349b5e25fSSatish Balay PetscFunctionReturn(0); 51449b5e25fSSatish Balay } 51549b5e25fSSatish Balay 51649b5e25fSSatish Balay 5174a2ae208SSatish Balay #undef __FUNCT__ 5184a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 51913f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 52049b5e25fSSatish Balay { 521045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 52213f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 52313f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 52413f74950SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2; 52549b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 52649b5e25fSSatish Balay 52749b5e25fSSatish Balay PetscFunctionBegin; 52849b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 52949b5e25fSSatish Balay row = im[k]; brow = row/bs; 53077431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 53177431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 53249b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 53349b5e25fSSatish Balay nrow = ailen[brow]; 53449b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 53577431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 53677431f27SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1); 53749b5e25fSSatish Balay col = in[l] ; 53849b5e25fSSatish Balay bcol = col/bs; 53949b5e25fSSatish Balay cidx = col%bs; 54049b5e25fSSatish Balay ridx = row%bs; 54149b5e25fSSatish Balay high = nrow; 54249b5e25fSSatish Balay low = 0; /* assume unsorted */ 54349b5e25fSSatish Balay while (high-low > 5) { 54449b5e25fSSatish Balay t = (low+high)/2; 54549b5e25fSSatish Balay if (rp[t] > bcol) high = t; 54649b5e25fSSatish Balay else low = t; 54749b5e25fSSatish Balay } 54849b5e25fSSatish Balay for (i=low; i<high; i++) { 54949b5e25fSSatish Balay if (rp[i] > bcol) break; 55049b5e25fSSatish Balay if (rp[i] == bcol) { 55149b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 55249b5e25fSSatish Balay goto finished; 55349b5e25fSSatish Balay } 55449b5e25fSSatish Balay } 55549b5e25fSSatish Balay *v++ = zero; 55649b5e25fSSatish Balay finished:; 55749b5e25fSSatish Balay } 55849b5e25fSSatish Balay } 55949b5e25fSSatish Balay PetscFunctionReturn(0); 56049b5e25fSSatish Balay } 56149b5e25fSSatish Balay 56249b5e25fSSatish Balay 5634a2ae208SSatish Balay #undef __FUNCT__ 5644a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 56513f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 56649b5e25fSSatish Balay { 5670880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5686849ba73SBarry Smith PetscErrorCode ierr; 569e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 57013f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 57113f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval; 5720880e062SHong Zhang PetscTruth roworiented=a->roworiented; 573f15d580aSBarry Smith const MatScalar *value = v; 574f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5750880e062SHong Zhang 57649b5e25fSSatish Balay PetscFunctionBegin; 5770880e062SHong Zhang if (roworiented) { 5780880e062SHong Zhang stepval = (n-1)*bs; 5790880e062SHong Zhang } else { 5800880e062SHong Zhang stepval = (m-1)*bs; 5810880e062SHong Zhang } 5820880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5830880e062SHong Zhang row = im[k]; 5840880e062SHong Zhang if (row < 0) continue; 5852515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 58677431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5870880e062SHong Zhang #endif 5880880e062SHong Zhang rp = aj + ai[row]; 5890880e062SHong Zhang ap = aa + bs2*ai[row]; 5900880e062SHong Zhang rmax = imax[row]; 5910880e062SHong Zhang nrow = ailen[row]; 5920880e062SHong Zhang low = 0; 593818f2c47SBarry Smith high = nrow; 5940880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5950880e062SHong Zhang if (in[l] < 0) continue; 5960880e062SHong Zhang col = in[l]; 5972515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 59877431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 599b1823623SSatish Balay #endif 6000880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6010880e062SHong Zhang if (roworiented) { 6020880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6030880e062SHong Zhang } else { 6040880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6050880e062SHong Zhang } 6067cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 607e2ee6c50SBarry Smith lastcol = col; 6080880e062SHong Zhang while (high-low > 7) { 6090880e062SHong Zhang t = (low+high)/2; 6100880e062SHong Zhang if (rp[t] > col) high = t; 6110880e062SHong Zhang else low = t; 6120880e062SHong Zhang } 6130880e062SHong Zhang for (i=low; i<high; i++) { 6140880e062SHong Zhang if (rp[i] > col) break; 6150880e062SHong Zhang if (rp[i] == col) { 6160880e062SHong Zhang bap = ap + bs2*i; 6170880e062SHong Zhang if (roworiented) { 6180880e062SHong Zhang if (is == ADD_VALUES) { 6190880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6200880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6210880e062SHong Zhang bap[jj] += *value++; 6220880e062SHong Zhang } 6230880e062SHong Zhang } 6240880e062SHong Zhang } else { 6250880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6260880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6270880e062SHong Zhang bap[jj] = *value++; 6280880e062SHong Zhang } 6290880e062SHong Zhang } 6300880e062SHong Zhang } 6310880e062SHong Zhang } else { 6320880e062SHong Zhang if (is == ADD_VALUES) { 6330880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6340880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6350880e062SHong Zhang *bap++ += *value++; 6360880e062SHong Zhang } 6370880e062SHong Zhang } 6380880e062SHong Zhang } else { 6390880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6400880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6410880e062SHong Zhang *bap++ = *value++; 6420880e062SHong Zhang } 6430880e062SHong Zhang } 6440880e062SHong Zhang } 6450880e062SHong Zhang } 6460880e062SHong Zhang goto noinsert2; 6470880e062SHong Zhang } 6480880e062SHong Zhang } 6490880e062SHong Zhang if (nonew == 1) goto noinsert2; 650085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 651ed1caa07SMatthew Knepley MatSeqXAIJReallocateAIJ(a,bs2,nrow,row,col,rmax,aa,ai,aj,a->mbs,rp,ap,imax,nonew); 652c03d1d03SSatish Balay N = nrow++ - 1; high++; 6530880e062SHong Zhang /* shift up all the later entries in this row */ 6540880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6550880e062SHong Zhang rp[ii+1] = rp[ii]; 6560880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6570880e062SHong Zhang } 6580880e062SHong Zhang if (N >= i) { 6590880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6600880e062SHong Zhang } 6610880e062SHong Zhang rp[i] = col; 6620880e062SHong Zhang bap = ap + bs2*i; 6630880e062SHong Zhang if (roworiented) { 6640880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6650880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6660880e062SHong Zhang bap[jj] = *value++; 6670880e062SHong Zhang } 6680880e062SHong Zhang } 6690880e062SHong Zhang } else { 6700880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6710880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6720880e062SHong Zhang *bap++ = *value++; 6730880e062SHong Zhang } 6740880e062SHong Zhang } 6750880e062SHong Zhang } 6760880e062SHong Zhang noinsert2:; 6770880e062SHong Zhang low = i; 6780880e062SHong Zhang } 6790880e062SHong Zhang ailen[row] = nrow; 6800880e062SHong Zhang } 6810880e062SHong Zhang PetscFunctionReturn(0); 68249b5e25fSSatish Balay } 68349b5e25fSSatish Balay 6844a2ae208SSatish Balay #undef __FUNCT__ 6854a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 686dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 68749b5e25fSSatish Balay { 68849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6896849ba73SBarry Smith PetscErrorCode ierr; 69013f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 69113f74950SBarry Smith PetscInt m = A->m,*ip,N,*ailen = a->ilen; 69213f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 69349b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 69449b5e25fSSatish Balay 69549b5e25fSSatish Balay PetscFunctionBegin; 69649b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 69749b5e25fSSatish Balay 69849b5e25fSSatish Balay if (m) rmax = ailen[0]; 69949b5e25fSSatish Balay for (i=1; i<mbs; i++) { 70049b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 70149b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 70249b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 70349b5e25fSSatish Balay if (fshift) { 70449b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 70549b5e25fSSatish Balay N = ailen[i]; 70649b5e25fSSatish Balay for (j=0; j<N; j++) { 70749b5e25fSSatish Balay ip[j-fshift] = ip[j]; 70849b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 70949b5e25fSSatish Balay } 71049b5e25fSSatish Balay } 71149b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 71249b5e25fSSatish Balay } 71349b5e25fSSatish Balay if (mbs) { 71449b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 71549b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 71649b5e25fSSatish Balay } 71749b5e25fSSatish Balay /* reset ilen and imax for each row */ 71849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 71949b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 72049b5e25fSSatish Balay } 7216c6c5352SBarry Smith a->nz = ai[mbs]; 72249b5e25fSSatish Balay 723b424e231SHong Zhang /* diagonals may have moved, reset it */ 724b424e231SHong Zhang if (a->diag) { 72513f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 72649b5e25fSSatish Balay } 727*ae15b995SBarry Smith ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->m,A->bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr); 728*ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 729*ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 73049b5e25fSSatish Balay a->reallocs = 0; 73149b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 73249b5e25fSSatish Balay PetscFunctionReturn(0); 73349b5e25fSSatish Balay } 73449b5e25fSSatish Balay 73549b5e25fSSatish Balay /* 73649b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 73749b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 73849b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 73949b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 74049b5e25fSSatish Balay */ 7414a2ae208SSatish Balay #undef __FUNCT__ 7424a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 74313f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 74449b5e25fSSatish Balay { 74513f74950SBarry Smith PetscInt i,j,k,row; 74649b5e25fSSatish Balay PetscTruth flg; 74749b5e25fSSatish Balay 74849b5e25fSSatish Balay PetscFunctionBegin; 74949b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 75049b5e25fSSatish Balay row = idx[i]; 75149b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 75249b5e25fSSatish Balay sizes[j] = 1; 75349b5e25fSSatish Balay i++; 75449b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 75549b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 75649b5e25fSSatish Balay i++; 75749b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 75849b5e25fSSatish Balay flg = PETSC_TRUE; 75949b5e25fSSatish Balay for (k=1; k<bs; k++) { 76049b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 76149b5e25fSSatish Balay flg = PETSC_FALSE; 76249b5e25fSSatish Balay break; 76349b5e25fSSatish Balay } 76449b5e25fSSatish Balay } 765abc0a331SBarry Smith if (flg) { /* No break in the bs */ 76649b5e25fSSatish Balay sizes[j] = bs; 76749b5e25fSSatish Balay i+= bs; 76849b5e25fSSatish Balay } else { 76949b5e25fSSatish Balay sizes[j] = 1; 77049b5e25fSSatish Balay i++; 77149b5e25fSSatish Balay } 77249b5e25fSSatish Balay } 77349b5e25fSSatish Balay } 77449b5e25fSSatish Balay *bs_max = j; 77549b5e25fSSatish Balay PetscFunctionReturn(0); 77649b5e25fSSatish Balay } 77749b5e25fSSatish Balay 77849b5e25fSSatish Balay 77949b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 78049b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 78149b5e25fSSatish Balay */ 78249b5e25fSSatish Balay 7834a2ae208SSatish Balay #undef __FUNCT__ 7844a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 78513f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 78649b5e25fSSatish Balay { 78749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7886849ba73SBarry Smith PetscErrorCode ierr; 789e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 79013f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 79113f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->bs,brow,bcol; 79213f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 79349b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 79449b5e25fSSatish Balay 79549b5e25fSSatish Balay PetscFunctionBegin; 79649b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 79749b5e25fSSatish Balay row = im[k]; /* row number */ 79849b5e25fSSatish Balay brow = row/bs; /* block row number */ 79949b5e25fSSatish Balay if (row < 0) continue; 8002515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 80177431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 80249b5e25fSSatish Balay #endif 80349b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 80449b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 80549b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 80649b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 80749b5e25fSSatish Balay low = 0; 80849b5e25fSSatish Balay 80949b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 81049b5e25fSSatish Balay if (in[l] < 0) continue; 8112515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 81277431f27SBarry Smith if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1); 81349b5e25fSSatish Balay #endif 81449b5e25fSSatish Balay col = in[l]; 81549b5e25fSSatish Balay bcol = col/bs; /* block col number */ 81649b5e25fSSatish Balay 817941593c8SHong Zhang if (brow > bcol) { 818941593c8SHong Zhang if (a->ignore_ltriangular){ 819941593c8SHong Zhang continue; /* ignore lower triangular values */ 820941593c8SHong Zhang } else { 821941593c8SHong 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)"); 822941593c8SHong Zhang } 823941593c8SHong Zhang } 824f4989cb3SHong Zhang 82549b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8268549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 82749b5e25fSSatish Balay /* element value a(k,l) */ 82849b5e25fSSatish Balay if (roworiented) { 82949b5e25fSSatish Balay value = v[l + k*n]; 83049b5e25fSSatish Balay } else { 83149b5e25fSSatish Balay value = v[k + l*m]; 83249b5e25fSSatish Balay } 83349b5e25fSSatish Balay 83449b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 8357cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 836e2ee6c50SBarry Smith lastcol = col; 83749b5e25fSSatish Balay while (high-low > 7) { 83849b5e25fSSatish Balay t = (low+high)/2; 83949b5e25fSSatish Balay if (rp[t] > bcol) high = t; 84049b5e25fSSatish Balay else low = t; 84149b5e25fSSatish Balay } 84249b5e25fSSatish Balay for (i=low; i<high; i++) { 84349b5e25fSSatish Balay if (rp[i] > bcol) break; 84449b5e25fSSatish Balay if (rp[i] == bcol) { 84549b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 84649b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 84749b5e25fSSatish Balay else *bap = value; 8488549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8498549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8508549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8518549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8528549e402SHong Zhang else *bap = value; 8538549e402SHong Zhang } 85449b5e25fSSatish Balay goto noinsert1; 85549b5e25fSSatish Balay } 85649b5e25fSSatish Balay } 85749b5e25fSSatish Balay 85849b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 859085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 860ed1caa07SMatthew Knepley MatSeqXAIJReallocateAIJ(a,bs2,nrow,brow,bcol,rmax,aa,ai,aj,a->mbs,rp,ap,imax,nonew); 86149b5e25fSSatish Balay 862c03d1d03SSatish Balay N = nrow++ - 1; high++; 86349b5e25fSSatish Balay /* shift up all the later entries in this row */ 86449b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 86549b5e25fSSatish Balay rp[ii+1] = rp[ii]; 86649b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 86749b5e25fSSatish Balay } 86849b5e25fSSatish Balay if (N>=i) { 86949b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 87049b5e25fSSatish Balay } 87149b5e25fSSatish Balay rp[i] = bcol; 87249b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 87349b5e25fSSatish Balay noinsert1:; 87449b5e25fSSatish Balay low = i; 8758549e402SHong Zhang } 87649b5e25fSSatish Balay } /* end of loop over added columns */ 87749b5e25fSSatish Balay ailen[brow] = nrow; 87849b5e25fSSatish Balay } /* end of loop over added rows */ 87949b5e25fSSatish Balay PetscFunctionReturn(0); 88049b5e25fSSatish Balay } 88149b5e25fSSatish Balay 8826849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 8836849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 88413f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt); 88513f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 88613f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]); 887f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar); 8886849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 8896849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 8906849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec); 8916849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 8926849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 8936849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat); 8946849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec); 89549b5e25fSSatish Balay 8966849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 8976849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 8986849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 8996849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 9006849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 9016849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 9026849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 9036849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 90449b5e25fSSatish Balay 9056849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 906d59c15a7SBarry Smith 9076849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 9086849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 9096849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 9106849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 9116849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 9126849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 9136849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 9146849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 91507f98182SSatish Balay 916af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,MatFactorInfo*,Mat*); 917af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,MatFactorInfo*,Mat*); 918af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,MatFactorInfo*,Mat*); 919af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,MatFactorInfo*,Mat*); 920af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,MatFactorInfo*,Mat*); 921af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,MatFactorInfo*,Mat*); 922af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,MatFactorInfo*,Mat*); 923af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,MatFactorInfo*,Mat*); 92413f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*); 92549b5e25fSSatish Balay 9266849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 9276849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 9286849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 9296849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 9306849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 9316849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 9326849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 9336849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 93449b5e25fSSatish Balay 9356849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 9366849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 9376849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 9386849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 9396849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 9406849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 9416849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 9426849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 94349b5e25fSSatish Balay 9444d101231SSatish Balay #ifdef HAVE_MatICCFactor 9456849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work! */ 9464a2ae208SSatish Balay #undef __FUNCT__ 9474d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 948dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 94949b5e25fSSatish Balay { 9504ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 95149b5e25fSSatish Balay Mat outA; 952dfbe8321SBarry Smith PetscErrorCode ierr; 95349b5e25fSSatish Balay PetscTruth row_identity,col_identity; 95449b5e25fSSatish Balay 95549b5e25fSSatish Balay PetscFunctionBegin; 95649b5e25fSSatish Balay outA = inA; 9571260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 95849b5e25fSSatish Balay 95949b5e25fSSatish Balay if (!a->diag) { 9601a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 96149b5e25fSSatish Balay } 96249b5e25fSSatish Balay /* 96349b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 96449b5e25fSSatish Balay for ILU(0) factorization with natural ordering 96549b5e25fSSatish Balay */ 96649b5e25fSSatish Balay switch (a->bs) { 96749b5e25fSSatish Balay case 1: 9680fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 96907f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 970d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 971*ae15b995SBarry Smith ierr = PetscInfo((inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr); 97249b5e25fSSatish Balay case 2: 9731a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9741a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9751a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 976*ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 97749b5e25fSSatish Balay break; 97849b5e25fSSatish Balay case 3: 9791a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9801a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9811a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 982*ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 98349b5e25fSSatish Balay break; 98449b5e25fSSatish Balay case 4: 9851a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 9861a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 9871a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 988*ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 98949b5e25fSSatish Balay break; 99049b5e25fSSatish Balay case 5: 9911a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 9921a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 9931a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 994*ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 99549b5e25fSSatish Balay break; 99649b5e25fSSatish Balay case 6: 9971a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 9981a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 9991a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 1000*ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 100149b5e25fSSatish Balay break; 100249b5e25fSSatish Balay case 7: 10031a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10041a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10051a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 1006*ae15b995SBarry Smith ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 100749b5e25fSSatish Balay break; 100849b5e25fSSatish Balay default: 100949b5e25fSSatish Balay a->row = row; 10101a3463dfSHong Zhang a->icol = col; 101149b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 101249b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 101349b5e25fSSatish Balay 101449b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 101549b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 101652e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 101749b5e25fSSatish Balay 101849b5e25fSSatish Balay if (!a->solve_work) { 101987828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 102052e6d16bSBarry Smith ierr = PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));CHKERRQ(ierr); 102149b5e25fSSatish Balay } 102249b5e25fSSatish Balay } 102349b5e25fSSatish Balay 1024af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 102549b5e25fSSatish Balay PetscFunctionReturn(0); 102649b5e25fSSatish Balay } 1027950f1e5bSHong Zhang #endif 1028950f1e5bSHong Zhang 10294a2ae208SSatish Balay #undef __FUNCT__ 10304a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 1031dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A) 103249b5e25fSSatish Balay { 103349b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 103449b5e25fSSatish Balay MPI_Comm comm = A->comm; 1035dfbe8321SBarry Smith PetscErrorCode ierr; 103649b5e25fSSatish Balay 103749b5e25fSSatish Balay PetscFunctionBegin; 103849b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 103949b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 104049b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 1041f5edf698SHong Zhang ierr = (*PetscHelpPrintf)(comm," -mat_ignore_lower_triangular: Ignore lower triangular values set by user\n");CHKERRQ(ierr); 1042f5edf698SHong Zhang ierr = (*PetscHelpPrintf)(comm," -mat_getrow_uppertriangular: Enable MatGetRow() for retrieving the upper triangular part of the row\n");CHKERRQ(ierr); 104349b5e25fSSatish Balay PetscFunctionReturn(0); 104449b5e25fSSatish Balay } 104549b5e25fSSatish Balay 104649b5e25fSSatish Balay EXTERN_C_BEGIN 10474a2ae208SSatish Balay #undef __FUNCT__ 10484a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1049be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 105049b5e25fSSatish Balay { 1051045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 105213f74950SBarry Smith PetscInt i,nz,n; 105349b5e25fSSatish Balay 105449b5e25fSSatish Balay PetscFunctionBegin; 10556c6c5352SBarry Smith nz = baij->maxnz; 1056c464158bSHong Zhang n = mat->n; 105749b5e25fSSatish Balay for (i=0; i<nz; i++) { 105849b5e25fSSatish Balay baij->j[i] = indices[i]; 105949b5e25fSSatish Balay } 10606c6c5352SBarry Smith baij->nz = nz; 106149b5e25fSSatish Balay for (i=0; i<n; i++) { 106249b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 106349b5e25fSSatish Balay } 106449b5e25fSSatish Balay PetscFunctionReturn(0); 106549b5e25fSSatish Balay } 106649b5e25fSSatish Balay EXTERN_C_END 106749b5e25fSSatish Balay 10684a2ae208SSatish Balay #undef __FUNCT__ 10694a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 107049b5e25fSSatish Balay /*@ 107119585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 107249b5e25fSSatish Balay in the matrix. 107349b5e25fSSatish Balay 107449b5e25fSSatish Balay Input Parameters: 107519585528SSatish Balay + mat - the SeqSBAIJ matrix 107649b5e25fSSatish Balay - indices - the column indices 107749b5e25fSSatish Balay 107849b5e25fSSatish Balay Level: advanced 107949b5e25fSSatish Balay 108049b5e25fSSatish Balay Notes: 108149b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 108249b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 108349b5e25fSSatish Balay of the MatSetValues() operation. 108449b5e25fSSatish Balay 108549b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1086d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 108749b5e25fSSatish Balay 1088ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 108949b5e25fSSatish Balay 1090ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 109149b5e25fSSatish Balay @*/ 1092be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 109349b5e25fSSatish Balay { 109413f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 109549b5e25fSSatish Balay 109649b5e25fSSatish Balay PetscFunctionBegin; 10974482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 10984482741eSBarry Smith PetscValidPointer(indices,2); 1099c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 110049b5e25fSSatish Balay if (f) { 110149b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 110249b5e25fSSatish Balay } else { 1103e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 110449b5e25fSSatish Balay } 110549b5e25fSSatish Balay PetscFunctionReturn(0); 110649b5e25fSSatish Balay } 110749b5e25fSSatish Balay 11084a2ae208SSatish Balay #undef __FUNCT__ 11093c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 11103c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 11113c896bc6SHong Zhang { 11123c896bc6SHong Zhang PetscErrorCode ierr; 11133c896bc6SHong Zhang 11143c896bc6SHong Zhang PetscFunctionBegin; 11153c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 11163c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 11173c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 11183c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 11193c896bc6SHong Zhang 11203c896bc6SHong Zhang if (a->i[A->m] != b->i[B->m]) { 11213c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 11223c896bc6SHong Zhang } 11233c896bc6SHong Zhang ierr = PetscMemcpy(b->a,a->a,(a->i[A->m])*sizeof(PetscScalar));CHKERRQ(ierr); 11243c896bc6SHong Zhang } else { 1125f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 11263c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1127f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 11283c896bc6SHong Zhang } 11293c896bc6SHong Zhang PetscFunctionReturn(0); 11303c896bc6SHong Zhang } 11313c896bc6SHong Zhang 11323c896bc6SHong Zhang #undef __FUNCT__ 11334a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1134dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1135273d9f13SBarry Smith { 1136dfbe8321SBarry Smith PetscErrorCode ierr; 1137273d9f13SBarry Smith 1138273d9f13SBarry Smith PetscFunctionBegin; 1139ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1140273d9f13SBarry Smith PetscFunctionReturn(0); 1141273d9f13SBarry Smith } 1142273d9f13SBarry Smith 1143a6ece127SHong Zhang #undef __FUNCT__ 1144a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1145dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1146a6ece127SHong Zhang { 1147a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1148a6ece127SHong Zhang PetscFunctionBegin; 1149a6ece127SHong Zhang *array = a->a; 1150a6ece127SHong Zhang PetscFunctionReturn(0); 1151a6ece127SHong Zhang } 1152a6ece127SHong Zhang 1153a6ece127SHong Zhang #undef __FUNCT__ 1154a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1155dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1156a6ece127SHong Zhang { 1157a6ece127SHong Zhang PetscFunctionBegin; 1158a6ece127SHong Zhang PetscFunctionReturn(0); 1159a6ece127SHong Zhang } 1160a6ece127SHong Zhang 116142ee4b1aSHong Zhang #include "petscblaslapack.h" 116242ee4b1aSHong Zhang #undef __FUNCT__ 116342ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1164f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 116542ee4b1aSHong Zhang { 116642ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1167dfbe8321SBarry Smith PetscErrorCode ierr; 1168521d7252SBarry Smith PetscInt i,bs=Y->bs,bs2,j; 11694ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 117042ee4b1aSHong Zhang 117142ee4b1aSHong Zhang PetscFunctionBegin; 117242ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1173f4df32b1SMatthew Knepley PetscScalar alpha = a; 1174f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1175c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1176c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1177c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1178c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1179c537a176SHong Zhang } 1180c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1181c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1182c4319e64SHong Zhang y->XtoY = X; 1183c537a176SHong Zhang } 1184c4319e64SHong Zhang bs2 = bs*bs; 11856c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1186c4319e64SHong Zhang j = 0; 1187c4319e64SHong Zhang while (j < bs2){ 1188f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1189c4319e64SHong Zhang j++; 1190c537a176SHong Zhang } 1191c4319e64SHong Zhang } 1192*ae15b995SBarry 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); 119342ee4b1aSHong Zhang } else { 1194f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1195f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1196f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 119742ee4b1aSHong Zhang } 119842ee4b1aSHong Zhang PetscFunctionReturn(0); 119942ee4b1aSHong Zhang } 120042ee4b1aSHong Zhang 1201efcf0fc3SBarry Smith #undef __FUNCT__ 1202efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1203dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1204efcf0fc3SBarry Smith { 1205efcf0fc3SBarry Smith PetscFunctionBegin; 1206efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1207efcf0fc3SBarry Smith PetscFunctionReturn(0); 1208efcf0fc3SBarry Smith } 1209efcf0fc3SBarry Smith 1210efcf0fc3SBarry Smith #undef __FUNCT__ 1211efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1212dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1213efcf0fc3SBarry Smith { 1214efcf0fc3SBarry Smith PetscFunctionBegin; 1215efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1216efcf0fc3SBarry Smith PetscFunctionReturn(0); 1217efcf0fc3SBarry Smith } 1218efcf0fc3SBarry Smith 1219efcf0fc3SBarry Smith #undef __FUNCT__ 1220efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1221dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1222efcf0fc3SBarry Smith { 1223efcf0fc3SBarry Smith PetscFunctionBegin; 1224efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1225efcf0fc3SBarry Smith PetscFunctionReturn(0); 1226efcf0fc3SBarry Smith } 1227efcf0fc3SBarry Smith 122899cafbc1SBarry Smith #undef __FUNCT__ 122999cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 123099cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 123199cafbc1SBarry Smith { 123299cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 123399cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 123499cafbc1SBarry Smith PetscScalar *aa = a->a; 123599cafbc1SBarry Smith 123699cafbc1SBarry Smith PetscFunctionBegin; 123799cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 123899cafbc1SBarry Smith PetscFunctionReturn(0); 123999cafbc1SBarry Smith } 124099cafbc1SBarry Smith 124199cafbc1SBarry Smith #undef __FUNCT__ 124299cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 124399cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 124499cafbc1SBarry Smith { 124599cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 124699cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 124799cafbc1SBarry Smith PetscScalar *aa = a->a; 124899cafbc1SBarry Smith 124999cafbc1SBarry Smith PetscFunctionBegin; 125099cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 125199cafbc1SBarry Smith PetscFunctionReturn(0); 125299cafbc1SBarry Smith } 125399cafbc1SBarry Smith 125449b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 125549b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 125649b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 125749b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 125849b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 125997304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1260431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1261e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 126249b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 126349b5e25fSSatish Balay 0, 126449b5e25fSSatish Balay 0, 126597304618SKris Buschelman /*10*/ 0, 126649b5e25fSSatish Balay 0, 12675f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1268d06b337dSHong Zhang MatRelax_SeqSBAIJ, 126949b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 127097304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 127149b5e25fSSatish Balay MatEqual_SeqSBAIJ, 127249b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 127349b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 127449b5e25fSSatish Balay MatNorm_SeqSBAIJ, 127597304618SKris Buschelman /*20*/ 0, 127649b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 127749b5e25fSSatish Balay 0, 127849b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 127949b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1280dcf5cc72SBarry Smith /*25*/ 0, 128149b5e25fSSatish Balay 0, 128249b5e25fSSatish Balay 0, 12835f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12845f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 128597304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1286c464158bSHong Zhang 0, 12874d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1288a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1289a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 129097304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 129149b5e25fSSatish Balay 0, 129249b5e25fSSatish Balay 0, 129349b5e25fSSatish Balay 0, 1294950f1e5bSHong Zhang 0, 129597304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 129649b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 129749b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 129849b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 12993c896bc6SHong Zhang MatCopy_SeqSBAIJ, 130097304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 130149b5e25fSSatish Balay MatScale_SeqSBAIJ, 130249b5e25fSSatish Balay 0, 130349b5e25fSSatish Balay 0, 130449b5e25fSSatish Balay 0, 1305521d7252SBarry Smith /*50*/ 0, 130649b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 130749b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 130849b5e25fSSatish Balay 0, 130949b5e25fSSatish Balay 0, 131097304618SKris Buschelman /*55*/ 0, 131149b5e25fSSatish Balay 0, 131249b5e25fSSatish Balay 0, 131349b5e25fSSatish Balay 0, 131449b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 131597304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 131649b5e25fSSatish Balay 0, 131749b5e25fSSatish Balay 0, 13188a124369SBarry Smith MatGetPetscMaps_Petsc, 1319d959ec07SHong Zhang 0, 132097304618SKris Buschelman /*65*/ 0, 1321d959ec07SHong Zhang 0, 1322d959ec07SHong Zhang 0, 1323d959ec07SHong Zhang 0, 1324d959ec07SHong Zhang 0, 132597304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13263e0d88b5SBarry Smith 0, 13273e0d88b5SBarry Smith 0, 13283e0d88b5SBarry Smith 0, 13293e0d88b5SBarry Smith 0, 133097304618SKris Buschelman /*75*/ 0, 13313e0d88b5SBarry Smith 0, 13323e0d88b5SBarry Smith 0, 13333e0d88b5SBarry Smith 0, 13343e0d88b5SBarry Smith 0, 133597304618SKris Buschelman /*80*/ 0, 13363e0d88b5SBarry Smith 0, 13373e0d88b5SBarry Smith 0, 13383e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 133997304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13403e0d88b5SBarry Smith #else 134197304618SKris Buschelman 0, 13423e0d88b5SBarry Smith #endif 1343865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1344865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1345865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1346efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1347865e5f61SKris Buschelman 0, 1348865e5f61SKris Buschelman 0, 1349865e5f61SKris Buschelman /*90*/ 0, 1350865e5f61SKris Buschelman 0, 1351865e5f61SKris Buschelman 0, 1352865e5f61SKris Buschelman 0, 1353865e5f61SKris Buschelman 0, 1354865e5f61SKris Buschelman /*95*/ 0, 1355865e5f61SKris Buschelman 0, 1356865e5f61SKris Buschelman 0, 135799cafbc1SBarry Smith 0, 135899cafbc1SBarry Smith 0, 135999cafbc1SBarry Smith /*100*/0, 136099cafbc1SBarry Smith 0, 136199cafbc1SBarry Smith 0, 136299cafbc1SBarry Smith 0, 136399cafbc1SBarry Smith 0, 136499cafbc1SBarry Smith /*105*/0, 136599cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1366f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1367f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 1368f5edf698SHong Zhang MatRestoreRowUpperTriangular_SeqSBAIJ 136999cafbc1SBarry Smith }; 1370be1d678aSKris Buschelman 137149b5e25fSSatish Balay EXTERN_C_BEGIN 13724a2ae208SSatish Balay #undef __FUNCT__ 13734a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1374be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 137549b5e25fSSatish Balay { 13764afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1377521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 1378dfbe8321SBarry Smith PetscErrorCode ierr; 137949b5e25fSSatish Balay 138049b5e25fSSatish Balay PetscFunctionBegin; 138149b5e25fSSatish Balay if (aij->nonew != 1) { 1382e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 138349b5e25fSSatish Balay } 138449b5e25fSSatish Balay 138549b5e25fSSatish Balay /* allocate space for values if not already there */ 138649b5e25fSSatish Balay if (!aij->saved_values) { 138787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 138849b5e25fSSatish Balay } 138949b5e25fSSatish Balay 139049b5e25fSSatish Balay /* copy values over */ 139187828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 139249b5e25fSSatish Balay PetscFunctionReturn(0); 139349b5e25fSSatish Balay } 139449b5e25fSSatish Balay EXTERN_C_END 139549b5e25fSSatish Balay 139649b5e25fSSatish Balay EXTERN_C_BEGIN 13974a2ae208SSatish Balay #undef __FUNCT__ 13984a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1399be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 140049b5e25fSSatish Balay { 14014afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 14026849ba73SBarry Smith PetscErrorCode ierr; 1403521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 140449b5e25fSSatish Balay 140549b5e25fSSatish Balay PetscFunctionBegin; 140649b5e25fSSatish Balay if (aij->nonew != 1) { 1407e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 140849b5e25fSSatish Balay } 140949b5e25fSSatish Balay if (!aij->saved_values) { 1410e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 141149b5e25fSSatish Balay } 141249b5e25fSSatish Balay 141349b5e25fSSatish Balay /* copy values over */ 141487828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 141549b5e25fSSatish Balay PetscFunctionReturn(0); 141649b5e25fSSatish Balay } 141749b5e25fSSatish Balay EXTERN_C_END 141849b5e25fSSatish Balay 14198549e402SHong Zhang EXTERN_C_BEGIN 14204a2ae208SSatish Balay #undef __FUNCT__ 1421a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1422be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 142349b5e25fSSatish Balay { 1424c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14256849ba73SBarry Smith PetscErrorCode ierr; 1426ab93d7beSBarry Smith PetscInt i,mbs,bs2; 1427ab93d7beSBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg; 142849b5e25fSSatish Balay 142949b5e25fSSatish Balay PetscFunctionBegin; 1430273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1431e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1432c464158bSHong Zhang mbs = B->m/bs; 143349b5e25fSSatish Balay bs2 = bs*bs; 143449b5e25fSSatish Balay 1435c464158bSHong Zhang if (mbs*bs != B->m) { 143629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 143749b5e25fSSatish Balay } 143849b5e25fSSatish Balay 1439ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1440ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1441ab93d7beSBarry Smith nz = 0; 1442ab93d7beSBarry Smith } 1443ab93d7beSBarry Smith 1444435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 144577431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 144649b5e25fSSatish Balay if (nnz) { 144749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 144877431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 144977431f27SBarry 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); 145049b5e25fSSatish Balay } 145149b5e25fSSatish Balay } 145249b5e25fSSatish Balay 1453e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 145449b5e25fSSatish Balay if (!flg) { 145549b5e25fSSatish Balay switch (bs) { 145649b5e25fSSatish Balay case 1: 14574ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 145849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1459d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1460e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 146149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 146249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1463431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1464431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 146549b5e25fSSatish Balay break; 146649b5e25fSSatish Balay case 2: 14674ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 146849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1469e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 147049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 147149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1472431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1473431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 147449b5e25fSSatish Balay break; 147549b5e25fSSatish Balay case 3: 14765f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 147749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1478e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 147949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 148049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1481431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1482431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 148349b5e25fSSatish Balay break; 148449b5e25fSSatish Balay case 4: 14855f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 148649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1487e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 148849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 148949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1490431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1491431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 149249b5e25fSSatish Balay break; 149349b5e25fSSatish Balay case 5: 14945f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 149549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1496e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 149749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 149849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1499431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1500431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 150149b5e25fSSatish Balay break; 150249b5e25fSSatish Balay case 6: 15035f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 150449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1505e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 150649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 150749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1508431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1509431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 151049b5e25fSSatish Balay break; 151149b5e25fSSatish Balay case 7: 1512de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 151349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1514e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1515de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 151649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1517431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1518431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 151949b5e25fSSatish Balay break; 152049b5e25fSSatish Balay } 152149b5e25fSSatish Balay } 152249b5e25fSSatish Balay 152349b5e25fSSatish Balay b->mbs = mbs; 15244afc71dfSHong Zhang b->nbs = mbs; 1525ab93d7beSBarry Smith if (!skipallocation) { 1526ab93d7beSBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1527ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1528ab93d7beSBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 152949b5e25fSSatish Balay if (!nnz) { 1530435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 153149b5e25fSSatish Balay else if (nz <= 0) nz = 1; 153249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15338cef66ccSHong Zhang b->imax[i] = nz; 153449b5e25fSSatish Balay } 1535153ea458SHong Zhang nz = nz*mbs; /* total nz */ 153649b5e25fSSatish Balay } else { 153749b5e25fSSatish Balay nz = 0; 15388cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 153949b5e25fSSatish Balay } 15406c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 154149b5e25fSSatish Balay 154249b5e25fSSatish Balay /* allocate the matrix space */ 1543a96a251dSBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->m+1,PetscInt,&b->i);CHKERRQ(ierr); 15446c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 154513f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 154649b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 154749b5e25fSSatish Balay 154849b5e25fSSatish Balay /* pointer to beginning of each row */ 154949b5e25fSSatish Balay b->i[0] = 0; 155049b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 155149b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 155249b5e25fSSatish Balay } 1553ab93d7beSBarry Smith } 155449b5e25fSSatish Balay 1555521d7252SBarry Smith B->bs = bs; 155649b5e25fSSatish Balay b->bs2 = bs2; 15576c6c5352SBarry Smith b->nz = 0; 15586c6c5352SBarry Smith b->maxnz = nz*bs2; 1559153ea458SHong Zhang 156016cdd363SHong Zhang b->inew = 0; 156116cdd363SHong Zhang b->jnew = 0; 156216cdd363SHong Zhang b->anew = 0; 156316cdd363SHong Zhang b->a2anew = 0; 15641a3463dfSHong Zhang b->permute = PETSC_FALSE; 1565c464158bSHong Zhang PetscFunctionReturn(0); 1566c464158bSHong Zhang } 1567a23d5eceSKris Buschelman EXTERN_C_END 1568153ea458SHong Zhang 1569d769727bSBarry Smith EXTERN_C_BEGIN 1570f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1571f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1572d769727bSBarry Smith EXTERN_C_END 1573d769727bSBarry Smith 15740bad9183SKris Buschelman /*MC 1575fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15760bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15770bad9183SKris Buschelman 15780bad9183SKris Buschelman Options Database Keys: 15790bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15800bad9183SKris Buschelman 15810bad9183SKris Buschelman Level: beginner 15820bad9183SKris Buschelman 15830bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15840bad9183SKris Buschelman M*/ 15850bad9183SKris Buschelman 1586a23d5eceSKris Buschelman EXTERN_C_BEGIN 1587a23d5eceSKris Buschelman #undef __FUNCT__ 1588a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1589be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1590a23d5eceSKris Buschelman { 1591a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1592dfbe8321SBarry Smith PetscErrorCode ierr; 159313f74950SBarry Smith PetscMPIInt size; 1594941593c8SHong Zhang PetscTruth flg; 1595a23d5eceSKris Buschelman 1596a23d5eceSKris Buschelman PetscFunctionBegin; 1597a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1598a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1599a23d5eceSKris Buschelman B->m = B->M = PetscMax(B->m,B->M); 1600a23d5eceSKris Buschelman B->n = B->N = PetscMax(B->n,B->N); 1601a23d5eceSKris Buschelman 1602a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1603a23d5eceSKris Buschelman B->data = (void*)b; 1604a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1605a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1606a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1607a23d5eceSKris Buschelman B->factor = 0; 1608a23d5eceSKris Buschelman B->mapping = 0; 1609a23d5eceSKris Buschelman b->row = 0; 1610a23d5eceSKris Buschelman b->icol = 0; 1611a23d5eceSKris Buschelman b->reallocs = 0; 1612a23d5eceSKris Buschelman b->saved_values = 0; 1613a23d5eceSKris Buschelman 1614a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1615a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 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); 1783f69a0ea3SMatthew Knepley ierr = MatSetSizes(C,A->m,A->n,A->m,A->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 179682327fa8SHong Zhang C->M = A->M; 179782327fa8SHong Zhang C->N = A->N; 1798521d7252SBarry Smith C->bs = A->bs; 179949b5e25fSSatish Balay c->bs2 = a->bs2; 180049b5e25fSSatish Balay c->mbs = a->mbs; 180149b5e25fSSatish Balay c->nbs = a->nbs; 180249b5e25fSSatish Balay 180313f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr); 180413f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr); 180549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 180649b5e25fSSatish Balay c->imax[i] = a->imax[i]; 180749b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 180849b5e25fSSatish Balay } 1809b40805acSSatish Balay ierr = PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 181049b5e25fSSatish Balay 181149b5e25fSSatish Balay /* allocate the matrix space */ 1812b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 181349b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 181413f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 1815b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 181649b5e25fSSatish Balay if (mbs > 0) { 181713f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 181849b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 181949b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 182049b5e25fSSatish Balay } else { 182149b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 182249b5e25fSSatish Balay } 182349b5e25fSSatish Balay } 182449b5e25fSSatish Balay 182549b5e25fSSatish Balay c->sorted = a->sorted; 182649b5e25fSSatish Balay c->roworiented = a->roworiented; 182749b5e25fSSatish Balay c->nonew = a->nonew; 182849b5e25fSSatish Balay 182949b5e25fSSatish Balay if (a->diag) { 183013f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 183152e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 183249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 183349b5e25fSSatish Balay c->diag[i] = a->diag[i]; 183449b5e25fSSatish Balay } 183549b5e25fSSatish Balay } else c->diag = 0; 18366c6c5352SBarry Smith c->nz = a->nz; 18376c6c5352SBarry Smith c->maxnz = a->maxnz; 183849b5e25fSSatish Balay c->solve_work = 0; 183949b5e25fSSatish Balay c->mult_work = 0; 184049b5e25fSSatish Balay *B = C; 1841b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 184249b5e25fSSatish Balay PetscFunctionReturn(0); 184349b5e25fSSatish Balay } 184449b5e25fSSatish Balay 18454a2ae208SSatish Balay #undef __FUNCT__ 18464a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1847f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A) 184849b5e25fSSatish Balay { 184949b5e25fSSatish Balay Mat_SeqSBAIJ *a; 185049b5e25fSSatish Balay Mat B; 18516849ba73SBarry Smith PetscErrorCode ierr; 185213f74950SBarry Smith int fd; 185313f74950SBarry Smith PetscMPIInt size; 185413f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 185513f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 185613f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 185713f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 185887828ca2SBarry Smith PetscScalar *aa; 185949b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 186049b5e25fSSatish Balay 186149b5e25fSSatish Balay PetscFunctionBegin; 1862b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 186349b5e25fSSatish Balay bs2 = bs*bs; 186449b5e25fSSatish Balay 186549b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 186629bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1867b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 186849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1869552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 187049b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 187149b5e25fSSatish Balay 187249b5e25fSSatish Balay if (header[3] < 0) { 187329bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 187449b5e25fSSatish Balay } 187549b5e25fSSatish Balay 187629bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 187749b5e25fSSatish Balay 187849b5e25fSSatish Balay /* 187949b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 188049b5e25fSSatish Balay divisible by the blocksize 188149b5e25fSSatish Balay */ 188249b5e25fSSatish Balay mbs = M/bs; 188349b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 188449b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 188549b5e25fSSatish Balay else mbs++; 188649b5e25fSSatish Balay if (extra_rows) { 1887*ae15b995SBarry Smith ierr = PetscInfo(0,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 188849b5e25fSSatish Balay } 188949b5e25fSSatish Balay 189049b5e25fSSatish Balay /* read in row lengths */ 189113f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 189249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 189349b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 189449b5e25fSSatish Balay 189549b5e25fSSatish Balay /* read in column indices */ 189613f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 189749b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 189849b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 189949b5e25fSSatish Balay 190049b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 190113f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 190213f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 190313f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 190413f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 190549b5e25fSSatish Balay masked = mask + mbs; 190649b5e25fSSatish Balay rowcount = 0; nzcount = 0; 190749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 190849b5e25fSSatish Balay nmask = 0; 190949b5e25fSSatish Balay for (j=0; j<bs; j++) { 191049b5e25fSSatish Balay kmax = rowlengths[rowcount]; 191149b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19122d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 191303630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 191449b5e25fSSatish Balay } 191549b5e25fSSatish Balay rowcount++; 191649b5e25fSSatish Balay } 1917574b2666SHong Zhang s_browlengths[i] += nmask; 1918574b2666SHong Zhang 191949b5e25fSSatish Balay /* zero out the mask elements we set */ 192049b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 192149b5e25fSSatish Balay } 192249b5e25fSSatish Balay 192349b5e25fSSatish Balay /* create our matrix */ 1924f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 1925f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 19269abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 1927ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 192849b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 192949b5e25fSSatish Balay 193049b5e25fSSatish Balay /* set matrix "i" values */ 193149b5e25fSSatish Balay a->i[0] = 0; 193249b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1933574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1934574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 193549b5e25fSSatish Balay } 19366c6c5352SBarry Smith a->nz = a->i[mbs]; 193749b5e25fSSatish Balay 193849b5e25fSSatish Balay /* read in nonzero values */ 193987828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 194049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 194149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 194249b5e25fSSatish Balay 194349b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 194449b5e25fSSatish Balay nzcount = 0; jcount = 0; 194549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 194649b5e25fSSatish Balay nzcountb = nzcount; 194749b5e25fSSatish Balay nmask = 0; 194849b5e25fSSatish Balay for (j=0; j<bs; j++) { 194949b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 195049b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19512d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 195203630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19532d703238SHong Zhang } 19542d703238SHong Zhang } 19552d703238SHong Zhang /* sort the masked values */ 19562d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19572d703238SHong Zhang 19582d703238SHong Zhang /* set "j" values into matrix */ 19592d703238SHong Zhang maskcount = 1; 19602d703238SHong Zhang for (j=0; j<nmask; j++) { 196149b5e25fSSatish Balay a->j[jcount++] = masked[j]; 196249b5e25fSSatish Balay mask[masked[j]] = maskcount++; 196349b5e25fSSatish Balay } 1964574b2666SHong Zhang 196549b5e25fSSatish Balay /* set "a" values into matrix */ 196649b5e25fSSatish Balay ishift = bs2*a->i[i]; 196749b5e25fSSatish Balay for (j=0; j<bs; j++) { 196849b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 196949b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1970574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1971574b2666SHong Zhang if (tmp >= i){ 197249b5e25fSSatish Balay block = mask[tmp] - 1; 197349b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 197449b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1975574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1976574b2666SHong Zhang } 1977574b2666SHong Zhang nzcountb++; 197849b5e25fSSatish Balay } 197949b5e25fSSatish Balay } 198049b5e25fSSatish Balay /* zero out the mask elements we set */ 198149b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 198249b5e25fSSatish Balay } 19836c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 198449b5e25fSSatish Balay 198549b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1986574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 198749b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 198849b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 198949b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 199049b5e25fSSatish Balay 19919abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19929abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 199349b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19949abb65ffSKris Buschelman *A = B; 199549b5e25fSSatish Balay PetscFunctionReturn(0); 199649b5e25fSSatish Balay } 1997574b2666SHong Zhang 1998d06b337dSHong Zhang #undef __FUNCT__ 1999d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 200013f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 2001d06b337dSHong Zhang { 2002d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2003d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 2004d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 20056849ba73SBarry Smith PetscErrorCode ierr; 2006521d7252SBarry Smith PetscInt m=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j; 2007521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 2008d06b337dSHong Zhang 2009d06b337dSHong Zhang PetscFunctionBegin; 2010b965ef7fSBarry Smith its = its*lits; 201177431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 2012d06b337dSHong Zhang 2013d06b337dSHong Zhang if (bs > 1) 2014d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 2015d06b337dSHong Zhang 20161ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2017d06b337dSHong Zhang if (xx != bb) { 20181ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 2019d06b337dSHong Zhang } else { 2020d06b337dSHong Zhang b = x; 2021d06b337dSHong Zhang } 2022d06b337dSHong Zhang 2023d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 2024d06b337dSHong Zhang 2025d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 20262798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2027d06b337dSHong Zhang for (i=0; i<m; i++) 2028d06b337dSHong Zhang t[i] = b[i]; 2029d06b337dSHong Zhang 2030d06b337dSHong Zhang for (i=0; i<m; i++){ 203144706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2032d06b337dSHong Zhang v = aa + ai[i] + 1; 2033d06b337dSHong Zhang vj = aj + ai[i] + 1; 2034d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2035d06b337dSHong Zhang x[i] = omega*t[i]/d; 2036d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 2037efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2038d06b337dSHong Zhang } 2039d06b337dSHong Zhang } 2040d06b337dSHong Zhang 20412798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 2042d06b337dSHong Zhang for (i=0; i<m; i++) 2043d06b337dSHong Zhang t[i] = b[i]; 2044d06b337dSHong Zhang 2045d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2046d06b337dSHong Zhang v = aa + ai[i] + 1; 2047d06b337dSHong Zhang vj = aj + ai[i] + 1; 2048d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2049d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2050efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2051d06b337dSHong Zhang } 2052d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2053d06b337dSHong Zhang d = *(aa + ai[i]); 2054d06b337dSHong Zhang v = aa + ai[i] + 1; 2055d06b337dSHong Zhang vj = aj + ai[i] + 1; 2056d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2057d06b337dSHong Zhang sum = t[i]; 2058d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2059efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2060d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2061d06b337dSHong Zhang } 2062d06b337dSHong Zhang } 2063d06b337dSHong Zhang its--; 2064d06b337dSHong Zhang } 2065d06b337dSHong Zhang 2066d06b337dSHong Zhang while (its--) { 206744706875SHong Zhang /* 206844706875SHong Zhang forward sweep: 206944706875SHong Zhang for i=0,...,m-1: 207044706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 207144706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 207244706875SHong Zhang b = b - x[i]*U^T(i,:); 2073d06b337dSHong Zhang 207444706875SHong Zhang */ 20752798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2076d06b337dSHong Zhang for (i=0; i<m; i++) 2077d06b337dSHong Zhang t[i] = b[i]; 2078d06b337dSHong Zhang 2079d06b337dSHong Zhang for (i=0; i<m; i++){ 208044706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2081d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2082d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2083d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2084d06b337dSHong Zhang sum = t[i]; 2085d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2086d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2087d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2088efee365bSSatish Balay ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2089d06b337dSHong Zhang } 2090d06b337dSHong Zhang } 2091d06b337dSHong Zhang 20922798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 209344706875SHong Zhang /* 209444706875SHong Zhang backward sweep: 209544706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 209644706875SHong Zhang for i=m-1,...,0: 209744706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 209844706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 209944706875SHong Zhang */ 2100d06b337dSHong Zhang for (i=0; i<m; i++) 2101d06b337dSHong Zhang t[i] = b[i]; 2102d06b337dSHong Zhang 2103d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2104d06b337dSHong Zhang v = aa + ai[i] + 1; 2105d06b337dSHong Zhang vj = aj + ai[i] + 1; 2106d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2107d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2108efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2109d06b337dSHong Zhang } 2110d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2111d06b337dSHong Zhang d = *(aa + ai[i]); 2112d06b337dSHong Zhang v = aa + ai[i] + 1; 2113d06b337dSHong Zhang vj = aj + ai[i] + 1; 2114d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2115d06b337dSHong Zhang sum = t[i]; 2116d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2117efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2118d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2119d06b337dSHong Zhang } 2120d06b337dSHong Zhang } 2121d06b337dSHong Zhang } 2122d06b337dSHong Zhang 2123d06b337dSHong Zhang ierr = PetscFree(t);CHKERRQ(ierr); 21241ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2125d06b337dSHong Zhang if (bb != xx) { 21261ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2127d06b337dSHong Zhang } 2128d06b337dSHong Zhang PetscFunctionReturn(0); 2129d06b337dSHong Zhang } 2130d06b337dSHong Zhang 2131d06b337dSHong Zhang 2132d06b337dSHong Zhang 2133d06b337dSHong Zhang 213449b5e25fSSatish Balay 213549b5e25fSSatish Balay 2136