1be1d678aSKris Buschelman #define PETSCMAT_DLL 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 4a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/inline/spops.h" 9aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1049b5e25fSSatish Balay 1149b5e25fSSatish Balay #define CHUNKSIZE 10 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay /* 1449b5e25fSSatish Balay Checks for missing diagonals 1549b5e25fSSatish Balay */ 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 18dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A) 1949b5e25fSSatish Balay { 20045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 216849ba73SBarry Smith PetscErrorCode ierr; 2213f74950SBarry Smith PetscInt *diag,*jj = a->j,i; 2349b5e25fSSatish Balay 2449b5e25fSSatish Balay PetscFunctionBegin; 25045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2649b5e25fSSatish Balay diag = a->diag; 2749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2877431f27SBarry Smith if (jj[diag[i]] != i) SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Matrix is missing diagonal number %D",i); 2949b5e25fSSatish Balay } 3049b5e25fSSatish Balay PetscFunctionReturn(0); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay 334a2ae208SSatish Balay #undef __FUNCT__ 344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 35dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 3649b5e25fSSatish Balay { 37045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 386849ba73SBarry Smith PetscErrorCode ierr; 3913f74950SBarry Smith PetscInt i,mbs = a->mbs; 4049b5e25fSSatish Balay 4149b5e25fSSatish Balay PetscFunctionBegin; 4249b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4349b5e25fSSatish Balay 4413f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 4552e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 46b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4749b5e25fSSatish Balay PetscFunctionReturn(0); 4849b5e25fSSatish Balay } 4949b5e25fSSatish Balay 504a2ae208SSatish Balay #undef __FUNCT__ 514a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 5213f74950SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 5349b5e25fSSatish Balay { 54a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5513f74950SBarry Smith PetscInt n = a->mbs,i; 5649b5e25fSSatish Balay 5749b5e25fSSatish Balay PetscFunctionBegin; 58d3e5a4abSHong Zhang *nn = n; 59a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 60a1373b80SHong Zhang 61a6ece127SHong Zhang if (oshift == 1) { 62a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 6313f74950SBarry Smith PetscInt nz = a->i[n]; 646c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 65a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 66a1373b80SHong Zhang *ia = a->i; *ja = a->j; 67a1373b80SHong Zhang } else { 68a1373b80SHong Zhang *ia = a->i; *ja = a->j; 69a6ece127SHong Zhang } 7049b5e25fSSatish Balay PetscFunctionReturn(0); 7149b5e25fSSatish Balay } 7249b5e25fSSatish Balay 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 7513f74950SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 7649b5e25fSSatish Balay { 77b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7813f74950SBarry Smith PetscInt i,n = a->mbs; 79a6ece127SHong Zhang 8049b5e25fSSatish Balay PetscFunctionBegin; 8149b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 82a6ece127SHong Zhang 83a6ece127SHong Zhang if (oshift == 1) { 8413f74950SBarry Smith PetscInt nz = a->i[n]-1; 856c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 86a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 87a6ece127SHong Zhang } 88a6ece127SHong Zhang PetscFunctionReturn(0); 8949b5e25fSSatish Balay } 9049b5e25fSSatish Balay 914a2ae208SSatish Balay #undef __FUNCT__ 924a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 93dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 9449b5e25fSSatish Balay { 9549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 96dfbe8321SBarry Smith PetscErrorCode ierr; 9749b5e25fSSatish Balay 9849b5e25fSSatish Balay PetscFunctionBegin; 9977431f27SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->m,a->nz); 10049b5e25fSSatish Balay if (!a->singlemalloc) { 101*a96a251dSBarry Smith ierr = PetscFree(a->a);CHKERRQ(ierr); 10249b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 10363c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 104*a96a251dSBarry Smith } else { 105*a96a251dSBarry Smith ierr = PetscFree3(a->a,a->i,a->j);CHKERRQ(ierr); 10649b5e25fSSatish Balay } 10749b5e25fSSatish Balay if (a->row) { 10849b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 10949b5e25fSSatish Balay } 11049b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 11149b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 11249b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 11349b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 114d59c15a7SBarry Smith if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 115d59c15a7SBarry Smith if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);} 116d59c15a7SBarry Smith if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 11749b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 118c4319e64SHong Zhang if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 1191a3463dfSHong Zhang 1201a3463dfSHong Zhang if (a->inew){ 1211a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1221a3463dfSHong Zhang a->inew = 0; 1231a3463dfSHong Zhang } 12449b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 125901853e0SKris Buschelman 126901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 127901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 128901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 129901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 130901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 131901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 13249b5e25fSSatish Balay PetscFunctionReturn(0); 13349b5e25fSSatish Balay } 13449b5e25fSSatish Balay 1354a2ae208SSatish Balay #undef __FUNCT__ 1364a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 137dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op) 13849b5e25fSSatish Balay { 139045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 14063ba0a88SBarry Smith PetscErrorCode ierr; 14149b5e25fSSatish Balay 14249b5e25fSSatish Balay PetscFunctionBegin; 1434d9d31abSKris Buschelman switch (op) { 1444d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1454d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1464d9d31abSKris Buschelman break; 1474d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1484d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1494d9d31abSKris Buschelman break; 1504d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1514d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1524d9d31abSKris Buschelman break; 1534d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1544d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1554d9d31abSKris Buschelman break; 1564d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1574d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1584d9d31abSKris Buschelman break; 1594d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1604d9d31abSKris Buschelman a->nonew = 1; 1614d9d31abSKris Buschelman break; 1624d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1634d9d31abSKris Buschelman a->nonew = -1; 1644d9d31abSKris Buschelman break; 1654d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1664d9d31abSKris Buschelman a->nonew = -2; 1674d9d31abSKris Buschelman break; 1684d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1694d9d31abSKris Buschelman a->nonew = 0; 1704d9d31abSKris Buschelman break; 1714d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1724d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1734d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1744d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1754d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 17663ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatSetOption_SeqSBAIJ:Option ignored\n"));CHKERRQ(ierr); 1774d9d31abSKris Buschelman break; 1784d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 17929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1809a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1819a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1829a4540c5SBarry Smith case MAT_HERMITIAN: 1839a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 18477e54ba9SKris Buschelman case MAT_SYMMETRIC: 18577e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1869a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1879a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1889a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 189941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 190941593c8SHong Zhang a->ignore_ltriangular = PETSC_TRUE; 191941593c8SHong Zhang break; 192941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 193941593c8SHong Zhang a->ignore_ltriangular = PETSC_FALSE; 19477e54ba9SKris Buschelman break; 1954d9d31abSKris Buschelman default: 19629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 19749b5e25fSSatish Balay } 19849b5e25fSSatish Balay PetscFunctionReturn(0); 19949b5e25fSSatish Balay } 20049b5e25fSSatish Balay 2014a2ae208SSatish Balay #undef __FUNCT__ 2024a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 20313f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 20449b5e25fSSatish Balay { 20549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2066849ba73SBarry Smith PetscErrorCode ierr; 20713f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 20849b5e25fSSatish Balay MatScalar *aa,*aa_i; 20987828ca2SBarry Smith PetscScalar *v_i; 21049b5e25fSSatish Balay 21149b5e25fSSatish Balay PetscFunctionBegin; 212521d7252SBarry Smith bs = A->bs; 21349b5e25fSSatish Balay ai = a->i; 21449b5e25fSSatish Balay aj = a->j; 21549b5e25fSSatish Balay aa = a->a; 21649b5e25fSSatish Balay bs2 = a->bs2; 21749b5e25fSSatish Balay 21877431f27SBarry Smith if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 21949b5e25fSSatish Balay 22049b5e25fSSatish Balay bn = row/bs; /* Block number */ 22149b5e25fSSatish Balay bp = row % bs; /* Block position */ 22249b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 22349b5e25fSSatish Balay *ncols = bs*M; 22449b5e25fSSatish Balay 22549b5e25fSSatish Balay if (v) { 22649b5e25fSSatish Balay *v = 0; 22749b5e25fSSatish Balay if (*ncols) { 22887828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 22949b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23049b5e25fSSatish Balay v_i = *v + i*bs; 23149b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 23249b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 23349b5e25fSSatish Balay } 23449b5e25fSSatish Balay } 23549b5e25fSSatish Balay } 23649b5e25fSSatish Balay 23749b5e25fSSatish Balay if (cols) { 23849b5e25fSSatish Balay *cols = 0; 23949b5e25fSSatish Balay if (*ncols) { 24013f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 24149b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 24249b5e25fSSatish Balay cols_i = *cols + i*bs; 24349b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 24449b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24549b5e25fSSatish Balay } 24649b5e25fSSatish Balay } 24749b5e25fSSatish Balay } 24849b5e25fSSatish Balay 24949b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2505ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2515ddb2528SHong Zhang #ifdef column_search 25249b5e25fSSatish Balay v_i = *v + M*bs; 25349b5e25fSSatish Balay cols_i = *cols + M*bs; 25449b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25549b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25649b5e25fSSatish Balay for (j=0; j<M; j++){ 25749b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 25849b5e25fSSatish Balay if (itmp == bn){ 25949b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 26049b5e25fSSatish Balay for (k=0; k<bs; k++) { 26149b5e25fSSatish Balay *cols_i++ = i*bs+k; 26249b5e25fSSatish Balay *v_i++ = aa_i[k]; 26349b5e25fSSatish Balay } 26449b5e25fSSatish Balay *ncols += bs; 26549b5e25fSSatish Balay break; 26649b5e25fSSatish Balay } 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay } 2695ddb2528SHong Zhang #endif 27049b5e25fSSatish Balay PetscFunctionReturn(0); 27149b5e25fSSatish Balay } 27249b5e25fSSatish Balay 2734a2ae208SSatish Balay #undef __FUNCT__ 2744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 27513f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 27649b5e25fSSatish Balay { 277dfbe8321SBarry Smith PetscErrorCode ierr; 27849b5e25fSSatish Balay 27949b5e25fSSatish Balay PetscFunctionBegin; 28049b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 28149b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 28249b5e25fSSatish Balay PetscFunctionReturn(0); 28349b5e25fSSatish Balay } 28449b5e25fSSatish Balay 2854a2ae208SSatish Balay #undef __FUNCT__ 2864a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 287dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 28849b5e25fSSatish Balay { 289dfbe8321SBarry Smith PetscErrorCode ierr; 29049b5e25fSSatish Balay PetscFunctionBegin; 291999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2928115998fSBarry Smith PetscFunctionReturn(0); 29349b5e25fSSatish Balay } 29449b5e25fSSatish Balay 2954a2ae208SSatish Balay #undef __FUNCT__ 2964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 2976849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 29849b5e25fSSatish Balay { 29949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 300dfbe8321SBarry Smith PetscErrorCode ierr; 301521d7252SBarry Smith PetscInt i,j,bs = A->bs,k,l,bs2=a->bs2; 302fb9695e5SSatish Balay char *name; 303f3ef73ceSBarry Smith PetscViewerFormat format; 30449b5e25fSSatish Balay 30549b5e25fSSatish Balay PetscFunctionBegin; 30680fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 307b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 308456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 30977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 310fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 31129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 312fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 313b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 31449b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31549b5e25fSSatish Balay for (j=0; j<bs; j++) { 31677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 31749b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31849b5e25fSSatish Balay for (l=0; l<bs; l++) { 31949b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 32049b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 32177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 32249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32349b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 32477431f27SBarry 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 (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 32777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32849b5e25fSSatish Balay } 32949b5e25fSSatish Balay #else 33049b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 33177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 33249b5e25fSSatish Balay } 33349b5e25fSSatish Balay #endif 33449b5e25fSSatish Balay } 33549b5e25fSSatish Balay } 336b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 33749b5e25fSSatish Balay } 33849b5e25fSSatish Balay } 339b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 34049b5e25fSSatish Balay } else { 341b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 34249b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 34349b5e25fSSatish Balay for (j=0; j<bs; j++) { 34477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 34549b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 34649b5e25fSSatish Balay for (l=0; l<bs; l++) { 34749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 34849b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 34977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 35049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35149b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 35277431f27SBarry 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 { 35577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35649b5e25fSSatish Balay } 35749b5e25fSSatish Balay #else 35877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35949b5e25fSSatish Balay #endif 36049b5e25fSSatish Balay } 36149b5e25fSSatish Balay } 362b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 36349b5e25fSSatish Balay } 36449b5e25fSSatish Balay } 365b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 36649b5e25fSSatish Balay } 367b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 36849b5e25fSSatish Balay PetscFunctionReturn(0); 36949b5e25fSSatish Balay } 37049b5e25fSSatish Balay 3714a2ae208SSatish Balay #undef __FUNCT__ 3724a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 3736849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 37449b5e25fSSatish Balay { 37549b5e25fSSatish Balay Mat A = (Mat) Aa; 37649b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 3776849ba73SBarry Smith PetscErrorCode ierr; 37813f74950SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2; 37913f74950SBarry Smith PetscMPIInt rank; 38049b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 38149b5e25fSSatish Balay MatScalar *aa; 38249b5e25fSSatish Balay MPI_Comm comm; 383b0a32e0cSBarry Smith PetscViewer viewer; 38449b5e25fSSatish Balay 38549b5e25fSSatish Balay PetscFunctionBegin; 38649b5e25fSSatish Balay /* 38749b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 38849b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 38949b5e25fSSatish Balay rest should return immediately. 39049b5e25fSSatish Balay */ 39149b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 39249b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 39349b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 39449b5e25fSSatish Balay 39549b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 39649b5e25fSSatish Balay 397b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 398b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 39949b5e25fSSatish Balay 40049b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 401b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 40249b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 40349b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 404c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 40549b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 40649b5e25fSSatish Balay aa = a->a + j*bs2; 40749b5e25fSSatish Balay for (k=0; k<bs; k++) { 40849b5e25fSSatish Balay for (l=0; l<bs; l++) { 40949b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 410b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 41149b5e25fSSatish Balay } 41249b5e25fSSatish Balay } 41349b5e25fSSatish Balay } 41449b5e25fSSatish Balay } 415b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 41649b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 41749b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 418c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 41949b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42049b5e25fSSatish Balay aa = a->a + j*bs2; 42149b5e25fSSatish Balay for (k=0; k<bs; k++) { 42249b5e25fSSatish Balay for (l=0; l<bs; l++) { 42349b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 424b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 42549b5e25fSSatish Balay } 42649b5e25fSSatish Balay } 42749b5e25fSSatish Balay } 42849b5e25fSSatish Balay } 42949b5e25fSSatish Balay 430b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 43149b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 43249b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 433c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 43449b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 43549b5e25fSSatish Balay aa = a->a + j*bs2; 43649b5e25fSSatish Balay for (k=0; k<bs; k++) { 43749b5e25fSSatish Balay for (l=0; l<bs; l++) { 43849b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 439b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44049b5e25fSSatish Balay } 44149b5e25fSSatish Balay } 44249b5e25fSSatish Balay } 44349b5e25fSSatish Balay } 44449b5e25fSSatish Balay PetscFunctionReturn(0); 44549b5e25fSSatish Balay } 44649b5e25fSSatish Balay 4474a2ae208SSatish Balay #undef __FUNCT__ 4484a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4496849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 45049b5e25fSSatish Balay { 451dfbe8321SBarry Smith PetscErrorCode ierr; 45249b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 453b0a32e0cSBarry Smith PetscDraw draw; 45449b5e25fSSatish Balay PetscTruth isnull; 45549b5e25fSSatish Balay 45649b5e25fSSatish Balay PetscFunctionBegin; 45749b5e25fSSatish Balay 458b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 459b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 46049b5e25fSSatish Balay 46149b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 462c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 46349b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 464b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 465b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 46649b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 46749b5e25fSSatish Balay PetscFunctionReturn(0); 46849b5e25fSSatish Balay } 46949b5e25fSSatish Balay 4704a2ae208SSatish Balay #undef __FUNCT__ 4714a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 472dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 47349b5e25fSSatish Balay { 474dfbe8321SBarry Smith PetscErrorCode ierr; 47532077d6dSBarry Smith PetscTruth iascii,isdraw; 47649b5e25fSSatish Balay 47749b5e25fSSatish Balay PetscFunctionBegin; 47832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 479fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 48032077d6dSBarry Smith if (iascii){ 48149b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 48249b5e25fSSatish Balay } else if (isdraw) { 48349b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 48449b5e25fSSatish Balay } else { 485a5e6ed63SBarry Smith Mat B; 486ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 487a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 488a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 48949b5e25fSSatish Balay } 49049b5e25fSSatish Balay PetscFunctionReturn(0); 49149b5e25fSSatish Balay } 49249b5e25fSSatish Balay 49349b5e25fSSatish Balay 4944a2ae208SSatish Balay #undef __FUNCT__ 4954a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 49613f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 49749b5e25fSSatish Balay { 498045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 49913f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 50013f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 50113f74950SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2; 50249b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 50349b5e25fSSatish Balay 50449b5e25fSSatish Balay PetscFunctionBegin; 50549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 50649b5e25fSSatish Balay row = im[k]; brow = row/bs; 50777431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 50877431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 50949b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 51049b5e25fSSatish Balay nrow = ailen[brow]; 51149b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 51277431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 51377431f27SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1); 51449b5e25fSSatish Balay col = in[l] ; 51549b5e25fSSatish Balay bcol = col/bs; 51649b5e25fSSatish Balay cidx = col%bs; 51749b5e25fSSatish Balay ridx = row%bs; 51849b5e25fSSatish Balay high = nrow; 51949b5e25fSSatish Balay low = 0; /* assume unsorted */ 52049b5e25fSSatish Balay while (high-low > 5) { 52149b5e25fSSatish Balay t = (low+high)/2; 52249b5e25fSSatish Balay if (rp[t] > bcol) high = t; 52349b5e25fSSatish Balay else low = t; 52449b5e25fSSatish Balay } 52549b5e25fSSatish Balay for (i=low; i<high; i++) { 52649b5e25fSSatish Balay if (rp[i] > bcol) break; 52749b5e25fSSatish Balay if (rp[i] == bcol) { 52849b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 52949b5e25fSSatish Balay goto finished; 53049b5e25fSSatish Balay } 53149b5e25fSSatish Balay } 53249b5e25fSSatish Balay *v++ = zero; 53349b5e25fSSatish Balay finished:; 53449b5e25fSSatish Balay } 53549b5e25fSSatish Balay } 53649b5e25fSSatish Balay PetscFunctionReturn(0); 53749b5e25fSSatish Balay } 53849b5e25fSSatish Balay 53949b5e25fSSatish Balay 5404a2ae208SSatish Balay #undef __FUNCT__ 5414a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 54213f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 54349b5e25fSSatish Balay { 5440880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5456849ba73SBarry Smith PetscErrorCode ierr; 546e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 54713f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 54813f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval; 5490880e062SHong Zhang PetscTruth roworiented=a->roworiented; 550f15d580aSBarry Smith const MatScalar *value = v; 551f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5520880e062SHong Zhang 55349b5e25fSSatish Balay PetscFunctionBegin; 5540880e062SHong Zhang if (roworiented) { 5550880e062SHong Zhang stepval = (n-1)*bs; 5560880e062SHong Zhang } else { 5570880e062SHong Zhang stepval = (m-1)*bs; 5580880e062SHong Zhang } 5590880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5600880e062SHong Zhang row = im[k]; 5610880e062SHong Zhang if (row < 0) continue; 5622515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 56377431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5640880e062SHong Zhang #endif 5650880e062SHong Zhang rp = aj + ai[row]; 5660880e062SHong Zhang ap = aa + bs2*ai[row]; 5670880e062SHong Zhang rmax = imax[row]; 5680880e062SHong Zhang nrow = ailen[row]; 5690880e062SHong Zhang low = 0; 570818f2c47SBarry Smith high = nrow; 5710880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5720880e062SHong Zhang if (in[l] < 0) continue; 5730880e062SHong Zhang col = in[l]; 5742515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 57577431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 576b1823623SSatish Balay #endif 5770880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 5780880e062SHong Zhang if (roworiented) { 5790880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 5800880e062SHong Zhang } else { 5810880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 5820880e062SHong Zhang } 583818f2c47SBarry Smith if (col < lastcol) low = 0; else high = nrow; 584e2ee6c50SBarry Smith lastcol = col; 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; 62777431f27SBarry 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 */ 63013f74950SBarry Smith PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6310880e062SHong Zhang MatScalar *new_a; 6320880e062SHong Zhang 63377431f27SBarry 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 */ 636*a96a251dSBarry Smith ierr = PetscMalloc3(bs2*new_nz,PetscScalar,&new_a,new_nz,PetscInt,&new_j,a->mbs+1,PetscInt,&new_i);CHKERRQ(ierr); 6370880e062SHong Zhang 6380880e062SHong Zhang /* copy over old data into new slots */ 6390880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6400880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 64113f74950SBarry Smith ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 6420880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 64313f74950SBarry Smith ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 6440880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 6450880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 6460880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 6470880e062SHong Zhang /* free up old matrix storage */ 6480880e062SHong Zhang if (!a->singlemalloc) { 649*a96a251dSBarry Smith ierr = PetscFree(a->a);CHKERRQ(ierr); 6500880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 6510880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 652*a96a251dSBarry Smith } else { 653*a96a251dSBarry Smith ierr = PetscFree3(a->a,a->i,a->j);CHKERRQ(ierr); 6540880e062SHong Zhang } 6550880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 6560880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 6570880e062SHong Zhang 6580880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 6590880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 66052e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar)));CHKERRQ(ierr); 6616c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 6620880e062SHong Zhang a->reallocs++; 6636c6c5352SBarry Smith a->nz++; 6640880e062SHong Zhang } 6650880e062SHong Zhang N = nrow++ - 1; 6660880e062SHong Zhang /* shift up all the later entries in this row */ 6670880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6680880e062SHong Zhang rp[ii+1] = rp[ii]; 6690880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6700880e062SHong Zhang } 6710880e062SHong Zhang if (N >= i) { 6720880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6730880e062SHong Zhang } 6740880e062SHong Zhang rp[i] = col; 6750880e062SHong Zhang bap = ap + bs2*i; 6760880e062SHong Zhang if (roworiented) { 6770880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6780880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6790880e062SHong Zhang bap[jj] = *value++; 6800880e062SHong Zhang } 6810880e062SHong Zhang } 6820880e062SHong Zhang } else { 6830880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6840880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6850880e062SHong Zhang *bap++ = *value++; 6860880e062SHong Zhang } 6870880e062SHong Zhang } 6880880e062SHong Zhang } 6890880e062SHong Zhang noinsert2:; 6900880e062SHong Zhang low = i; 6910880e062SHong Zhang } 6920880e062SHong Zhang ailen[row] = nrow; 6930880e062SHong Zhang } 6940880e062SHong Zhang PetscFunctionReturn(0); 69549b5e25fSSatish Balay } 69649b5e25fSSatish Balay 6974a2ae208SSatish Balay #undef __FUNCT__ 6984a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 699dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 70049b5e25fSSatish Balay { 70149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7026849ba73SBarry Smith PetscErrorCode ierr; 70313f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 70413f74950SBarry Smith PetscInt m = A->m,*ip,N,*ailen = a->ilen; 70513f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 70649b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 70749b5e25fSSatish Balay 70849b5e25fSSatish Balay PetscFunctionBegin; 70949b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 71049b5e25fSSatish Balay 71149b5e25fSSatish Balay if (m) rmax = ailen[0]; 71249b5e25fSSatish Balay for (i=1; i<mbs; i++) { 71349b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 71449b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 71549b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 71649b5e25fSSatish Balay if (fshift) { 71749b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 71849b5e25fSSatish Balay N = ailen[i]; 71949b5e25fSSatish Balay for (j=0; j<N; j++) { 72049b5e25fSSatish Balay ip[j-fshift] = ip[j]; 72149b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 72249b5e25fSSatish Balay } 72349b5e25fSSatish Balay } 72449b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 72549b5e25fSSatish Balay } 72649b5e25fSSatish Balay if (mbs) { 72749b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 72849b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 72949b5e25fSSatish Balay } 73049b5e25fSSatish Balay /* reset ilen and imax for each row */ 73149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 73249b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 73349b5e25fSSatish Balay } 7346c6c5352SBarry Smith a->nz = ai[mbs]; 73549b5e25fSSatish Balay 736b424e231SHong Zhang /* diagonals may have moved, reset it */ 737b424e231SHong Zhang if (a->diag) { 73813f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 73949b5e25fSSatish Balay } 74063ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n", 74163ba0a88SBarry Smith m,A->m,A->bs,fshift*bs2,a->nz*bs2));CHKERRQ(ierr); 74263ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %D\n",a->reallocs));CHKERRQ(ierr); 74363ba0a88SBarry Smith ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %D\n",rmax));CHKERRQ(ierr); 74449b5e25fSSatish Balay a->reallocs = 0; 74549b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 74649b5e25fSSatish Balay PetscFunctionReturn(0); 74749b5e25fSSatish Balay } 74849b5e25fSSatish Balay 74949b5e25fSSatish Balay /* 75049b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 75149b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 75249b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 75349b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 75449b5e25fSSatish Balay */ 7554a2ae208SSatish Balay #undef __FUNCT__ 7564a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 75713f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 75849b5e25fSSatish Balay { 75913f74950SBarry Smith PetscInt i,j,k,row; 76049b5e25fSSatish Balay PetscTruth flg; 76149b5e25fSSatish Balay 76249b5e25fSSatish Balay PetscFunctionBegin; 76349b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 76449b5e25fSSatish Balay row = idx[i]; 76549b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 76649b5e25fSSatish Balay sizes[j] = 1; 76749b5e25fSSatish Balay i++; 76849b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 76949b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 77049b5e25fSSatish Balay i++; 77149b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 77249b5e25fSSatish Balay flg = PETSC_TRUE; 77349b5e25fSSatish Balay for (k=1; k<bs; k++) { 77449b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 77549b5e25fSSatish Balay flg = PETSC_FALSE; 77649b5e25fSSatish Balay break; 77749b5e25fSSatish Balay } 77849b5e25fSSatish Balay } 779abc0a331SBarry Smith if (flg) { /* No break in the bs */ 78049b5e25fSSatish Balay sizes[j] = bs; 78149b5e25fSSatish Balay i+= bs; 78249b5e25fSSatish Balay } else { 78349b5e25fSSatish Balay sizes[j] = 1; 78449b5e25fSSatish Balay i++; 78549b5e25fSSatish Balay } 78649b5e25fSSatish Balay } 78749b5e25fSSatish Balay } 78849b5e25fSSatish Balay *bs_max = j; 78949b5e25fSSatish Balay PetscFunctionReturn(0); 79049b5e25fSSatish Balay } 79149b5e25fSSatish Balay 79249b5e25fSSatish Balay 79349b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 79449b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 79549b5e25fSSatish Balay */ 79649b5e25fSSatish Balay 7974a2ae208SSatish Balay #undef __FUNCT__ 7984a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 79913f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 80049b5e25fSSatish Balay { 80149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8026849ba73SBarry Smith PetscErrorCode ierr; 803e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 80413f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 80513f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->bs,brow,bcol; 80613f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 80749b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 80849b5e25fSSatish Balay 80949b5e25fSSatish Balay PetscFunctionBegin; 81049b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 81149b5e25fSSatish Balay row = im[k]; /* row number */ 81249b5e25fSSatish Balay brow = row/bs; /* block row number */ 81349b5e25fSSatish Balay if (row < 0) continue; 8142515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 81577431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 81649b5e25fSSatish Balay #endif 81749b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 81849b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 81949b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 82049b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 82149b5e25fSSatish Balay low = 0; 82249b5e25fSSatish Balay 82349b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 82449b5e25fSSatish Balay if (in[l] < 0) continue; 8252515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 82677431f27SBarry Smith if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1); 82749b5e25fSSatish Balay #endif 82849b5e25fSSatish Balay col = in[l]; 82949b5e25fSSatish Balay bcol = col/bs; /* block col number */ 83049b5e25fSSatish Balay 831941593c8SHong Zhang if (brow > bcol) { 832941593c8SHong Zhang if (a->ignore_ltriangular){ 833941593c8SHong Zhang continue; /* ignore lower triangular values */ 834941593c8SHong Zhang } else { 835941593c8SHong Zhang SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR)"); 836941593c8SHong Zhang } 837941593c8SHong Zhang } 838f4989cb3SHong Zhang 83949b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8408549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 84149b5e25fSSatish Balay /* element value a(k,l) */ 84249b5e25fSSatish Balay if (roworiented) { 84349b5e25fSSatish Balay value = v[l + k*n]; 84449b5e25fSSatish Balay } else { 84549b5e25fSSatish Balay value = v[k + l*m]; 84649b5e25fSSatish Balay } 84749b5e25fSSatish Balay 84849b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 849e2ee6c50SBarry Smith if (col < lastcol) low = 0; high = nrow; 850e2ee6c50SBarry Smith lastcol = col; 85149b5e25fSSatish Balay while (high-low > 7) { 85249b5e25fSSatish Balay t = (low+high)/2; 85349b5e25fSSatish Balay if (rp[t] > bcol) high = t; 85449b5e25fSSatish Balay else low = t; 85549b5e25fSSatish Balay } 85649b5e25fSSatish Balay for (i=low; i<high; i++) { 85749b5e25fSSatish Balay if (rp[i] > bcol) break; 85849b5e25fSSatish Balay if (rp[i] == bcol) { 85949b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 86049b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 86149b5e25fSSatish Balay else *bap = value; 8628549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8638549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8648549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8658549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8668549e402SHong Zhang else *bap = value; 8678549e402SHong Zhang } 86849b5e25fSSatish Balay goto noinsert1; 86949b5e25fSSatish Balay } 87049b5e25fSSatish Balay } 87149b5e25fSSatish Balay 87249b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 87377431f27SBarry Smith else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 87449b5e25fSSatish Balay if (nrow >= rmax) { 87549b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 87613f74950SBarry Smith PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 87749b5e25fSSatish Balay MatScalar *new_a; 87849b5e25fSSatish Balay 87977431f27SBarry Smith if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 88049b5e25fSSatish Balay 88149b5e25fSSatish Balay /* Malloc new storage space */ 882*a96a251dSBarry Smith ierr = PetscMalloc3(bs2*new_nz,PetscScalar,&new_a,new_nz,PetscInt,&new_j,a->mbs+1,PetscInt,&new_i);CHKERRQ(ierr); 88349b5e25fSSatish Balay 88449b5e25fSSatish Balay /* copy over old data into new slots */ 88549b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 88649b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 88713f74950SBarry Smith ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 88849b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 88913f74950SBarry Smith ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 89049b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 89149b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 89249b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 89349b5e25fSSatish Balay /* free up old matrix storage */ 89449b5e25fSSatish Balay if (!a->singlemalloc) { 895*a96a251dSBarry Smith ierr = PetscFree(a->a);CHKERRQ(ierr); 89649b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 89749b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 898*a96a251dSBarry Smith } else { 899*a96a251dSBarry Smith ierr = PetscFree3(a->a,a->i,a->j);CHKERRQ(ierr); 90049b5e25fSSatish Balay } 90149b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 90249b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 90349b5e25fSSatish Balay 90449b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 90549b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 90652e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar)));CHKERRQ(ierr); 9076c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 90849b5e25fSSatish Balay a->reallocs++; 9096c6c5352SBarry Smith a->nz++; 91049b5e25fSSatish Balay } 91149b5e25fSSatish Balay 91249b5e25fSSatish Balay N = nrow++ - 1; 91349b5e25fSSatish Balay /* shift up all the later entries in this row */ 91449b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 91549b5e25fSSatish Balay rp[ii+1] = rp[ii]; 91649b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 91749b5e25fSSatish Balay } 91849b5e25fSSatish Balay if (N>=i) { 91949b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 92049b5e25fSSatish Balay } 92149b5e25fSSatish Balay rp[i] = bcol; 92249b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 92349b5e25fSSatish Balay noinsert1:; 92449b5e25fSSatish Balay low = i; 9258549e402SHong Zhang } 92649b5e25fSSatish Balay } /* end of loop over added columns */ 92749b5e25fSSatish Balay ailen[brow] = nrow; 92849b5e25fSSatish Balay } /* end of loop over added rows */ 92949b5e25fSSatish Balay PetscFunctionReturn(0); 93049b5e25fSSatish Balay } 93149b5e25fSSatish Balay 9326849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 9336849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 93413f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt); 93513f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 93613f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]); 9376849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat); 9386849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 9396849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 9406849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec); 9416849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 9426849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 9436849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat); 9446849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec); 94549b5e25fSSatish Balay 9466849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 9476849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 9486849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 9496849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 9506849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 9516849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 9526849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 9536849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 95449b5e25fSSatish Balay 9556849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 956d59c15a7SBarry Smith 9576849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 9586849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 9596849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 9606849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 9616849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 9626849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 9636849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 9646849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 96507f98182SSatish Balay 966af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,MatFactorInfo*,Mat*); 967af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,MatFactorInfo*,Mat*); 968af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,MatFactorInfo*,Mat*); 969af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,MatFactorInfo*,Mat*); 970af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,MatFactorInfo*,Mat*); 971af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,MatFactorInfo*,Mat*); 972af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,MatFactorInfo*,Mat*); 973af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,MatFactorInfo*,Mat*); 97413f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*); 97549b5e25fSSatish Balay 9766849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 9776849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 9786849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 9796849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 9806849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 9816849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 9826849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 9836849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 98449b5e25fSSatish Balay 9856849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 9866849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 9876849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 9886849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 9896849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 9906849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 9916849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 9926849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 99349b5e25fSSatish Balay 9944d101231SSatish Balay #ifdef HAVE_MatICCFactor 9956849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work! */ 9964a2ae208SSatish Balay #undef __FUNCT__ 9974d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 998dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 99949b5e25fSSatish Balay { 10004ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 100149b5e25fSSatish Balay Mat outA; 1002dfbe8321SBarry Smith PetscErrorCode ierr; 100349b5e25fSSatish Balay PetscTruth row_identity,col_identity; 100449b5e25fSSatish Balay 100549b5e25fSSatish Balay PetscFunctionBegin; 100649b5e25fSSatish Balay outA = inA; 10071260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 100849b5e25fSSatish Balay 100949b5e25fSSatish Balay if (!a->diag) { 10101a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 101149b5e25fSSatish Balay } 101249b5e25fSSatish Balay /* 101349b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 101449b5e25fSSatish Balay for ILU(0) factorization with natural ordering 101549b5e25fSSatish Balay */ 101649b5e25fSSatish Balay switch (a->bs) { 101749b5e25fSSatish Balay case 1: 10180fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 101907f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1020d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 102163ba0a88SBarry Smith ierr = PetscLoginfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"));CHKERRQ(ierr); 102249b5e25fSSatish Balay case 2: 10231a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10241a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10251a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 102663ba0a88SBarry Smith ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"));CHKERRQ(ierr); 102749b5e25fSSatish Balay break; 102849b5e25fSSatish Balay case 3: 10291a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10301a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10311a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 103263ba0a88SBarry Smith ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"));CHKERRQ(ierr); 103349b5e25fSSatish Balay break; 103449b5e25fSSatish Balay case 4: 10351a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10361a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10371a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 103863ba0a88SBarry Smith ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"));CHKERRQ(ierr); 103949b5e25fSSatish Balay break; 104049b5e25fSSatish Balay case 5: 10411a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10421a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10431a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 104463ba0a88SBarry Smith ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"));CHKERRQ(ierr); 104549b5e25fSSatish Balay break; 104649b5e25fSSatish Balay case 6: 10471a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10481a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10491a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 105063ba0a88SBarry Smith ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"));CHKERRQ(ierr); 105149b5e25fSSatish Balay break; 105249b5e25fSSatish Balay case 7: 10531a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10541a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10551a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 105663ba0a88SBarry Smith ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"));CHKERRQ(ierr); 105749b5e25fSSatish Balay break; 105849b5e25fSSatish Balay default: 105949b5e25fSSatish Balay a->row = row; 10601a3463dfSHong Zhang a->icol = col; 106149b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 106249b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 106349b5e25fSSatish Balay 106449b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 106549b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 106652e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 106749b5e25fSSatish Balay 106849b5e25fSSatish Balay if (!a->solve_work) { 106987828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 107052e6d16bSBarry Smith ierr = PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));CHKERRQ(ierr); 107149b5e25fSSatish Balay } 107249b5e25fSSatish Balay } 107349b5e25fSSatish Balay 1074af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 107549b5e25fSSatish Balay PetscFunctionReturn(0); 107649b5e25fSSatish Balay } 1077950f1e5bSHong Zhang #endif 1078950f1e5bSHong Zhang 10794a2ae208SSatish Balay #undef __FUNCT__ 10804a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 1081dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A) 108249b5e25fSSatish Balay { 108349b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 108449b5e25fSSatish Balay MPI_Comm comm = A->comm; 1085dfbe8321SBarry Smith PetscErrorCode ierr; 108649b5e25fSSatish Balay 108749b5e25fSSatish Balay PetscFunctionBegin; 108849b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 108949b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 109049b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 1091941593c8SHong Zhang ierr = (*PetscHelpPrintf)(comm," -mat_ignore_ltriangular: Ignore lower triangular values set by user\n");CHKERRQ(ierr); 109249b5e25fSSatish Balay PetscFunctionReturn(0); 109349b5e25fSSatish Balay } 109449b5e25fSSatish Balay 109549b5e25fSSatish Balay EXTERN_C_BEGIN 10964a2ae208SSatish Balay #undef __FUNCT__ 10974a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1098be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 109949b5e25fSSatish Balay { 1100045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 110113f74950SBarry Smith PetscInt i,nz,n; 110249b5e25fSSatish Balay 110349b5e25fSSatish Balay PetscFunctionBegin; 11046c6c5352SBarry Smith nz = baij->maxnz; 1105c464158bSHong Zhang n = mat->n; 110649b5e25fSSatish Balay for (i=0; i<nz; i++) { 110749b5e25fSSatish Balay baij->j[i] = indices[i]; 110849b5e25fSSatish Balay } 11096c6c5352SBarry Smith baij->nz = nz; 111049b5e25fSSatish Balay for (i=0; i<n; i++) { 111149b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 111249b5e25fSSatish Balay } 111349b5e25fSSatish Balay PetscFunctionReturn(0); 111449b5e25fSSatish Balay } 111549b5e25fSSatish Balay EXTERN_C_END 111649b5e25fSSatish Balay 11174a2ae208SSatish Balay #undef __FUNCT__ 11184a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 111949b5e25fSSatish Balay /*@ 112019585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 112149b5e25fSSatish Balay in the matrix. 112249b5e25fSSatish Balay 112349b5e25fSSatish Balay Input Parameters: 112419585528SSatish Balay + mat - the SeqSBAIJ matrix 112549b5e25fSSatish Balay - indices - the column indices 112649b5e25fSSatish Balay 112749b5e25fSSatish Balay Level: advanced 112849b5e25fSSatish Balay 112949b5e25fSSatish Balay Notes: 113049b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 113149b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 113249b5e25fSSatish Balay of the MatSetValues() operation. 113349b5e25fSSatish Balay 113449b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1135d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 113649b5e25fSSatish Balay 1137ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 113849b5e25fSSatish Balay 1139ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 114049b5e25fSSatish Balay @*/ 1141be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 114249b5e25fSSatish Balay { 114313f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 114449b5e25fSSatish Balay 114549b5e25fSSatish Balay PetscFunctionBegin; 11464482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 11474482741eSBarry Smith PetscValidPointer(indices,2); 1148c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 114949b5e25fSSatish Balay if (f) { 115049b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 115149b5e25fSSatish Balay } else { 1152e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 115349b5e25fSSatish Balay } 115449b5e25fSSatish Balay PetscFunctionReturn(0); 115549b5e25fSSatish Balay } 115649b5e25fSSatish Balay 11574a2ae208SSatish Balay #undef __FUNCT__ 11584a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1159dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1160273d9f13SBarry Smith { 1161dfbe8321SBarry Smith PetscErrorCode ierr; 1162273d9f13SBarry Smith 1163273d9f13SBarry Smith PetscFunctionBegin; 1164273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1165273d9f13SBarry Smith PetscFunctionReturn(0); 1166273d9f13SBarry Smith } 1167273d9f13SBarry Smith 1168a6ece127SHong Zhang #undef __FUNCT__ 1169a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1170dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1171a6ece127SHong Zhang { 1172a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1173a6ece127SHong Zhang PetscFunctionBegin; 1174a6ece127SHong Zhang *array = a->a; 1175a6ece127SHong Zhang PetscFunctionReturn(0); 1176a6ece127SHong Zhang } 1177a6ece127SHong Zhang 1178a6ece127SHong Zhang #undef __FUNCT__ 1179a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1180dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1181a6ece127SHong Zhang { 1182a6ece127SHong Zhang PetscFunctionBegin; 1183a6ece127SHong Zhang PetscFunctionReturn(0); 1184a6ece127SHong Zhang } 1185a6ece127SHong Zhang 118642ee4b1aSHong Zhang #include "petscblaslapack.h" 118742ee4b1aSHong Zhang #undef __FUNCT__ 118842ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1189dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 119042ee4b1aSHong Zhang { 119142ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1192dfbe8321SBarry Smith PetscErrorCode ierr; 1193521d7252SBarry Smith PetscInt i,bs=Y->bs,bs2,j; 11944ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 119542ee4b1aSHong Zhang 119642ee4b1aSHong Zhang PetscFunctionBegin; 119742ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 119871044d3cSBarry Smith BLASaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one); 1199c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1200c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1201c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1202c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1203c537a176SHong Zhang } 1204c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1205c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1206c4319e64SHong Zhang y->XtoY = X; 1207c537a176SHong Zhang } 1208c4319e64SHong Zhang bs2 = bs*bs; 12096c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1210c4319e64SHong Zhang j = 0; 1211c4319e64SHong Zhang while (j < bs2){ 12126550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1213c4319e64SHong Zhang j++; 1214c537a176SHong Zhang } 1215c4319e64SHong Zhang } 121663ba0a88SBarry Smith ierr = 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)));CHKERRQ(ierr); 121742ee4b1aSHong Zhang } else { 121842ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 121942ee4b1aSHong Zhang } 122042ee4b1aSHong Zhang PetscFunctionReturn(0); 122142ee4b1aSHong Zhang } 122242ee4b1aSHong Zhang 1223efcf0fc3SBarry Smith #undef __FUNCT__ 1224efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1225dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1226efcf0fc3SBarry Smith { 1227efcf0fc3SBarry Smith PetscFunctionBegin; 1228efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1229efcf0fc3SBarry Smith PetscFunctionReturn(0); 1230efcf0fc3SBarry Smith } 1231efcf0fc3SBarry Smith 1232efcf0fc3SBarry Smith #undef __FUNCT__ 1233efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1234dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,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__ "MatIsHermitian_SeqSBAIJ" 1243dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1244efcf0fc3SBarry Smith { 1245efcf0fc3SBarry Smith PetscFunctionBegin; 1246efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1247efcf0fc3SBarry Smith PetscFunctionReturn(0); 1248efcf0fc3SBarry Smith } 1249efcf0fc3SBarry Smith 125049b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 125149b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 125249b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 125349b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 125449b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 125597304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1256e005ede5SBarry Smith MatMult_SeqSBAIJ_N, 1257e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 125849b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 125949b5e25fSSatish Balay 0, 126049b5e25fSSatish Balay 0, 126197304618SKris Buschelman /*10*/ 0, 126249b5e25fSSatish Balay 0, 12635f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1264d06b337dSHong Zhang MatRelax_SeqSBAIJ, 126549b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 126697304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 126749b5e25fSSatish Balay MatEqual_SeqSBAIJ, 126849b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 126949b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 127049b5e25fSSatish Balay MatNorm_SeqSBAIJ, 127197304618SKris Buschelman /*20*/ 0, 127249b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 127349b5e25fSSatish Balay 0, 127449b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 127549b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1276dcf5cc72SBarry Smith /*25*/ 0, 127749b5e25fSSatish Balay 0, 127849b5e25fSSatish Balay 0, 12795f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12805f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 128197304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1282c464158bSHong Zhang 0, 12834d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1284a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1285a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 128697304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 128749b5e25fSSatish Balay 0, 128849b5e25fSSatish Balay 0, 128949b5e25fSSatish Balay 0, 1290950f1e5bSHong Zhang 0, 129197304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 129249b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 129349b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 129449b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 129549b5e25fSSatish Balay 0, 129697304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 129749b5e25fSSatish Balay MatScale_SeqSBAIJ, 129849b5e25fSSatish Balay 0, 129949b5e25fSSatish Balay 0, 130049b5e25fSSatish Balay 0, 1301521d7252SBarry Smith /*50*/ 0, 130249b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 130349b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 130449b5e25fSSatish Balay 0, 130549b5e25fSSatish Balay 0, 130697304618SKris Buschelman /*55*/ 0, 130749b5e25fSSatish Balay 0, 130849b5e25fSSatish Balay 0, 130949b5e25fSSatish Balay 0, 131049b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 131197304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 131249b5e25fSSatish Balay 0, 131349b5e25fSSatish Balay 0, 13148a124369SBarry Smith MatGetPetscMaps_Petsc, 1315d959ec07SHong Zhang 0, 131697304618SKris Buschelman /*65*/ 0, 1317d959ec07SHong Zhang 0, 1318d959ec07SHong Zhang 0, 1319d959ec07SHong Zhang 0, 1320d959ec07SHong Zhang 0, 132197304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13223e0d88b5SBarry Smith 0, 13233e0d88b5SBarry Smith 0, 13243e0d88b5SBarry Smith 0, 13253e0d88b5SBarry Smith 0, 132697304618SKris Buschelman /*75*/ 0, 13273e0d88b5SBarry Smith 0, 13283e0d88b5SBarry Smith 0, 13293e0d88b5SBarry Smith 0, 13303e0d88b5SBarry Smith 0, 133197304618SKris Buschelman /*80*/ 0, 13323e0d88b5SBarry Smith 0, 13333e0d88b5SBarry Smith 0, 13343e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 133597304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13363e0d88b5SBarry Smith #else 133797304618SKris Buschelman 0, 13383e0d88b5SBarry Smith #endif 1339865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1340865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1341865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1342efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1343865e5f61SKris Buschelman 0, 1344865e5f61SKris Buschelman 0, 1345865e5f61SKris Buschelman /*90*/ 0, 1346865e5f61SKris Buschelman 0, 1347865e5f61SKris Buschelman 0, 1348865e5f61SKris Buschelman 0, 1349865e5f61SKris Buschelman 0, 1350865e5f61SKris Buschelman /*95*/ 0, 1351865e5f61SKris Buschelman 0, 1352865e5f61SKris Buschelman 0, 1353865e5f61SKris Buschelman 0}; 1354be1d678aSKris Buschelman 135549b5e25fSSatish Balay EXTERN_C_BEGIN 13564a2ae208SSatish Balay #undef __FUNCT__ 13574a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1358be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 135949b5e25fSSatish Balay { 13604afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1361521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 1362dfbe8321SBarry Smith PetscErrorCode ierr; 136349b5e25fSSatish Balay 136449b5e25fSSatish Balay PetscFunctionBegin; 136549b5e25fSSatish Balay if (aij->nonew != 1) { 1366e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 136749b5e25fSSatish Balay } 136849b5e25fSSatish Balay 136949b5e25fSSatish Balay /* allocate space for values if not already there */ 137049b5e25fSSatish Balay if (!aij->saved_values) { 137187828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 137249b5e25fSSatish Balay } 137349b5e25fSSatish Balay 137449b5e25fSSatish Balay /* copy values over */ 137587828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 137649b5e25fSSatish Balay PetscFunctionReturn(0); 137749b5e25fSSatish Balay } 137849b5e25fSSatish Balay EXTERN_C_END 137949b5e25fSSatish Balay 138049b5e25fSSatish Balay EXTERN_C_BEGIN 13814a2ae208SSatish Balay #undef __FUNCT__ 13824a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1383be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 138449b5e25fSSatish Balay { 13854afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13866849ba73SBarry Smith PetscErrorCode ierr; 1387521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 138849b5e25fSSatish Balay 138949b5e25fSSatish Balay PetscFunctionBegin; 139049b5e25fSSatish Balay if (aij->nonew != 1) { 1391e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 139249b5e25fSSatish Balay } 139349b5e25fSSatish Balay if (!aij->saved_values) { 1394e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 139549b5e25fSSatish Balay } 139649b5e25fSSatish Balay 139749b5e25fSSatish Balay /* copy values over */ 139887828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 139949b5e25fSSatish Balay PetscFunctionReturn(0); 140049b5e25fSSatish Balay } 140149b5e25fSSatish Balay EXTERN_C_END 140249b5e25fSSatish Balay 14038549e402SHong Zhang EXTERN_C_BEGIN 14044a2ae208SSatish Balay #undef __FUNCT__ 1405a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1406be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 140749b5e25fSSatish Balay { 1408c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14096849ba73SBarry Smith PetscErrorCode ierr; 141013f74950SBarry Smith PetscInt i,len,mbs,bs2; 141149b5e25fSSatish Balay PetscTruth flg; 141249b5e25fSSatish Balay 141349b5e25fSSatish Balay PetscFunctionBegin; 1414273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1415e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1416c464158bSHong Zhang mbs = B->m/bs; 141749b5e25fSSatish Balay bs2 = bs*bs; 141849b5e25fSSatish Balay 1419c464158bSHong Zhang if (mbs*bs != B->m) { 142029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 142149b5e25fSSatish Balay } 142249b5e25fSSatish Balay 1423435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 142477431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 142549b5e25fSSatish Balay if (nnz) { 142649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 142777431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 142877431f27SBarry 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); 142949b5e25fSSatish Balay } 143049b5e25fSSatish Balay } 143149b5e25fSSatish Balay 1432e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 143349b5e25fSSatish Balay if (!flg) { 143449b5e25fSSatish Balay switch (bs) { 143549b5e25fSSatish Balay case 1: 14364ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 143749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1438d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1439e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 144049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 144149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 144249b5e25fSSatish Balay break; 144349b5e25fSSatish Balay case 2: 14444ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 144549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1446e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 144749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 144849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 144949b5e25fSSatish Balay break; 145049b5e25fSSatish Balay case 3: 14515f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 145249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1453e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 145449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 145549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 145649b5e25fSSatish Balay break; 145749b5e25fSSatish Balay case 4: 14585f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 145949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1460e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 146149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 146249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 146349b5e25fSSatish Balay break; 146449b5e25fSSatish Balay case 5: 14655f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 146649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1467e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 146849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 146949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 147049b5e25fSSatish Balay break; 147149b5e25fSSatish Balay case 6: 14725f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 147349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1474e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 147549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 147649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 147749b5e25fSSatish Balay break; 147849b5e25fSSatish Balay case 7: 1479de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 148049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1481e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1482de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 148349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 148449b5e25fSSatish Balay break; 148549b5e25fSSatish Balay } 148649b5e25fSSatish Balay } 148749b5e25fSSatish Balay 148849b5e25fSSatish Balay b->mbs = mbs; 14894afc71dfSHong Zhang b->nbs = mbs; 149013f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->imax);CHKERRQ(ierr); 149149b5e25fSSatish Balay if (!nnz) { 1492435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 149349b5e25fSSatish Balay else if (nz <= 0) nz = 1; 149449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14958cef66ccSHong Zhang b->imax[i] = nz; 149649b5e25fSSatish Balay } 1497153ea458SHong Zhang nz = nz*mbs; /* total nz */ 149849b5e25fSSatish Balay } else { 149949b5e25fSSatish Balay nz = 0; 15008cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 150149b5e25fSSatish Balay } 15026c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 150349b5e25fSSatish Balay 150449b5e25fSSatish Balay /* allocate the matrix space */ 1505*a96a251dSBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->m+1,PetscInt,&b->i);CHKERRQ(ierr); 15066c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 150713f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 150849b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 150949b5e25fSSatish Balay 151049b5e25fSSatish Balay /* pointer to beginning of each row */ 151149b5e25fSSatish Balay b->i[0] = 0; 151249b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 151349b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 151449b5e25fSSatish Balay } 151549b5e25fSSatish Balay 151649b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 151713f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->ilen);CHKERRQ(ierr); 151852e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 151949b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 152049b5e25fSSatish Balay 1521521d7252SBarry Smith B->bs = bs; 152249b5e25fSSatish Balay b->bs2 = bs2; 15236c6c5352SBarry Smith b->nz = 0; 15246c6c5352SBarry Smith b->maxnz = nz*bs2; 1525153ea458SHong Zhang 152616cdd363SHong Zhang b->inew = 0; 152716cdd363SHong Zhang b->jnew = 0; 152816cdd363SHong Zhang b->anew = 0; 152916cdd363SHong Zhang b->a2anew = 0; 15301a3463dfSHong Zhang b->permute = PETSC_FALSE; 1531c464158bSHong Zhang PetscFunctionReturn(0); 1532c464158bSHong Zhang } 1533a23d5eceSKris Buschelman EXTERN_C_END 1534153ea458SHong Zhang 1535d769727bSBarry Smith EXTERN_C_BEGIN 1536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,MatReuse,Mat*); 1537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 1538d769727bSBarry Smith EXTERN_C_END 1539d769727bSBarry Smith 15400bad9183SKris Buschelman /*MC 1541fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15420bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15430bad9183SKris Buschelman 15440bad9183SKris Buschelman Options Database Keys: 15450bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15460bad9183SKris Buschelman 15470bad9183SKris Buschelman Level: beginner 15480bad9183SKris Buschelman 15490bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15500bad9183SKris Buschelman M*/ 15510bad9183SKris Buschelman 1552a23d5eceSKris Buschelman EXTERN_C_BEGIN 1553a23d5eceSKris Buschelman #undef __FUNCT__ 1554a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1555be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1556a23d5eceSKris Buschelman { 1557a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1558dfbe8321SBarry Smith PetscErrorCode ierr; 155913f74950SBarry Smith PetscMPIInt size; 1560941593c8SHong Zhang PetscTruth flg; 1561a23d5eceSKris Buschelman 1562a23d5eceSKris Buschelman PetscFunctionBegin; 1563a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1564a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1565a23d5eceSKris Buschelman B->m = B->M = PetscMax(B->m,B->M); 1566a23d5eceSKris Buschelman B->n = B->N = PetscMax(B->n,B->N); 1567a23d5eceSKris Buschelman 1568a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1569a23d5eceSKris Buschelman B->data = (void*)b; 1570a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1571a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1572a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1573a23d5eceSKris Buschelman B->factor = 0; 1574a23d5eceSKris Buschelman B->mapping = 0; 1575a23d5eceSKris Buschelman b->row = 0; 1576a23d5eceSKris Buschelman b->icol = 0; 1577a23d5eceSKris Buschelman b->reallocs = 0; 1578a23d5eceSKris Buschelman b->saved_values = 0; 1579a23d5eceSKris Buschelman 1580a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1581a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1582a23d5eceSKris Buschelman 1583a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1584a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1585a23d5eceSKris Buschelman b->nonew = 0; 1586a23d5eceSKris Buschelman b->diag = 0; 1587a23d5eceSKris Buschelman b->solve_work = 0; 1588a23d5eceSKris Buschelman b->mult_work = 0; 1589a23d5eceSKris Buschelman B->spptr = 0; 1590a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1591a23d5eceSKris Buschelman b->xtoy = 0; 1592a23d5eceSKris Buschelman b->XtoY = 0; 1593a23d5eceSKris Buschelman 1594a23d5eceSKris Buschelman b->inew = 0; 1595a23d5eceSKris Buschelman b->jnew = 0; 1596a23d5eceSKris Buschelman b->anew = 0; 1597a23d5eceSKris Buschelman b->a2anew = 0; 1598a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1599a23d5eceSKris Buschelman 1600941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 1601941593c8SHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr); 1602941593c8SHong Zhang if (flg) b->ignore_ltriangular = PETSC_TRUE; 1603941593c8SHong Zhang 1604a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1605a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1606a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1607a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1608a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1609a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1610a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1611a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1612a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16134e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16144e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16154e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1616a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1617a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1618a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1619a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1620a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1621a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 162223ce1328SBarry Smith 162323ce1328SBarry Smith B->symmetric = PETSC_TRUE; 162423ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 162523ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 162623ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 1627a23d5eceSKris Buschelman PetscFunctionReturn(0); 1628a23d5eceSKris Buschelman } 1629a23d5eceSKris Buschelman EXTERN_C_END 1630a23d5eceSKris Buschelman 1631a23d5eceSKris Buschelman #undef __FUNCT__ 1632a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1633a23d5eceSKris Buschelman /*@C 1634a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1635a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1636a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1637a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1638a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1639a23d5eceSKris Buschelman 1640a23d5eceSKris Buschelman Collective on Mat 1641a23d5eceSKris Buschelman 1642a23d5eceSKris Buschelman Input Parameters: 1643a23d5eceSKris Buschelman + A - the symmetric matrix 1644a23d5eceSKris Buschelman . bs - size of block 1645a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1646a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1647a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1648a23d5eceSKris Buschelman 1649a23d5eceSKris Buschelman Options Database Keys: 1650a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1651a23d5eceSKris Buschelman block calculations (much slower) 1652a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1653a23d5eceSKris Buschelman 1654a23d5eceSKris Buschelman Level: intermediate 1655a23d5eceSKris Buschelman 1656a23d5eceSKris Buschelman Notes: 1657a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1658a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1659a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1660a23d5eceSKris Buschelman matrices. 1661a23d5eceSKris Buschelman 166249a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 166349a6f317SBarry Smith 166449a6f317SBarry Smith 1665a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1666a23d5eceSKris Buschelman @*/ 1667be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 166813f74950SBarry Smith { 166913f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1670a23d5eceSKris Buschelman 1671a23d5eceSKris Buschelman PetscFunctionBegin; 1672a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1673a23d5eceSKris Buschelman if (f) { 1674a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1675a23d5eceSKris Buschelman } 1676a23d5eceSKris Buschelman PetscFunctionReturn(0); 1677a23d5eceSKris Buschelman } 167849b5e25fSSatish Balay 16794a2ae208SSatish Balay #undef __FUNCT__ 16804a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1681c464158bSHong Zhang /*@C 1682c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1683c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1684c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1685c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1686c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 168749b5e25fSSatish Balay 1688c464158bSHong Zhang Collective on MPI_Comm 1689c464158bSHong Zhang 1690c464158bSHong Zhang Input Parameters: 1691c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1692c464158bSHong Zhang . bs - size of block 1693c464158bSHong Zhang . m - number of rows, or number of columns 1694c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1695744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1696744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1697c464158bSHong Zhang 1698c464158bSHong Zhang Output Parameter: 1699c464158bSHong Zhang . A - the symmetric matrix 1700c464158bSHong Zhang 1701c464158bSHong Zhang Options Database Keys: 1702c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1703c464158bSHong Zhang block calculations (much slower) 1704c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1705c464158bSHong Zhang 1706c464158bSHong Zhang Level: intermediate 1707c464158bSHong Zhang 1708c464158bSHong Zhang Notes: 1709c464158bSHong Zhang 1710c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1711c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1712c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1713c464158bSHong Zhang matrices. 1714c464158bSHong Zhang 171549a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 171649a6f317SBarry Smith 1717c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1718c464158bSHong Zhang @*/ 1719be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1720c464158bSHong Zhang { 1721dfbe8321SBarry Smith PetscErrorCode ierr; 1722c464158bSHong Zhang 1723c464158bSHong Zhang PetscFunctionBegin; 1724c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1725c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1726273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 172749b5e25fSSatish Balay PetscFunctionReturn(0); 172849b5e25fSSatish Balay } 172949b5e25fSSatish Balay 17304a2ae208SSatish Balay #undef __FUNCT__ 17314a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1732dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 173349b5e25fSSatish Balay { 173449b5e25fSSatish Balay Mat C; 173549b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17366849ba73SBarry Smith PetscErrorCode ierr; 173713f74950SBarry Smith PetscInt i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 173849b5e25fSSatish Balay 173949b5e25fSSatish Balay PetscFunctionBegin; 174029bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 174149b5e25fSSatish Balay 174249b5e25fSSatish Balay *B = 0; 1743692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1744be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17451d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1746692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1747692f9cbeSHong Zhang 1748273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 174949b5e25fSSatish Balay C->factor = A->factor; 175049b5e25fSSatish Balay c->row = 0; 175149b5e25fSSatish Balay c->icol = 0; 175249b5e25fSSatish Balay c->saved_values = 0; 175349b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 175449b5e25fSSatish Balay C->assembled = PETSC_TRUE; 175549b5e25fSSatish Balay 175682327fa8SHong Zhang C->M = A->M; 175782327fa8SHong Zhang C->N = A->N; 1758521d7252SBarry Smith C->bs = A->bs; 175949b5e25fSSatish Balay c->bs2 = a->bs2; 176049b5e25fSSatish Balay c->mbs = a->mbs; 176149b5e25fSSatish Balay c->nbs = a->nbs; 176249b5e25fSSatish Balay 176313f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr); 176413f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr); 176549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 176649b5e25fSSatish Balay c->imax[i] = a->imax[i]; 176749b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 176849b5e25fSSatish Balay } 176949b5e25fSSatish Balay 177049b5e25fSSatish Balay /* allocate the matrix space */ 1771*a96a251dSBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 177249b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 177313f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 177449b5e25fSSatish Balay if (mbs > 0) { 177513f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 177649b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 177749b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 177849b5e25fSSatish Balay } else { 177949b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 178049b5e25fSSatish Balay } 178149b5e25fSSatish Balay } 178249b5e25fSSatish Balay 178352e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 178449b5e25fSSatish Balay c->sorted = a->sorted; 178549b5e25fSSatish Balay c->roworiented = a->roworiented; 178649b5e25fSSatish Balay c->nonew = a->nonew; 178749b5e25fSSatish Balay 178849b5e25fSSatish Balay if (a->diag) { 178913f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 179052e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 179149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 179249b5e25fSSatish Balay c->diag[i] = a->diag[i]; 179349b5e25fSSatish Balay } 179449b5e25fSSatish Balay } else c->diag = 0; 17956c6c5352SBarry Smith c->nz = a->nz; 17966c6c5352SBarry Smith c->maxnz = a->maxnz; 179749b5e25fSSatish Balay c->solve_work = 0; 179849b5e25fSSatish Balay c->mult_work = 0; 179949b5e25fSSatish Balay *B = C; 1800b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 180149b5e25fSSatish Balay PetscFunctionReturn(0); 180249b5e25fSSatish Balay } 180349b5e25fSSatish Balay 18044a2ae208SSatish Balay #undef __FUNCT__ 18054a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1806dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A) 180749b5e25fSSatish Balay { 180849b5e25fSSatish Balay Mat_SeqSBAIJ *a; 180949b5e25fSSatish Balay Mat B; 18106849ba73SBarry Smith PetscErrorCode ierr; 181113f74950SBarry Smith int fd; 181213f74950SBarry Smith PetscMPIInt size; 181313f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 181413f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 181513f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 181613f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 181787828ca2SBarry Smith PetscScalar *aa; 181849b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 181949b5e25fSSatish Balay 182049b5e25fSSatish Balay PetscFunctionBegin; 1821b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 182249b5e25fSSatish Balay bs2 = bs*bs; 182349b5e25fSSatish Balay 182449b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 182529bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1826b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 182749b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1828552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 182949b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 183049b5e25fSSatish Balay 183149b5e25fSSatish Balay if (header[3] < 0) { 183229bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 183349b5e25fSSatish Balay } 183449b5e25fSSatish Balay 183529bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 183649b5e25fSSatish Balay 183749b5e25fSSatish Balay /* 183849b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 183949b5e25fSSatish Balay divisible by the blocksize 184049b5e25fSSatish Balay */ 184149b5e25fSSatish Balay mbs = M/bs; 184249b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 184349b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 184449b5e25fSSatish Balay else mbs++; 184549b5e25fSSatish Balay if (extra_rows) { 184663ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"));CHKERRQ(ierr); 184749b5e25fSSatish Balay } 184849b5e25fSSatish Balay 184949b5e25fSSatish Balay /* read in row lengths */ 185013f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 185149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 185249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 185349b5e25fSSatish Balay 185449b5e25fSSatish Balay /* read in column indices */ 185513f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 185649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 185749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 185849b5e25fSSatish Balay 185949b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 186013f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 186113f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 186213f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 186313f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 186449b5e25fSSatish Balay masked = mask + mbs; 186549b5e25fSSatish Balay rowcount = 0; nzcount = 0; 186649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 186749b5e25fSSatish Balay nmask = 0; 186849b5e25fSSatish Balay for (j=0; j<bs; j++) { 186949b5e25fSSatish Balay kmax = rowlengths[rowcount]; 187049b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18712d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 187203630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 187349b5e25fSSatish Balay } 187449b5e25fSSatish Balay rowcount++; 187549b5e25fSSatish Balay } 1876574b2666SHong Zhang s_browlengths[i] += nmask; 1877574b2666SHong Zhang 187849b5e25fSSatish Balay /* zero out the mask elements we set */ 187949b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 188049b5e25fSSatish Balay } 188149b5e25fSSatish Balay 188249b5e25fSSatish Balay /* create our matrix */ 18839abb65ffSKris Buschelman ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr); 18849abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 188578473fd7SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr); 188649b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 188749b5e25fSSatish Balay 188849b5e25fSSatish Balay /* set matrix "i" values */ 188949b5e25fSSatish Balay a->i[0] = 0; 189049b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1891574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1892574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 189349b5e25fSSatish Balay } 18946c6c5352SBarry Smith a->nz = a->i[mbs]; 189549b5e25fSSatish Balay 189649b5e25fSSatish Balay /* read in nonzero values */ 189787828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 189849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 189949b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 190049b5e25fSSatish Balay 190149b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 190249b5e25fSSatish Balay nzcount = 0; jcount = 0; 190349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 190449b5e25fSSatish Balay nzcountb = nzcount; 190549b5e25fSSatish Balay nmask = 0; 190649b5e25fSSatish Balay for (j=0; j<bs; j++) { 190749b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 190849b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19092d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 191003630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19112d703238SHong Zhang } 19122d703238SHong Zhang } 19132d703238SHong Zhang /* sort the masked values */ 19142d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19152d703238SHong Zhang 19162d703238SHong Zhang /* set "j" values into matrix */ 19172d703238SHong Zhang maskcount = 1; 19182d703238SHong Zhang for (j=0; j<nmask; j++) { 191949b5e25fSSatish Balay a->j[jcount++] = masked[j]; 192049b5e25fSSatish Balay mask[masked[j]] = maskcount++; 192149b5e25fSSatish Balay } 1922574b2666SHong Zhang 192349b5e25fSSatish Balay /* set "a" values into matrix */ 192449b5e25fSSatish Balay ishift = bs2*a->i[i]; 192549b5e25fSSatish Balay for (j=0; j<bs; j++) { 192649b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 192749b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1928574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1929574b2666SHong Zhang if (tmp >= i){ 193049b5e25fSSatish Balay block = mask[tmp] - 1; 193149b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 193249b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1933574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1934574b2666SHong Zhang } 1935574b2666SHong Zhang nzcountb++; 193649b5e25fSSatish Balay } 193749b5e25fSSatish Balay } 193849b5e25fSSatish Balay /* zero out the mask elements we set */ 193949b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 194049b5e25fSSatish Balay } 19416c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 194249b5e25fSSatish Balay 194349b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1944574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 194549b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 194649b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 194749b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 194849b5e25fSSatish Balay 19499abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19509abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195149b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19529abb65ffSKris Buschelman *A = B; 195349b5e25fSSatish Balay PetscFunctionReturn(0); 195449b5e25fSSatish Balay } 1955574b2666SHong Zhang 1956d06b337dSHong Zhang #undef __FUNCT__ 1957d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 195813f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1959d06b337dSHong Zhang { 1960d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1961d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1962d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19636849ba73SBarry Smith PetscErrorCode ierr; 1964521d7252SBarry Smith PetscInt m=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j; 1965521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 1966d06b337dSHong Zhang 1967d06b337dSHong Zhang PetscFunctionBegin; 1968b965ef7fSBarry Smith its = its*lits; 196977431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1970d06b337dSHong Zhang 1971d06b337dSHong Zhang if (bs > 1) 1972d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1973d06b337dSHong Zhang 19741ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1975d06b337dSHong Zhang if (xx != bb) { 19761ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1977d06b337dSHong Zhang } else { 1978d06b337dSHong Zhang b = x; 1979d06b337dSHong Zhang } 1980d06b337dSHong Zhang 1981d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1982d06b337dSHong Zhang 1983d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19842798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1985d06b337dSHong Zhang for (i=0; i<m; i++) 1986d06b337dSHong Zhang t[i] = b[i]; 1987d06b337dSHong Zhang 1988d06b337dSHong Zhang for (i=0; i<m; i++){ 198944706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1990d06b337dSHong Zhang v = aa + ai[i] + 1; 1991d06b337dSHong Zhang vj = aj + ai[i] + 1; 1992d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1993d06b337dSHong Zhang x[i] = omega*t[i]/d; 1994d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 1995efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 1996d06b337dSHong Zhang } 1997d06b337dSHong Zhang } 1998d06b337dSHong Zhang 19992798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 2000d06b337dSHong Zhang for (i=0; i<m; i++) 2001d06b337dSHong Zhang t[i] = b[i]; 2002d06b337dSHong Zhang 2003d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2004d06b337dSHong Zhang v = aa + ai[i] + 1; 2005d06b337dSHong Zhang vj = aj + ai[i] + 1; 2006d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2007d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2008efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2009d06b337dSHong Zhang } 2010d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2011d06b337dSHong Zhang d = *(aa + ai[i]); 2012d06b337dSHong Zhang v = aa + ai[i] + 1; 2013d06b337dSHong Zhang vj = aj + ai[i] + 1; 2014d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2015d06b337dSHong Zhang sum = t[i]; 2016d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2017efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2018d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2019d06b337dSHong Zhang } 2020d06b337dSHong Zhang } 2021d06b337dSHong Zhang its--; 2022d06b337dSHong Zhang } 2023d06b337dSHong Zhang 2024d06b337dSHong Zhang while (its--) { 202544706875SHong Zhang /* 202644706875SHong Zhang forward sweep: 202744706875SHong Zhang for i=0,...,m-1: 202844706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 202944706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 203044706875SHong Zhang b = b - x[i]*U^T(i,:); 2031d06b337dSHong Zhang 203244706875SHong Zhang */ 20332798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2034d06b337dSHong Zhang for (i=0; i<m; i++) 2035d06b337dSHong Zhang t[i] = b[i]; 2036d06b337dSHong Zhang 2037d06b337dSHong Zhang for (i=0; i<m; i++){ 203844706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2039d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2040d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2041d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2042d06b337dSHong Zhang sum = t[i]; 2043d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2044d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2045d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2046efee365bSSatish Balay ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr); 2047d06b337dSHong Zhang } 2048d06b337dSHong Zhang } 2049d06b337dSHong Zhang 20502798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 205144706875SHong Zhang /* 205244706875SHong Zhang backward sweep: 205344706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 205444706875SHong Zhang for i=m-1,...,0: 205544706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 205644706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 205744706875SHong Zhang */ 2058d06b337dSHong Zhang for (i=0; i<m; i++) 2059d06b337dSHong Zhang t[i] = b[i]; 2060d06b337dSHong Zhang 2061d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2062d06b337dSHong Zhang v = aa + ai[i] + 1; 2063d06b337dSHong Zhang vj = aj + ai[i] + 1; 2064d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2065d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 2066efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2067d06b337dSHong Zhang } 2068d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2069d06b337dSHong Zhang d = *(aa + ai[i]); 2070d06b337dSHong Zhang v = aa + ai[i] + 1; 2071d06b337dSHong Zhang vj = aj + ai[i] + 1; 2072d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2073d06b337dSHong Zhang sum = t[i]; 2074d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 2075efee365bSSatish Balay ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr); 2076d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2077d06b337dSHong Zhang } 2078d06b337dSHong Zhang } 2079d06b337dSHong Zhang } 2080d06b337dSHong Zhang 2081d06b337dSHong Zhang ierr = PetscFree(t);CHKERRQ(ierr); 20821ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2083d06b337dSHong Zhang if (bb != xx) { 20841ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2085d06b337dSHong Zhang } 2086d06b337dSHong Zhang PetscFunctionReturn(0); 2087d06b337dSHong Zhang } 2088d06b337dSHong Zhang 2089d06b337dSHong Zhang 2090d06b337dSHong Zhang 2091d06b337dSHong Zhang 209249b5e25fSSatish Balay 209349b5e25fSSatish Balay 2094