149b5e25fSSatish Balay 249b5e25fSSatish Balay /* 3a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 449b5e25fSSatish Balay matrix storage format. 549b5e25fSSatish Balay */ 69e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 749b5e25fSSatish Balay #include "src/inline/spops.h" 8aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 949b5e25fSSatish Balay 1049b5e25fSSatish Balay #define CHUNKSIZE 10 1149b5e25fSSatish Balay 1249b5e25fSSatish Balay /* 1349b5e25fSSatish Balay Checks for missing diagonals 1449b5e25fSSatish Balay */ 154a2ae208SSatish Balay #undef __FUNCT__ 164a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 17dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A) 1849b5e25fSSatish Balay { 19045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 206849ba73SBarry Smith PetscErrorCode ierr; 216849ba73SBarry Smith int *diag,*jj = a->j,i; 2249b5e25fSSatish Balay 2349b5e25fSSatish Balay PetscFunctionBegin; 24045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2549b5e25fSSatish Balay diag = a->diag; 2649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 27*77431f27SBarry Smith if (jj[diag[i]] != i) SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Matrix is missing diagonal number %D",i); 2849b5e25fSSatish Balay } 2949b5e25fSSatish Balay PetscFunctionReturn(0); 3049b5e25fSSatish Balay } 3149b5e25fSSatish Balay 324a2ae208SSatish Balay #undef __FUNCT__ 334a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 34dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 3549b5e25fSSatish Balay { 36045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 376849ba73SBarry Smith PetscErrorCode ierr; 386849ba73SBarry Smith int i,mbs = a->mbs; 3949b5e25fSSatish Balay 4049b5e25fSSatish Balay PetscFunctionBegin; 4149b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4249b5e25fSSatish Balay 43b424e231SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr); 44b424e231SHong Zhang PetscLogObjectMemory(A,(mbs+1)*sizeof(int)); 45b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4649b5e25fSSatish Balay PetscFunctionReturn(0); 4749b5e25fSSatish Balay } 4849b5e25fSSatish Balay 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 516849ba73SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done) 5249b5e25fSSatish Balay { 53a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 549e37e433SSatish Balay int n = a->mbs,i; 5549b5e25fSSatish Balay 5649b5e25fSSatish Balay PetscFunctionBegin; 57d3e5a4abSHong Zhang *nn = n; 58a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 59a1373b80SHong Zhang 60a6ece127SHong Zhang if (oshift == 1) { 61a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 626c6c5352SBarry Smith int nz = a->i[n]; 636c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 64a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 65a1373b80SHong Zhang *ia = a->i; *ja = a->j; 66a1373b80SHong Zhang } else { 67a1373b80SHong Zhang *ia = a->i; *ja = a->j; 68a6ece127SHong Zhang } 6949b5e25fSSatish Balay PetscFunctionReturn(0); 7049b5e25fSSatish Balay } 7149b5e25fSSatish Balay 724a2ae208SSatish Balay #undef __FUNCT__ 734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 746849ba73SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done) 7549b5e25fSSatish Balay { 76b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 77a15ac5e4SHong Zhang int i,n = a->mbs; 78a6ece127SHong Zhang 7949b5e25fSSatish Balay PetscFunctionBegin; 8049b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 81a6ece127SHong Zhang 82a6ece127SHong Zhang if (oshift == 1) { 836c6c5352SBarry Smith int nz = a->i[n]-1; 846c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 85a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 86a6ece127SHong Zhang } 87a6ece127SHong Zhang PetscFunctionReturn(0); 8849b5e25fSSatish Balay } 8949b5e25fSSatish Balay 904a2ae208SSatish Balay #undef __FUNCT__ 914a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 92dfbe8321SBarry Smith PetscErrorCode MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9349b5e25fSSatish Balay { 9449b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 9549b5e25fSSatish Balay 9649b5e25fSSatish Balay PetscFunctionBegin; 9749b5e25fSSatish Balay *bs = sbaij->bs; 9849b5e25fSSatish Balay PetscFunctionReturn(0); 9949b5e25fSSatish Balay } 10049b5e25fSSatish Balay 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 103dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 10449b5e25fSSatish Balay { 10549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 106dfbe8321SBarry Smith PetscErrorCode ierr; 10749b5e25fSSatish Balay 10849b5e25fSSatish Balay PetscFunctionBegin; 109*77431f27SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->m,a->nz); 11049b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 11149b5e25fSSatish Balay if (!a->singlemalloc) { 11249b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 11363c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 11449b5e25fSSatish Balay } 11549b5e25fSSatish Balay if (a->row) { 11649b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 11749b5e25fSSatish Balay } 11849b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 11949b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 12049b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 12149b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 122d59c15a7SBarry Smith if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 123d59c15a7SBarry Smith if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);} 124d59c15a7SBarry Smith if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 12549b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 126c4319e64SHong Zhang if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 1271a3463dfSHong Zhang 1281a3463dfSHong Zhang if (a->inew){ 1291a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1301a3463dfSHong Zhang a->inew = 0; 1311a3463dfSHong Zhang } 13249b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 133901853e0SKris Buschelman 134901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 135901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 136901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 137901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 138901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 139901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 14049b5e25fSSatish Balay PetscFunctionReturn(0); 14149b5e25fSSatish Balay } 14249b5e25fSSatish Balay 1434a2ae208SSatish Balay #undef __FUNCT__ 1444a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 145dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op) 14649b5e25fSSatish Balay { 147045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 14849b5e25fSSatish Balay 14949b5e25fSSatish Balay PetscFunctionBegin; 1504d9d31abSKris Buschelman switch (op) { 1514d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1524d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1554d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1584d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1614d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1644d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1674d9d31abSKris Buschelman a->nonew = 1; 1684d9d31abSKris Buschelman break; 1694d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1704d9d31abSKris Buschelman a->nonew = -1; 1714d9d31abSKris Buschelman break; 1724d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1734d9d31abSKris Buschelman a->nonew = -2; 1744d9d31abSKris Buschelman break; 1754d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1764d9d31abSKris Buschelman a->nonew = 0; 1774d9d31abSKris Buschelman break; 1784d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1794d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1804d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1814d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1824d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 183b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1844d9d31abSKris Buschelman break; 1854d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 18629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1879a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1889a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1899a4540c5SBarry Smith case MAT_HERMITIAN: 1909a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 19177e54ba9SKris Buschelman case MAT_SYMMETRIC: 19277e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1939a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1949a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1959a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 19677e54ba9SKris Buschelman break; 1974d9d31abSKris Buschelman default: 19829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 19949b5e25fSSatish Balay } 20049b5e25fSSatish Balay PetscFunctionReturn(0); 20149b5e25fSSatish Balay } 20249b5e25fSSatish Balay 2034a2ae208SSatish Balay #undef __FUNCT__ 2044a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 205dfbe8321SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v) 20649b5e25fSSatish Balay { 20749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2086849ba73SBarry Smith PetscErrorCode ierr; 2096849ba73SBarry Smith int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 21049b5e25fSSatish Balay MatScalar *aa,*aa_i; 21187828ca2SBarry Smith PetscScalar *v_i; 21249b5e25fSSatish Balay 21349b5e25fSSatish Balay PetscFunctionBegin; 21449b5e25fSSatish Balay bs = a->bs; 21549b5e25fSSatish Balay ai = a->i; 21649b5e25fSSatish Balay aj = a->j; 21749b5e25fSSatish Balay aa = a->a; 21849b5e25fSSatish Balay bs2 = a->bs2; 21949b5e25fSSatish Balay 220*77431f27SBarry Smith if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 22149b5e25fSSatish Balay 22249b5e25fSSatish Balay bn = row/bs; /* Block number */ 22349b5e25fSSatish Balay bp = row % bs; /* Block position */ 22449b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 22549b5e25fSSatish Balay *ncols = bs*M; 22649b5e25fSSatish Balay 22749b5e25fSSatish Balay if (v) { 22849b5e25fSSatish Balay *v = 0; 22949b5e25fSSatish Balay if (*ncols) { 23087828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 23149b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23249b5e25fSSatish Balay v_i = *v + i*bs; 23349b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 23449b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 23549b5e25fSSatish Balay } 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay } 23849b5e25fSSatish Balay 23949b5e25fSSatish Balay if (cols) { 24049b5e25fSSatish Balay *cols = 0; 24149b5e25fSSatish Balay if (*ncols) { 24282502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 24349b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 24449b5e25fSSatish Balay cols_i = *cols + i*bs; 24549b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 24649b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24749b5e25fSSatish Balay } 24849b5e25fSSatish Balay } 24949b5e25fSSatish Balay } 25049b5e25fSSatish Balay 25149b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2525ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2535ddb2528SHong Zhang #ifdef column_search 25449b5e25fSSatish Balay v_i = *v + M*bs; 25549b5e25fSSatish Balay cols_i = *cols + M*bs; 25649b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25749b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25849b5e25fSSatish Balay for (j=0; j<M; j++){ 25949b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 26049b5e25fSSatish Balay if (itmp == bn){ 26149b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 26249b5e25fSSatish Balay for (k=0; k<bs; k++) { 26349b5e25fSSatish Balay *cols_i++ = i*bs+k; 26449b5e25fSSatish Balay *v_i++ = aa_i[k]; 26549b5e25fSSatish Balay } 26649b5e25fSSatish Balay *ncols += bs; 26749b5e25fSSatish Balay break; 26849b5e25fSSatish Balay } 26949b5e25fSSatish Balay } 27049b5e25fSSatish Balay } 2715ddb2528SHong Zhang #endif 27249b5e25fSSatish Balay 27349b5e25fSSatish Balay PetscFunctionReturn(0); 27449b5e25fSSatish Balay } 27549b5e25fSSatish Balay 2764a2ae208SSatish Balay #undef __FUNCT__ 2774a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 278dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 27949b5e25fSSatish Balay { 280dfbe8321SBarry Smith PetscErrorCode ierr; 28149b5e25fSSatish Balay 28249b5e25fSSatish Balay PetscFunctionBegin; 28349b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 28449b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 28549b5e25fSSatish Balay PetscFunctionReturn(0); 28649b5e25fSSatish Balay } 28749b5e25fSSatish Balay 2884a2ae208SSatish Balay #undef __FUNCT__ 2894a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 290dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 29149b5e25fSSatish Balay { 292dfbe8321SBarry Smith PetscErrorCode ierr; 29349b5e25fSSatish Balay PetscFunctionBegin; 294999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2958115998fSBarry Smith PetscFunctionReturn(0); 29649b5e25fSSatish Balay } 29749b5e25fSSatish Balay 2984a2ae208SSatish Balay #undef __FUNCT__ 2994a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3006849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 30149b5e25fSSatish Balay { 30249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 303dfbe8321SBarry Smith PetscErrorCode ierr; 304dfbe8321SBarry Smith int i,j,bs = a->bs,k,l,bs2=a->bs2; 305fb9695e5SSatish Balay char *name; 306f3ef73ceSBarry Smith PetscViewerFormat format; 30749b5e25fSSatish Balay 30849b5e25fSSatish Balay PetscFunctionBegin; 30980fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 310b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 311456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 312*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 313fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 31429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 315fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 316b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 31749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31849b5e25fSSatish Balay for (j=0; j<bs; j++) { 319*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 32049b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 32149b5e25fSSatish Balay for (l=0; l<bs; l++) { 32249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 32349b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 324*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 32549b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32649b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 327*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 32849b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32949b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 330*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 33149b5e25fSSatish Balay } 33249b5e25fSSatish Balay #else 33349b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 334*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 33549b5e25fSSatish Balay } 33649b5e25fSSatish Balay #endif 33749b5e25fSSatish Balay } 33849b5e25fSSatish Balay } 339b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 34049b5e25fSSatish Balay } 34149b5e25fSSatish Balay } 342b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 34349b5e25fSSatish Balay } else { 344b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 34549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 34649b5e25fSSatish Balay for (j=0; j<bs; j++) { 347*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 34849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 34949b5e25fSSatish Balay for (l=0; l<bs; l++) { 35049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 35149b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 352*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 35349b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35449b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 355*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 35649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35749b5e25fSSatish Balay } else { 358*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35949b5e25fSSatish Balay } 36049b5e25fSSatish Balay #else 361*77431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 36249b5e25fSSatish Balay #endif 36349b5e25fSSatish Balay } 36449b5e25fSSatish Balay } 365b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 36649b5e25fSSatish Balay } 36749b5e25fSSatish Balay } 368b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 36949b5e25fSSatish Balay } 370b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 37149b5e25fSSatish Balay PetscFunctionReturn(0); 37249b5e25fSSatish Balay } 37349b5e25fSSatish Balay 3744a2ae208SSatish Balay #undef __FUNCT__ 3754a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 3766849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 37749b5e25fSSatish Balay { 37849b5e25fSSatish Balay Mat A = (Mat) Aa; 37949b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 3806849ba73SBarry Smith PetscErrorCode ierr; 3816849ba73SBarry Smith int row,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 38249b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 38349b5e25fSSatish Balay MatScalar *aa; 38449b5e25fSSatish Balay MPI_Comm comm; 385b0a32e0cSBarry Smith PetscViewer viewer; 38649b5e25fSSatish Balay 38749b5e25fSSatish Balay PetscFunctionBegin; 38849b5e25fSSatish Balay /* 38949b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 39049b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 39149b5e25fSSatish Balay rest should return immediately. 39249b5e25fSSatish Balay */ 39349b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 39449b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 39549b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 39649b5e25fSSatish Balay 39749b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 39849b5e25fSSatish Balay 399b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 400b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 40149b5e25fSSatish Balay 40249b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 403b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 40449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 40549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 406c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 40749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 40849b5e25fSSatish Balay aa = a->a + j*bs2; 40949b5e25fSSatish Balay for (k=0; k<bs; k++) { 41049b5e25fSSatish Balay for (l=0; l<bs; l++) { 41149b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 412b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 41349b5e25fSSatish Balay } 41449b5e25fSSatish Balay } 41549b5e25fSSatish Balay } 41649b5e25fSSatish Balay } 417b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 41849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 41949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 420c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 42149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42249b5e25fSSatish Balay aa = a->a + j*bs2; 42349b5e25fSSatish Balay for (k=0; k<bs; k++) { 42449b5e25fSSatish Balay for (l=0; l<bs; l++) { 42549b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 426b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 42749b5e25fSSatish Balay } 42849b5e25fSSatish Balay } 42949b5e25fSSatish Balay } 43049b5e25fSSatish Balay } 43149b5e25fSSatish Balay 432b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 43349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 43449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 435c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 43649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 43749b5e25fSSatish Balay aa = a->a + j*bs2; 43849b5e25fSSatish Balay for (k=0; k<bs; k++) { 43949b5e25fSSatish Balay for (l=0; l<bs; l++) { 44049b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 441b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44249b5e25fSSatish Balay } 44349b5e25fSSatish Balay } 44449b5e25fSSatish Balay } 44549b5e25fSSatish Balay } 44649b5e25fSSatish Balay PetscFunctionReturn(0); 44749b5e25fSSatish Balay } 44849b5e25fSSatish Balay 4494a2ae208SSatish Balay #undef __FUNCT__ 4504a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4516849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 45249b5e25fSSatish Balay { 453dfbe8321SBarry Smith PetscErrorCode ierr; 45449b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 455b0a32e0cSBarry Smith PetscDraw draw; 45649b5e25fSSatish Balay PetscTruth isnull; 45749b5e25fSSatish Balay 45849b5e25fSSatish Balay PetscFunctionBegin; 45949b5e25fSSatish Balay 460b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 461b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 46249b5e25fSSatish Balay 46349b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 464c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 46549b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 466b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 467b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 46849b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 46949b5e25fSSatish Balay PetscFunctionReturn(0); 47049b5e25fSSatish Balay } 47149b5e25fSSatish Balay 4724a2ae208SSatish Balay #undef __FUNCT__ 4734a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 474dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 47549b5e25fSSatish Balay { 476dfbe8321SBarry Smith PetscErrorCode ierr; 47732077d6dSBarry Smith PetscTruth iascii,isdraw; 47849b5e25fSSatish Balay 47949b5e25fSSatish Balay PetscFunctionBegin; 48032077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 481fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 48232077d6dSBarry Smith if (iascii){ 48349b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 48449b5e25fSSatish Balay } else if (isdraw) { 48549b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 48649b5e25fSSatish Balay } else { 487a5e6ed63SBarry Smith Mat B; 488a5e6ed63SBarry Smith ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr); 489a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 490a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 49149b5e25fSSatish Balay } 49249b5e25fSSatish Balay PetscFunctionReturn(0); 49349b5e25fSSatish Balay } 49449b5e25fSSatish Balay 49549b5e25fSSatish Balay 4964a2ae208SSatish Balay #undef __FUNCT__ 4974a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 498dfbe8321SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[]) 49949b5e25fSSatish Balay { 500045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 50149b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 50249b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 50349b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 50449b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 50549b5e25fSSatish Balay 50649b5e25fSSatish Balay PetscFunctionBegin; 50749b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 50849b5e25fSSatish Balay row = im[k]; brow = row/bs; 509*77431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 510*77431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 51149b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 51249b5e25fSSatish Balay nrow = ailen[brow]; 51349b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 514*77431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 515*77431f27SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1); 51649b5e25fSSatish Balay col = in[l] ; 51749b5e25fSSatish Balay bcol = col/bs; 51849b5e25fSSatish Balay cidx = col%bs; 51949b5e25fSSatish Balay ridx = row%bs; 52049b5e25fSSatish Balay high = nrow; 52149b5e25fSSatish Balay low = 0; /* assume unsorted */ 52249b5e25fSSatish Balay while (high-low > 5) { 52349b5e25fSSatish Balay t = (low+high)/2; 52449b5e25fSSatish Balay if (rp[t] > bcol) high = t; 52549b5e25fSSatish Balay else low = t; 52649b5e25fSSatish Balay } 52749b5e25fSSatish Balay for (i=low; i<high; i++) { 52849b5e25fSSatish Balay if (rp[i] > bcol) break; 52949b5e25fSSatish Balay if (rp[i] == bcol) { 53049b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 53149b5e25fSSatish Balay goto finished; 53249b5e25fSSatish Balay } 53349b5e25fSSatish Balay } 53449b5e25fSSatish Balay *v++ = zero; 53549b5e25fSSatish Balay finished:; 53649b5e25fSSatish Balay } 53749b5e25fSSatish Balay } 53849b5e25fSSatish Balay PetscFunctionReturn(0); 53949b5e25fSSatish Balay } 54049b5e25fSSatish Balay 54149b5e25fSSatish Balay 5424a2ae208SSatish Balay #undef __FUNCT__ 5434a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 544dfbe8321SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is) 54549b5e25fSSatish Balay { 5460880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5476849ba73SBarry Smith PetscErrorCode ierr; 5480880e062SHong Zhang int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 5490880e062SHong Zhang int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 5506849ba73SBarry Smith int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval; 5510880e062SHong Zhang PetscTruth roworiented=a->roworiented; 552f15d580aSBarry Smith const MatScalar *value = v; 553f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5540880e062SHong Zhang 55549b5e25fSSatish Balay PetscFunctionBegin; 5560880e062SHong Zhang if (roworiented) { 5570880e062SHong Zhang stepval = (n-1)*bs; 5580880e062SHong Zhang } else { 5590880e062SHong Zhang stepval = (m-1)*bs; 5600880e062SHong Zhang } 5610880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5620880e062SHong Zhang row = im[k]; 5630880e062SHong Zhang if (row < 0) continue; 5640880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 565*77431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5660880e062SHong Zhang #endif 5670880e062SHong Zhang rp = aj + ai[row]; 5680880e062SHong Zhang ap = aa + bs2*ai[row]; 5690880e062SHong Zhang rmax = imax[row]; 5700880e062SHong Zhang nrow = ailen[row]; 5710880e062SHong Zhang low = 0; 5720880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5730880e062SHong Zhang if (in[l] < 0) continue; 5740880e062SHong Zhang col = in[l]; 575b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g) 576*77431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 577b1823623SSatish Balay #endif 5780880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 5790880e062SHong Zhang if (roworiented) { 5800880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 5810880e062SHong Zhang } else { 5820880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 5830880e062SHong Zhang } 5840880e062SHong Zhang if (!sorted) low = 0; high = nrow; 5850880e062SHong Zhang while (high-low > 7) { 5860880e062SHong Zhang t = (low+high)/2; 5870880e062SHong Zhang if (rp[t] > col) high = t; 5880880e062SHong Zhang else low = t; 5890880e062SHong Zhang } 5900880e062SHong Zhang for (i=low; i<high; i++) { 5910880e062SHong Zhang if (rp[i] > col) break; 5920880e062SHong Zhang if (rp[i] == col) { 5930880e062SHong Zhang bap = ap + bs2*i; 5940880e062SHong Zhang if (roworiented) { 5950880e062SHong Zhang if (is == ADD_VALUES) { 5960880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5970880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5980880e062SHong Zhang bap[jj] += *value++; 5990880e062SHong Zhang } 6000880e062SHong Zhang } 6010880e062SHong Zhang } else { 6020880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6030880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6040880e062SHong Zhang bap[jj] = *value++; 6050880e062SHong Zhang } 6060880e062SHong Zhang } 6070880e062SHong Zhang } 6080880e062SHong Zhang } else { 6090880e062SHong Zhang if (is == ADD_VALUES) { 6100880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6110880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6120880e062SHong Zhang *bap++ += *value++; 6130880e062SHong Zhang } 6140880e062SHong Zhang } 6150880e062SHong Zhang } else { 6160880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6170880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6180880e062SHong Zhang *bap++ = *value++; 6190880e062SHong Zhang } 6200880e062SHong Zhang } 6210880e062SHong Zhang } 6220880e062SHong Zhang } 6230880e062SHong Zhang goto noinsert2; 6240880e062SHong Zhang } 6250880e062SHong Zhang } 6260880e062SHong Zhang if (nonew == 1) goto noinsert2; 627*77431f27SBarry Smith else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 6280880e062SHong Zhang if (nrow >= rmax) { 6290880e062SHong Zhang /* there is no extra room in row, therefore enlarge */ 6300880e062SHong Zhang int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6310880e062SHong Zhang MatScalar *new_a; 6320880e062SHong Zhang 633*77431f27SBarry Smith if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 6340880e062SHong Zhang 6350880e062SHong Zhang /* malloc new storage space */ 6360880e062SHong Zhang len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 6370880e062SHong Zhang ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 6380880e062SHong Zhang new_j = (int*)(new_a + bs2*new_nz); 6390880e062SHong Zhang new_i = new_j + new_nz; 6400880e062SHong Zhang 6410880e062SHong Zhang /* copy over old data into new slots */ 6420880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6430880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 6440880e062SHong Zhang ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 6450880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 6460880e062SHong Zhang ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 6470880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 6480880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 6490880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 6500880e062SHong Zhang /* free up old matrix storage */ 6510880e062SHong Zhang ierr = PetscFree(a->a);CHKERRQ(ierr); 6520880e062SHong Zhang if (!a->singlemalloc) { 6530880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 6540880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 6550880e062SHong Zhang } 6560880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 6570880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 6580880e062SHong Zhang 6590880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 6600880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 6610880e062SHong Zhang PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 6626c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 6630880e062SHong Zhang a->reallocs++; 6646c6c5352SBarry Smith a->nz++; 6650880e062SHong Zhang } 6660880e062SHong Zhang N = nrow++ - 1; 6670880e062SHong Zhang /* shift up all the later entries in this row */ 6680880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6690880e062SHong Zhang rp[ii+1] = rp[ii]; 6700880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6710880e062SHong Zhang } 6720880e062SHong Zhang if (N >= i) { 6730880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6740880e062SHong Zhang } 6750880e062SHong Zhang rp[i] = col; 6760880e062SHong Zhang bap = ap + bs2*i; 6770880e062SHong Zhang if (roworiented) { 6780880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6790880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6800880e062SHong Zhang bap[jj] = *value++; 6810880e062SHong Zhang } 6820880e062SHong Zhang } 6830880e062SHong Zhang } else { 6840880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6850880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6860880e062SHong Zhang *bap++ = *value++; 6870880e062SHong Zhang } 6880880e062SHong Zhang } 6890880e062SHong Zhang } 6900880e062SHong Zhang noinsert2:; 6910880e062SHong Zhang low = i; 6920880e062SHong Zhang } 6930880e062SHong Zhang ailen[row] = nrow; 6940880e062SHong Zhang } 6950880e062SHong Zhang PetscFunctionReturn(0); 69649b5e25fSSatish Balay } 69749b5e25fSSatish Balay 6984a2ae208SSatish Balay #undef __FUNCT__ 6994a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 700dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 70149b5e25fSSatish Balay { 70249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7036849ba73SBarry Smith PetscErrorCode ierr; 70449b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 705c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 7066849ba73SBarry Smith int mbs = a->mbs,bs2 = a->bs2,rmax = 0; 70749b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 70849b5e25fSSatish Balay 70949b5e25fSSatish Balay PetscFunctionBegin; 71049b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 71149b5e25fSSatish Balay 71249b5e25fSSatish Balay if (m) rmax = ailen[0]; 71349b5e25fSSatish Balay for (i=1; i<mbs; i++) { 71449b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 71549b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 71649b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 71749b5e25fSSatish Balay if (fshift) { 71849b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 71949b5e25fSSatish Balay N = ailen[i]; 72049b5e25fSSatish Balay for (j=0; j<N; j++) { 72149b5e25fSSatish Balay ip[j-fshift] = ip[j]; 72249b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 72349b5e25fSSatish Balay } 72449b5e25fSSatish Balay } 72549b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 72649b5e25fSSatish Balay } 72749b5e25fSSatish Balay if (mbs) { 72849b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 72949b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 73049b5e25fSSatish Balay } 73149b5e25fSSatish Balay /* reset ilen and imax for each row */ 73249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 73349b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 73449b5e25fSSatish Balay } 7356c6c5352SBarry Smith a->nz = ai[mbs]; 73649b5e25fSSatish Balay 737b424e231SHong Zhang /* diagonals may have moved, reset it */ 738b424e231SHong Zhang if (a->diag) { 739b424e231SHong Zhang ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr); 74049b5e25fSSatish Balay } 741*77431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n", 7426c6c5352SBarry Smith m,A->m,a->bs,fshift*bs2,a->nz*bs2); 743*77431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %D\n", 74449b5e25fSSatish Balay a->reallocs); 745*77431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %D\n",rmax); 74649b5e25fSSatish Balay a->reallocs = 0; 74749b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 74849b5e25fSSatish Balay 74949b5e25fSSatish Balay PetscFunctionReturn(0); 75049b5e25fSSatish Balay } 75149b5e25fSSatish Balay 75249b5e25fSSatish Balay /* 75349b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 75449b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 75549b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 75649b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 75749b5e25fSSatish Balay */ 7584a2ae208SSatish Balay #undef __FUNCT__ 7594a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 760dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 76149b5e25fSSatish Balay { 76249b5e25fSSatish Balay int i,j,k,row; 76349b5e25fSSatish Balay PetscTruth flg; 76449b5e25fSSatish Balay 76549b5e25fSSatish Balay PetscFunctionBegin; 76649b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 76749b5e25fSSatish Balay row = idx[i]; 76849b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 76949b5e25fSSatish Balay sizes[j] = 1; 77049b5e25fSSatish Balay i++; 77149b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 77249b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 77349b5e25fSSatish Balay i++; 77449b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 77549b5e25fSSatish Balay flg = PETSC_TRUE; 77649b5e25fSSatish Balay for (k=1; k<bs; k++) { 77749b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 77849b5e25fSSatish Balay flg = PETSC_FALSE; 77949b5e25fSSatish Balay break; 78049b5e25fSSatish Balay } 78149b5e25fSSatish Balay } 78249b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 78349b5e25fSSatish Balay sizes[j] = bs; 78449b5e25fSSatish Balay i+= bs; 78549b5e25fSSatish Balay } else { 78649b5e25fSSatish Balay sizes[j] = 1; 78749b5e25fSSatish Balay i++; 78849b5e25fSSatish Balay } 78949b5e25fSSatish Balay } 79049b5e25fSSatish Balay } 79149b5e25fSSatish Balay *bs_max = j; 79249b5e25fSSatish Balay PetscFunctionReturn(0); 79349b5e25fSSatish Balay } 79449b5e25fSSatish Balay 7954a2ae208SSatish Balay #undef __FUNCT__ 7964a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 797dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag) 79849b5e25fSSatish Balay { 79949b5e25fSSatish Balay PetscFunctionBegin; 800c0f24835SHong Zhang SETERRQ(PETSC_ERR_SUP,"No support for this function yet"); 80149b5e25fSSatish Balay } 80249b5e25fSSatish Balay 80349b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 80449b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 80549b5e25fSSatish Balay */ 80649b5e25fSSatish Balay 8074a2ae208SSatish Balay #undef __FUNCT__ 8084a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 809dfbe8321SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is) 81049b5e25fSSatish Balay { 81149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8126849ba73SBarry Smith PetscErrorCode ierr; 81349b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 81449b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 81549b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 8166849ba73SBarry Smith int ridx,cidx,bs2=a->bs2; 81749b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 81849b5e25fSSatish Balay 81949b5e25fSSatish Balay PetscFunctionBegin; 82049b5e25fSSatish Balay 82149b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 82249b5e25fSSatish Balay row = im[k]; /* row number */ 82349b5e25fSSatish Balay brow = row/bs; /* block row number */ 82449b5e25fSSatish Balay if (row < 0) continue; 82549b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 826*77431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 82749b5e25fSSatish Balay #endif 82849b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 82949b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 83049b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 83149b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 83249b5e25fSSatish Balay low = 0; 83349b5e25fSSatish Balay 83449b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 83549b5e25fSSatish Balay if (in[l] < 0) continue; 83649b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 837*77431f27SBarry Smith if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1); 83849b5e25fSSatish Balay #endif 83949b5e25fSSatish Balay col = in[l]; 84049b5e25fSSatish Balay bcol = col/bs; /* block col number */ 84149b5e25fSSatish Balay 84249b5e25fSSatish Balay if (brow <= bcol){ 84349b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8448549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 84549b5e25fSSatish Balay /* element value a(k,l) */ 84649b5e25fSSatish Balay if (roworiented) { 84749b5e25fSSatish Balay value = v[l + k*n]; 84849b5e25fSSatish Balay } else { 84949b5e25fSSatish Balay value = v[k + l*m]; 85049b5e25fSSatish Balay } 85149b5e25fSSatish Balay 85249b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 85349b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 85449b5e25fSSatish Balay while (high-low > 7) { 85549b5e25fSSatish Balay t = (low+high)/2; 85649b5e25fSSatish Balay if (rp[t] > bcol) high = t; 85749b5e25fSSatish Balay else low = t; 85849b5e25fSSatish Balay } 85949b5e25fSSatish Balay for (i=low; i<high; i++) { 860*77431f27SBarry Smith /* printf("The loop of i=low.., rp[%D]=%D\n",i,rp[i]); */ 86149b5e25fSSatish Balay if (rp[i] > bcol) break; 86249b5e25fSSatish Balay if (rp[i] == bcol) { 86349b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 86449b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 86549b5e25fSSatish Balay else *bap = value; 8668549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8678549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8688549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8698549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8708549e402SHong Zhang else *bap = value; 8718549e402SHong Zhang } 87249b5e25fSSatish Balay goto noinsert1; 87349b5e25fSSatish Balay } 87449b5e25fSSatish Balay } 87549b5e25fSSatish Balay 87649b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 877*77431f27SBarry Smith else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 87849b5e25fSSatish Balay if (nrow >= rmax) { 87949b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 88049b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 88149b5e25fSSatish Balay MatScalar *new_a; 88249b5e25fSSatish Balay 883*77431f27SBarry Smith if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 88449b5e25fSSatish Balay 88549b5e25fSSatish Balay /* Malloc new storage space */ 88649b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 88782502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 88849b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 88949b5e25fSSatish Balay new_i = new_j + new_nz; 89049b5e25fSSatish Balay 89149b5e25fSSatish Balay /* copy over old data into new slots */ 89249b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 89349b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 89449b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 89549b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 89649b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 89749b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 89849b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 89949b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 90049b5e25fSSatish Balay /* free up old matrix storage */ 90149b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 90249b5e25fSSatish Balay if (!a->singlemalloc) { 90349b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 90449b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 90549b5e25fSSatish Balay } 90649b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 90749b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 90849b5e25fSSatish Balay 90949b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 91049b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 911b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 9126c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 91349b5e25fSSatish Balay a->reallocs++; 9146c6c5352SBarry Smith a->nz++; 91549b5e25fSSatish Balay } 91649b5e25fSSatish Balay 91749b5e25fSSatish Balay N = nrow++ - 1; 91849b5e25fSSatish Balay /* shift up all the later entries in this row */ 91949b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 92049b5e25fSSatish Balay rp[ii+1] = rp[ii]; 92149b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 92249b5e25fSSatish Balay } 92349b5e25fSSatish Balay if (N>=i) { 92449b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 92549b5e25fSSatish Balay } 92649b5e25fSSatish Balay rp[i] = bcol; 92749b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 92849b5e25fSSatish Balay noinsert1:; 92949b5e25fSSatish Balay low = i; 93049b5e25fSSatish Balay /* } */ 9318549e402SHong Zhang } 93249b5e25fSSatish Balay } /* end of if .. if.. */ 93349b5e25fSSatish Balay } /* end of loop over added columns */ 93449b5e25fSSatish Balay ailen[brow] = nrow; 93549b5e25fSSatish Balay } /* end of loop over added rows */ 93649b5e25fSSatish Balay 93749b5e25fSSatish Balay PetscFunctionReturn(0); 93849b5e25fSSatish Balay } 93949b5e25fSSatish Balay 9406849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 9416849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 9426849ba73SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int); 9436849ba73SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 9446849ba73SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]); 9456849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat); 9466849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 9476849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 9486849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec); 9496849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 9506849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 9516849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat); 9526849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec); 95349b5e25fSSatish Balay 9546849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 9556849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 9566849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 9576849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 9586849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 9596849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 9606849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 9616849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 96249b5e25fSSatish Balay 9636849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 964d59c15a7SBarry Smith 9656849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 9666849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 9676849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 9686849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 9696849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 9706849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 9716849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 9726849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 97307f98182SSatish Balay 9746849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9756849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9766849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9776849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9786849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9796849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9806849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 9816849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 9826849ba73SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*); 98349b5e25fSSatish Balay 9846849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 9856849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 9866849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 9876849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 9886849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 9896849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 9906849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 9916849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 99249b5e25fSSatish Balay 9936849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 9946849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 9956849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 9966849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 9976849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 9986849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 9996849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 10006849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 100149b5e25fSSatish Balay 10024d101231SSatish Balay #ifdef HAVE_MatICCFactor 10036849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work! */ 10044a2ae208SSatish Balay #undef __FUNCT__ 10054d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 1006dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 100749b5e25fSSatish Balay { 10084ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 100949b5e25fSSatish Balay Mat outA; 1010dfbe8321SBarry Smith PetscErrorCode ierr; 101149b5e25fSSatish Balay PetscTruth row_identity,col_identity; 101249b5e25fSSatish Balay 101349b5e25fSSatish Balay PetscFunctionBegin; 101449b5e25fSSatish Balay outA = inA; 10151260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 101649b5e25fSSatish Balay 101749b5e25fSSatish Balay if (!a->diag) { 10181a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 101949b5e25fSSatish Balay } 102049b5e25fSSatish Balay /* 102149b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 102249b5e25fSSatish Balay for ILU(0) factorization with natural ordering 102349b5e25fSSatish Balay */ 102449b5e25fSSatish Balay switch (a->bs) { 102549b5e25fSSatish Balay case 1: 10260fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 102707f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1028d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 10294d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 103049b5e25fSSatish Balay case 2: 10311a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10321a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10331a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10344d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 103549b5e25fSSatish Balay break; 103649b5e25fSSatish Balay case 3: 10371a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10381a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10391a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10404d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 104149b5e25fSSatish Balay break; 104249b5e25fSSatish Balay case 4: 10431a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10441a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10451a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10464d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 104749b5e25fSSatish Balay break; 104849b5e25fSSatish Balay case 5: 10491a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10501a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10511a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10524d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 105349b5e25fSSatish Balay break; 105449b5e25fSSatish Balay case 6: 10551a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10561a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10571a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10584d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 105949b5e25fSSatish Balay break; 106049b5e25fSSatish Balay case 7: 10611a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10621a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10631a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10644d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 106549b5e25fSSatish Balay break; 106649b5e25fSSatish Balay default: 106749b5e25fSSatish Balay a->row = row; 10681a3463dfSHong Zhang a->icol = col; 106949b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 107049b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 107149b5e25fSSatish Balay 107249b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 107349b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1074b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 107549b5e25fSSatish Balay 107649b5e25fSSatish Balay if (!a->solve_work) { 107787828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 107887828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 107949b5e25fSSatish Balay } 108049b5e25fSSatish Balay } 108149b5e25fSSatish Balay 10821a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 108349b5e25fSSatish Balay 108449b5e25fSSatish Balay PetscFunctionReturn(0); 108549b5e25fSSatish Balay } 1086950f1e5bSHong Zhang #endif 1087950f1e5bSHong Zhang 10884a2ae208SSatish Balay #undef __FUNCT__ 10894a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 1090dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A) 109149b5e25fSSatish Balay { 109249b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 109349b5e25fSSatish Balay MPI_Comm comm = A->comm; 1094dfbe8321SBarry Smith PetscErrorCode ierr; 109549b5e25fSSatish Balay 109649b5e25fSSatish Balay PetscFunctionBegin; 109749b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 109849b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 109949b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 110049b5e25fSSatish Balay PetscFunctionReturn(0); 110149b5e25fSSatish Balay } 110249b5e25fSSatish Balay 110349b5e25fSSatish Balay EXTERN_C_BEGIN 11044a2ae208SSatish Balay #undef __FUNCT__ 11054a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1106dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 110749b5e25fSSatish Balay { 1108045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 110949b5e25fSSatish Balay int i,nz,n; 111049b5e25fSSatish Balay 111149b5e25fSSatish Balay PetscFunctionBegin; 11126c6c5352SBarry Smith nz = baij->maxnz; 1113c464158bSHong Zhang n = mat->n; 111449b5e25fSSatish Balay for (i=0; i<nz; i++) { 111549b5e25fSSatish Balay baij->j[i] = indices[i]; 111649b5e25fSSatish Balay } 11176c6c5352SBarry Smith baij->nz = nz; 111849b5e25fSSatish Balay for (i=0; i<n; i++) { 111949b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 112049b5e25fSSatish Balay } 112149b5e25fSSatish Balay 112249b5e25fSSatish Balay PetscFunctionReturn(0); 112349b5e25fSSatish Balay } 112449b5e25fSSatish Balay EXTERN_C_END 112549b5e25fSSatish Balay 11264a2ae208SSatish Balay #undef __FUNCT__ 11274a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 112849b5e25fSSatish Balay /*@ 112919585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 113049b5e25fSSatish Balay in the matrix. 113149b5e25fSSatish Balay 113249b5e25fSSatish Balay Input Parameters: 113319585528SSatish Balay + mat - the SeqSBAIJ matrix 113449b5e25fSSatish Balay - indices - the column indices 113549b5e25fSSatish Balay 113649b5e25fSSatish Balay Level: advanced 113749b5e25fSSatish Balay 113849b5e25fSSatish Balay Notes: 113949b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 114049b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 114149b5e25fSSatish Balay of the MatSetValues() operation. 114249b5e25fSSatish Balay 114349b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 114419585528SSatish Balay MatCreateSeqSBAIJ(). 114549b5e25fSSatish Balay 1146ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 114749b5e25fSSatish Balay 1148ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 114949b5e25fSSatish Balay @*/ 1150dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 115149b5e25fSSatish Balay { 1152dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,int *); 115349b5e25fSSatish Balay 115449b5e25fSSatish Balay PetscFunctionBegin; 11554482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 11564482741eSBarry Smith PetscValidPointer(indices,2); 1157c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 115849b5e25fSSatish Balay if (f) { 115949b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 116049b5e25fSSatish Balay } else { 1161e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 116249b5e25fSSatish Balay } 116349b5e25fSSatish Balay PetscFunctionReturn(0); 116449b5e25fSSatish Balay } 116549b5e25fSSatish Balay 11664a2ae208SSatish Balay #undef __FUNCT__ 11674a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1168dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1169273d9f13SBarry Smith { 1170dfbe8321SBarry Smith PetscErrorCode ierr; 1171273d9f13SBarry Smith 1172273d9f13SBarry Smith PetscFunctionBegin; 1173273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1174273d9f13SBarry Smith PetscFunctionReturn(0); 1175273d9f13SBarry Smith } 1176273d9f13SBarry Smith 1177a6ece127SHong Zhang #undef __FUNCT__ 1178a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1179dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1180a6ece127SHong Zhang { 1181a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1182a6ece127SHong Zhang PetscFunctionBegin; 1183a6ece127SHong Zhang *array = a->a; 1184a6ece127SHong Zhang PetscFunctionReturn(0); 1185a6ece127SHong Zhang } 1186a6ece127SHong Zhang 1187a6ece127SHong Zhang #undef __FUNCT__ 1188a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1189dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1190a6ece127SHong Zhang { 1191a6ece127SHong Zhang PetscFunctionBegin; 1192a6ece127SHong Zhang PetscFunctionReturn(0); 1193a6ece127SHong Zhang } 1194a6ece127SHong Zhang 119542ee4b1aSHong Zhang #include "petscblaslapack.h" 119642ee4b1aSHong Zhang #undef __FUNCT__ 119742ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1198dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 119942ee4b1aSHong Zhang { 120042ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1201dfbe8321SBarry Smith PetscErrorCode ierr; 12024ce68768SBarry Smith int i,bs=y->bs,bs2,j; 12034ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 120442ee4b1aSHong Zhang 120542ee4b1aSHong Zhang PetscFunctionBegin; 120642ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 12074ce68768SBarry Smith BLaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one); 1208c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1209c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1210c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1211c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1212c537a176SHong Zhang } 1213c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1214c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1215c4319e64SHong Zhang y->XtoY = X; 1216c537a176SHong Zhang } 1217c4319e64SHong Zhang bs2 = bs*bs; 12186c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1219c4319e64SHong Zhang j = 0; 1220c4319e64SHong Zhang while (j < bs2){ 12216550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1222c4319e64SHong Zhang j++; 1223c537a176SHong Zhang } 1224c4319e64SHong Zhang } 1225*77431f27SBarry Smith PetscLogInfo(0,"MatAXPY_SeqSBAIJ: 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)); 122642ee4b1aSHong Zhang } else { 122742ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 122842ee4b1aSHong Zhang } 122942ee4b1aSHong Zhang PetscFunctionReturn(0); 123042ee4b1aSHong Zhang } 123142ee4b1aSHong Zhang 1232efcf0fc3SBarry Smith #undef __FUNCT__ 1233efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1234dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1235efcf0fc3SBarry Smith { 1236efcf0fc3SBarry Smith PetscFunctionBegin; 1237efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1238efcf0fc3SBarry Smith PetscFunctionReturn(0); 1239efcf0fc3SBarry Smith } 1240efcf0fc3SBarry Smith 1241efcf0fc3SBarry Smith #undef __FUNCT__ 1242efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1243dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1244efcf0fc3SBarry Smith { 1245efcf0fc3SBarry Smith PetscFunctionBegin; 1246efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1247efcf0fc3SBarry Smith PetscFunctionReturn(0); 1248efcf0fc3SBarry Smith } 1249efcf0fc3SBarry Smith 1250efcf0fc3SBarry Smith #undef __FUNCT__ 1251efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1252dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1253efcf0fc3SBarry Smith { 1254efcf0fc3SBarry Smith PetscFunctionBegin; 1255efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1256efcf0fc3SBarry Smith PetscFunctionReturn(0); 1257efcf0fc3SBarry Smith } 1258efcf0fc3SBarry Smith 125949b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 126049b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 126149b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 126249b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 126349b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 126497304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1265e005ede5SBarry Smith MatMult_SeqSBAIJ_N, 1266e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 126749b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 126849b5e25fSSatish Balay 0, 126949b5e25fSSatish Balay 0, 127097304618SKris Buschelman /*10*/ 0, 127149b5e25fSSatish Balay 0, 12725f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1273d06b337dSHong Zhang MatRelax_SeqSBAIJ, 127449b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 127597304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 127649b5e25fSSatish Balay MatEqual_SeqSBAIJ, 127749b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 127849b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 127949b5e25fSSatish Balay MatNorm_SeqSBAIJ, 128097304618SKris Buschelman /*20*/ 0, 128149b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 128249b5e25fSSatish Balay 0, 128349b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 128449b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 128597304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ, 128649b5e25fSSatish Balay 0, 128749b5e25fSSatish Balay 0, 12885f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12895f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 129097304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1291c464158bSHong Zhang 0, 12924d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1293a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1294a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 129597304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 129649b5e25fSSatish Balay 0, 129749b5e25fSSatish Balay 0, 129849b5e25fSSatish Balay 0, 1299950f1e5bSHong Zhang 0, 130097304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 130149b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 130249b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 130349b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 130449b5e25fSSatish Balay 0, 130597304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 130649b5e25fSSatish Balay MatScale_SeqSBAIJ, 130749b5e25fSSatish Balay 0, 130849b5e25fSSatish Balay 0, 130949b5e25fSSatish Balay 0, 131097304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ, 131149b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 131249b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 131349b5e25fSSatish Balay 0, 131449b5e25fSSatish Balay 0, 131597304618SKris Buschelman /*55*/ 0, 131649b5e25fSSatish Balay 0, 131749b5e25fSSatish Balay 0, 131849b5e25fSSatish Balay 0, 131949b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 132097304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 132149b5e25fSSatish Balay 0, 132249b5e25fSSatish Balay 0, 13238a124369SBarry Smith MatGetPetscMaps_Petsc, 1324d959ec07SHong Zhang 0, 132597304618SKris Buschelman /*65*/ 0, 1326d959ec07SHong Zhang 0, 1327d959ec07SHong Zhang 0, 1328d959ec07SHong Zhang 0, 1329d959ec07SHong Zhang 0, 133097304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13313e0d88b5SBarry Smith 0, 13323e0d88b5SBarry Smith 0, 13333e0d88b5SBarry Smith 0, 13343e0d88b5SBarry Smith 0, 133597304618SKris Buschelman /*75*/ 0, 13363e0d88b5SBarry Smith 0, 13373e0d88b5SBarry Smith 0, 13383e0d88b5SBarry Smith 0, 13393e0d88b5SBarry Smith 0, 134097304618SKris Buschelman /*80*/ 0, 13413e0d88b5SBarry Smith 0, 13423e0d88b5SBarry Smith 0, 13433e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 134497304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13453e0d88b5SBarry Smith #else 134697304618SKris Buschelman 0, 13473e0d88b5SBarry Smith #endif 1348865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1349865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1350865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1351efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1352865e5f61SKris Buschelman 0, 1353865e5f61SKris Buschelman 0, 1354865e5f61SKris Buschelman /*90*/ 0, 1355865e5f61SKris Buschelman 0, 1356865e5f61SKris Buschelman 0, 1357865e5f61SKris Buschelman 0, 1358865e5f61SKris Buschelman 0, 1359865e5f61SKris Buschelman /*95*/ 0, 1360865e5f61SKris Buschelman 0, 1361865e5f61SKris Buschelman 0, 1362865e5f61SKris Buschelman 0}; 136349b5e25fSSatish Balay 136449b5e25fSSatish Balay EXTERN_C_BEGIN 13654a2ae208SSatish Balay #undef __FUNCT__ 13664a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1367dfbe8321SBarry Smith PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat) 136849b5e25fSSatish Balay { 13694afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1370c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 1371dfbe8321SBarry Smith PetscErrorCode ierr; 137249b5e25fSSatish Balay 137349b5e25fSSatish Balay PetscFunctionBegin; 137449b5e25fSSatish Balay if (aij->nonew != 1) { 1375e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 137649b5e25fSSatish Balay } 137749b5e25fSSatish Balay 137849b5e25fSSatish Balay /* allocate space for values if not already there */ 137949b5e25fSSatish Balay if (!aij->saved_values) { 138087828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 138149b5e25fSSatish Balay } 138249b5e25fSSatish Balay 138349b5e25fSSatish Balay /* copy values over */ 138487828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 138549b5e25fSSatish Balay PetscFunctionReturn(0); 138649b5e25fSSatish Balay } 138749b5e25fSSatish Balay EXTERN_C_END 138849b5e25fSSatish Balay 138949b5e25fSSatish Balay EXTERN_C_BEGIN 13904a2ae208SSatish Balay #undef __FUNCT__ 13914a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1392dfbe8321SBarry Smith PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat) 139349b5e25fSSatish Balay { 13944afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13956849ba73SBarry Smith PetscErrorCode ierr; 13966849ba73SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2; 139749b5e25fSSatish Balay 139849b5e25fSSatish Balay PetscFunctionBegin; 139949b5e25fSSatish Balay if (aij->nonew != 1) { 1400e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 140149b5e25fSSatish Balay } 140249b5e25fSSatish Balay if (!aij->saved_values) { 1403e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 140449b5e25fSSatish Balay } 140549b5e25fSSatish Balay 140649b5e25fSSatish Balay /* copy values over */ 140787828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 140849b5e25fSSatish Balay PetscFunctionReturn(0); 140949b5e25fSSatish Balay } 141049b5e25fSSatish Balay EXTERN_C_END 141149b5e25fSSatish Balay 14128549e402SHong Zhang EXTERN_C_BEGIN 14134a2ae208SSatish Balay #undef __FUNCT__ 1414a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1415dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz) 141649b5e25fSSatish Balay { 1417c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14186849ba73SBarry Smith PetscErrorCode ierr; 14196849ba73SBarry Smith int i,len,mbs,bs2; 142049b5e25fSSatish Balay PetscTruth flg; 142149b5e25fSSatish Balay 142249b5e25fSSatish Balay PetscFunctionBegin; 1423273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1424e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1425c464158bSHong Zhang mbs = B->m/bs; 142649b5e25fSSatish Balay bs2 = bs*bs; 142749b5e25fSSatish Balay 1428c464158bSHong Zhang if (mbs*bs != B->m) { 142929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 143049b5e25fSSatish Balay } 143149b5e25fSSatish Balay 1432435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1433*77431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 143449b5e25fSSatish Balay if (nnz) { 143549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 1436*77431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 1437*77431f27SBarry 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); 143849b5e25fSSatish Balay } 143949b5e25fSSatish Balay } 144049b5e25fSSatish Balay 1441e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 144249b5e25fSSatish Balay if (!flg) { 144349b5e25fSSatish Balay switch (bs) { 144449b5e25fSSatish Balay case 1: 14454ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 144649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1447d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1448e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 144949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 145049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 145149b5e25fSSatish Balay break; 145249b5e25fSSatish Balay case 2: 14534ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 145449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1455e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 145649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 145749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 145849b5e25fSSatish Balay break; 145949b5e25fSSatish Balay case 3: 14605f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 146149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1462e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 146349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 146449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 146549b5e25fSSatish Balay break; 146649b5e25fSSatish Balay case 4: 14675f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 146849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1469e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 147049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 147149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 147249b5e25fSSatish Balay break; 147349b5e25fSSatish Balay case 5: 14745f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 147549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1476e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 147749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 147849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 147949b5e25fSSatish Balay break; 148049b5e25fSSatish Balay case 6: 14815f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 148249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1483e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 148449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 148549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 148649b5e25fSSatish Balay break; 148749b5e25fSSatish Balay case 7: 1488de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 148949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1490e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1491de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 149249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 149349b5e25fSSatish Balay break; 149449b5e25fSSatish Balay } 149549b5e25fSSatish Balay } 149649b5e25fSSatish Balay 149749b5e25fSSatish Balay b->mbs = mbs; 14984afc71dfSHong Zhang b->nbs = mbs; 1499b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 150049b5e25fSSatish Balay if (!nnz) { 1501435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 150249b5e25fSSatish Balay else if (nz <= 0) nz = 1; 150349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15048cef66ccSHong Zhang b->imax[i] = nz; 150549b5e25fSSatish Balay } 1506153ea458SHong Zhang nz = nz*mbs; /* total nz */ 150749b5e25fSSatish Balay } else { 150849b5e25fSSatish Balay nz = 0; 15098cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 151049b5e25fSSatish Balay } 15116c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 151249b5e25fSSatish Balay 151349b5e25fSSatish Balay /* allocate the matrix space */ 15146c6c5352SBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 151582502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 15166c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15176c6c5352SBarry Smith b->j = (int*)(b->a + nz*bs2); 15186c6c5352SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15196c6c5352SBarry Smith b->i = b->j + nz; 152049b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 152149b5e25fSSatish Balay 152249b5e25fSSatish Balay /* pointer to beginning of each row */ 152349b5e25fSSatish Balay b->i[0] = 0; 152449b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 152549b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 152649b5e25fSSatish Balay } 152749b5e25fSSatish Balay 152849b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 15295033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1530b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 153149b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 153249b5e25fSSatish Balay 153349b5e25fSSatish Balay b->bs = bs; 153449b5e25fSSatish Balay b->bs2 = bs2; 15356c6c5352SBarry Smith b->nz = 0; 15366c6c5352SBarry Smith b->maxnz = nz*bs2; 1537153ea458SHong Zhang 153816cdd363SHong Zhang b->inew = 0; 153916cdd363SHong Zhang b->jnew = 0; 154016cdd363SHong Zhang b->anew = 0; 154116cdd363SHong Zhang b->a2anew = 0; 15421a3463dfSHong Zhang b->permute = PETSC_FALSE; 1543c464158bSHong Zhang PetscFunctionReturn(0); 1544c464158bSHong Zhang } 1545a23d5eceSKris Buschelman EXTERN_C_END 1546153ea458SHong Zhang 1547d769727bSBarry Smith EXTERN_C_BEGIN 1548dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,Mat*); 1549dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,Mat*); 1550d769727bSBarry Smith EXTERN_C_END 1551d769727bSBarry Smith 15520bad9183SKris Buschelman /*MC 1553fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15540bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15550bad9183SKris Buschelman 15560bad9183SKris Buschelman Options Database Keys: 15570bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15580bad9183SKris Buschelman 15590bad9183SKris Buschelman Level: beginner 15600bad9183SKris Buschelman 15610bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15620bad9183SKris Buschelman M*/ 15630bad9183SKris Buschelman 1564a23d5eceSKris Buschelman EXTERN_C_BEGIN 1565a23d5eceSKris Buschelman #undef __FUNCT__ 1566a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1567dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqSBAIJ(Mat B) 1568a23d5eceSKris Buschelman { 1569a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1570dfbe8321SBarry Smith PetscErrorCode ierr; 1571dfbe8321SBarry Smith int size; 1572a23d5eceSKris Buschelman 1573a23d5eceSKris Buschelman PetscFunctionBegin; 1574a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1575a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1576a23d5eceSKris Buschelman B->m = B->M = PetscMax(B->m,B->M); 1577a23d5eceSKris Buschelman B->n = B->N = PetscMax(B->n,B->N); 1578a23d5eceSKris Buschelman 1579a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1580a23d5eceSKris Buschelman B->data = (void*)b; 1581a23d5eceSKris Buschelman ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1582a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1583a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1584a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1585a23d5eceSKris Buschelman B->factor = 0; 1586a23d5eceSKris Buschelman B->lupivotthreshold = 1.0; 1587a23d5eceSKris Buschelman B->mapping = 0; 1588a23d5eceSKris Buschelman b->row = 0; 1589a23d5eceSKris Buschelman b->icol = 0; 1590a23d5eceSKris Buschelman b->reallocs = 0; 1591a23d5eceSKris Buschelman b->saved_values = 0; 1592a23d5eceSKris Buschelman 1593a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1594a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1595a23d5eceSKris Buschelman 1596a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1597a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1598a23d5eceSKris Buschelman b->nonew = 0; 1599a23d5eceSKris Buschelman b->diag = 0; 1600a23d5eceSKris Buschelman b->solve_work = 0; 1601a23d5eceSKris Buschelman b->mult_work = 0; 1602a23d5eceSKris Buschelman B->spptr = 0; 1603a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1604a23d5eceSKris Buschelman b->xtoy = 0; 1605a23d5eceSKris Buschelman b->XtoY = 0; 1606a23d5eceSKris Buschelman 1607a23d5eceSKris Buschelman b->inew = 0; 1608a23d5eceSKris Buschelman b->jnew = 0; 1609a23d5eceSKris Buschelman b->anew = 0; 1610a23d5eceSKris Buschelman b->a2anew = 0; 1611a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1612a23d5eceSKris Buschelman 1613a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1614a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1615a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1616a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1617a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1618a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1619a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1620a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1621a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16224e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16234e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16244e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1625a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1626a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1627a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1628a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1629a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1630a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 163123ce1328SBarry Smith 163223ce1328SBarry Smith B->symmetric = PETSC_TRUE; 163323ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 163423ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 163523ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 1636a23d5eceSKris Buschelman PetscFunctionReturn(0); 1637a23d5eceSKris Buschelman } 1638a23d5eceSKris Buschelman EXTERN_C_END 1639a23d5eceSKris Buschelman 1640a23d5eceSKris Buschelman #undef __FUNCT__ 1641a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1642a23d5eceSKris Buschelman /*@C 1643a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1644a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1645a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1646a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1647a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1648a23d5eceSKris Buschelman 1649a23d5eceSKris Buschelman Collective on Mat 1650a23d5eceSKris Buschelman 1651a23d5eceSKris Buschelman Input Parameters: 1652a23d5eceSKris Buschelman + A - the symmetric matrix 1653a23d5eceSKris Buschelman . bs - size of block 1654a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1655a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1656a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1657a23d5eceSKris Buschelman 1658a23d5eceSKris Buschelman Options Database Keys: 1659a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1660a23d5eceSKris Buschelman block calculations (much slower) 1661a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1662a23d5eceSKris Buschelman 1663a23d5eceSKris Buschelman Level: intermediate 1664a23d5eceSKris Buschelman 1665a23d5eceSKris Buschelman Notes: 1666a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1667a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1668a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1669a23d5eceSKris Buschelman matrices. 1670a23d5eceSKris Buschelman 1671a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1672a23d5eceSKris Buschelman @*/ 1673dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) { 1674dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,int,int,const int[]); 1675a23d5eceSKris Buschelman 1676a23d5eceSKris Buschelman PetscFunctionBegin; 1677a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1678a23d5eceSKris Buschelman if (f) { 1679a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1680a23d5eceSKris Buschelman } 1681a23d5eceSKris Buschelman PetscFunctionReturn(0); 1682a23d5eceSKris Buschelman } 168349b5e25fSSatish Balay 16844a2ae208SSatish Balay #undef __FUNCT__ 16854a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1686c464158bSHong Zhang /*@C 1687c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1688c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1689c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1690c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1691c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 169249b5e25fSSatish Balay 1693c464158bSHong Zhang Collective on MPI_Comm 1694c464158bSHong Zhang 1695c464158bSHong Zhang Input Parameters: 1696c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1697c464158bSHong Zhang . bs - size of block 1698c464158bSHong Zhang . m - number of rows, or number of columns 1699c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1700744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1701744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1702c464158bSHong Zhang 1703c464158bSHong Zhang Output Parameter: 1704c464158bSHong Zhang . A - the symmetric matrix 1705c464158bSHong Zhang 1706c464158bSHong Zhang Options Database Keys: 1707c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1708c464158bSHong Zhang block calculations (much slower) 1709c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1710c464158bSHong Zhang 1711c464158bSHong Zhang Level: intermediate 1712c464158bSHong Zhang 1713c464158bSHong Zhang Notes: 1714c464158bSHong Zhang 1715c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1716c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1717c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1718c464158bSHong Zhang matrices. 1719c464158bSHong Zhang 1720c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1721c464158bSHong Zhang @*/ 1722dfbe8321SBarry Smith PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A) 1723c464158bSHong Zhang { 1724dfbe8321SBarry Smith PetscErrorCode ierr; 1725c464158bSHong Zhang 1726c464158bSHong Zhang PetscFunctionBegin; 1727c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1728c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1729273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 173049b5e25fSSatish Balay PetscFunctionReturn(0); 173149b5e25fSSatish Balay } 173249b5e25fSSatish Balay 17334a2ae208SSatish Balay #undef __FUNCT__ 17344a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1735dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 173649b5e25fSSatish Balay { 173749b5e25fSSatish Balay Mat C; 173849b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17396849ba73SBarry Smith PetscErrorCode ierr; 17406849ba73SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 174149b5e25fSSatish Balay 174249b5e25fSSatish Balay PetscFunctionBegin; 174329bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 174449b5e25fSSatish Balay 174549b5e25fSSatish Balay *B = 0; 1746692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1747be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17481d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1749692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1750692f9cbeSHong Zhang 1751273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 175249b5e25fSSatish Balay C->factor = A->factor; 175349b5e25fSSatish Balay c->row = 0; 175449b5e25fSSatish Balay c->icol = 0; 175549b5e25fSSatish Balay c->saved_values = 0; 175649b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 175749b5e25fSSatish Balay C->assembled = PETSC_TRUE; 175849b5e25fSSatish Balay 175982327fa8SHong Zhang C->M = A->M; 176082327fa8SHong Zhang C->N = A->N; 176149b5e25fSSatish Balay c->bs = a->bs; 176249b5e25fSSatish Balay c->bs2 = a->bs2; 176349b5e25fSSatish Balay c->mbs = a->mbs; 176449b5e25fSSatish Balay c->nbs = a->nbs; 176549b5e25fSSatish Balay 1766b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1767b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 176849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 176949b5e25fSSatish Balay c->imax[i] = a->imax[i]; 177049b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 177149b5e25fSSatish Balay } 177249b5e25fSSatish Balay 177349b5e25fSSatish Balay /* allocate the matrix space */ 177449b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 177549b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 177682502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 177749b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 177849b5e25fSSatish Balay c->i = c->j + nz; 177949b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 178049b5e25fSSatish Balay if (mbs > 0) { 178149b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 178249b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 178349b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 178449b5e25fSSatish Balay } else { 178549b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 178649b5e25fSSatish Balay } 178749b5e25fSSatish Balay } 178849b5e25fSSatish Balay 1789b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 179049b5e25fSSatish Balay c->sorted = a->sorted; 179149b5e25fSSatish Balay c->roworiented = a->roworiented; 179249b5e25fSSatish Balay c->nonew = a->nonew; 179349b5e25fSSatish Balay 179449b5e25fSSatish Balay if (a->diag) { 1795b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1796b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 179749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 179849b5e25fSSatish Balay c->diag[i] = a->diag[i]; 179949b5e25fSSatish Balay } 180049b5e25fSSatish Balay } else c->diag = 0; 18016c6c5352SBarry Smith c->nz = a->nz; 18026c6c5352SBarry Smith c->maxnz = a->maxnz; 180349b5e25fSSatish Balay c->solve_work = 0; 180449b5e25fSSatish Balay c->mult_work = 0; 180549b5e25fSSatish Balay *B = C; 1806b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 180749b5e25fSSatish Balay PetscFunctionReturn(0); 180849b5e25fSSatish Balay } 180949b5e25fSSatish Balay 18104a2ae208SSatish Balay #undef __FUNCT__ 18114a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1812dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A) 181349b5e25fSSatish Balay { 181449b5e25fSSatish Balay Mat_SeqSBAIJ *a; 181549b5e25fSSatish Balay Mat B; 18166849ba73SBarry Smith PetscErrorCode ierr; 18176849ba73SBarry Smith int i,nz,fd,header[4],size,*rowlengths=0,M,N,bs=1; 18183f7326a9SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 181949b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 182049b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 182187828ca2SBarry Smith PetscScalar *aa; 182249b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 182349b5e25fSSatish Balay 182449b5e25fSSatish Balay PetscFunctionBegin; 1825b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 182649b5e25fSSatish Balay bs2 = bs*bs; 182749b5e25fSSatish Balay 182849b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 182929bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1830b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 183149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1832552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 183349b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 183449b5e25fSSatish Balay 183549b5e25fSSatish Balay if (header[3] < 0) { 183629bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 183749b5e25fSSatish Balay } 183849b5e25fSSatish Balay 183929bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 184049b5e25fSSatish Balay 184149b5e25fSSatish Balay /* 184249b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 184349b5e25fSSatish Balay divisible by the blocksize 184449b5e25fSSatish Balay */ 184549b5e25fSSatish Balay mbs = M/bs; 184649b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 184749b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 184849b5e25fSSatish Balay else mbs++; 184949b5e25fSSatish Balay if (extra_rows) { 1850b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 185149b5e25fSSatish Balay } 185249b5e25fSSatish Balay 185349b5e25fSSatish Balay /* read in row lengths */ 1854b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 185549b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 185649b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 185749b5e25fSSatish Balay 185849b5e25fSSatish Balay /* read in column indices */ 1859b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 186049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 186149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 186249b5e25fSSatish Balay 186349b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 186482502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 1865a91a614fSBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 186682502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 186749b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 186849b5e25fSSatish Balay masked = mask + mbs; 186949b5e25fSSatish Balay rowcount = 0; nzcount = 0; 187049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 187149b5e25fSSatish Balay nmask = 0; 187249b5e25fSSatish Balay for (j=0; j<bs; j++) { 187349b5e25fSSatish Balay kmax = rowlengths[rowcount]; 187449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18752d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 187603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 187749b5e25fSSatish Balay } 187849b5e25fSSatish Balay rowcount++; 187949b5e25fSSatish Balay } 1880574b2666SHong Zhang s_browlengths[i] += nmask; 1881574b2666SHong Zhang 188249b5e25fSSatish Balay /* zero out the mask elements we set */ 188349b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 188449b5e25fSSatish Balay } 188549b5e25fSSatish Balay 188649b5e25fSSatish Balay /* create our matrix */ 18879abb65ffSKris Buschelman ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr); 18889abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 188978473fd7SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr); 189049b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 189149b5e25fSSatish Balay 189249b5e25fSSatish Balay /* set matrix "i" values */ 189349b5e25fSSatish Balay a->i[0] = 0; 189449b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1895574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1896574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 189749b5e25fSSatish Balay } 18986c6c5352SBarry Smith a->nz = a->i[mbs]; 189949b5e25fSSatish Balay 190049b5e25fSSatish Balay /* read in nonzero values */ 190187828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 190249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 190349b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 190449b5e25fSSatish Balay 190549b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 190649b5e25fSSatish Balay nzcount = 0; jcount = 0; 190749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 190849b5e25fSSatish Balay nzcountb = nzcount; 190949b5e25fSSatish Balay nmask = 0; 191049b5e25fSSatish Balay for (j=0; j<bs; j++) { 191149b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 191249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19132d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 191403630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19152d703238SHong Zhang } 19162d703238SHong Zhang } 19172d703238SHong Zhang /* sort the masked values */ 19182d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19192d703238SHong Zhang 19202d703238SHong Zhang /* set "j" values into matrix */ 19212d703238SHong Zhang maskcount = 1; 19222d703238SHong Zhang for (j=0; j<nmask; j++) { 192349b5e25fSSatish Balay a->j[jcount++] = masked[j]; 192449b5e25fSSatish Balay mask[masked[j]] = maskcount++; 192549b5e25fSSatish Balay } 1926574b2666SHong Zhang 192749b5e25fSSatish Balay /* set "a" values into matrix */ 192849b5e25fSSatish Balay ishift = bs2*a->i[i]; 192949b5e25fSSatish Balay for (j=0; j<bs; j++) { 193049b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 193149b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1932574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1933574b2666SHong Zhang if (tmp >= i){ 193449b5e25fSSatish Balay block = mask[tmp] - 1; 193549b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 193649b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1937574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1938574b2666SHong Zhang } 1939574b2666SHong Zhang nzcountb++; 194049b5e25fSSatish Balay } 194149b5e25fSSatish Balay } 194249b5e25fSSatish Balay /* zero out the mask elements we set */ 194349b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 194449b5e25fSSatish Balay } 19456c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 194649b5e25fSSatish Balay 194749b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1948574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 194949b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 195049b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 195149b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 195249b5e25fSSatish Balay 19539abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19549abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195549b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19569abb65ffSKris Buschelman *A = B; 195749b5e25fSSatish Balay PetscFunctionReturn(0); 195849b5e25fSSatish Balay } 1959574b2666SHong Zhang 1960d06b337dSHong Zhang #undef __FUNCT__ 1961d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1962dfbe8321SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx) 1963d06b337dSHong Zhang { 1964d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1965d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1966d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19676849ba73SBarry Smith PetscErrorCode ierr; 19686849ba73SBarry Smith int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j; 1969d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1970d06b337dSHong Zhang 1971d06b337dSHong Zhang PetscFunctionBegin; 1972b965ef7fSBarry Smith its = its*lits; 1973*77431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1974d06b337dSHong Zhang 1975d06b337dSHong Zhang if (bs > 1) 1976d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1977d06b337dSHong Zhang 19781ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1979d06b337dSHong Zhang if (xx != bb) { 19801ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1981d06b337dSHong Zhang } else { 1982d06b337dSHong Zhang b = x; 1983d06b337dSHong Zhang } 1984d06b337dSHong Zhang 1985d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1986d06b337dSHong Zhang 1987d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19882798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1989d06b337dSHong Zhang for (i=0; i<m; i++) 1990d06b337dSHong Zhang t[i] = b[i]; 1991d06b337dSHong Zhang 1992d06b337dSHong Zhang for (i=0; i<m; i++){ 199344706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1994d06b337dSHong Zhang v = aa + ai[i] + 1; 1995d06b337dSHong Zhang vj = aj + ai[i] + 1; 1996d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1997d06b337dSHong Zhang x[i] = omega*t[i]/d; 1998d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 199944706875SHong Zhang PetscLogFlops(2*nz-1); 2000d06b337dSHong Zhang } 2001d06b337dSHong Zhang } 2002d06b337dSHong Zhang 20032798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 2004d06b337dSHong Zhang for (i=0; i<m; i++) 2005d06b337dSHong Zhang t[i] = b[i]; 2006d06b337dSHong Zhang 2007d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2008d06b337dSHong Zhang v = aa + ai[i] + 1; 2009d06b337dSHong Zhang vj = aj + ai[i] + 1; 2010d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2011d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 201244706875SHong Zhang PetscLogFlops(2*nz-1); 2013d06b337dSHong Zhang } 2014d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2015d06b337dSHong Zhang d = *(aa + ai[i]); 2016d06b337dSHong Zhang v = aa + ai[i] + 1; 2017d06b337dSHong Zhang vj = aj + ai[i] + 1; 2018d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2019d06b337dSHong Zhang sum = t[i]; 2020d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 202144706875SHong Zhang PetscLogFlops(2*nz-1); 2022d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2023d06b337dSHong Zhang } 2024d06b337dSHong Zhang } 2025d06b337dSHong Zhang its--; 2026d06b337dSHong Zhang } 2027d06b337dSHong Zhang 2028d06b337dSHong Zhang while (its--) { 202944706875SHong Zhang /* 203044706875SHong Zhang forward sweep: 203144706875SHong Zhang for i=0,...,m-1: 203244706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 203344706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 203444706875SHong Zhang b = b - x[i]*U^T(i,:); 2035d06b337dSHong Zhang 203644706875SHong Zhang */ 20372798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2038d06b337dSHong Zhang for (i=0; i<m; i++) 2039d06b337dSHong Zhang t[i] = b[i]; 2040d06b337dSHong Zhang 2041d06b337dSHong Zhang for (i=0; i<m; i++){ 204244706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2043d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2044d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2045d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2046d06b337dSHong Zhang sum = t[i]; 2047d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2048d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2049d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 205044706875SHong Zhang PetscLogFlops(4*nz-2); 2051d06b337dSHong Zhang } 2052d06b337dSHong Zhang } 2053d06b337dSHong Zhang 20542798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 205544706875SHong Zhang /* 205644706875SHong Zhang backward sweep: 205744706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 205844706875SHong Zhang for i=m-1,...,0: 205944706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 206044706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 206144706875SHong Zhang */ 2062d06b337dSHong Zhang for (i=0; i<m; i++) 2063d06b337dSHong Zhang t[i] = b[i]; 2064d06b337dSHong Zhang 2065d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2066d06b337dSHong Zhang v = aa + ai[i] + 1; 2067d06b337dSHong Zhang vj = aj + ai[i] + 1; 2068d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2069d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 207044706875SHong Zhang PetscLogFlops(2*nz-1); 2071d06b337dSHong Zhang } 2072d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2073d06b337dSHong Zhang d = *(aa + ai[i]); 2074d06b337dSHong Zhang v = aa + ai[i] + 1; 2075d06b337dSHong Zhang vj = aj + ai[i] + 1; 2076d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2077d06b337dSHong Zhang sum = t[i]; 2078d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 207944706875SHong Zhang PetscLogFlops(2*nz-1); 2080d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2081d06b337dSHong Zhang } 2082d06b337dSHong Zhang } 2083d06b337dSHong Zhang } 2084d06b337dSHong Zhang 2085d06b337dSHong Zhang ierr = PetscFree(t);CHKERRQ(ierr); 20861ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2087d06b337dSHong Zhang if (bb != xx) { 20881ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2089d06b337dSHong Zhang } 20902798e883SHong Zhang 2091d06b337dSHong Zhang PetscFunctionReturn(0); 2092d06b337dSHong Zhang } 2093d06b337dSHong Zhang 2094d06b337dSHong Zhang 2095d06b337dSHong Zhang 2096d06b337dSHong Zhang 209749b5e25fSSatish Balay 209849b5e25fSSatish Balay 2099