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 */ 77c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 87c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 949b5e25fSSatish Balay 10b5b17502SBarry Smith #include "../src/mat/impls/sbaij/seq/relax.h" 1170dcbbb9SBarry Smith #define USESHORT 12b5b17502SBarry Smith #include "../src/mat/impls/sbaij/seq/relax.h" 1370dcbbb9SBarry Smith 14db4efbfdSBarry Smith extern PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat,PetscTruth); 1549b5e25fSSatish Balay #define CHUNKSIZE 10 1649b5e25fSSatish Balay 17b5b17502SBarry Smith 1849b5e25fSSatish Balay /* 1949b5e25fSSatish Balay Checks for missing diagonals 2049b5e25fSSatish Balay */ 214a2ae208SSatish Balay #undef __FUNCT__ 224a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 232af78befSBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscTruth *missing,PetscInt *dd) 2449b5e25fSSatish Balay { 25045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 266849ba73SBarry Smith PetscErrorCode ierr; 2713f74950SBarry Smith PetscInt *diag,*jj = a->j,i; 2849b5e25fSSatish Balay 2949b5e25fSSatish Balay PetscFunctionBegin; 30045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 3149b5e25fSSatish Balay diag = a->diag; 322af78befSBarry Smith *missing = PETSC_FALSE; 3349b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 342af78befSBarry Smith if (jj[diag[i]] != i) { 352af78befSBarry Smith *missing = PETSC_TRUE; 362af78befSBarry Smith if (dd) *dd = i; 372af78befSBarry Smith break; 382af78befSBarry Smith } 3949b5e25fSSatish Balay } 4049b5e25fSSatish Balay PetscFunctionReturn(0); 4149b5e25fSSatish Balay } 4249b5e25fSSatish Balay 434a2ae208SSatish Balay #undef __FUNCT__ 444a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 45dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 4649b5e25fSSatish Balay { 47045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 486849ba73SBarry Smith PetscErrorCode ierr; 4909f38230SBarry Smith PetscInt i; 5049b5e25fSSatish Balay 5149b5e25fSSatish Balay PetscFunctionBegin; 5209f38230SBarry Smith if (!a->diag) { 5309f38230SBarry Smith ierr = PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 54c760cd28SBarry Smith ierr = PetscLogObjectMemory(A,a->mbs*sizeof(PetscInt));CHKERRQ(ierr); 55c760cd28SBarry Smith a->free_diag = PETSC_TRUE; 5609f38230SBarry Smith } 5709f38230SBarry Smith for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i]; 5849b5e25fSSatish Balay PetscFunctionReturn(0); 5949b5e25fSSatish Balay } 6049b5e25fSSatish Balay 614a2ae208SSatish Balay #undef __FUNCT__ 624a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 638f7157efSSatish Balay static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 6449b5e25fSSatish Balay { 65a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 66d0f46423SBarry Smith PetscInt i,j,n = a->mbs,nz = a->i[n],bs = A->rmap->bs; 678f7157efSSatish Balay PetscErrorCode ierr; 6849b5e25fSSatish Balay 6949b5e25fSSatish Balay PetscFunctionBegin; 70d3e5a4abSHong Zhang *nn = n; 71a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 728f7157efSSatish Balay if (!blockcompressed) { 738f7157efSSatish Balay /* malloc & create the natural set of indices */ 74f1d0d59dSSatish Balay ierr = PetscMalloc2((n+1)*bs,PetscInt,ia,nz*bs,PetscInt,ja);CHKERRQ(ierr); 758f7157efSSatish Balay for (i=0; i<n+1; i++) { 768f7157efSSatish Balay for (j=0; j<bs; j++) { 778f7157efSSatish Balay *ia[i*bs+j] = a->i[i]*bs+j+oshift; 788f7157efSSatish Balay } 798f7157efSSatish Balay } 808f7157efSSatish Balay for (i=0; i<nz; i++) { 818f7157efSSatish Balay for (j=0; j<bs; j++) { 828f7157efSSatish Balay *ja[i*bs+j] = a->j[i]*bs+j+oshift; 838f7157efSSatish Balay } 848f7157efSSatish Balay } 858f7157efSSatish Balay } else { /* blockcompressed */ 86a6ece127SHong Zhang if (oshift == 1) { 87a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 886c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 89a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 908f7157efSSatish Balay } 91a1373b80SHong Zhang *ia = a->i; *ja = a->j; 92a6ece127SHong Zhang } 938f7157efSSatish Balay 9449b5e25fSSatish Balay PetscFunctionReturn(0); 9549b5e25fSSatish Balay } 9649b5e25fSSatish Balay 974a2ae208SSatish Balay #undef __FUNCT__ 984a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 998f7157efSSatish Balay static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 10049b5e25fSSatish Balay { 101b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1028f7157efSSatish Balay PetscInt i,n = a->mbs,nz = a->i[n]; 1038f7157efSSatish Balay PetscErrorCode ierr; 104a6ece127SHong Zhang 10549b5e25fSSatish Balay PetscFunctionBegin; 10649b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 107a6ece127SHong Zhang 1088f7157efSSatish Balay if (!blockcompressed) { 1098f7157efSSatish Balay ierr = PetscFree2(*ia,*ja);CHKERRQ(ierr); 1108f7157efSSatish Balay } else if (oshift == 1) { /* blockcompressed */ 1116c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 112a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 113a6ece127SHong Zhang } 1148f7157efSSatish Balay 115a6ece127SHong Zhang PetscFunctionReturn(0); 11649b5e25fSSatish Balay } 11749b5e25fSSatish Balay 1184a2ae208SSatish Balay #undef __FUNCT__ 1194a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 120dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 12149b5e25fSSatish Balay { 12249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 123dfbe8321SBarry Smith PetscErrorCode ierr; 12449b5e25fSSatish Balay 12549b5e25fSSatish Balay PetscFunctionBegin; 126a9f03627SSatish Balay #if defined(PETSC_USE_LOG) 127d0f46423SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz); 128a9f03627SSatish Balay #endif 129e6b907acSBarry Smith ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr); 1307f53bb6cSHong Zhang if (a->free_diag){ierr = PetscFree(a->diag);CHKERRQ(ierr);} 1319bfd6278SHong Zhang if (a->row) {ierr = ISDestroy(a->row);CHKERRQ(ierr);} 1329bfd6278SHong Zhang if (a->col){ierr = ISDestroy(a->col);CHKERRQ(ierr);} 1339bfd6278SHong Zhang if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 134061b2667SBarry Smith if (a->idiag) {ierr = PetscFree(a->idiag);CHKERRQ(ierr);} 1350def2e27SBarry Smith if (a->inode.size) {ierr = PetscFree(a->inode.size);CHKERRQ(ierr);} 136c760cd28SBarry Smith if (a->free_diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 137c760cd28SBarry Smith if (a->free_imax_ilen) {ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);} 13805b42c5fSBarry Smith ierr = PetscFree(a->solve_work);CHKERRQ(ierr); 13941f059aeSBarry Smith ierr = PetscFree(a->sor_work);CHKERRQ(ierr); 14005b42c5fSBarry Smith ierr = PetscFree(a->solves_work);CHKERRQ(ierr); 14105b42c5fSBarry Smith ierr = PetscFree(a->mult_work);CHKERRQ(ierr); 14205b42c5fSBarry Smith ierr = PetscFree(a->saved_values);CHKERRQ(ierr); 14305b42c5fSBarry Smith ierr = PetscFree(a->xtoy);CHKERRQ(ierr); 1444da8f245SBarry Smith if (a->free_jshort) {ierr = PetscFree(a->jshort);CHKERRQ(ierr);} 1451a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1464da8f245SBarry Smith if (a->parent) {ierr = MatDestroy(a->parent);CHKERRQ(ierr);} 14749b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 148901853e0SKris Buschelman 149dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 150901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 151901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 152901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 153901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 154901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 155901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 15649b5e25fSSatish Balay PetscFunctionReturn(0); 15749b5e25fSSatish Balay } 15849b5e25fSSatish Balay 1594a2ae208SSatish Balay #undef __FUNCT__ 1604a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 1614e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscTruth flg) 16249b5e25fSSatish Balay { 163045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 16463ba0a88SBarry Smith PetscErrorCode ierr; 16549b5e25fSSatish Balay 16649b5e25fSSatish Balay PetscFunctionBegin; 1674d9d31abSKris Buschelman switch (op) { 1684d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1694e0d8c25SBarry Smith a->roworiented = flg; 1704d9d31abSKris Buschelman break; 171a9817697SBarry Smith case MAT_KEEP_NONZERO_PATTERN: 172a9817697SBarry Smith a->keepnonzeropattern = flg; 1734d9d31abSKris Buschelman break; 174512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 175512a5fc5SBarry Smith a->nonew = (flg ? 0 : 1); 1764d9d31abSKris Buschelman break; 1774d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1784e0d8c25SBarry Smith a->nonew = (flg ? -1 : 0); 1794d9d31abSKris Buschelman break; 1804d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1814e0d8c25SBarry Smith a->nonew = (flg ? -2 : 0); 1824d9d31abSKris Buschelman break; 18328b2fa4aSMatthew Knepley case MAT_UNUSED_NONZERO_LOCATION_ERR: 18428b2fa4aSMatthew Knepley a->nounused = (flg ? -1 : 0); 18528b2fa4aSMatthew Knepley break; 1864e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1874d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1884d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 189290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 1904d9d31abSKris Buschelman break; 1919a4540c5SBarry Smith case MAT_HERMITIAN: 192eeffb40dSHong Zhang if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first"); 1930bc54ff2SBarry Smith if (A->cmap->n < 65536 && A->cmap->bs == 1) { 194eeffb40dSHong Zhang A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort; 1950bc54ff2SBarry Smith } else if (A->cmap->bs == 1) { 196eeffb40dSHong Zhang A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian; 1970bc54ff2SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1"); 198eeffb40dSHong Zhang break; 19977e54ba9SKris Buschelman case MAT_SYMMETRIC: 20077e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 2019a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 2024e0d8c25SBarry Smith if (!flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 203290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);CHKERRQ(ierr); 204290bbb0aSBarry Smith break; 205941593c8SHong Zhang case MAT_IGNORE_LOWER_TRIANGULAR: 2064e0d8c25SBarry Smith a->ignore_ltriangular = flg; 207941593c8SHong Zhang break; 208941593c8SHong Zhang case MAT_ERROR_LOWER_TRIANGULAR: 2094e0d8c25SBarry Smith a->ignore_ltriangular = flg; 21077e54ba9SKris Buschelman break; 211f5edf698SHong Zhang case MAT_GETROW_UPPERTRIANGULAR: 2124e0d8c25SBarry Smith a->getrow_utriangular = flg; 213f5edf698SHong Zhang break; 2144d9d31abSKris Buschelman default: 215ad86a440SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 21649b5e25fSSatish Balay } 21749b5e25fSSatish Balay PetscFunctionReturn(0); 21849b5e25fSSatish Balay } 21949b5e25fSSatish Balay 2204a2ae208SSatish Balay #undef __FUNCT__ 2214a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 22213f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 22349b5e25fSSatish Balay { 22449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2256849ba73SBarry Smith PetscErrorCode ierr; 22613f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 22749b5e25fSSatish Balay MatScalar *aa,*aa_i; 22887828ca2SBarry Smith PetscScalar *v_i; 22949b5e25fSSatish Balay 23049b5e25fSSatish Balay PetscFunctionBegin; 2314e0d8c25SBarry Smith if (A && !a->getrow_utriangular) SETERRQ(PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()"); 232f5edf698SHong Zhang /* Get the upper triangular part of the row */ 233d0f46423SBarry Smith bs = A->rmap->bs; 23449b5e25fSSatish Balay ai = a->i; 23549b5e25fSSatish Balay aj = a->j; 23649b5e25fSSatish Balay aa = a->a; 23749b5e25fSSatish Balay bs2 = a->bs2; 23849b5e25fSSatish Balay 239d0f46423SBarry Smith if (row < 0 || row >= A->rmap->N) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 24049b5e25fSSatish Balay 24149b5e25fSSatish Balay bn = row/bs; /* Block number */ 24249b5e25fSSatish Balay bp = row % bs; /* Block position */ 24349b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 24449b5e25fSSatish Balay *ncols = bs*M; 24549b5e25fSSatish Balay 24649b5e25fSSatish Balay if (v) { 24749b5e25fSSatish Balay *v = 0; 24849b5e25fSSatish Balay if (*ncols) { 24987828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 25049b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 25149b5e25fSSatish Balay v_i = *v + i*bs; 25249b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 25349b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 25449b5e25fSSatish Balay } 25549b5e25fSSatish Balay } 25649b5e25fSSatish Balay } 25749b5e25fSSatish Balay 25849b5e25fSSatish Balay if (cols) { 25949b5e25fSSatish Balay *cols = 0; 26049b5e25fSSatish Balay if (*ncols) { 26113f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 26249b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 26349b5e25fSSatish Balay cols_i = *cols + i*bs; 26449b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 26549b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 26649b5e25fSSatish Balay } 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay } 26949b5e25fSSatish Balay 27049b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2715ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2725ddb2528SHong Zhang #ifdef column_search 27349b5e25fSSatish Balay v_i = *v + M*bs; 27449b5e25fSSatish Balay cols_i = *cols + M*bs; 27549b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 27649b5e25fSSatish Balay M = ai[i+1] - ai[i]; 27749b5e25fSSatish Balay for (j=0; j<M; j++){ 27849b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 27949b5e25fSSatish Balay if (itmp == bn){ 28049b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 28149b5e25fSSatish Balay for (k=0; k<bs; k++) { 28249b5e25fSSatish Balay *cols_i++ = i*bs+k; 28349b5e25fSSatish Balay *v_i++ = aa_i[k]; 28449b5e25fSSatish Balay } 28549b5e25fSSatish Balay *ncols += bs; 28649b5e25fSSatish Balay break; 28749b5e25fSSatish Balay } 28849b5e25fSSatish Balay } 28949b5e25fSSatish Balay } 2905ddb2528SHong Zhang #endif 29149b5e25fSSatish Balay PetscFunctionReturn(0); 29249b5e25fSSatish Balay } 29349b5e25fSSatish Balay 2944a2ae208SSatish Balay #undef __FUNCT__ 2954a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 29613f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 29749b5e25fSSatish Balay { 298dfbe8321SBarry Smith PetscErrorCode ierr; 29949b5e25fSSatish Balay 30049b5e25fSSatish Balay PetscFunctionBegin; 30105b42c5fSBarry Smith if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);} 30205b42c5fSBarry Smith if (v) {ierr = PetscFree(*v);CHKERRQ(ierr);} 30349b5e25fSSatish Balay PetscFunctionReturn(0); 30449b5e25fSSatish Balay } 30549b5e25fSSatish Balay 3064a2ae208SSatish Balay #undef __FUNCT__ 307f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ" 308f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A) 309f5edf698SHong Zhang { 310f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 311f5edf698SHong Zhang 312f5edf698SHong Zhang PetscFunctionBegin; 313f5edf698SHong Zhang a->getrow_utriangular = PETSC_TRUE; 314f5edf698SHong Zhang PetscFunctionReturn(0); 315f5edf698SHong Zhang } 316f5edf698SHong Zhang #undef __FUNCT__ 317f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ" 318f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A) 319f5edf698SHong Zhang { 320f5edf698SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 321f5edf698SHong Zhang 322f5edf698SHong Zhang PetscFunctionBegin; 323f5edf698SHong Zhang a->getrow_utriangular = PETSC_FALSE; 324f5edf698SHong Zhang PetscFunctionReturn(0); 325f5edf698SHong Zhang } 326f5edf698SHong Zhang 327f5edf698SHong Zhang #undef __FUNCT__ 3284a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 329fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B) 33049b5e25fSSatish Balay { 331dfbe8321SBarry Smith PetscErrorCode ierr; 33249b5e25fSSatish Balay PetscFunctionBegin; 333815cbec1SBarry Smith if (reuse == MAT_INITIAL_MATRIX || *B != A) { 334999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 335fc4dec0aSBarry Smith } 3368115998fSBarry Smith PetscFunctionReturn(0); 33749b5e25fSSatish Balay } 33849b5e25fSSatish Balay 3394a2ae208SSatish Balay #undef __FUNCT__ 3404a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 3416849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 34249b5e25fSSatish Balay { 34349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 344dfbe8321SBarry Smith PetscErrorCode ierr; 345d0f46423SBarry Smith PetscInt i,j,bs = A->rmap->bs,k,l,bs2=a->bs2; 3462dcb1b2aSMatthew Knepley const char *name; 347f3ef73ceSBarry Smith PetscViewerFormat format; 34849b5e25fSSatish Balay 34949b5e25fSSatish Balay PetscFunctionBegin; 35080fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 351b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 352456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 35377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 354fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 355d2507d54SMatthew Knepley Mat aij; 356d2507d54SMatthew Knepley 35770d5e725SHong Zhang if (A->factor && bs>1){ 35870d5e725SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr); 35970d5e725SHong Zhang PetscFunctionReturn(0); 36070d5e725SHong Zhang } 361c9f458caSMatthew Knepley ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr); 362c9f458caSMatthew Knepley ierr = MatView(aij,viewer);CHKERRQ(ierr); 363c9f458caSMatthew Knepley ierr = MatDestroy(aij);CHKERRQ(ierr); 364fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 365b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 36749b5e25fSSatish Balay for (j=0; j<bs; j++) { 36877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 36949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 37049b5e25fSSatish Balay for (l=0; l<bs; l++) { 37149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37249b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 373a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 37449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37549b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 376a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 37749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37849b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 379a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38049b5e25fSSatish Balay } 38149b5e25fSSatish Balay #else 38249b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 383a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38449b5e25fSSatish Balay } 38549b5e25fSSatish Balay #endif 38649b5e25fSSatish Balay } 38749b5e25fSSatish Balay } 388b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 38949b5e25fSSatish Balay } 39049b5e25fSSatish Balay } 391b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 392c1490034SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 393c1490034SHong Zhang PetscFunctionReturn(0); 39449b5e25fSSatish Balay } else { 3958608aa04SHong Zhang if (A->factor && bs>1){ 3968608aa04SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored. MatView_SeqSBAIJ_ASCII() may not display complete or logically correct entries!\n");CHKERRQ(ierr); 3978608aa04SHong Zhang } 398b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 39949b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 40049b5e25fSSatish Balay for (j=0; j<bs; j++) { 40177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 40249b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 40349b5e25fSSatish Balay for (l=0; l<bs; l++) { 40449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 40549b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 406a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l, 40749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40849b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 409a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l, 41049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 41149b5e25fSSatish Balay } else { 412a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 41349b5e25fSSatish Balay } 41449b5e25fSSatish Balay #else 415e9f7bc9eSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 41649b5e25fSSatish Balay #endif 41749b5e25fSSatish Balay } 41849b5e25fSSatish Balay } 419b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 42049b5e25fSSatish Balay } 42149b5e25fSSatish Balay } 422b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 42349b5e25fSSatish Balay } 424b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 42549b5e25fSSatish Balay PetscFunctionReturn(0); 42649b5e25fSSatish Balay } 42749b5e25fSSatish Balay 4284a2ae208SSatish Balay #undef __FUNCT__ 4294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 4306849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 43149b5e25fSSatish Balay { 43249b5e25fSSatish Balay Mat A = (Mat) Aa; 43349b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 4346849ba73SBarry Smith PetscErrorCode ierr; 435d0f46423SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2; 43613f74950SBarry Smith PetscMPIInt rank; 43749b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 43849b5e25fSSatish Balay MatScalar *aa; 43949b5e25fSSatish Balay MPI_Comm comm; 440b0a32e0cSBarry Smith PetscViewer viewer; 44149b5e25fSSatish Balay 44249b5e25fSSatish Balay PetscFunctionBegin; 44349b5e25fSSatish Balay /* 44449b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 44549b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 44649b5e25fSSatish Balay rest should return immediately. 44749b5e25fSSatish Balay */ 44849b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 44949b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 45049b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 45149b5e25fSSatish Balay 45249b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 45349b5e25fSSatish Balay 454b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 455b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 45649b5e25fSSatish Balay 45749b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 458b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 45949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 461d0f46423SBarry Smith y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 46249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46349b5e25fSSatish Balay aa = a->a + j*bs2; 46449b5e25fSSatish Balay for (k=0; k<bs; k++) { 46549b5e25fSSatish Balay for (l=0; l<bs; l++) { 46649b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 467b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46849b5e25fSSatish Balay } 46949b5e25fSSatish Balay } 47049b5e25fSSatish Balay } 47149b5e25fSSatish Balay } 472b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 47349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 47449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 475d0f46423SBarry Smith y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 47649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47749b5e25fSSatish Balay aa = a->a + j*bs2; 47849b5e25fSSatish Balay for (k=0; k<bs; k++) { 47949b5e25fSSatish Balay for (l=0; l<bs; l++) { 48049b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 481b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 48249b5e25fSSatish Balay } 48349b5e25fSSatish Balay } 48449b5e25fSSatish Balay } 48549b5e25fSSatish Balay } 48649b5e25fSSatish Balay 487b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 48849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 48949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 490d0f46423SBarry Smith y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0; 49149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 49249b5e25fSSatish Balay aa = a->a + j*bs2; 49349b5e25fSSatish Balay for (k=0; k<bs; k++) { 49449b5e25fSSatish Balay for (l=0; l<bs; l++) { 49549b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 496b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 49749b5e25fSSatish Balay } 49849b5e25fSSatish Balay } 49949b5e25fSSatish Balay } 50049b5e25fSSatish Balay } 50149b5e25fSSatish Balay PetscFunctionReturn(0); 50249b5e25fSSatish Balay } 50349b5e25fSSatish Balay 5044a2ae208SSatish Balay #undef __FUNCT__ 5054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 5066849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 50749b5e25fSSatish Balay { 508dfbe8321SBarry Smith PetscErrorCode ierr; 50949b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 510b0a32e0cSBarry Smith PetscDraw draw; 51149b5e25fSSatish Balay PetscTruth isnull; 51249b5e25fSSatish Balay 51349b5e25fSSatish Balay PetscFunctionBegin; 514b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 515b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 51649b5e25fSSatish Balay 51749b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 518d0f46423SBarry Smith xr = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0; 51949b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 520b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 521b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 52249b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 52349b5e25fSSatish Balay PetscFunctionReturn(0); 52449b5e25fSSatish Balay } 52549b5e25fSSatish Balay 5264a2ae208SSatish Balay #undef __FUNCT__ 5274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 528dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 52949b5e25fSSatish Balay { 530dfbe8321SBarry Smith PetscErrorCode ierr; 53132077d6dSBarry Smith PetscTruth iascii,isdraw; 53249b5e25fSSatish Balay 53349b5e25fSSatish Balay PetscFunctionBegin; 53432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 535fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 53632077d6dSBarry Smith if (iascii){ 53749b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 53849b5e25fSSatish Balay } else if (isdraw) { 53949b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 54049b5e25fSSatish Balay } else { 541a5e6ed63SBarry Smith Mat B; 542ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 543a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 544a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 54549b5e25fSSatish Balay } 54649b5e25fSSatish Balay PetscFunctionReturn(0); 54749b5e25fSSatish Balay } 54849b5e25fSSatish Balay 54949b5e25fSSatish Balay 5504a2ae208SSatish Balay #undef __FUNCT__ 5514a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 55213f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 55349b5e25fSSatish Balay { 554045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 55513f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 55613f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 557d0f46423SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2; 55897e567efSBarry Smith MatScalar *ap,*aa = a->a; 55949b5e25fSSatish Balay 56049b5e25fSSatish Balay PetscFunctionBegin; 56149b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 56249b5e25fSSatish Balay row = im[k]; brow = row/bs; 56397e567efSBarry Smith if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */ 564d0f46423SBarry Smith if (row >= A->rmap->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 56549b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 56649b5e25fSSatish Balay nrow = ailen[brow]; 56749b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 56897e567efSBarry Smith if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */ 569d0f46423SBarry Smith if (in[l] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1); 57049b5e25fSSatish Balay col = in[l] ; 57149b5e25fSSatish Balay bcol = col/bs; 57249b5e25fSSatish Balay cidx = col%bs; 57349b5e25fSSatish Balay ridx = row%bs; 57449b5e25fSSatish Balay high = nrow; 57549b5e25fSSatish Balay low = 0; /* assume unsorted */ 57649b5e25fSSatish Balay while (high-low > 5) { 57749b5e25fSSatish Balay t = (low+high)/2; 57849b5e25fSSatish Balay if (rp[t] > bcol) high = t; 57949b5e25fSSatish Balay else low = t; 58049b5e25fSSatish Balay } 58149b5e25fSSatish Balay for (i=low; i<high; i++) { 58249b5e25fSSatish Balay if (rp[i] > bcol) break; 58349b5e25fSSatish Balay if (rp[i] == bcol) { 58449b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 58549b5e25fSSatish Balay goto finished; 58649b5e25fSSatish Balay } 58749b5e25fSSatish Balay } 58897e567efSBarry Smith *v++ = 0.0; 58949b5e25fSSatish Balay finished:; 59049b5e25fSSatish Balay } 59149b5e25fSSatish Balay } 59249b5e25fSSatish Balay PetscFunctionReturn(0); 59349b5e25fSSatish Balay } 59449b5e25fSSatish Balay 59549b5e25fSSatish Balay 5964a2ae208SSatish Balay #undef __FUNCT__ 5974a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 59813f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 59949b5e25fSSatish Balay { 6000880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6016849ba73SBarry Smith PetscErrorCode ierr; 602e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 60313f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 604d0f46423SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval; 6050880e062SHong Zhang PetscTruth roworiented=a->roworiented; 606dd6ea824SBarry Smith const PetscScalar *value = v; 607f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 6080880e062SHong Zhang 60949b5e25fSSatish Balay PetscFunctionBegin; 6100880e062SHong Zhang if (roworiented) { 6110880e062SHong Zhang stepval = (n-1)*bs; 6120880e062SHong Zhang } else { 6130880e062SHong Zhang stepval = (m-1)*bs; 6140880e062SHong Zhang } 6150880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 6160880e062SHong Zhang row = im[k]; 6170880e062SHong Zhang if (row < 0) continue; 6182515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 61977431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 6200880e062SHong Zhang #endif 6210880e062SHong Zhang rp = aj + ai[row]; 6220880e062SHong Zhang ap = aa + bs2*ai[row]; 6230880e062SHong Zhang rmax = imax[row]; 6240880e062SHong Zhang nrow = ailen[row]; 6250880e062SHong Zhang low = 0; 626818f2c47SBarry Smith high = nrow; 6270880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6280880e062SHong Zhang if (in[l] < 0) continue; 6290880e062SHong Zhang col = in[l]; 6302515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 63177431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 632b1823623SSatish Balay #endif 633b98bf0e1SJed Brown if (col < row) { 634b98bf0e1SJed Brown if (a->ignore_ltriangular) { 635b98bf0e1SJed Brown continue; /* ignore lower triangular block */ 636b98bf0e1SJed Brown } else { 637b98bf0e1SJed Brown SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)"); 638b98bf0e1SJed Brown } 639b98bf0e1SJed Brown } 6400880e062SHong Zhang if (roworiented) { 6410880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6420880e062SHong Zhang } else { 6430880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6440880e062SHong Zhang } 6457cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 646e2ee6c50SBarry Smith lastcol = col; 6470880e062SHong Zhang while (high-low > 7) { 6480880e062SHong Zhang t = (low+high)/2; 6490880e062SHong Zhang if (rp[t] > col) high = t; 6500880e062SHong Zhang else low = t; 6510880e062SHong Zhang } 6520880e062SHong Zhang for (i=low; i<high; i++) { 6530880e062SHong Zhang if (rp[i] > col) break; 6540880e062SHong Zhang if (rp[i] == col) { 6550880e062SHong Zhang bap = ap + bs2*i; 6560880e062SHong Zhang if (roworiented) { 6570880e062SHong Zhang if (is == ADD_VALUES) { 6580880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6590880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6600880e062SHong Zhang bap[jj] += *value++; 6610880e062SHong Zhang } 6620880e062SHong Zhang } 6630880e062SHong Zhang } else { 6640880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6650880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6660880e062SHong Zhang bap[jj] = *value++; 6670880e062SHong Zhang } 6680880e062SHong Zhang } 6690880e062SHong Zhang } 6700880e062SHong Zhang } else { 6710880e062SHong Zhang if (is == ADD_VALUES) { 6720880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6730880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6740880e062SHong Zhang *bap++ += *value++; 6750880e062SHong Zhang } 6760880e062SHong Zhang } 6770880e062SHong Zhang } else { 6780880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6790880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6800880e062SHong Zhang *bap++ = *value++; 6810880e062SHong Zhang } 6820880e062SHong Zhang } 6830880e062SHong Zhang } 6840880e062SHong Zhang } 6850880e062SHong Zhang goto noinsert2; 6860880e062SHong Zhang } 6870880e062SHong Zhang } 6880880e062SHong Zhang if (nonew == 1) goto noinsert2; 689085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 690421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 691c03d1d03SSatish Balay N = nrow++ - 1; high++; 6920880e062SHong Zhang /* shift up all the later entries in this row */ 6930880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6940880e062SHong Zhang rp[ii+1] = rp[ii]; 6950880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6960880e062SHong Zhang } 6970880e062SHong Zhang if (N >= i) { 6980880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6990880e062SHong Zhang } 7000880e062SHong Zhang rp[i] = col; 7010880e062SHong Zhang bap = ap + bs2*i; 7020880e062SHong Zhang if (roworiented) { 7030880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7040880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 7050880e062SHong Zhang bap[jj] = *value++; 7060880e062SHong Zhang } 7070880e062SHong Zhang } 7080880e062SHong Zhang } else { 7090880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7100880e062SHong Zhang for (jj=0; jj<bs; jj++) { 7110880e062SHong Zhang *bap++ = *value++; 7120880e062SHong Zhang } 7130880e062SHong Zhang } 7140880e062SHong Zhang } 7150880e062SHong Zhang noinsert2:; 7160880e062SHong Zhang low = i; 7170880e062SHong Zhang } 7180880e062SHong Zhang ailen[row] = nrow; 7190880e062SHong Zhang } 7200880e062SHong Zhang PetscFunctionReturn(0); 72149b5e25fSSatish Balay } 72249b5e25fSSatish Balay 72364831d72SBarry Smith /* 72464831d72SBarry Smith This is not yet used 72564831d72SBarry Smith */ 7264a2ae208SSatish Balay #undef __FUNCT__ 7274108e4d5SBarry Smith #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode" 7284108e4d5SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A) 7290def2e27SBarry Smith { 7300def2e27SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7310def2e27SBarry Smith PetscErrorCode ierr; 7320def2e27SBarry Smith const PetscInt *ai = a->i, *aj = a->j,*cols; 7330def2e27SBarry Smith PetscInt i = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts; 7340def2e27SBarry Smith PetscTruth flag; 7350def2e27SBarry Smith 7360def2e27SBarry Smith PetscFunctionBegin; 7370def2e27SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&ns);CHKERRQ(ierr); 7380def2e27SBarry Smith while (i < m){ 7390def2e27SBarry Smith nzx = ai[i+1] - ai[i]; /* Number of nonzeros */ 7400def2e27SBarry Smith /* Limits the number of elements in a node to 'a->inode.limit' */ 7410def2e27SBarry Smith for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) { 7420def2e27SBarry Smith nzy = ai[j+1] - ai[j]; 7430def2e27SBarry Smith if (nzy != (nzx - j + i)) break; 7440def2e27SBarry Smith ierr = PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);CHKERRQ(ierr); 7450def2e27SBarry Smith if (!flag) break; 7460def2e27SBarry Smith } 7470def2e27SBarry Smith ns[node_count++] = blk_size; 7480def2e27SBarry Smith i = j; 7490def2e27SBarry Smith } 7500def2e27SBarry Smith if (!a->inode.size && m && node_count > .9*m) { 7510def2e27SBarry Smith ierr = PetscFree(ns);CHKERRQ(ierr); 7520def2e27SBarry Smith ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr); 7530def2e27SBarry Smith } else { 7540def2e27SBarry Smith a->inode.node_count = node_count; 7550def2e27SBarry Smith ierr = PetscMalloc(node_count*sizeof(PetscInt),&a->inode.size);CHKERRQ(ierr); 756c760cd28SBarry Smith ierr = PetscLogObjectMemory(A,node_count*sizeof(PetscInt));CHKERRQ(ierr); 7570def2e27SBarry Smith ierr = PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt)); 7580def2e27SBarry Smith ierr = PetscFree(ns);CHKERRQ(ierr); 7590def2e27SBarry Smith ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr); 7600def2e27SBarry Smith 7610def2e27SBarry Smith /* count collections of adjacent columns in each inode */ 7620def2e27SBarry Smith row = 0; 7630def2e27SBarry Smith cnt = 0; 7640def2e27SBarry Smith for (i=0; i<node_count; i++) { 7650def2e27SBarry Smith cols = aj + ai[row] + a->inode.size[i]; 7660def2e27SBarry Smith nz = ai[row+1] - ai[row] - a->inode.size[i]; 7670def2e27SBarry Smith for (j=1; j<nz; j++) { 7680def2e27SBarry Smith if (cols[j] != cols[j-1]+1) { 7690def2e27SBarry Smith cnt++; 7700def2e27SBarry Smith } 7710def2e27SBarry Smith } 7720def2e27SBarry Smith cnt++; 7730def2e27SBarry Smith row += a->inode.size[i]; 7740def2e27SBarry Smith } 7750def2e27SBarry Smith ierr = PetscMalloc(2*cnt*sizeof(PetscInt),&counts);CHKERRQ(ierr); 7760def2e27SBarry Smith cnt = 0; 7770def2e27SBarry Smith row = 0; 7780def2e27SBarry Smith for (i=0; i<node_count; i++) { 7790def2e27SBarry Smith cols = aj + ai[row] + a->inode.size[i]; 7800def2e27SBarry Smith CHKMEMQ; 7810def2e27SBarry Smith counts[2*cnt] = cols[0]; 7820def2e27SBarry Smith CHKMEMQ; 7830def2e27SBarry Smith nz = ai[row+1] - ai[row] - a->inode.size[i]; 7840def2e27SBarry Smith cnt2 = 1; 7850def2e27SBarry Smith for (j=1; j<nz; j++) { 7860def2e27SBarry Smith if (cols[j] != cols[j-1]+1) { 7870def2e27SBarry Smith CHKMEMQ; 7880def2e27SBarry Smith counts[2*(cnt++)+1] = cnt2; 7890def2e27SBarry Smith counts[2*cnt] = cols[j]; 7900def2e27SBarry Smith CHKMEMQ; 7910def2e27SBarry Smith cnt2 = 1; 7920def2e27SBarry Smith } else cnt2++; 7930def2e27SBarry Smith } 7940def2e27SBarry Smith CHKMEMQ; 7950def2e27SBarry Smith counts[2*(cnt++)+1] = cnt2; 7960def2e27SBarry Smith CHKMEMQ; 7970def2e27SBarry Smith row += a->inode.size[i]; 7980def2e27SBarry Smith } 7990def2e27SBarry Smith ierr = PetscIntView(2*cnt,counts,0); 8000def2e27SBarry Smith } 80138702af4SBarry Smith PetscFunctionReturn(0); 80238702af4SBarry Smith } 80338702af4SBarry Smith 80438702af4SBarry Smith #undef __FUNCT__ 8054a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 806dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 80749b5e25fSSatish Balay { 80849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8096849ba73SBarry Smith PetscErrorCode ierr; 81013f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 811d0f46423SBarry Smith PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen; 81213f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 81349b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 81449b5e25fSSatish Balay 81549b5e25fSSatish Balay PetscFunctionBegin; 81649b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 81749b5e25fSSatish Balay 81849b5e25fSSatish Balay if (m) rmax = ailen[0]; 81949b5e25fSSatish Balay for (i=1; i<mbs; i++) { 82049b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 82149b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 82249b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 82349b5e25fSSatish Balay if (fshift) { 82449b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 82549b5e25fSSatish Balay N = ailen[i]; 82649b5e25fSSatish Balay for (j=0; j<N; j++) { 82749b5e25fSSatish Balay ip[j-fshift] = ip[j]; 82849b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 82949b5e25fSSatish Balay } 83049b5e25fSSatish Balay } 83149b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 83249b5e25fSSatish Balay } 83349b5e25fSSatish Balay if (mbs) { 83449b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 83549b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 83649b5e25fSSatish Balay } 83749b5e25fSSatish Balay /* reset ilen and imax for each row */ 83849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 83949b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 84049b5e25fSSatish Balay } 8416c6c5352SBarry Smith a->nz = ai[mbs]; 84249b5e25fSSatish Balay 843b424e231SHong Zhang /* diagonals may have moved, reset it */ 844b424e231SHong Zhang if (a->diag) { 8452ed38d0bSJed Brown ierr = PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));CHKERRQ(ierr); 84649b5e25fSSatish Balay } 84728b2fa4aSMatthew Knepley if (fshift && a->nounused == -1) { 84828b2fa4aSMatthew Knepley SETERRQ4(PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D block size %D, %D unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2); 84928b2fa4aSMatthew Knepley } 850d0f46423SBarry Smith ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap->N,A->rmap->bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr); 851ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 852ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 85349b5e25fSSatish Balay a->reallocs = 0; 85449b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 855061b2667SBarry Smith a->idiagvalid = PETSC_FALSE; 85638702af4SBarry Smith 85738702af4SBarry Smith if (A->cmap->n < 65536 && A->cmap->bs == 1) { 85838702af4SBarry Smith if (!a->jshort) { 85938702af4SBarry Smith ierr = PetscMalloc(a->i[A->rmap->n]*sizeof(unsigned short),&a->jshort);CHKERRQ(ierr); 860c760cd28SBarry Smith ierr = PetscLogObjectMemory(A,a->i[A->rmap->n]*sizeof(unsigned short));CHKERRQ(ierr); 86138702af4SBarry Smith for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i]; 86238702af4SBarry Smith A->ops->mult = MatMult_SeqSBAIJ_1_ushort; 86341f059aeSBarry Smith A->ops->sor = MatSOR_SeqSBAIJ_ushort; 8644da8f245SBarry Smith a->free_jshort = PETSC_TRUE; 86538702af4SBarry Smith } 86638702af4SBarry Smith } 86749b5e25fSSatish Balay PetscFunctionReturn(0); 86849b5e25fSSatish Balay } 86949b5e25fSSatish Balay 87049b5e25fSSatish Balay /* 87149b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 87249b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 87349b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 87449b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 87549b5e25fSSatish Balay */ 8764a2ae208SSatish Balay #undef __FUNCT__ 8774a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 87813f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 87949b5e25fSSatish Balay { 88013f74950SBarry Smith PetscInt i,j,k,row; 88149b5e25fSSatish Balay PetscTruth flg; 88249b5e25fSSatish Balay 88349b5e25fSSatish Balay PetscFunctionBegin; 88449b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 88549b5e25fSSatish Balay row = idx[i]; 88649b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 88749b5e25fSSatish Balay sizes[j] = 1; 88849b5e25fSSatish Balay i++; 88949b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 89049b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 89149b5e25fSSatish Balay i++; 89249b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 89349b5e25fSSatish Balay flg = PETSC_TRUE; 89449b5e25fSSatish Balay for (k=1; k<bs; k++) { 89549b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 89649b5e25fSSatish Balay flg = PETSC_FALSE; 89749b5e25fSSatish Balay break; 89849b5e25fSSatish Balay } 89949b5e25fSSatish Balay } 900abc0a331SBarry Smith if (flg) { /* No break in the bs */ 90149b5e25fSSatish Balay sizes[j] = bs; 90249b5e25fSSatish Balay i+= bs; 90349b5e25fSSatish Balay } else { 90449b5e25fSSatish Balay sizes[j] = 1; 90549b5e25fSSatish Balay i++; 90649b5e25fSSatish Balay } 90749b5e25fSSatish Balay } 90849b5e25fSSatish Balay } 90949b5e25fSSatish Balay *bs_max = j; 91049b5e25fSSatish Balay PetscFunctionReturn(0); 91149b5e25fSSatish Balay } 91249b5e25fSSatish Balay 91349b5e25fSSatish Balay 91449b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 91549b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 91649b5e25fSSatish Balay */ 91749b5e25fSSatish Balay 9184a2ae208SSatish Balay #undef __FUNCT__ 9194a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 92013f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 92149b5e25fSSatish Balay { 92249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 9236849ba73SBarry Smith PetscErrorCode ierr; 924e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 92513f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 926d0f46423SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 92713f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 92849b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 92949b5e25fSSatish Balay 93049b5e25fSSatish Balay PetscFunctionBegin; 93171fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 93249b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 93349b5e25fSSatish Balay row = im[k]; /* row number */ 93449b5e25fSSatish Balay brow = row/bs; /* block row number */ 93549b5e25fSSatish Balay if (row < 0) continue; 9362515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 937d0f46423SBarry Smith if (row >= A->rmap->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 93849b5e25fSSatish Balay #endif 93949b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 94049b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 94149b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 94249b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 94349b5e25fSSatish Balay low = 0; 94449b5e25fSSatish Balay 94549b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 94649b5e25fSSatish Balay if (in[l] < 0) continue; 9472515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 948d0f46423SBarry Smith if (in[l] >= A->rmap->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap->N-1); 94949b5e25fSSatish Balay #endif 95049b5e25fSSatish Balay col = in[l]; 95149b5e25fSSatish Balay bcol = col/bs; /* block col number */ 95249b5e25fSSatish Balay 953941593c8SHong Zhang if (brow > bcol) { 954941593c8SHong Zhang if (a->ignore_ltriangular){ 955941593c8SHong Zhang continue; /* ignore lower triangular values */ 956941593c8SHong Zhang } else { 9574e0d8c25SBarry Smith SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)"); 958941593c8SHong Zhang } 959941593c8SHong Zhang } 960f4989cb3SHong Zhang 96149b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 9628549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 96349b5e25fSSatish Balay /* element value a(k,l) */ 96449b5e25fSSatish Balay if (roworiented) { 96549b5e25fSSatish Balay value = v[l + k*n]; 96649b5e25fSSatish Balay } else { 96749b5e25fSSatish Balay value = v[k + l*m]; 96849b5e25fSSatish Balay } 96949b5e25fSSatish Balay 97049b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 9717cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 972e2ee6c50SBarry Smith lastcol = col; 97349b5e25fSSatish Balay while (high-low > 7) { 97449b5e25fSSatish Balay t = (low+high)/2; 97549b5e25fSSatish Balay if (rp[t] > bcol) high = t; 97649b5e25fSSatish Balay else low = t; 97749b5e25fSSatish Balay } 97849b5e25fSSatish Balay for (i=low; i<high; i++) { 97949b5e25fSSatish Balay if (rp[i] > bcol) break; 98049b5e25fSSatish Balay if (rp[i] == bcol) { 98149b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 98249b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 98349b5e25fSSatish Balay else *bap = value; 9848549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 9858549e402SHong Zhang if (brow == bcol && ridx < cidx){ 9868549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 9878549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 9888549e402SHong Zhang else *bap = value; 9898549e402SHong Zhang } 99049b5e25fSSatish Balay goto noinsert1; 99149b5e25fSSatish Balay } 99249b5e25fSSatish Balay } 99349b5e25fSSatish Balay 99449b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 995085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 996421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 99749b5e25fSSatish Balay 998c03d1d03SSatish Balay N = nrow++ - 1; high++; 99949b5e25fSSatish Balay /* shift up all the later entries in this row */ 100049b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 100149b5e25fSSatish Balay rp[ii+1] = rp[ii]; 100249b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 100349b5e25fSSatish Balay } 100449b5e25fSSatish Balay if (N>=i) { 100549b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 100649b5e25fSSatish Balay } 100749b5e25fSSatish Balay rp[i] = bcol; 100849b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 100949b5e25fSSatish Balay noinsert1:; 101049b5e25fSSatish Balay low = i; 10118549e402SHong Zhang } 101249b5e25fSSatish Balay } /* end of loop over added columns */ 101349b5e25fSSatish Balay ailen[brow] = nrow; 101449b5e25fSSatish Balay } /* end of loop over added rows */ 101549b5e25fSSatish Balay PetscFunctionReturn(0); 101649b5e25fSSatish Balay } 101749b5e25fSSatish Balay 10184a2ae208SSatish Balay #undef __FUNCT__ 10194d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 10200481f469SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info) 102149b5e25fSSatish Balay { 10224ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 102349b5e25fSSatish Balay Mat outA; 1024dfbe8321SBarry Smith PetscErrorCode ierr; 1025c84f5b01SHong Zhang PetscTruth row_identity; 102649b5e25fSSatish Balay 102749b5e25fSSatish Balay PetscFunctionBegin; 1028c84f5b01SHong Zhang if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 1029c84f5b01SHong Zhang ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1030c84f5b01SHong Zhang if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 1031d0f46423SBarry Smith if (inA->rmap->bs != 1) SETERRQ1(PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */ 1032c84f5b01SHong Zhang 103349b5e25fSSatish Balay outA = inA; 1034c078aec8SLisandro Dalcin inA->factor = MAT_FACTOR_ICC; 103549b5e25fSSatish Balay 10361a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 1037db4efbfdSBarry Smith ierr = MatSeqSBAIJSetNumericFactorization(inA,row_identity);CHKERRQ(ierr); 103849b5e25fSSatish Balay 1039c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1040c3122656SLisandro Dalcin if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } 1041c84f5b01SHong Zhang a->row = row; 1042c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1043c3122656SLisandro Dalcin if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } 1044c84f5b01SHong Zhang a->col = row; 1045c84f5b01SHong Zhang 1046c84f5b01SHong Zhang /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 1047c84f5b01SHong Zhang if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 104852e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 104949b5e25fSSatish Balay 105049b5e25fSSatish Balay if (!a->solve_work) { 1051d0f46423SBarry Smith ierr = PetscMalloc((inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 1052d0f46423SBarry Smith ierr = PetscLogObjectMemory(inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));CHKERRQ(ierr); 105349b5e25fSSatish Balay } 105449b5e25fSSatish Balay 1055719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(outA,inA,info);CHKERRQ(ierr); 105649b5e25fSSatish Balay PetscFunctionReturn(0); 105749b5e25fSSatish Balay } 1058950f1e5bSHong Zhang 105949b5e25fSSatish Balay EXTERN_C_BEGIN 10604a2ae208SSatish Balay #undef __FUNCT__ 10614a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1062be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 106349b5e25fSSatish Balay { 1064045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 106513f74950SBarry Smith PetscInt i,nz,n; 106649b5e25fSSatish Balay 106749b5e25fSSatish Balay PetscFunctionBegin; 10686c6c5352SBarry Smith nz = baij->maxnz; 1069d0f46423SBarry Smith n = mat->cmap->n; 107049b5e25fSSatish Balay for (i=0; i<nz; i++) { 107149b5e25fSSatish Balay baij->j[i] = indices[i]; 107249b5e25fSSatish Balay } 10736c6c5352SBarry Smith baij->nz = nz; 107449b5e25fSSatish Balay for (i=0; i<n; i++) { 107549b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 107649b5e25fSSatish Balay } 107749b5e25fSSatish Balay PetscFunctionReturn(0); 107849b5e25fSSatish Balay } 107949b5e25fSSatish Balay EXTERN_C_END 108049b5e25fSSatish Balay 10814a2ae208SSatish Balay #undef __FUNCT__ 10824a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 108349b5e25fSSatish Balay /*@ 108419585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 108549b5e25fSSatish Balay in the matrix. 108649b5e25fSSatish Balay 108749b5e25fSSatish Balay Input Parameters: 108819585528SSatish Balay + mat - the SeqSBAIJ matrix 108949b5e25fSSatish Balay - indices - the column indices 109049b5e25fSSatish Balay 109149b5e25fSSatish Balay Level: advanced 109249b5e25fSSatish Balay 109349b5e25fSSatish Balay Notes: 109449b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 109549b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 109649b5e25fSSatish Balay of the MatSetValues() operation. 109749b5e25fSSatish Balay 109849b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1099d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 110049b5e25fSSatish Balay 1101ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 110249b5e25fSSatish Balay 1103ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 110449b5e25fSSatish Balay @*/ 1105be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 110649b5e25fSSatish Balay { 110713f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 110849b5e25fSSatish Balay 110949b5e25fSSatish Balay PetscFunctionBegin; 11104482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 11114482741eSBarry Smith PetscValidPointer(indices,2); 1112c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 111349b5e25fSSatish Balay if (f) { 111449b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 111549b5e25fSSatish Balay } else { 1116e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 111749b5e25fSSatish Balay } 111849b5e25fSSatish Balay PetscFunctionReturn(0); 111949b5e25fSSatish Balay } 112049b5e25fSSatish Balay 11214a2ae208SSatish Balay #undef __FUNCT__ 11223c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 11233c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 11243c896bc6SHong Zhang { 11253c896bc6SHong Zhang PetscErrorCode ierr; 11263c896bc6SHong Zhang 11273c896bc6SHong Zhang PetscFunctionBegin; 11283c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 11293c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 11303c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 11313c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 11323c896bc6SHong Zhang 1133d0f46423SBarry Smith if (a->i[A->rmap->N] != b->i[B->rmap->N]) { 11343c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 11353c896bc6SHong Zhang } 1136d0f46423SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));CHKERRQ(ierr); 11373c896bc6SHong Zhang } else { 1138f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 11393c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1140f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 11413c896bc6SHong Zhang } 11423c896bc6SHong Zhang PetscFunctionReturn(0); 11433c896bc6SHong Zhang } 11443c896bc6SHong Zhang 11453c896bc6SHong Zhang #undef __FUNCT__ 11464a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1147dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1148273d9f13SBarry Smith { 1149dfbe8321SBarry Smith PetscErrorCode ierr; 1150273d9f13SBarry Smith 1151273d9f13SBarry Smith PetscFunctionBegin; 1152db4efbfdSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,-PetscMax(A->rmap->bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1153273d9f13SBarry Smith PetscFunctionReturn(0); 1154273d9f13SBarry Smith } 1155273d9f13SBarry Smith 1156a6ece127SHong Zhang #undef __FUNCT__ 1157a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1158dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1159a6ece127SHong Zhang { 1160a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1161a6ece127SHong Zhang PetscFunctionBegin; 1162a6ece127SHong Zhang *array = a->a; 1163a6ece127SHong Zhang PetscFunctionReturn(0); 1164a6ece127SHong Zhang } 1165a6ece127SHong Zhang 1166a6ece127SHong Zhang #undef __FUNCT__ 1167a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1168dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1169a6ece127SHong Zhang { 1170a6ece127SHong Zhang PetscFunctionBegin; 1171a6ece127SHong Zhang PetscFunctionReturn(0); 1172a6ece127SHong Zhang } 1173a6ece127SHong Zhang 117442ee4b1aSHong Zhang #include "petscblaslapack.h" 117542ee4b1aSHong Zhang #undef __FUNCT__ 117642ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1177f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 117842ee4b1aSHong Zhang { 117942ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1180dfbe8321SBarry Smith PetscErrorCode ierr; 1181d0f46423SBarry Smith PetscInt i,bs=Y->rmap->bs,bs2,j; 11820805154bSBarry Smith PetscBLASInt one = 1,bnz = PetscBLASIntCast(x->nz); 118342ee4b1aSHong Zhang 118442ee4b1aSHong Zhang PetscFunctionBegin; 118542ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1186f4df32b1SMatthew Knepley PetscScalar alpha = a; 1187f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1188c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1189c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1190c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1191c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1192c537a176SHong Zhang } 1193c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1194c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1195c4319e64SHong Zhang y->XtoY = X; 1196c537a176SHong Zhang } 1197c4319e64SHong Zhang bs2 = bs*bs; 11986c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1199c4319e64SHong Zhang j = 0; 1200c4319e64SHong Zhang while (j < bs2){ 1201f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1202c4319e64SHong Zhang j++; 1203c537a176SHong Zhang } 1204c4319e64SHong Zhang } 12051e2582c4SBarry Smith ierr = PetscInfo3(Y,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %G\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));CHKERRQ(ierr); 120642ee4b1aSHong Zhang } else { 1207f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1208f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1209f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 121042ee4b1aSHong Zhang } 121142ee4b1aSHong Zhang PetscFunctionReturn(0); 121242ee4b1aSHong Zhang } 121342ee4b1aSHong Zhang 1214efcf0fc3SBarry Smith #undef __FUNCT__ 12156363de48SJed Brown #define __FUNCT__ "MatSetBlockSize_SeqSBAIJ" 12166363de48SJed Brown PetscErrorCode MatSetBlockSize_SeqSBAIJ(Mat A,PetscInt bs) 12176363de48SJed Brown { 12186363de48SJed Brown PetscInt rbs,cbs; 12196363de48SJed Brown PetscErrorCode ierr; 12206363de48SJed Brown 12216363de48SJed Brown PetscFunctionBegin; 12226363de48SJed Brown ierr = PetscLayoutGetBlockSize(A->rmap,&rbs);CHKERRQ(ierr); 12236363de48SJed Brown ierr = PetscLayoutGetBlockSize(A->cmap,&cbs);CHKERRQ(ierr); 12246363de48SJed Brown if (rbs != bs) SETERRQ2(PETSC_ERR_ARG_SIZ,"Attempt to set block size %d with SBAIJ %d",bs,rbs); 12256363de48SJed Brown if (cbs != bs) SETERRQ2(PETSC_ERR_ARG_SIZ,"Attempt to set block size %d with SBAIJ %d",bs,cbs); 12266363de48SJed Brown PetscFunctionReturn(0); 12276363de48SJed Brown } 12286363de48SJed Brown 12296363de48SJed Brown #undef __FUNCT__ 1230efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1231dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1232efcf0fc3SBarry Smith { 1233efcf0fc3SBarry Smith PetscFunctionBegin; 1234efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1235efcf0fc3SBarry Smith PetscFunctionReturn(0); 1236efcf0fc3SBarry Smith } 1237efcf0fc3SBarry Smith 1238efcf0fc3SBarry Smith #undef __FUNCT__ 1239efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1240dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1241efcf0fc3SBarry Smith { 1242efcf0fc3SBarry Smith PetscFunctionBegin; 1243efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1244efcf0fc3SBarry Smith PetscFunctionReturn(0); 1245efcf0fc3SBarry Smith } 1246efcf0fc3SBarry Smith 1247efcf0fc3SBarry Smith #undef __FUNCT__ 1248efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1249ab5e4463SMatthew Knepley PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1250efcf0fc3SBarry Smith { 1251efcf0fc3SBarry Smith PetscFunctionBegin; 1252efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1253efcf0fc3SBarry Smith PetscFunctionReturn(0); 1254efcf0fc3SBarry Smith } 1255efcf0fc3SBarry Smith 125699cafbc1SBarry Smith #undef __FUNCT__ 125799cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 125899cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 125999cafbc1SBarry Smith { 126099cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 126199cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 1262dd6ea824SBarry Smith MatScalar *aa = a->a; 126399cafbc1SBarry Smith 126499cafbc1SBarry Smith PetscFunctionBegin; 126599cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 126699cafbc1SBarry Smith PetscFunctionReturn(0); 126799cafbc1SBarry Smith } 126899cafbc1SBarry Smith 126999cafbc1SBarry Smith #undef __FUNCT__ 127099cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 127199cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 127299cafbc1SBarry Smith { 127399cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 127499cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 1275dd6ea824SBarry Smith MatScalar *aa = a->a; 127699cafbc1SBarry Smith 127799cafbc1SBarry Smith PetscFunctionBegin; 127899cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 127999cafbc1SBarry Smith PetscFunctionReturn(0); 128099cafbc1SBarry Smith } 128199cafbc1SBarry Smith 128249b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 128349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 128449b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 128549b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 128649b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 128797304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1288431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1289e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 1290db4efbfdSBarry Smith 0, 129149b5e25fSSatish Balay 0, 129249b5e25fSSatish Balay 0, 129397304618SKris Buschelman /*10*/ 0, 129449b5e25fSSatish Balay 0, 1295c078aec8SLisandro Dalcin MatCholeskyFactor_SeqSBAIJ, 129641f059aeSBarry Smith MatSOR_SeqSBAIJ, 129749b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 129897304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 129949b5e25fSSatish Balay MatEqual_SeqSBAIJ, 130049b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 130149b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 130249b5e25fSSatish Balay MatNorm_SeqSBAIJ, 130397304618SKris Buschelman /*20*/ 0, 130449b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 130549b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 130649b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1307d519adbfSMatthew Knepley /*24*/ 0, 130849b5e25fSSatish Balay 0, 130949b5e25fSSatish Balay 0, 1310db4efbfdSBarry Smith 0, 1311db4efbfdSBarry Smith 0, 1312d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqSBAIJ, 1313c464158bSHong Zhang 0, 1314db4efbfdSBarry Smith 0, 1315a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1316a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 1317d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqSBAIJ, 1318719d5645SBarry Smith 0, 1319719d5645SBarry Smith 0, 132049b5e25fSSatish Balay 0, 1321c84f5b01SHong Zhang MatICCFactor_SeqSBAIJ, 1322d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqSBAIJ, 132349b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 132449b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 132549b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 13263c896bc6SHong Zhang MatCopy_SeqSBAIJ, 1327d519adbfSMatthew Knepley /*44*/ 0, 132849b5e25fSSatish Balay MatScale_SeqSBAIJ, 132949b5e25fSSatish Balay 0, 133049b5e25fSSatish Balay 0, 133149b5e25fSSatish Balay 0, 13326363de48SJed Brown /*49*/ MatSetBlockSize_SeqSBAIJ, 133349b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 133449b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 133549b5e25fSSatish Balay 0, 133649b5e25fSSatish Balay 0, 1337d519adbfSMatthew Knepley /*54*/ 0, 133849b5e25fSSatish Balay 0, 133949b5e25fSSatish Balay 0, 134049b5e25fSSatish Balay 0, 134149b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 1342d519adbfSMatthew Knepley /*59*/ MatGetSubMatrix_SeqSBAIJ, 134349b5e25fSSatish Balay 0, 134449b5e25fSSatish Balay 0, 1345357abbc8SBarry Smith 0, 1346d959ec07SHong Zhang 0, 1347d519adbfSMatthew Knepley /*64*/ 0, 1348d959ec07SHong Zhang 0, 1349d959ec07SHong Zhang 0, 1350d959ec07SHong Zhang 0, 1351d959ec07SHong Zhang 0, 1352d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqSBAIJ, 13533e0d88b5SBarry Smith 0, 13543e0d88b5SBarry Smith 0, 13553e0d88b5SBarry Smith 0, 13563e0d88b5SBarry Smith 0, 1357d519adbfSMatthew Knepley /*74*/ 0, 13583e0d88b5SBarry Smith 0, 13593e0d88b5SBarry Smith 0, 13603e0d88b5SBarry Smith 0, 13613e0d88b5SBarry Smith 0, 1362d519adbfSMatthew Knepley /*79*/ 0, 13633e0d88b5SBarry Smith 0, 13643e0d88b5SBarry Smith 0, 13653e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 136697304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13673e0d88b5SBarry Smith #else 136897304618SKris Buschelman 0, 13693e0d88b5SBarry Smith #endif 1370865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1371d519adbfSMatthew Knepley /*84*/ MatIsSymmetric_SeqSBAIJ, 1372865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1373efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1374865e5f61SKris Buschelman 0, 1375865e5f61SKris Buschelman 0, 1376d519adbfSMatthew Knepley /*89*/ 0, 1377865e5f61SKris Buschelman 0, 1378865e5f61SKris Buschelman 0, 1379865e5f61SKris Buschelman 0, 1380865e5f61SKris Buschelman 0, 1381d519adbfSMatthew Knepley /*94*/ 0, 1382865e5f61SKris Buschelman 0, 1383865e5f61SKris Buschelman 0, 138499cafbc1SBarry Smith 0, 138599cafbc1SBarry Smith 0, 1386d519adbfSMatthew Knepley /*99*/ 0, 138799cafbc1SBarry Smith 0, 138899cafbc1SBarry Smith 0, 138999cafbc1SBarry Smith 0, 139099cafbc1SBarry Smith 0, 1391d519adbfSMatthew Knepley /*104*/0, 139299cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1393f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1394f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 13952af78befSBarry Smith MatRestoreRowUpperTriangular_SeqSBAIJ, 1396d519adbfSMatthew Knepley /*109*/0, 13972af78befSBarry Smith 0, 13982af78befSBarry Smith 0, 13992af78befSBarry Smith 0, 1400547795f9SHong Zhang MatMissingDiagonal_SeqSBAIJ, 1401547795f9SHong Zhang /*114*/0, 1402547795f9SHong Zhang 0, 1403547795f9SHong Zhang 0, 1404547795f9SHong Zhang 0, 1405547795f9SHong Zhang 0, 1406547795f9SHong Zhang /*119*/0, 1407547795f9SHong Zhang 0, 1408547795f9SHong Zhang 0 140999cafbc1SBarry Smith }; 1410be1d678aSKris Buschelman 141149b5e25fSSatish Balay EXTERN_C_BEGIN 14124a2ae208SSatish Balay #undef __FUNCT__ 14134a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1414be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 141549b5e25fSSatish Balay { 14164afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1417d0f46423SBarry Smith PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2; 1418dfbe8321SBarry Smith PetscErrorCode ierr; 141949b5e25fSSatish Balay 142049b5e25fSSatish Balay PetscFunctionBegin; 142149b5e25fSSatish Balay if (aij->nonew != 1) { 1422512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 142349b5e25fSSatish Balay } 142449b5e25fSSatish Balay 142549b5e25fSSatish Balay /* allocate space for values if not already there */ 142649b5e25fSSatish Balay if (!aij->saved_values) { 142787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 142849b5e25fSSatish Balay } 142949b5e25fSSatish Balay 143049b5e25fSSatish Balay /* copy values over */ 143187828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 143249b5e25fSSatish Balay PetscFunctionReturn(0); 143349b5e25fSSatish Balay } 143449b5e25fSSatish Balay EXTERN_C_END 143549b5e25fSSatish Balay 143649b5e25fSSatish Balay EXTERN_C_BEGIN 14374a2ae208SSatish Balay #undef __FUNCT__ 14384a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1439be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 144049b5e25fSSatish Balay { 14414afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 14426849ba73SBarry Smith PetscErrorCode ierr; 1443d0f46423SBarry Smith PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2; 144449b5e25fSSatish Balay 144549b5e25fSSatish Balay PetscFunctionBegin; 144649b5e25fSSatish Balay if (aij->nonew != 1) { 1447512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 144849b5e25fSSatish Balay } 144949b5e25fSSatish Balay if (!aij->saved_values) { 1450e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 145149b5e25fSSatish Balay } 145249b5e25fSSatish Balay 145349b5e25fSSatish Balay /* copy values over */ 145487828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 145549b5e25fSSatish Balay PetscFunctionReturn(0); 145649b5e25fSSatish Balay } 145749b5e25fSSatish Balay EXTERN_C_END 145849b5e25fSSatish Balay 14598549e402SHong Zhang EXTERN_C_BEGIN 14604a2ae208SSatish Balay #undef __FUNCT__ 1461a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1462be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 146349b5e25fSSatish Balay { 1464c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14656849ba73SBarry Smith PetscErrorCode ierr; 1466db4efbfdSBarry Smith PetscInt i,mbs,bs2, newbs = PetscAbs(bs); 146790d69ab7SBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg = PETSC_FALSE; 146849b5e25fSSatish Balay 146949b5e25fSSatish Balay PetscFunctionBegin; 1470273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1471db4efbfdSBarry Smith if (bs < 0) { 1472db4efbfdSBarry Smith ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for MPISBAIJ matrix","Mat");CHKERRQ(ierr); 1473db4efbfdSBarry Smith ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSeqSBAIJSetPreallocation",newbs,&newbs,PETSC_NULL);CHKERRQ(ierr); 1474db4efbfdSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 1475db4efbfdSBarry Smith bs = PetscAbs(bs); 1476db4efbfdSBarry Smith } 1477db4efbfdSBarry Smith if (nnz && newbs != bs) { 1478db4efbfdSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot change blocksize from command line if setting nnz"); 1479db4efbfdSBarry Smith } 1480db4efbfdSBarry Smith bs = newbs; 1481db4efbfdSBarry Smith 148226283091SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr); 148326283091SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr); 148426283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 148526283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 1486899cda47SBarry Smith 1487d0f46423SBarry Smith mbs = B->rmap->N/bs; 148849b5e25fSSatish Balay bs2 = bs*bs; 148949b5e25fSSatish Balay 1490d0f46423SBarry Smith if (mbs*bs != B->rmap->N) { 149129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 149249b5e25fSSatish Balay } 149349b5e25fSSatish Balay 1494ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1495ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1496ab93d7beSBarry Smith nz = 0; 1497ab93d7beSBarry Smith } 1498ab93d7beSBarry Smith 1499435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 150077431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 150149b5e25fSSatish Balay if (nnz) { 150249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 150377431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 150477431f27SBarry 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); 150549b5e25fSSatish Balay } 150649b5e25fSSatish Balay } 150749b5e25fSSatish Balay 1508db4efbfdSBarry Smith B->ops->mult = MatMult_SeqSBAIJ_N; 1509db4efbfdSBarry Smith B->ops->multadd = MatMultAdd_SeqSBAIJ_N; 1510db4efbfdSBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_N; 1511db4efbfdSBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N; 151290d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,PETSC_NULL);CHKERRQ(ierr); 151349b5e25fSSatish Balay if (!flg) { 151449b5e25fSSatish Balay switch (bs) { 151549b5e25fSSatish Balay case 1: 151649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 151749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1518431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1519431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 152049b5e25fSSatish Balay break; 152149b5e25fSSatish Balay case 2: 152249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 152349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1524431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1525431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 152649b5e25fSSatish Balay break; 152749b5e25fSSatish Balay case 3: 152849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 152949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1530431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1531431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 153249b5e25fSSatish Balay break; 153349b5e25fSSatish Balay case 4: 153449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 153549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1536431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1537431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 153849b5e25fSSatish Balay break; 153949b5e25fSSatish Balay case 5: 154049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 154149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1542431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1543431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 154449b5e25fSSatish Balay break; 154549b5e25fSSatish Balay case 6: 154649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 154749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1548431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1549431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 155049b5e25fSSatish Balay break; 155149b5e25fSSatish Balay case 7: 1552de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 155349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1554431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1555431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 155649b5e25fSSatish Balay break; 155749b5e25fSSatish Balay } 155849b5e25fSSatish Balay } 155949b5e25fSSatish Balay 156049b5e25fSSatish Balay b->mbs = mbs; 15614afc71dfSHong Zhang b->nbs = mbs; 1562ab93d7beSBarry Smith if (!skipallocation) { 15632ee49352SLisandro Dalcin if (!b->imax) { 1564ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1565c760cd28SBarry Smith b->free_imax_ilen = PETSC_TRUE; 15662ee49352SLisandro Dalcin ierr = PetscLogObjectMemory(B,2*mbs*sizeof(PetscInt));CHKERRQ(ierr); 15672ee49352SLisandro Dalcin } 156849b5e25fSSatish Balay if (!nnz) { 1569435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 157049b5e25fSSatish Balay else if (nz <= 0) nz = 1; 157149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15728cef66ccSHong Zhang b->imax[i] = nz; 157349b5e25fSSatish Balay } 1574153ea458SHong Zhang nz = nz*mbs; /* total nz */ 157549b5e25fSSatish Balay } else { 157649b5e25fSSatish Balay nz = 0; 15778cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 157849b5e25fSSatish Balay } 15792ee49352SLisandro Dalcin /* b->ilen will count nonzeros in each block row so far. */ 15802ee49352SLisandro Dalcin for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 15816c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 158249b5e25fSSatish Balay 158349b5e25fSSatish Balay /* allocate the matrix space */ 15842ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);CHKERRQ(ierr); 1585d0f46423SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap->N+1,PetscInt,&b->i);CHKERRQ(ierr); 1586d0f46423SBarry Smith ierr = PetscLogObjectMemory(B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 15876c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 158813f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 158949b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 159049b5e25fSSatish Balay 159149b5e25fSSatish Balay /* pointer to beginning of each row */ 1592*e60cf9a0SBarry Smith b->i[0] = 0; 159349b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 159449b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 159549b5e25fSSatish Balay } 1596e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1597e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1598e811da20SHong Zhang } else { 1599e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1600e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1601ab93d7beSBarry Smith } 160249b5e25fSSatish Balay 1603d0f46423SBarry Smith B->rmap->bs = bs; 160449b5e25fSSatish Balay b->bs2 = bs2; 16056c6c5352SBarry Smith b->nz = 0; 16066c6c5352SBarry Smith b->maxnz = nz*bs2; 1607153ea458SHong Zhang 160816cdd363SHong Zhang b->inew = 0; 160916cdd363SHong Zhang b->jnew = 0; 161016cdd363SHong Zhang b->anew = 0; 161116cdd363SHong Zhang b->a2anew = 0; 16121a3463dfSHong Zhang b->permute = PETSC_FALSE; 1613c464158bSHong Zhang PetscFunctionReturn(0); 1614c464158bSHong Zhang } 1615a23d5eceSKris Buschelman EXTERN_C_END 1616153ea458SHong Zhang 1617db4efbfdSBarry Smith #undef __FUNCT__ 1618db4efbfdSBarry Smith #define __FUNCT__ "MatSeqBAIJSetNumericFactorization" 1619db4efbfdSBarry Smith /* 1620db4efbfdSBarry Smith This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization 1621db4efbfdSBarry Smith */ 1622db4efbfdSBarry Smith PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat B,PetscTruth natural) 1623db4efbfdSBarry Smith { 1624db4efbfdSBarry Smith PetscErrorCode ierr; 162590d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE; 1626db4efbfdSBarry Smith PetscInt bs = B->rmap->bs; 1627db4efbfdSBarry Smith 1628db4efbfdSBarry Smith PetscFunctionBegin; 162990d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,PETSC_NULL);CHKERRQ(ierr); 1630db4efbfdSBarry Smith if (flg) bs = 8; 1631db4efbfdSBarry Smith 1632db4efbfdSBarry Smith if (!natural) { 1633db4efbfdSBarry Smith switch (bs) { 1634db4efbfdSBarry Smith case 1: 1635db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 1636db4efbfdSBarry Smith break; 1637db4efbfdSBarry Smith case 2: 1638db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 1639db4efbfdSBarry Smith break; 1640db4efbfdSBarry Smith case 3: 1641db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 1642db4efbfdSBarry Smith break; 1643db4efbfdSBarry Smith case 4: 1644db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 1645db4efbfdSBarry Smith break; 1646db4efbfdSBarry Smith case 5: 1647db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 1648db4efbfdSBarry Smith break; 1649db4efbfdSBarry Smith case 6: 1650db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 1651db4efbfdSBarry Smith break; 1652db4efbfdSBarry Smith case 7: 1653db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 1654db4efbfdSBarry Smith break; 1655db4efbfdSBarry Smith default: 1656db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N; 1657db4efbfdSBarry Smith break; 1658db4efbfdSBarry Smith } 1659db4efbfdSBarry Smith } else { 1660db4efbfdSBarry Smith switch (bs) { 1661db4efbfdSBarry Smith case 1: 1662db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 1663db4efbfdSBarry Smith break; 1664db4efbfdSBarry Smith case 2: 1665db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 1666db4efbfdSBarry Smith break; 1667db4efbfdSBarry Smith case 3: 1668db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 1669db4efbfdSBarry Smith break; 1670db4efbfdSBarry Smith case 4: 1671db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 1672db4efbfdSBarry Smith break; 1673db4efbfdSBarry Smith case 5: 1674db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 1675db4efbfdSBarry Smith break; 1676db4efbfdSBarry Smith case 6: 1677db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 1678db4efbfdSBarry Smith break; 1679db4efbfdSBarry Smith case 7: 1680db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 1681db4efbfdSBarry Smith break; 1682db4efbfdSBarry Smith default: 1683db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 1684db4efbfdSBarry Smith break; 1685db4efbfdSBarry Smith } 1686db4efbfdSBarry Smith } 1687db4efbfdSBarry Smith PetscFunctionReturn(0); 1688db4efbfdSBarry Smith } 1689db4efbfdSBarry Smith 1690d769727bSBarry Smith EXTERN_C_BEGIN 1691f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1692f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1693d769727bSBarry Smith EXTERN_C_END 1694d769727bSBarry Smith 1695e631078cSBarry Smith 1696e631078cSBarry Smith EXTERN_C_BEGIN 16975c9eb25fSBarry Smith #undef __FUNCT__ 16985c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_seqsbaij_petsc" 16995c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 17005c9eb25fSBarry Smith { 1701d0f46423SBarry Smith PetscInt n = A->rmap->n; 17025c9eb25fSBarry Smith PetscErrorCode ierr; 17035c9eb25fSBarry Smith 17045c9eb25fSBarry Smith PetscFunctionBegin; 17055c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 17065c9eb25fSBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 17075c9eb25fSBarry Smith if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 17085c9eb25fSBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 17095c9eb25fSBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1710db4efbfdSBarry Smith (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; 1711db4efbfdSBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; 17125c9eb25fSBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 1713719d5645SBarry Smith (*B)->factor = ftype; 17145c9eb25fSBarry Smith PetscFunctionReturn(0); 17155c9eb25fSBarry Smith } 1716e631078cSBarry Smith EXTERN_C_END 17175c9eb25fSBarry Smith 17185c9eb25fSBarry Smith EXTERN_C_BEGIN 1719db4efbfdSBarry Smith #undef __FUNCT__ 1720db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqsbaij_petsc" 1721db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 1722db4efbfdSBarry Smith { 1723db4efbfdSBarry Smith PetscFunctionBegin; 1724db4efbfdSBarry Smith if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 1725db4efbfdSBarry Smith *flg = PETSC_TRUE; 1726db4efbfdSBarry Smith } else { 1727db4efbfdSBarry Smith *flg = PETSC_FALSE; 1728db4efbfdSBarry Smith } 1729db4efbfdSBarry Smith PetscFunctionReturn(0); 1730db4efbfdSBarry Smith } 1731db4efbfdSBarry Smith EXTERN_C_END 1732db4efbfdSBarry Smith 1733db4efbfdSBarry Smith EXTERN_C_BEGIN 1734611f576cSBarry Smith #if defined(PETSC_HAVE_MUMPS) 17355c9eb25fSBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat,MatFactorType,Mat*); 1736611f576cSBarry Smith #endif 1737611f576cSBarry Smith #if defined(PETSC_HAVE_SPOOLES) 17385c9eb25fSBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_spooles(Mat,MatFactorType,Mat*); 1739611f576cSBarry Smith #endif 1740b5e56a35SBarry Smith #if defined(PETSC_HAVE_PASTIX) 1741b5e56a35SBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*); 1742b5e56a35SBarry Smith #endif 17435c9eb25fSBarry Smith EXTERN_C_END 17445c9eb25fSBarry Smith 17450bad9183SKris Buschelman /*MC 1746fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 17470bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 17480bad9183SKris Buschelman 17490bad9183SKris Buschelman Options Database Keys: 17500bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 17510bad9183SKris Buschelman 17520bad9183SKris Buschelman Level: beginner 17530bad9183SKris Buschelman 17540bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 17550bad9183SKris Buschelman M*/ 17560bad9183SKris Buschelman 1757a23d5eceSKris Buschelman EXTERN_C_BEGIN 1758a23d5eceSKris Buschelman #undef __FUNCT__ 1759a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1760be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1761a23d5eceSKris Buschelman { 1762a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1763dfbe8321SBarry Smith PetscErrorCode ierr; 176413f74950SBarry Smith PetscMPIInt size; 17650def2e27SBarry Smith PetscTruth no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE; 1766a23d5eceSKris Buschelman 1767a23d5eceSKris Buschelman PetscFunctionBegin; 17687adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 1769a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1770a23d5eceSKris Buschelman 177138f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1772a23d5eceSKris Buschelman B->data = (void*)b; 1773a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1774a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1775a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1776a23d5eceSKris Buschelman B->mapping = 0; 1777a23d5eceSKris Buschelman b->row = 0; 1778a23d5eceSKris Buschelman b->icol = 0; 1779a23d5eceSKris Buschelman b->reallocs = 0; 1780a23d5eceSKris Buschelman b->saved_values = 0; 17810def2e27SBarry Smith b->inode.limit = 5; 17820def2e27SBarry Smith b->inode.max_limit = 5; 1783a23d5eceSKris Buschelman 1784a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1785a23d5eceSKris Buschelman b->nonew = 0; 1786a23d5eceSKris Buschelman b->diag = 0; 1787a23d5eceSKris Buschelman b->solve_work = 0; 1788a23d5eceSKris Buschelman b->mult_work = 0; 1789a23d5eceSKris Buschelman B->spptr = 0; 1790a9817697SBarry Smith b->keepnonzeropattern = PETSC_FALSE; 1791a23d5eceSKris Buschelman b->xtoy = 0; 1792a23d5eceSKris Buschelman b->XtoY = 0; 1793a23d5eceSKris Buschelman 1794a23d5eceSKris Buschelman b->inew = 0; 1795a23d5eceSKris Buschelman b->jnew = 0; 1796a23d5eceSKris Buschelman b->anew = 0; 1797a23d5eceSKris Buschelman b->a2anew = 0; 1798a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1799a23d5eceSKris Buschelman 1800941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 180190d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,PETSC_NULL);CHKERRQ(ierr); 1802941593c8SHong Zhang 1803f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 180490d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,PETSC_NULL);CHKERRQ(ierr); 1805f5edf698SHong Zhang 1806b5e56a35SBarry Smith #if defined(PETSC_HAVE_PASTIX) 1807ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_pastix_C", 1808b5e56a35SBarry Smith "MatGetFactor_seqsbaij_pastix", 1809b5e56a35SBarry Smith MatGetFactor_seqsbaij_pastix);CHKERRQ(ierr); 1810b5e56a35SBarry Smith #endif 1811611f576cSBarry Smith #if defined(PETSC_HAVE_SPOOLES) 1812ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_spooles_C", 18135c9eb25fSBarry Smith "MatGetFactor_seqsbaij_spooles", 18145c9eb25fSBarry Smith MatGetFactor_seqsbaij_spooles);CHKERRQ(ierr); 1815611f576cSBarry Smith #endif 1816611f576cSBarry Smith #if defined(PETSC_HAVE_MUMPS) 1817ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_mumps_C", 18185c9eb25fSBarry Smith "MatGetFactor_seqsbaij_mumps", 18195c9eb25fSBarry Smith MatGetFactor_seqsbaij_mumps);CHKERRQ(ierr); 1820611f576cSBarry Smith #endif 1821ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactorAvailable_petsc_C", 1822db4efbfdSBarry Smith "MatGetFactorAvailable_seqsbaij_petsc", 1823db4efbfdSBarry Smith MatGetFactorAvailable_seqsbaij_petsc);CHKERRQ(ierr); 1824ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 18255c9eb25fSBarry Smith "MatGetFactor_seqsbaij_petsc", 18265c9eb25fSBarry Smith MatGetFactor_seqsbaij_petsc);CHKERRQ(ierr); 1827a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1828a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1829a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1830a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1831a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1832a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1833a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1834a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1835a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 18364e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 18374e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 18384e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1839a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1840a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1841a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1842a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1843a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1844a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 184523ce1328SBarry Smith 184623ce1328SBarry Smith B->symmetric = PETSC_TRUE; 184723ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 184823ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 184923ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 185017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 18510def2e27SBarry Smith 18520def2e27SBarry Smith ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");CHKERRQ(ierr); 18530def2e27SBarry Smith ierr = PetscOptionsTruth("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);CHKERRQ(ierr); 18540def2e27SBarry Smith if (no_unroll) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);} 18550def2e27SBarry Smith ierr = PetscOptionsTruth("-mat_no_inode","Do not optimize for inodes (slower)",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);CHKERRQ(ierr); 18560def2e27SBarry Smith if (no_inode) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);} 18570def2e27SBarry Smith ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",PETSC_NULL,b->inode.limit,&b->inode.limit,PETSC_NULL);CHKERRQ(ierr); 18580def2e27SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 18590def2e27SBarry Smith b->inode.use = (PetscTruth)(!(no_unroll || no_inode)); 18600def2e27SBarry Smith if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit; 18610def2e27SBarry Smith 1862a23d5eceSKris Buschelman PetscFunctionReturn(0); 1863a23d5eceSKris Buschelman } 1864a23d5eceSKris Buschelman EXTERN_C_END 1865a23d5eceSKris Buschelman 1866a23d5eceSKris Buschelman #undef __FUNCT__ 1867a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1868a23d5eceSKris Buschelman /*@C 1869a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1870a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1871a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1872a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1873a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1874a23d5eceSKris Buschelman 1875a23d5eceSKris Buschelman Collective on Mat 1876a23d5eceSKris Buschelman 1877a23d5eceSKris Buschelman Input Parameters: 1878a23d5eceSKris Buschelman + A - the symmetric matrix 1879a23d5eceSKris Buschelman . bs - size of block 1880a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1881a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1882a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1883a23d5eceSKris Buschelman 1884a23d5eceSKris Buschelman Options Database Keys: 1885a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1886a23d5eceSKris Buschelman block calculations (much slower) 1887db4efbfdSBarry Smith . -mat_block_size - size of the blocks to use (only works if a negative bs is passed in 1888a23d5eceSKris Buschelman 1889a23d5eceSKris Buschelman Level: intermediate 1890a23d5eceSKris Buschelman 1891a23d5eceSKris Buschelman Notes: 1892a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1893a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1894a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1895a23d5eceSKris Buschelman matrices. 1896a23d5eceSKris Buschelman 1897aa95bbe8SBarry Smith You can call MatGetInfo() to get information on how effective the preallocation was; 1898aa95bbe8SBarry Smith for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 1899aa95bbe8SBarry Smith You can also run with the option -info and look for messages with the string 1900aa95bbe8SBarry Smith malloc in them to see if additional memory allocation was needed. 1901aa95bbe8SBarry Smith 190249a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 190349a6f317SBarry Smith 190449a6f317SBarry Smith 1905a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1906a23d5eceSKris Buschelman @*/ 1907be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 190813f74950SBarry Smith { 190913f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1910a23d5eceSKris Buschelman 1911a23d5eceSKris Buschelman PetscFunctionBegin; 1912a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1913a23d5eceSKris Buschelman if (f) { 1914a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1915a23d5eceSKris Buschelman } 1916a23d5eceSKris Buschelman PetscFunctionReturn(0); 1917a23d5eceSKris Buschelman } 191849b5e25fSSatish Balay 19194a2ae208SSatish Balay #undef __FUNCT__ 19204a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1921c464158bSHong Zhang /*@C 1922c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1923c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1924c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1925c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1926c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 192749b5e25fSSatish Balay 1928c464158bSHong Zhang Collective on MPI_Comm 1929c464158bSHong Zhang 1930c464158bSHong Zhang Input Parameters: 1931c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1932c464158bSHong Zhang . bs - size of block 1933c464158bSHong Zhang . m - number of rows, or number of columns 1934c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1935744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1936744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1937c464158bSHong Zhang 1938c464158bSHong Zhang Output Parameter: 1939c464158bSHong Zhang . A - the symmetric matrix 1940c464158bSHong Zhang 1941c464158bSHong Zhang Options Database Keys: 1942c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1943c464158bSHong Zhang block calculations (much slower) 1944c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1945c464158bSHong Zhang 1946c464158bSHong Zhang Level: intermediate 1947c464158bSHong Zhang 1948175b88e8SBarry Smith It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 1949ae1d86c5SBarry Smith MatXXXXSetPreallocation() paradgm instead of this routine directly. 1950175b88e8SBarry Smith [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 1951175b88e8SBarry Smith 1952c464158bSHong Zhang Notes: 19536d6d819aSHong Zhang The number of rows and columns must be divisible by blocksize. 19546d6d819aSHong Zhang This matrix type does not support complex Hermitian operation. 1955c464158bSHong Zhang 1956c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1957c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1958c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1959c464158bSHong Zhang matrices. 1960c464158bSHong Zhang 196149a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 196249a6f317SBarry Smith 1963c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1964c464158bSHong Zhang @*/ 1965be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1966c464158bSHong Zhang { 1967dfbe8321SBarry Smith PetscErrorCode ierr; 1968c464158bSHong Zhang 1969c464158bSHong Zhang PetscFunctionBegin; 1970f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1971f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1972c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1973ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 197449b5e25fSSatish Balay PetscFunctionReturn(0); 197549b5e25fSSatish Balay } 197649b5e25fSSatish Balay 19774a2ae208SSatish Balay #undef __FUNCT__ 19784a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1979dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 198049b5e25fSSatish Balay { 198149b5e25fSSatish Balay Mat C; 198249b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 19836849ba73SBarry Smith PetscErrorCode ierr; 1984b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 198549b5e25fSSatish Balay 198649b5e25fSSatish Balay PetscFunctionBegin; 198729bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 198849b5e25fSSatish Balay 198949b5e25fSSatish Balay *B = 0; 19907adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 1991d0f46423SBarry Smith ierr = MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);CHKERRQ(ierr); 19927adad957SLisandro Dalcin ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 19931d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1994692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1995692f9cbeSHong Zhang 1996273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 199749b5e25fSSatish Balay C->factor = A->factor; 199849b5e25fSSatish Balay c->row = 0; 199949b5e25fSSatish Balay c->icol = 0; 200049b5e25fSSatish Balay c->saved_values = 0; 2001a9817697SBarry Smith c->keepnonzeropattern = a->keepnonzeropattern; 200249b5e25fSSatish Balay C->assembled = PETSC_TRUE; 200349b5e25fSSatish Balay 200426283091SBarry Smith ierr = PetscLayoutCopy(A->rmap,&C->rmap);CHKERRQ(ierr); 200526283091SBarry Smith ierr = PetscLayoutCopy(A->cmap,&C->cmap);CHKERRQ(ierr); 200649b5e25fSSatish Balay c->bs2 = a->bs2; 200749b5e25fSSatish Balay c->mbs = a->mbs; 200849b5e25fSSatish Balay c->nbs = a->nbs; 200949b5e25fSSatish Balay 2010c760cd28SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 2011c760cd28SBarry Smith c->imax = a->imax; 2012c760cd28SBarry Smith c->ilen = a->ilen; 2013c760cd28SBarry Smith c->free_imax_ilen = PETSC_FALSE; 2014c760cd28SBarry Smith } else { 20158777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 2016c760cd28SBarry Smith ierr = PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 201749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 201849b5e25fSSatish Balay c->imax[i] = a->imax[i]; 201949b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 202049b5e25fSSatish Balay } 2021c760cd28SBarry Smith c->free_imax_ilen = PETSC_TRUE; 2022c760cd28SBarry Smith } 202349b5e25fSSatish Balay 202449b5e25fSSatish Balay /* allocate the matrix space */ 20254da8f245SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 20264da8f245SBarry Smith ierr = PetscMalloc(bs2*nz*sizeof(MatScalar),&c->a);CHKERRQ(ierr); 20274da8f245SBarry Smith ierr = PetscLogObjectMemory(C,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 20284da8f245SBarry Smith c->singlemalloc = PETSC_FALSE; 20294da8f245SBarry Smith c->free_ij = PETSC_FALSE; 20304da8f245SBarry Smith c->parent = A; 20314da8f245SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 20324da8f245SBarry Smith ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 20334da8f245SBarry Smith ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 20344da8f245SBarry Smith } else { 2035b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 203613f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 2037b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 20384da8f245SBarry Smith c->singlemalloc = PETSC_TRUE; 20394da8f245SBarry Smith c->free_ij = PETSC_TRUE; 20404da8f245SBarry Smith } 204149b5e25fSSatish Balay if (mbs > 0) { 20424da8f245SBarry Smith if (cpvalues != MAT_SHARE_NONZERO_PATTERN) { 204313f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 20444da8f245SBarry Smith } 204549b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 204649b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 204749b5e25fSSatish Balay } else { 204849b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 204949b5e25fSSatish Balay } 2050a1c3900fSBarry Smith if (a->jshort) { 20514da8f245SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 20524da8f245SBarry Smith c->jshort = a->jshort; 20534da8f245SBarry Smith c->free_jshort = PETSC_FALSE; 20544da8f245SBarry Smith } else { 2055a1c3900fSBarry Smith ierr = PetscMalloc(nz*sizeof(unsigned short),&c->jshort);CHKERRQ(ierr); 2056c760cd28SBarry Smith ierr = PetscLogObjectMemory(C,nz*sizeof(unsigned short));CHKERRQ(ierr); 2057a1c3900fSBarry Smith ierr = PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));CHKERRQ(ierr); 20584da8f245SBarry Smith c->free_jshort = PETSC_TRUE; 20594da8f245SBarry Smith } 2060a1c3900fSBarry Smith } 206149b5e25fSSatish Balay } 206249b5e25fSSatish Balay 206349b5e25fSSatish Balay c->roworiented = a->roworiented; 206449b5e25fSSatish Balay c->nonew = a->nonew; 206549b5e25fSSatish Balay 206649b5e25fSSatish Balay if (a->diag) { 2067c760cd28SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 2068c760cd28SBarry Smith c->diag = a->diag; 2069c760cd28SBarry Smith c->free_diag = PETSC_FALSE; 2070c760cd28SBarry Smith } else { 20712ed38d0bSJed Brown ierr = PetscMalloc(mbs*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 20722ed38d0bSJed Brown ierr = PetscLogObjectMemory(C,mbs*sizeof(PetscInt));CHKERRQ(ierr); 207349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 207449b5e25fSSatish Balay c->diag[i] = a->diag[i]; 207549b5e25fSSatish Balay } 2076c760cd28SBarry Smith c->free_diag = PETSC_TRUE; 2077c760cd28SBarry Smith } 207849b5e25fSSatish Balay } else c->diag = 0; 20796c6c5352SBarry Smith c->nz = a->nz; 20806c6c5352SBarry Smith c->maxnz = a->maxnz; 208149b5e25fSSatish Balay c->solve_work = 0; 208249b5e25fSSatish Balay c->mult_work = 0; 2083e6b907acSBarry Smith c->free_a = PETSC_TRUE; 208449b5e25fSSatish Balay *B = C; 20857adad957SLisandro Dalcin ierr = PetscFListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);CHKERRQ(ierr); 208649b5e25fSSatish Balay PetscFunctionReturn(0); 208749b5e25fSSatish Balay } 208849b5e25fSSatish Balay 20894a2ae208SSatish Balay #undef __FUNCT__ 20904a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 2091a313700dSBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, const MatType type,Mat *A) 209249b5e25fSSatish Balay { 209349b5e25fSSatish Balay Mat_SeqSBAIJ *a; 209449b5e25fSSatish Balay Mat B; 20956849ba73SBarry Smith PetscErrorCode ierr; 209613f74950SBarry Smith int fd; 209713f74950SBarry Smith PetscMPIInt size; 209813f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 209913f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 210013f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 210113f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 210287828ca2SBarry Smith PetscScalar *aa; 210349b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 210449b5e25fSSatish Balay 210549b5e25fSSatish Balay PetscFunctionBegin; 2106b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 210749b5e25fSSatish Balay bs2 = bs*bs; 210849b5e25fSSatish Balay 210949b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 211029bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 2111b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 211249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 2113552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 211449b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 211549b5e25fSSatish Balay 211649b5e25fSSatish Balay if (header[3] < 0) { 211729bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 211849b5e25fSSatish Balay } 211949b5e25fSSatish Balay 212029bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 212149b5e25fSSatish Balay 212249b5e25fSSatish Balay /* 212349b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 212449b5e25fSSatish Balay divisible by the blocksize 212549b5e25fSSatish Balay */ 212649b5e25fSSatish Balay mbs = M/bs; 212749b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 212849b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 212949b5e25fSSatish Balay else mbs++; 213049b5e25fSSatish Balay if (extra_rows) { 21311e2582c4SBarry Smith ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 213249b5e25fSSatish Balay } 213349b5e25fSSatish Balay 213449b5e25fSSatish Balay /* read in row lengths */ 213513f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 213649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 213749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 213849b5e25fSSatish Balay 213949b5e25fSSatish Balay /* read in column indices */ 214013f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 214149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 214249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 214349b5e25fSSatish Balay 214449b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 214513f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 214613f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 214774ed9c26SBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&mask,mbs,PetscInt,&masked);CHKERRQ(ierr); 214813f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 214974ed9c26SBarry Smith rowcount = 0; 215074ed9c26SBarry Smith nzcount = 0; 215149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 215249b5e25fSSatish Balay nmask = 0; 215349b5e25fSSatish Balay for (j=0; j<bs; j++) { 215449b5e25fSSatish Balay kmax = rowlengths[rowcount]; 215549b5e25fSSatish Balay for (k=0; k<kmax; k++) { 21562d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 215703630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 215849b5e25fSSatish Balay } 215949b5e25fSSatish Balay rowcount++; 216049b5e25fSSatish Balay } 2161574b2666SHong Zhang s_browlengths[i] += nmask; 2162574b2666SHong Zhang 216349b5e25fSSatish Balay /* zero out the mask elements we set */ 216449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 216549b5e25fSSatish Balay } 216649b5e25fSSatish Balay 216749b5e25fSSatish Balay /* create our matrix */ 2168f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 2169f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 21709abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 2171ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 217249b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 217349b5e25fSSatish Balay 217449b5e25fSSatish Balay /* set matrix "i" values */ 2175*e60cf9a0SBarry Smith a->i[0] = 0; 217649b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 2177574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 2178574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 217949b5e25fSSatish Balay } 21806c6c5352SBarry Smith a->nz = a->i[mbs]; 218149b5e25fSSatish Balay 218249b5e25fSSatish Balay /* read in nonzero values */ 218387828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 218449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 218549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 218649b5e25fSSatish Balay 218749b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 218849b5e25fSSatish Balay nzcount = 0; jcount = 0; 218949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 219049b5e25fSSatish Balay nzcountb = nzcount; 219149b5e25fSSatish Balay nmask = 0; 219249b5e25fSSatish Balay for (j=0; j<bs; j++) { 219349b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 219449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 21952d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 219603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 21972d703238SHong Zhang } 21982d703238SHong Zhang } 21992d703238SHong Zhang /* sort the masked values */ 22002d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 22012d703238SHong Zhang 22022d703238SHong Zhang /* set "j" values into matrix */ 22032d703238SHong Zhang maskcount = 1; 22042d703238SHong Zhang for (j=0; j<nmask; j++) { 220549b5e25fSSatish Balay a->j[jcount++] = masked[j]; 220649b5e25fSSatish Balay mask[masked[j]] = maskcount++; 220749b5e25fSSatish Balay } 2208574b2666SHong Zhang 220949b5e25fSSatish Balay /* set "a" values into matrix */ 221049b5e25fSSatish Balay ishift = bs2*a->i[i]; 221149b5e25fSSatish Balay for (j=0; j<bs; j++) { 221249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 221349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 2214574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 2215574b2666SHong Zhang if (tmp >= i){ 221649b5e25fSSatish Balay block = mask[tmp] - 1; 221749b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 221849b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 2219574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 2220574b2666SHong Zhang } 2221574b2666SHong Zhang nzcountb++; 222249b5e25fSSatish Balay } 222349b5e25fSSatish Balay } 222449b5e25fSSatish Balay /* zero out the mask elements we set */ 222549b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 222649b5e25fSSatish Balay } 22276c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 222849b5e25fSSatish Balay 222949b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2230574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 223149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 223249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 223374ed9c26SBarry Smith ierr = PetscFree2(mask,masked);CHKERRQ(ierr); 223449b5e25fSSatish Balay 22359abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 22369abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223749b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 22389abb65ffSKris Buschelman *A = B; 223949b5e25fSSatish Balay PetscFunctionReturn(0); 224049b5e25fSSatish Balay } 2241574b2666SHong Zhang 2242d06b337dSHong Zhang #undef __FUNCT__ 2243c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2244c75a6043SHong Zhang /*@ 2245c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2246c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2247c75a6043SHong Zhang 2248c75a6043SHong Zhang Collective on MPI_Comm 2249c75a6043SHong Zhang 2250c75a6043SHong Zhang Input Parameters: 2251c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2252c75a6043SHong Zhang . bs - size of block 2253c75a6043SHong Zhang . m - number of rows 2254c75a6043SHong Zhang . n - number of columns 2255c75a6043SHong Zhang . i - row indices 2256c75a6043SHong Zhang . j - column indices 2257c75a6043SHong Zhang - a - matrix values 2258c75a6043SHong Zhang 2259c75a6043SHong Zhang Output Parameter: 2260c75a6043SHong Zhang . mat - the matrix 2261c75a6043SHong Zhang 2262c75a6043SHong Zhang Level: intermediate 2263c75a6043SHong Zhang 2264c75a6043SHong Zhang Notes: 2265c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2266c75a6043SHong Zhang once the matrix is destroyed 2267c75a6043SHong Zhang 2268c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2269c75a6043SHong Zhang 2270c75a6043SHong Zhang The i and j indices are 0 based 2271c75a6043SHong Zhang 2272c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2273c75a6043SHong Zhang 2274c75a6043SHong Zhang @*/ 2275c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2276c75a6043SHong Zhang { 2277c75a6043SHong Zhang PetscErrorCode ierr; 2278c75a6043SHong Zhang PetscInt ii; 2279c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2280c75a6043SHong Zhang 2281c75a6043SHong Zhang PetscFunctionBegin; 2282c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2283c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2284c75a6043SHong Zhang 2285c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2286c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2287c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2288c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2289c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2290c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2291c760cd28SBarry Smith ierr = PetscLogObjectMemory(*mat,2*m*sizeof(PetscInt));CHKERRQ(ierr); 2292c75a6043SHong Zhang 2293c75a6043SHong Zhang sbaij->i = i; 2294c75a6043SHong Zhang sbaij->j = j; 2295c75a6043SHong Zhang sbaij->a = a; 2296c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2297c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2298e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2299e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2300c75a6043SHong Zhang 2301c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2302c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2303c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2304c75a6043SHong Zhang if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]); 2305c75a6043SHong Zhang #endif 2306c75a6043SHong Zhang } 2307c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2308c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2309c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2310c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2311c75a6043SHong Zhang } 2312c75a6043SHong Zhang #endif 2313c75a6043SHong Zhang 2314c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2315c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2316c75a6043SHong Zhang PetscFunctionReturn(0); 2317c75a6043SHong Zhang } 2318d06b337dSHong Zhang 2319d06b337dSHong Zhang 2320d06b337dSHong Zhang 232149b5e25fSSatish Balay 232249b5e25fSSatish Balay 2323