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: 192*eeffb40dSHong Zhang if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first"); 193*eeffb40dSHong Zhang #if defined(USESHORT) 194*eeffb40dSHong Zhang A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort; 195*eeffb40dSHong Zhang #else 196*eeffb40dSHong Zhang A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian; 197*eeffb40dSHong Zhang #endif 198*eeffb40dSHong 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 6330880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6340880e062SHong Zhang if (roworiented) { 6350880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6360880e062SHong Zhang } else { 6370880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6380880e062SHong Zhang } 6397cd84e04SBarry Smith if (col <= lastcol) low = 0; else high = nrow; 640e2ee6c50SBarry Smith lastcol = col; 6410880e062SHong Zhang while (high-low > 7) { 6420880e062SHong Zhang t = (low+high)/2; 6430880e062SHong Zhang if (rp[t] > col) high = t; 6440880e062SHong Zhang else low = t; 6450880e062SHong Zhang } 6460880e062SHong Zhang for (i=low; i<high; i++) { 6470880e062SHong Zhang if (rp[i] > col) break; 6480880e062SHong Zhang if (rp[i] == col) { 6490880e062SHong Zhang bap = ap + bs2*i; 6500880e062SHong Zhang if (roworiented) { 6510880e062SHong Zhang if (is == ADD_VALUES) { 6520880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6530880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6540880e062SHong Zhang bap[jj] += *value++; 6550880e062SHong Zhang } 6560880e062SHong Zhang } 6570880e062SHong Zhang } else { 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 } 6640880e062SHong Zhang } else { 6650880e062SHong Zhang if (is == ADD_VALUES) { 6660880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6670880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6680880e062SHong Zhang *bap++ += *value++; 6690880e062SHong Zhang } 6700880e062SHong Zhang } 6710880e062SHong Zhang } else { 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 } 6780880e062SHong Zhang } 6790880e062SHong Zhang goto noinsert2; 6800880e062SHong Zhang } 6810880e062SHong Zhang } 6820880e062SHong Zhang if (nonew == 1) goto noinsert2; 683085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 684421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 685c03d1d03SSatish Balay N = nrow++ - 1; high++; 6860880e062SHong Zhang /* shift up all the later entries in this row */ 6870880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6880880e062SHong Zhang rp[ii+1] = rp[ii]; 6890880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6900880e062SHong Zhang } 6910880e062SHong Zhang if (N >= i) { 6920880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6930880e062SHong Zhang } 6940880e062SHong Zhang rp[i] = col; 6950880e062SHong Zhang bap = ap + bs2*i; 6960880e062SHong Zhang if (roworiented) { 6970880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6980880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6990880e062SHong Zhang bap[jj] = *value++; 7000880e062SHong Zhang } 7010880e062SHong Zhang } 7020880e062SHong Zhang } else { 7030880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7040880e062SHong Zhang for (jj=0; jj<bs; jj++) { 7050880e062SHong Zhang *bap++ = *value++; 7060880e062SHong Zhang } 7070880e062SHong Zhang } 7080880e062SHong Zhang } 7090880e062SHong Zhang noinsert2:; 7100880e062SHong Zhang low = i; 7110880e062SHong Zhang } 7120880e062SHong Zhang ailen[row] = nrow; 7130880e062SHong Zhang } 7140880e062SHong Zhang PetscFunctionReturn(0); 71549b5e25fSSatish Balay } 71649b5e25fSSatish Balay 71764831d72SBarry Smith /* 71864831d72SBarry Smith This is not yet used 71964831d72SBarry Smith */ 7204a2ae208SSatish Balay #undef __FUNCT__ 7214108e4d5SBarry Smith #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode" 7224108e4d5SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A) 7230def2e27SBarry Smith { 7240def2e27SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7250def2e27SBarry Smith PetscErrorCode ierr; 7260def2e27SBarry Smith const PetscInt *ai = a->i, *aj = a->j,*cols; 7270def2e27SBarry Smith PetscInt i = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts; 7280def2e27SBarry Smith PetscTruth flag; 7290def2e27SBarry Smith 7300def2e27SBarry Smith PetscFunctionBegin; 7310def2e27SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&ns);CHKERRQ(ierr); 7320def2e27SBarry Smith while (i < m){ 7330def2e27SBarry Smith nzx = ai[i+1] - ai[i]; /* Number of nonzeros */ 7340def2e27SBarry Smith /* Limits the number of elements in a node to 'a->inode.limit' */ 7350def2e27SBarry Smith for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) { 7360def2e27SBarry Smith nzy = ai[j+1] - ai[j]; 7370def2e27SBarry Smith if (nzy != (nzx - j + i)) break; 7380def2e27SBarry Smith ierr = PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);CHKERRQ(ierr); 7390def2e27SBarry Smith if (!flag) break; 7400def2e27SBarry Smith } 7410def2e27SBarry Smith ns[node_count++] = blk_size; 7420def2e27SBarry Smith i = j; 7430def2e27SBarry Smith } 7440def2e27SBarry Smith if (!a->inode.size && m && node_count > .9*m) { 7450def2e27SBarry Smith ierr = PetscFree(ns);CHKERRQ(ierr); 7460def2e27SBarry Smith ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr); 7470def2e27SBarry Smith } else { 7480def2e27SBarry Smith a->inode.node_count = node_count; 7490def2e27SBarry Smith ierr = PetscMalloc(node_count*sizeof(PetscInt),&a->inode.size);CHKERRQ(ierr); 750c760cd28SBarry Smith ierr = PetscLogObjectMemory(A,node_count*sizeof(PetscInt));CHKERRQ(ierr); 7510def2e27SBarry Smith ierr = PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt)); 7520def2e27SBarry Smith ierr = PetscFree(ns);CHKERRQ(ierr); 7530def2e27SBarry Smith ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr); 7540def2e27SBarry Smith 7550def2e27SBarry Smith /* count collections of adjacent columns in each inode */ 7560def2e27SBarry Smith row = 0; 7570def2e27SBarry Smith cnt = 0; 7580def2e27SBarry Smith for (i=0; i<node_count; i++) { 7590def2e27SBarry Smith cols = aj + ai[row] + a->inode.size[i]; 7600def2e27SBarry Smith nz = ai[row+1] - ai[row] - a->inode.size[i]; 7610def2e27SBarry Smith for (j=1; j<nz; j++) { 7620def2e27SBarry Smith if (cols[j] != cols[j-1]+1) { 7630def2e27SBarry Smith cnt++; 7640def2e27SBarry Smith } 7650def2e27SBarry Smith } 7660def2e27SBarry Smith cnt++; 7670def2e27SBarry Smith row += a->inode.size[i]; 7680def2e27SBarry Smith } 7690def2e27SBarry Smith ierr = PetscMalloc(2*cnt*sizeof(PetscInt),&counts);CHKERRQ(ierr); 7700def2e27SBarry Smith cnt = 0; 7710def2e27SBarry Smith row = 0; 7720def2e27SBarry Smith for (i=0; i<node_count; i++) { 7730def2e27SBarry Smith cols = aj + ai[row] + a->inode.size[i]; 7740def2e27SBarry Smith CHKMEMQ; 7750def2e27SBarry Smith counts[2*cnt] = cols[0]; 7760def2e27SBarry Smith CHKMEMQ; 7770def2e27SBarry Smith nz = ai[row+1] - ai[row] - a->inode.size[i]; 7780def2e27SBarry Smith cnt2 = 1; 7790def2e27SBarry Smith for (j=1; j<nz; j++) { 7800def2e27SBarry Smith if (cols[j] != cols[j-1]+1) { 7810def2e27SBarry Smith CHKMEMQ; 7820def2e27SBarry Smith counts[2*(cnt++)+1] = cnt2; 7830def2e27SBarry Smith counts[2*cnt] = cols[j]; 7840def2e27SBarry Smith CHKMEMQ; 7850def2e27SBarry Smith cnt2 = 1; 7860def2e27SBarry Smith } else cnt2++; 7870def2e27SBarry Smith } 7880def2e27SBarry Smith CHKMEMQ; 7890def2e27SBarry Smith counts[2*(cnt++)+1] = cnt2; 7900def2e27SBarry Smith CHKMEMQ; 7910def2e27SBarry Smith row += a->inode.size[i]; 7920def2e27SBarry Smith } 7930def2e27SBarry Smith ierr = PetscIntView(2*cnt,counts,0); 7940def2e27SBarry Smith } 79538702af4SBarry Smith PetscFunctionReturn(0); 79638702af4SBarry Smith } 79738702af4SBarry Smith 79838702af4SBarry Smith #undef __FUNCT__ 7994a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 800dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 80149b5e25fSSatish Balay { 80249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 8036849ba73SBarry Smith PetscErrorCode ierr; 80413f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 805d0f46423SBarry Smith PetscInt m = A->rmap->N,*ip,N,*ailen = a->ilen; 80613f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 80749b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 80849b5e25fSSatish Balay 80949b5e25fSSatish Balay PetscFunctionBegin; 81049b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 81149b5e25fSSatish Balay 81249b5e25fSSatish Balay if (m) rmax = ailen[0]; 81349b5e25fSSatish Balay for (i=1; i<mbs; i++) { 81449b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 81549b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 81649b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 81749b5e25fSSatish Balay if (fshift) { 81849b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 81949b5e25fSSatish Balay N = ailen[i]; 82049b5e25fSSatish Balay for (j=0; j<N; j++) { 82149b5e25fSSatish Balay ip[j-fshift] = ip[j]; 82249b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 82349b5e25fSSatish Balay } 82449b5e25fSSatish Balay } 82549b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 82649b5e25fSSatish Balay } 82749b5e25fSSatish Balay if (mbs) { 82849b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 82949b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 83049b5e25fSSatish Balay } 83149b5e25fSSatish Balay /* reset ilen and imax for each row */ 83249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 83349b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 83449b5e25fSSatish Balay } 8356c6c5352SBarry Smith a->nz = ai[mbs]; 83649b5e25fSSatish Balay 837b424e231SHong Zhang /* diagonals may have moved, reset it */ 838b424e231SHong Zhang if (a->diag) { 83913f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 84049b5e25fSSatish Balay } 84128b2fa4aSMatthew Knepley if (fshift && a->nounused == -1) { 84228b2fa4aSMatthew 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); 84328b2fa4aSMatthew Knepley } 844d0f46423SBarry 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); 845ae15b995SBarry Smith ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr); 846ae15b995SBarry Smith ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr); 84749b5e25fSSatish Balay a->reallocs = 0; 84849b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 849061b2667SBarry Smith a->idiagvalid = PETSC_FALSE; 85038702af4SBarry Smith 85138702af4SBarry Smith if (A->cmap->n < 65536 && A->cmap->bs == 1) { 85238702af4SBarry Smith if (!a->jshort) { 85338702af4SBarry Smith ierr = PetscMalloc(a->i[A->rmap->n]*sizeof(unsigned short),&a->jshort);CHKERRQ(ierr); 854c760cd28SBarry Smith ierr = PetscLogObjectMemory(A,a->i[A->rmap->n]*sizeof(unsigned short));CHKERRQ(ierr); 85538702af4SBarry Smith for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i]; 85638702af4SBarry Smith A->ops->mult = MatMult_SeqSBAIJ_1_ushort; 85741f059aeSBarry Smith A->ops->sor = MatSOR_SeqSBAIJ_ushort; 8584da8f245SBarry Smith a->free_jshort = PETSC_TRUE; 85938702af4SBarry Smith } 86038702af4SBarry Smith } 86149b5e25fSSatish Balay PetscFunctionReturn(0); 86249b5e25fSSatish Balay } 86349b5e25fSSatish Balay 86449b5e25fSSatish Balay /* 86549b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 86649b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 86749b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 86849b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 86949b5e25fSSatish Balay */ 8704a2ae208SSatish Balay #undef __FUNCT__ 8714a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 87213f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 87349b5e25fSSatish Balay { 87413f74950SBarry Smith PetscInt i,j,k,row; 87549b5e25fSSatish Balay PetscTruth flg; 87649b5e25fSSatish Balay 87749b5e25fSSatish Balay PetscFunctionBegin; 87849b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 87949b5e25fSSatish Balay row = idx[i]; 88049b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 88149b5e25fSSatish Balay sizes[j] = 1; 88249b5e25fSSatish Balay i++; 88349b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 88449b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 88549b5e25fSSatish Balay i++; 88649b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 88749b5e25fSSatish Balay flg = PETSC_TRUE; 88849b5e25fSSatish Balay for (k=1; k<bs; k++) { 88949b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 89049b5e25fSSatish Balay flg = PETSC_FALSE; 89149b5e25fSSatish Balay break; 89249b5e25fSSatish Balay } 89349b5e25fSSatish Balay } 894abc0a331SBarry Smith if (flg) { /* No break in the bs */ 89549b5e25fSSatish Balay sizes[j] = bs; 89649b5e25fSSatish Balay i+= bs; 89749b5e25fSSatish Balay } else { 89849b5e25fSSatish Balay sizes[j] = 1; 89949b5e25fSSatish Balay i++; 90049b5e25fSSatish Balay } 90149b5e25fSSatish Balay } 90249b5e25fSSatish Balay } 90349b5e25fSSatish Balay *bs_max = j; 90449b5e25fSSatish Balay PetscFunctionReturn(0); 90549b5e25fSSatish Balay } 90649b5e25fSSatish Balay 90749b5e25fSSatish Balay 90849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 90949b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 91049b5e25fSSatish Balay */ 91149b5e25fSSatish Balay 9124a2ae208SSatish Balay #undef __FUNCT__ 9134a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 91413f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 91549b5e25fSSatish Balay { 91649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 9176849ba73SBarry Smith PetscErrorCode ierr; 918e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 91913f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 920d0f46423SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol; 92113f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 92249b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 92349b5e25fSSatish Balay 92449b5e25fSSatish Balay PetscFunctionBegin; 92549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 92649b5e25fSSatish Balay row = im[k]; /* row number */ 92749b5e25fSSatish Balay brow = row/bs; /* block row number */ 92849b5e25fSSatish Balay if (row < 0) continue; 9292515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 930d0f46423SBarry Smith if (row >= A->rmap->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1); 93149b5e25fSSatish Balay #endif 93249b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 93349b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 93449b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 93549b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 93649b5e25fSSatish Balay low = 0; 93749b5e25fSSatish Balay 93849b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 93949b5e25fSSatish Balay if (in[l] < 0) continue; 9402515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 941d0f46423SBarry 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); 94249b5e25fSSatish Balay #endif 94349b5e25fSSatish Balay col = in[l]; 94449b5e25fSSatish Balay bcol = col/bs; /* block col number */ 94549b5e25fSSatish Balay 946941593c8SHong Zhang if (brow > bcol) { 947941593c8SHong Zhang if (a->ignore_ltriangular){ 948941593c8SHong Zhang continue; /* ignore lower triangular values */ 949941593c8SHong Zhang } else { 9504e0d8c25SBarry 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)"); 951941593c8SHong Zhang } 952941593c8SHong Zhang } 953f4989cb3SHong Zhang 95449b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 9558549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 95649b5e25fSSatish Balay /* element value a(k,l) */ 95749b5e25fSSatish Balay if (roworiented) { 95849b5e25fSSatish Balay value = v[l + k*n]; 95949b5e25fSSatish Balay } else { 96049b5e25fSSatish Balay value = v[k + l*m]; 96149b5e25fSSatish Balay } 96249b5e25fSSatish Balay 96349b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 9647cd84e04SBarry Smith if (col <= lastcol) low = 0; high = nrow; 965e2ee6c50SBarry Smith lastcol = col; 96649b5e25fSSatish Balay while (high-low > 7) { 96749b5e25fSSatish Balay t = (low+high)/2; 96849b5e25fSSatish Balay if (rp[t] > bcol) high = t; 96949b5e25fSSatish Balay else low = t; 97049b5e25fSSatish Balay } 97149b5e25fSSatish Balay for (i=low; i<high; i++) { 97249b5e25fSSatish Balay if (rp[i] > bcol) break; 97349b5e25fSSatish Balay if (rp[i] == bcol) { 97449b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 97549b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 97649b5e25fSSatish Balay else *bap = value; 9778549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 9788549e402SHong Zhang if (brow == bcol && ridx < cidx){ 9798549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 9808549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 9818549e402SHong Zhang else *bap = value; 9828549e402SHong Zhang } 98349b5e25fSSatish Balay goto noinsert1; 98449b5e25fSSatish Balay } 98549b5e25fSSatish Balay } 98649b5e25fSSatish Balay 98749b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 988085a36d4SBarry Smith if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 989421e10b8SBarry Smith MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); 99049b5e25fSSatish Balay 991c03d1d03SSatish Balay N = nrow++ - 1; high++; 99249b5e25fSSatish Balay /* shift up all the later entries in this row */ 99349b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 99449b5e25fSSatish Balay rp[ii+1] = rp[ii]; 99549b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 99649b5e25fSSatish Balay } 99749b5e25fSSatish Balay if (N>=i) { 99849b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 99949b5e25fSSatish Balay } 100049b5e25fSSatish Balay rp[i] = bcol; 100149b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 100249b5e25fSSatish Balay noinsert1:; 100349b5e25fSSatish Balay low = i; 10048549e402SHong Zhang } 100549b5e25fSSatish Balay } /* end of loop over added columns */ 100649b5e25fSSatish Balay ailen[brow] = nrow; 100749b5e25fSSatish Balay } /* end of loop over added rows */ 100849b5e25fSSatish Balay PetscFunctionReturn(0); 100949b5e25fSSatish Balay } 101049b5e25fSSatish Balay 10114a2ae208SSatish Balay #undef __FUNCT__ 10124d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 10130481f469SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info) 101449b5e25fSSatish Balay { 10154ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 101649b5e25fSSatish Balay Mat outA; 1017dfbe8321SBarry Smith PetscErrorCode ierr; 1018c84f5b01SHong Zhang PetscTruth row_identity; 101949b5e25fSSatish Balay 102049b5e25fSSatish Balay PetscFunctionBegin; 1021c84f5b01SHong Zhang if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc"); 1022c84f5b01SHong Zhang ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1023c84f5b01SHong Zhang if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 1024d0f46423SBarry 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()! */ 1025c84f5b01SHong Zhang 102649b5e25fSSatish Balay outA = inA; 1027c078aec8SLisandro Dalcin inA->factor = MAT_FACTOR_ICC; 102849b5e25fSSatish Balay 10291a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 1030db4efbfdSBarry Smith ierr = MatSeqSBAIJSetNumericFactorization(inA,row_identity);CHKERRQ(ierr); 103149b5e25fSSatish Balay 1032c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1033c3122656SLisandro Dalcin if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } 1034c84f5b01SHong Zhang a->row = row; 1035c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1036c3122656SLisandro Dalcin if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } 1037c84f5b01SHong Zhang a->col = row; 1038c84f5b01SHong Zhang 1039c84f5b01SHong Zhang /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */ 1040c84f5b01SHong Zhang if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);} 104152e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 104249b5e25fSSatish Balay 104349b5e25fSSatish Balay if (!a->solve_work) { 1044d0f46423SBarry Smith ierr = PetscMalloc((inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 1045d0f46423SBarry Smith ierr = PetscLogObjectMemory(inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));CHKERRQ(ierr); 104649b5e25fSSatish Balay } 104749b5e25fSSatish Balay 1048719d5645SBarry Smith ierr = MatCholeskyFactorNumeric(outA,inA,info);CHKERRQ(ierr); 104949b5e25fSSatish Balay PetscFunctionReturn(0); 105049b5e25fSSatish Balay } 1051950f1e5bSHong Zhang 105249b5e25fSSatish Balay EXTERN_C_BEGIN 10534a2ae208SSatish Balay #undef __FUNCT__ 10544a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 1055be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 105649b5e25fSSatish Balay { 1057045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 105813f74950SBarry Smith PetscInt i,nz,n; 105949b5e25fSSatish Balay 106049b5e25fSSatish Balay PetscFunctionBegin; 10616c6c5352SBarry Smith nz = baij->maxnz; 1062d0f46423SBarry Smith n = mat->cmap->n; 106349b5e25fSSatish Balay for (i=0; i<nz; i++) { 106449b5e25fSSatish Balay baij->j[i] = indices[i]; 106549b5e25fSSatish Balay } 10666c6c5352SBarry Smith baij->nz = nz; 106749b5e25fSSatish Balay for (i=0; i<n; i++) { 106849b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 106949b5e25fSSatish Balay } 107049b5e25fSSatish Balay PetscFunctionReturn(0); 107149b5e25fSSatish Balay } 107249b5e25fSSatish Balay EXTERN_C_END 107349b5e25fSSatish Balay 10744a2ae208SSatish Balay #undef __FUNCT__ 10754a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 107649b5e25fSSatish Balay /*@ 107719585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 107849b5e25fSSatish Balay in the matrix. 107949b5e25fSSatish Balay 108049b5e25fSSatish Balay Input Parameters: 108119585528SSatish Balay + mat - the SeqSBAIJ matrix 108249b5e25fSSatish Balay - indices - the column indices 108349b5e25fSSatish Balay 108449b5e25fSSatish Balay Level: advanced 108549b5e25fSSatish Balay 108649b5e25fSSatish Balay Notes: 108749b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 108849b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 108949b5e25fSSatish Balay of the MatSetValues() operation. 109049b5e25fSSatish Balay 109149b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1092d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 109349b5e25fSSatish Balay 1094ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 109549b5e25fSSatish Balay 1096ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 109749b5e25fSSatish Balay @*/ 1098be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 109949b5e25fSSatish Balay { 110013f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 110149b5e25fSSatish Balay 110249b5e25fSSatish Balay PetscFunctionBegin; 11034482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 11044482741eSBarry Smith PetscValidPointer(indices,2); 1105c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 110649b5e25fSSatish Balay if (f) { 110749b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 110849b5e25fSSatish Balay } else { 1109e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 111049b5e25fSSatish Balay } 111149b5e25fSSatish Balay PetscFunctionReturn(0); 111249b5e25fSSatish Balay } 111349b5e25fSSatish Balay 11144a2ae208SSatish Balay #undef __FUNCT__ 11153c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ" 11163c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str) 11173c896bc6SHong Zhang { 11183c896bc6SHong Zhang PetscErrorCode ierr; 11193c896bc6SHong Zhang 11203c896bc6SHong Zhang PetscFunctionBegin; 11213c896bc6SHong Zhang /* If the two matrices have the same copy implementation, use fast copy. */ 11223c896bc6SHong Zhang if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { 11233c896bc6SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 11243c896bc6SHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 11253c896bc6SHong Zhang 1126d0f46423SBarry Smith if (a->i[A->rmap->N] != b->i[B->rmap->N]) { 11273c896bc6SHong Zhang SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); 11283c896bc6SHong Zhang } 1129d0f46423SBarry Smith ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));CHKERRQ(ierr); 11303c896bc6SHong Zhang } else { 1131f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr); 11323c896bc6SHong Zhang ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1133f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr); 11343c896bc6SHong Zhang } 11353c896bc6SHong Zhang PetscFunctionReturn(0); 11363c896bc6SHong Zhang } 11373c896bc6SHong Zhang 11383c896bc6SHong Zhang #undef __FUNCT__ 11394a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1140dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1141273d9f13SBarry Smith { 1142dfbe8321SBarry Smith PetscErrorCode ierr; 1143273d9f13SBarry Smith 1144273d9f13SBarry Smith PetscFunctionBegin; 1145db4efbfdSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,-PetscMax(A->rmap->bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr); 1146273d9f13SBarry Smith PetscFunctionReturn(0); 1147273d9f13SBarry Smith } 1148273d9f13SBarry Smith 1149a6ece127SHong Zhang #undef __FUNCT__ 1150a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1151dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1152a6ece127SHong Zhang { 1153a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1154a6ece127SHong Zhang PetscFunctionBegin; 1155a6ece127SHong Zhang *array = a->a; 1156a6ece127SHong Zhang PetscFunctionReturn(0); 1157a6ece127SHong Zhang } 1158a6ece127SHong Zhang 1159a6ece127SHong Zhang #undef __FUNCT__ 1160a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1161dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1162a6ece127SHong Zhang { 1163a6ece127SHong Zhang PetscFunctionBegin; 1164a6ece127SHong Zhang PetscFunctionReturn(0); 1165a6ece127SHong Zhang } 1166a6ece127SHong Zhang 116742ee4b1aSHong Zhang #include "petscblaslapack.h" 116842ee4b1aSHong Zhang #undef __FUNCT__ 116942ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1170f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) 117142ee4b1aSHong Zhang { 117242ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1173dfbe8321SBarry Smith PetscErrorCode ierr; 1174d0f46423SBarry Smith PetscInt i,bs=Y->rmap->bs,bs2,j; 11750805154bSBarry Smith PetscBLASInt one = 1,bnz = PetscBLASIntCast(x->nz); 117642ee4b1aSHong Zhang 117742ee4b1aSHong Zhang PetscFunctionBegin; 117842ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1179f4df32b1SMatthew Knepley PetscScalar alpha = a; 1180f4df32b1SMatthew Knepley BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); 1181c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1182c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1183c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1184c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1185c537a176SHong Zhang } 1186c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1187c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1188c4319e64SHong Zhang y->XtoY = X; 1189c537a176SHong Zhang } 1190c4319e64SHong Zhang bs2 = bs*bs; 11916c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1192c4319e64SHong Zhang j = 0; 1193c4319e64SHong Zhang while (j < bs2){ 1194f4df32b1SMatthew Knepley y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]); 1195c4319e64SHong Zhang j++; 1196c537a176SHong Zhang } 1197c4319e64SHong Zhang } 11981e2582c4SBarry 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); 119942ee4b1aSHong Zhang } else { 1200f5edf698SHong Zhang ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr); 1201f4df32b1SMatthew Knepley ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1202f5edf698SHong Zhang ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr); 120342ee4b1aSHong Zhang } 120442ee4b1aSHong Zhang PetscFunctionReturn(0); 120542ee4b1aSHong Zhang } 120642ee4b1aSHong Zhang 1207efcf0fc3SBarry Smith #undef __FUNCT__ 1208efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1209dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1210efcf0fc3SBarry Smith { 1211efcf0fc3SBarry Smith PetscFunctionBegin; 1212efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1213efcf0fc3SBarry Smith PetscFunctionReturn(0); 1214efcf0fc3SBarry Smith } 1215efcf0fc3SBarry Smith 1216efcf0fc3SBarry Smith #undef __FUNCT__ 1217efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1218dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1219efcf0fc3SBarry Smith { 1220efcf0fc3SBarry Smith PetscFunctionBegin; 1221efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1222efcf0fc3SBarry Smith PetscFunctionReturn(0); 1223efcf0fc3SBarry Smith } 1224efcf0fc3SBarry Smith 1225efcf0fc3SBarry Smith #undef __FUNCT__ 1226efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1227ab5e4463SMatthew Knepley PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1228efcf0fc3SBarry Smith { 1229efcf0fc3SBarry Smith PetscFunctionBegin; 1230efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1231efcf0fc3SBarry Smith PetscFunctionReturn(0); 1232efcf0fc3SBarry Smith } 1233efcf0fc3SBarry Smith 123499cafbc1SBarry Smith #undef __FUNCT__ 123599cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ" 123699cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A) 123799cafbc1SBarry Smith { 123899cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 123999cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 1240dd6ea824SBarry Smith MatScalar *aa = a->a; 124199cafbc1SBarry Smith 124299cafbc1SBarry Smith PetscFunctionBegin; 124399cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 124499cafbc1SBarry Smith PetscFunctionReturn(0); 124599cafbc1SBarry Smith } 124699cafbc1SBarry Smith 124799cafbc1SBarry Smith #undef __FUNCT__ 124899cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ" 124999cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A) 125099cafbc1SBarry Smith { 125199cafbc1SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 125299cafbc1SBarry Smith PetscInt i,nz = a->bs2*a->i[a->mbs]; 1253dd6ea824SBarry Smith MatScalar *aa = a->a; 125499cafbc1SBarry Smith 125599cafbc1SBarry Smith PetscFunctionBegin; 125699cafbc1SBarry Smith for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 125799cafbc1SBarry Smith PetscFunctionReturn(0); 125899cafbc1SBarry Smith } 125999cafbc1SBarry Smith 126049b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 126149b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 126249b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 126349b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 126449b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 126597304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1266431c96f7SBarry Smith MatMult_SeqSBAIJ_N, /* transpose versions are same as non-transpose versions */ 1267e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 1268db4efbfdSBarry Smith 0, 126949b5e25fSSatish Balay 0, 127049b5e25fSSatish Balay 0, 127197304618SKris Buschelman /*10*/ 0, 127249b5e25fSSatish Balay 0, 1273c078aec8SLisandro Dalcin MatCholeskyFactor_SeqSBAIJ, 127441f059aeSBarry Smith MatSOR_SeqSBAIJ, 127549b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 127697304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 127749b5e25fSSatish Balay MatEqual_SeqSBAIJ, 127849b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 127949b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 128049b5e25fSSatish Balay MatNorm_SeqSBAIJ, 128197304618SKris Buschelman /*20*/ 0, 128249b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 128349b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 128449b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1285d519adbfSMatthew Knepley /*24*/ 0, 128649b5e25fSSatish Balay 0, 128749b5e25fSSatish Balay 0, 1288db4efbfdSBarry Smith 0, 1289db4efbfdSBarry Smith 0, 1290d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqSBAIJ, 1291c464158bSHong Zhang 0, 1292db4efbfdSBarry Smith 0, 1293a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1294a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 1295d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqSBAIJ, 1296719d5645SBarry Smith 0, 1297719d5645SBarry Smith 0, 129849b5e25fSSatish Balay 0, 1299c84f5b01SHong Zhang MatICCFactor_SeqSBAIJ, 1300d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqSBAIJ, 130149b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 130249b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 130349b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 13043c896bc6SHong Zhang MatCopy_SeqSBAIJ, 1305d519adbfSMatthew Knepley /*44*/ 0, 130649b5e25fSSatish Balay MatScale_SeqSBAIJ, 130749b5e25fSSatish Balay 0, 130849b5e25fSSatish Balay 0, 130949b5e25fSSatish Balay 0, 1310d519adbfSMatthew Knepley /*49*/ 0, 131149b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 131249b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 131349b5e25fSSatish Balay 0, 131449b5e25fSSatish Balay 0, 1315d519adbfSMatthew Knepley /*54*/ 0, 131649b5e25fSSatish Balay 0, 131749b5e25fSSatish Balay 0, 131849b5e25fSSatish Balay 0, 131949b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 1320d519adbfSMatthew Knepley /*59*/ MatGetSubMatrix_SeqSBAIJ, 132149b5e25fSSatish Balay 0, 132249b5e25fSSatish Balay 0, 1323357abbc8SBarry Smith 0, 1324d959ec07SHong Zhang 0, 1325d519adbfSMatthew Knepley /*64*/ 0, 1326d959ec07SHong Zhang 0, 1327d959ec07SHong Zhang 0, 1328d959ec07SHong Zhang 0, 1329d959ec07SHong Zhang 0, 1330d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqSBAIJ, 13313e0d88b5SBarry Smith 0, 13323e0d88b5SBarry Smith 0, 13333e0d88b5SBarry Smith 0, 13343e0d88b5SBarry Smith 0, 1335d519adbfSMatthew Knepley /*74*/ 0, 13363e0d88b5SBarry Smith 0, 13373e0d88b5SBarry Smith 0, 13383e0d88b5SBarry Smith 0, 13393e0d88b5SBarry Smith 0, 1340d519adbfSMatthew Knepley /*79*/ 0, 13413e0d88b5SBarry Smith 0, 13423e0d88b5SBarry Smith 0, 13433e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 134497304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13453e0d88b5SBarry Smith #else 134697304618SKris Buschelman 0, 13473e0d88b5SBarry Smith #endif 1348865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1349d519adbfSMatthew Knepley /*84*/ MatIsSymmetric_SeqSBAIJ, 1350865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1351efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1352865e5f61SKris Buschelman 0, 1353865e5f61SKris Buschelman 0, 1354d519adbfSMatthew Knepley /*89*/ 0, 1355865e5f61SKris Buschelman 0, 1356865e5f61SKris Buschelman 0, 1357865e5f61SKris Buschelman 0, 1358865e5f61SKris Buschelman 0, 1359d519adbfSMatthew Knepley /*94*/ 0, 1360865e5f61SKris Buschelman 0, 1361865e5f61SKris Buschelman 0, 136299cafbc1SBarry Smith 0, 136399cafbc1SBarry Smith 0, 1364d519adbfSMatthew Knepley /*99*/ 0, 136599cafbc1SBarry Smith 0, 136699cafbc1SBarry Smith 0, 136799cafbc1SBarry Smith 0, 136899cafbc1SBarry Smith 0, 1369d519adbfSMatthew Knepley /*104*/0, 137099cafbc1SBarry Smith MatRealPart_SeqSBAIJ, 1371f5edf698SHong Zhang MatImaginaryPart_SeqSBAIJ, 1372f5edf698SHong Zhang MatGetRowUpperTriangular_SeqSBAIJ, 13732af78befSBarry Smith MatRestoreRowUpperTriangular_SeqSBAIJ, 1374d519adbfSMatthew Knepley /*109*/0, 13752af78befSBarry Smith 0, 13762af78befSBarry Smith 0, 13772af78befSBarry Smith 0, 13782af78befSBarry Smith MatMissingDiagonal_SeqSBAIJ 137999cafbc1SBarry Smith }; 1380be1d678aSKris Buschelman 138149b5e25fSSatish Balay EXTERN_C_BEGIN 13824a2ae208SSatish Balay #undef __FUNCT__ 13834a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1384be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat) 138549b5e25fSSatish Balay { 13864afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1387d0f46423SBarry Smith PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2; 1388dfbe8321SBarry Smith PetscErrorCode ierr; 138949b5e25fSSatish Balay 139049b5e25fSSatish Balay PetscFunctionBegin; 139149b5e25fSSatish Balay if (aij->nonew != 1) { 1392512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 139349b5e25fSSatish Balay } 139449b5e25fSSatish Balay 139549b5e25fSSatish Balay /* allocate space for values if not already there */ 139649b5e25fSSatish Balay if (!aij->saved_values) { 139787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 139849b5e25fSSatish Balay } 139949b5e25fSSatish Balay 140049b5e25fSSatish Balay /* copy values over */ 140187828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 140249b5e25fSSatish Balay PetscFunctionReturn(0); 140349b5e25fSSatish Balay } 140449b5e25fSSatish Balay EXTERN_C_END 140549b5e25fSSatish Balay 140649b5e25fSSatish Balay EXTERN_C_BEGIN 14074a2ae208SSatish Balay #undef __FUNCT__ 14084a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1409be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat) 141049b5e25fSSatish Balay { 14114afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 14126849ba73SBarry Smith PetscErrorCode ierr; 1413d0f46423SBarry Smith PetscInt nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2; 141449b5e25fSSatish Balay 141549b5e25fSSatish Balay PetscFunctionBegin; 141649b5e25fSSatish Balay if (aij->nonew != 1) { 1417512a5fc5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first"); 141849b5e25fSSatish Balay } 141949b5e25fSSatish Balay if (!aij->saved_values) { 1420e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 142149b5e25fSSatish Balay } 142249b5e25fSSatish Balay 142349b5e25fSSatish Balay /* copy values over */ 142487828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 142549b5e25fSSatish Balay PetscFunctionReturn(0); 142649b5e25fSSatish Balay } 142749b5e25fSSatish Balay EXTERN_C_END 142849b5e25fSSatish Balay 14298549e402SHong Zhang EXTERN_C_BEGIN 14304a2ae208SSatish Balay #undef __FUNCT__ 1431a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1432be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 143349b5e25fSSatish Balay { 1434c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14356849ba73SBarry Smith PetscErrorCode ierr; 1436db4efbfdSBarry Smith PetscInt i,mbs,bs2, newbs = PetscAbs(bs); 143790d69ab7SBarry Smith PetscTruth skipallocation = PETSC_FALSE,flg = PETSC_FALSE; 143849b5e25fSSatish Balay 143949b5e25fSSatish Balay PetscFunctionBegin; 1440273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1441db4efbfdSBarry Smith if (bs < 0) { 1442db4efbfdSBarry Smith ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for MPISBAIJ matrix","Mat");CHKERRQ(ierr); 1443db4efbfdSBarry Smith ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSeqSBAIJSetPreallocation",newbs,&newbs,PETSC_NULL);CHKERRQ(ierr); 1444db4efbfdSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 1445db4efbfdSBarry Smith bs = PetscAbs(bs); 1446db4efbfdSBarry Smith } 1447db4efbfdSBarry Smith if (nnz && newbs != bs) { 1448db4efbfdSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot change blocksize from command line if setting nnz"); 1449db4efbfdSBarry Smith } 1450db4efbfdSBarry Smith bs = newbs; 1451db4efbfdSBarry Smith 145226283091SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr); 145326283091SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr); 145426283091SBarry Smith ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 145526283091SBarry Smith ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 1456899cda47SBarry Smith 1457d0f46423SBarry Smith mbs = B->rmap->N/bs; 145849b5e25fSSatish Balay bs2 = bs*bs; 145949b5e25fSSatish Balay 1460d0f46423SBarry Smith if (mbs*bs != B->rmap->N) { 146129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 146249b5e25fSSatish Balay } 146349b5e25fSSatish Balay 1464ab93d7beSBarry Smith if (nz == MAT_SKIP_ALLOCATION) { 1465ab93d7beSBarry Smith skipallocation = PETSC_TRUE; 1466ab93d7beSBarry Smith nz = 0; 1467ab93d7beSBarry Smith } 1468ab93d7beSBarry Smith 1469435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 147077431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 147149b5e25fSSatish Balay if (nnz) { 147249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 147377431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 147477431f27SBarry 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); 147549b5e25fSSatish Balay } 147649b5e25fSSatish Balay } 147749b5e25fSSatish Balay 1478db4efbfdSBarry Smith B->ops->mult = MatMult_SeqSBAIJ_N; 1479db4efbfdSBarry Smith B->ops->multadd = MatMultAdd_SeqSBAIJ_N; 1480db4efbfdSBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_N; 1481db4efbfdSBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N; 148290d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,PETSC_NULL);CHKERRQ(ierr); 148349b5e25fSSatish Balay if (!flg) { 148449b5e25fSSatish Balay switch (bs) { 148549b5e25fSSatish Balay case 1: 148649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 148749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 1488431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_1; 1489431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1; 149049b5e25fSSatish Balay break; 149149b5e25fSSatish Balay case 2: 149249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 149349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 1494431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_2; 1495431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2; 149649b5e25fSSatish Balay break; 149749b5e25fSSatish Balay case 3: 149849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 149949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 1500431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_3; 1501431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3; 150249b5e25fSSatish Balay break; 150349b5e25fSSatish Balay case 4: 150449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 150549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 1506431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_4; 1507431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4; 150849b5e25fSSatish Balay break; 150949b5e25fSSatish Balay case 5: 151049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 151149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 1512431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_5; 1513431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5; 151449b5e25fSSatish Balay break; 151549b5e25fSSatish Balay case 6: 151649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 151749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 1518431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_6; 1519431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6; 152049b5e25fSSatish Balay break; 152149b5e25fSSatish Balay case 7: 1522de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 152349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 1524431c96f7SBarry Smith B->ops->multtranspose = MatMult_SeqSBAIJ_7; 1525431c96f7SBarry Smith B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7; 152649b5e25fSSatish Balay break; 152749b5e25fSSatish Balay } 152849b5e25fSSatish Balay } 152949b5e25fSSatish Balay 153049b5e25fSSatish Balay b->mbs = mbs; 15314afc71dfSHong Zhang b->nbs = mbs; 1532ab93d7beSBarry Smith if (!skipallocation) { 15332ee49352SLisandro Dalcin if (!b->imax) { 1534ab93d7beSBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr); 1535c760cd28SBarry Smith b->free_imax_ilen = PETSC_TRUE; 15362ee49352SLisandro Dalcin ierr = PetscLogObjectMemory(B,2*mbs*sizeof(PetscInt));CHKERRQ(ierr); 15372ee49352SLisandro Dalcin } 153849b5e25fSSatish Balay if (!nnz) { 1539435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 154049b5e25fSSatish Balay else if (nz <= 0) nz = 1; 154149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15428cef66ccSHong Zhang b->imax[i] = nz; 154349b5e25fSSatish Balay } 1544153ea458SHong Zhang nz = nz*mbs; /* total nz */ 154549b5e25fSSatish Balay } else { 154649b5e25fSSatish Balay nz = 0; 15478cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 154849b5e25fSSatish Balay } 15492ee49352SLisandro Dalcin /* b->ilen will count nonzeros in each block row so far. */ 15502ee49352SLisandro Dalcin for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 15516c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 155249b5e25fSSatish Balay 155349b5e25fSSatish Balay /* allocate the matrix space */ 15542ee49352SLisandro Dalcin ierr = MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);CHKERRQ(ierr); 1555d0f46423SBarry Smith ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap->N+1,PetscInt,&b->i);CHKERRQ(ierr); 1556d0f46423SBarry Smith ierr = PetscLogObjectMemory(B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 15576c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 155813f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 155949b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 156049b5e25fSSatish Balay 156149b5e25fSSatish Balay /* pointer to beginning of each row */ 156249b5e25fSSatish Balay b->i[0] = 0; 156349b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 156449b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 156549b5e25fSSatish Balay } 1566e6b907acSBarry Smith b->free_a = PETSC_TRUE; 1567e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 1568e811da20SHong Zhang } else { 1569e6b907acSBarry Smith b->free_a = PETSC_FALSE; 1570e6b907acSBarry Smith b->free_ij = PETSC_FALSE; 1571ab93d7beSBarry Smith } 157249b5e25fSSatish Balay 1573d0f46423SBarry Smith B->rmap->bs = bs; 157449b5e25fSSatish Balay b->bs2 = bs2; 15756c6c5352SBarry Smith b->nz = 0; 15766c6c5352SBarry Smith b->maxnz = nz*bs2; 1577153ea458SHong Zhang 157816cdd363SHong Zhang b->inew = 0; 157916cdd363SHong Zhang b->jnew = 0; 158016cdd363SHong Zhang b->anew = 0; 158116cdd363SHong Zhang b->a2anew = 0; 15821a3463dfSHong Zhang b->permute = PETSC_FALSE; 1583c464158bSHong Zhang PetscFunctionReturn(0); 1584c464158bSHong Zhang } 1585a23d5eceSKris Buschelman EXTERN_C_END 1586153ea458SHong Zhang 1587db4efbfdSBarry Smith #undef __FUNCT__ 1588db4efbfdSBarry Smith #define __FUNCT__ "MatSeqBAIJSetNumericFactorization" 1589db4efbfdSBarry Smith /* 1590db4efbfdSBarry Smith This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization 1591db4efbfdSBarry Smith */ 1592db4efbfdSBarry Smith PetscErrorCode MatSeqSBAIJSetNumericFactorization(Mat B,PetscTruth natural) 1593db4efbfdSBarry Smith { 1594db4efbfdSBarry Smith PetscErrorCode ierr; 159590d69ab7SBarry Smith PetscTruth flg = PETSC_FALSE; 1596db4efbfdSBarry Smith PetscInt bs = B->rmap->bs; 1597db4efbfdSBarry Smith 1598db4efbfdSBarry Smith PetscFunctionBegin; 159990d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,PETSC_NULL);CHKERRQ(ierr); 1600db4efbfdSBarry Smith if (flg) bs = 8; 1601db4efbfdSBarry Smith 1602db4efbfdSBarry Smith if (!natural) { 1603db4efbfdSBarry Smith switch (bs) { 1604db4efbfdSBarry Smith case 1: 1605db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 1606db4efbfdSBarry Smith break; 1607db4efbfdSBarry Smith case 2: 1608db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 1609db4efbfdSBarry Smith break; 1610db4efbfdSBarry Smith case 3: 1611db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 1612db4efbfdSBarry Smith break; 1613db4efbfdSBarry Smith case 4: 1614db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 1615db4efbfdSBarry Smith break; 1616db4efbfdSBarry Smith case 5: 1617db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 1618db4efbfdSBarry Smith break; 1619db4efbfdSBarry Smith case 6: 1620db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 1621db4efbfdSBarry Smith break; 1622db4efbfdSBarry Smith case 7: 1623db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 1624db4efbfdSBarry Smith break; 1625db4efbfdSBarry Smith default: 1626db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N; 1627db4efbfdSBarry Smith break; 1628db4efbfdSBarry Smith } 1629db4efbfdSBarry Smith } else { 1630db4efbfdSBarry Smith switch (bs) { 1631db4efbfdSBarry Smith case 1: 1632db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 1633db4efbfdSBarry Smith break; 1634db4efbfdSBarry Smith case 2: 1635db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 1636db4efbfdSBarry Smith break; 1637db4efbfdSBarry Smith case 3: 1638db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 1639db4efbfdSBarry Smith break; 1640db4efbfdSBarry Smith case 4: 1641db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 1642db4efbfdSBarry Smith break; 1643db4efbfdSBarry Smith case 5: 1644db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 1645db4efbfdSBarry Smith break; 1646db4efbfdSBarry Smith case 6: 1647db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 1648db4efbfdSBarry Smith break; 1649db4efbfdSBarry Smith case 7: 1650db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 1651db4efbfdSBarry Smith break; 1652db4efbfdSBarry Smith default: 1653db4efbfdSBarry Smith B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 1654db4efbfdSBarry Smith break; 1655db4efbfdSBarry Smith } 1656db4efbfdSBarry Smith } 1657db4efbfdSBarry Smith PetscFunctionReturn(0); 1658db4efbfdSBarry Smith } 1659db4efbfdSBarry Smith 1660d769727bSBarry Smith EXTERN_C_BEGIN 1661f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*); 1662f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*); 1663d769727bSBarry Smith EXTERN_C_END 1664d769727bSBarry Smith 1665e631078cSBarry Smith 1666e631078cSBarry Smith EXTERN_C_BEGIN 16675c9eb25fSBarry Smith #undef __FUNCT__ 16685c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_seqsbaij_petsc" 16695c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 16705c9eb25fSBarry Smith { 1671d0f46423SBarry Smith PetscInt n = A->rmap->n; 16725c9eb25fSBarry Smith PetscErrorCode ierr; 16735c9eb25fSBarry Smith 16745c9eb25fSBarry Smith PetscFunctionBegin; 16755c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 16765c9eb25fSBarry Smith ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 16775c9eb25fSBarry Smith if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 16785c9eb25fSBarry Smith ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 16795c9eb25fSBarry Smith ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1680db4efbfdSBarry Smith (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; 1681db4efbfdSBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; 16825c9eb25fSBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 1683719d5645SBarry Smith (*B)->factor = ftype; 16845c9eb25fSBarry Smith PetscFunctionReturn(0); 16855c9eb25fSBarry Smith } 1686e631078cSBarry Smith EXTERN_C_END 16875c9eb25fSBarry Smith 16885c9eb25fSBarry Smith EXTERN_C_BEGIN 1689db4efbfdSBarry Smith #undef __FUNCT__ 1690db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactorAvailable_seqsbaij_petsc" 1691db4efbfdSBarry Smith PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 1692db4efbfdSBarry Smith { 1693db4efbfdSBarry Smith PetscFunctionBegin; 1694db4efbfdSBarry Smith if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 1695db4efbfdSBarry Smith *flg = PETSC_TRUE; 1696db4efbfdSBarry Smith } else { 1697db4efbfdSBarry Smith *flg = PETSC_FALSE; 1698db4efbfdSBarry Smith } 1699db4efbfdSBarry Smith PetscFunctionReturn(0); 1700db4efbfdSBarry Smith } 1701db4efbfdSBarry Smith EXTERN_C_END 1702db4efbfdSBarry Smith 1703db4efbfdSBarry Smith EXTERN_C_BEGIN 1704611f576cSBarry Smith #if defined(PETSC_HAVE_MUMPS) 17055c9eb25fSBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat,MatFactorType,Mat*); 1706611f576cSBarry Smith #endif 1707611f576cSBarry Smith #if defined(PETSC_HAVE_SPOOLES) 17085c9eb25fSBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_spooles(Mat,MatFactorType,Mat*); 1709611f576cSBarry Smith #endif 1710b5e56a35SBarry Smith #if defined(PETSC_HAVE_PASTIX) 1711b5e56a35SBarry Smith extern PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*); 1712b5e56a35SBarry Smith #endif 17135c9eb25fSBarry Smith EXTERN_C_END 17145c9eb25fSBarry Smith 17150bad9183SKris Buschelman /*MC 1716fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 17170bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 17180bad9183SKris Buschelman 17190bad9183SKris Buschelman Options Database Keys: 17200bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 17210bad9183SKris Buschelman 17220bad9183SKris Buschelman Level: beginner 17230bad9183SKris Buschelman 17240bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 17250bad9183SKris Buschelman M*/ 17260bad9183SKris Buschelman 1727a23d5eceSKris Buschelman EXTERN_C_BEGIN 1728a23d5eceSKris Buschelman #undef __FUNCT__ 1729a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1730be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B) 1731a23d5eceSKris Buschelman { 1732a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1733dfbe8321SBarry Smith PetscErrorCode ierr; 173413f74950SBarry Smith PetscMPIInt size; 17350def2e27SBarry Smith PetscTruth no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE; 1736a23d5eceSKris Buschelman 1737a23d5eceSKris Buschelman PetscFunctionBegin; 17387adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 1739a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1740a23d5eceSKris Buschelman 174138f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1742a23d5eceSKris Buschelman B->data = (void*)b; 1743a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1744a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1745a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1746a23d5eceSKris Buschelman B->mapping = 0; 1747a23d5eceSKris Buschelman b->row = 0; 1748a23d5eceSKris Buschelman b->icol = 0; 1749a23d5eceSKris Buschelman b->reallocs = 0; 1750a23d5eceSKris Buschelman b->saved_values = 0; 17510def2e27SBarry Smith b->inode.limit = 5; 17520def2e27SBarry Smith b->inode.max_limit = 5; 1753a23d5eceSKris Buschelman 1754a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1755a23d5eceSKris Buschelman b->nonew = 0; 1756a23d5eceSKris Buschelman b->diag = 0; 1757a23d5eceSKris Buschelman b->solve_work = 0; 1758a23d5eceSKris Buschelman b->mult_work = 0; 1759a23d5eceSKris Buschelman B->spptr = 0; 1760a9817697SBarry Smith b->keepnonzeropattern = PETSC_FALSE; 1761a23d5eceSKris Buschelman b->xtoy = 0; 1762a23d5eceSKris Buschelman b->XtoY = 0; 1763a23d5eceSKris Buschelman 1764a23d5eceSKris Buschelman b->inew = 0; 1765a23d5eceSKris Buschelman b->jnew = 0; 1766a23d5eceSKris Buschelman b->anew = 0; 1767a23d5eceSKris Buschelman b->a2anew = 0; 1768a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1769a23d5eceSKris Buschelman 1770941593c8SHong Zhang b->ignore_ltriangular = PETSC_FALSE; 177190d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,PETSC_NULL);CHKERRQ(ierr); 1772941593c8SHong Zhang 1773f5edf698SHong Zhang b->getrow_utriangular = PETSC_FALSE; 177490d69ab7SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,PETSC_NULL);CHKERRQ(ierr); 1775f5edf698SHong Zhang 1776b5e56a35SBarry Smith #if defined(PETSC_HAVE_PASTIX) 1777ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_pastix_C", 1778b5e56a35SBarry Smith "MatGetFactor_seqsbaij_pastix", 1779b5e56a35SBarry Smith MatGetFactor_seqsbaij_pastix);CHKERRQ(ierr); 1780b5e56a35SBarry Smith #endif 1781611f576cSBarry Smith #if defined(PETSC_HAVE_SPOOLES) 1782ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_spooles_C", 17835c9eb25fSBarry Smith "MatGetFactor_seqsbaij_spooles", 17845c9eb25fSBarry Smith MatGetFactor_seqsbaij_spooles);CHKERRQ(ierr); 1785611f576cSBarry Smith #endif 1786611f576cSBarry Smith #if defined(PETSC_HAVE_MUMPS) 1787ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_mumps_C", 17885c9eb25fSBarry Smith "MatGetFactor_seqsbaij_mumps", 17895c9eb25fSBarry Smith MatGetFactor_seqsbaij_mumps);CHKERRQ(ierr); 1790611f576cSBarry Smith #endif 1791ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactorAvailable_petsc_C", 1792db4efbfdSBarry Smith "MatGetFactorAvailable_seqsbaij_petsc", 1793db4efbfdSBarry Smith MatGetFactorAvailable_seqsbaij_petsc);CHKERRQ(ierr); 1794ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 17955c9eb25fSBarry Smith "MatGetFactor_seqsbaij_petsc", 17965c9eb25fSBarry Smith MatGetFactor_seqsbaij_petsc);CHKERRQ(ierr); 1797a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1798a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1799a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1800a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1801a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1802a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1803a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1804a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1805a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 18064e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 18074e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 18084e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1809a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1810a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1811a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1812a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1813a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1814a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 181523ce1328SBarry Smith 181623ce1328SBarry Smith B->symmetric = PETSC_TRUE; 181723ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 181823ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 181923ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 182017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr); 18210def2e27SBarry Smith 18220def2e27SBarry Smith ierr = PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");CHKERRQ(ierr); 18230def2e27SBarry Smith ierr = PetscOptionsTruth("-mat_no_unroll","Do not optimize for inodes (slower)",PETSC_NULL,no_unroll,&no_unroll,PETSC_NULL);CHKERRQ(ierr); 18240def2e27SBarry Smith if (no_unroll) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);} 18250def2e27SBarry Smith ierr = PetscOptionsTruth("-mat_no_inode","Do not optimize for inodes (slower)",PETSC_NULL,no_inode,&no_inode,PETSC_NULL);CHKERRQ(ierr); 18260def2e27SBarry Smith if (no_inode) {ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);} 18270def2e27SBarry 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); 18280def2e27SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 18290def2e27SBarry Smith b->inode.use = (PetscTruth)(!(no_unroll || no_inode)); 18300def2e27SBarry Smith if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit; 18310def2e27SBarry Smith 1832a23d5eceSKris Buschelman PetscFunctionReturn(0); 1833a23d5eceSKris Buschelman } 1834a23d5eceSKris Buschelman EXTERN_C_END 1835a23d5eceSKris Buschelman 1836a23d5eceSKris Buschelman #undef __FUNCT__ 1837a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1838a23d5eceSKris Buschelman /*@C 1839a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1840a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1841a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1842a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1843a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1844a23d5eceSKris Buschelman 1845a23d5eceSKris Buschelman Collective on Mat 1846a23d5eceSKris Buschelman 1847a23d5eceSKris Buschelman Input Parameters: 1848a23d5eceSKris Buschelman + A - the symmetric matrix 1849a23d5eceSKris Buschelman . bs - size of block 1850a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1851a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1852a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1853a23d5eceSKris Buschelman 1854a23d5eceSKris Buschelman Options Database Keys: 1855a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1856a23d5eceSKris Buschelman block calculations (much slower) 1857db4efbfdSBarry Smith . -mat_block_size - size of the blocks to use (only works if a negative bs is passed in 1858a23d5eceSKris Buschelman 1859a23d5eceSKris Buschelman Level: intermediate 1860a23d5eceSKris Buschelman 1861a23d5eceSKris Buschelman Notes: 1862a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1863a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1864a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1865a23d5eceSKris Buschelman matrices. 1866a23d5eceSKris Buschelman 1867aa95bbe8SBarry Smith You can call MatGetInfo() to get information on how effective the preallocation was; 1868aa95bbe8SBarry Smith for example the fields mallocs,nz_allocated,nz_used,nz_unneeded; 1869aa95bbe8SBarry Smith You can also run with the option -info and look for messages with the string 1870aa95bbe8SBarry Smith malloc in them to see if additional memory allocation was needed. 1871aa95bbe8SBarry Smith 187249a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 187349a6f317SBarry Smith 187449a6f317SBarry Smith 1875a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1876a23d5eceSKris Buschelman @*/ 1877be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 187813f74950SBarry Smith { 187913f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1880a23d5eceSKris Buschelman 1881a23d5eceSKris Buschelman PetscFunctionBegin; 1882a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1883a23d5eceSKris Buschelman if (f) { 1884a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1885a23d5eceSKris Buschelman } 1886a23d5eceSKris Buschelman PetscFunctionReturn(0); 1887a23d5eceSKris Buschelman } 188849b5e25fSSatish Balay 18894a2ae208SSatish Balay #undef __FUNCT__ 18904a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1891c464158bSHong Zhang /*@C 1892c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1893c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1894c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1895c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1896c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 189749b5e25fSSatish Balay 1898c464158bSHong Zhang Collective on MPI_Comm 1899c464158bSHong Zhang 1900c464158bSHong Zhang Input Parameters: 1901c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1902c464158bSHong Zhang . bs - size of block 1903c464158bSHong Zhang . m - number of rows, or number of columns 1904c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1905744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1906744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1907c464158bSHong Zhang 1908c464158bSHong Zhang Output Parameter: 1909c464158bSHong Zhang . A - the symmetric matrix 1910c464158bSHong Zhang 1911c464158bSHong Zhang Options Database Keys: 1912c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1913c464158bSHong Zhang block calculations (much slower) 1914c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1915c464158bSHong Zhang 1916c464158bSHong Zhang Level: intermediate 1917c464158bSHong Zhang 1918175b88e8SBarry Smith It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(), 1919ae1d86c5SBarry Smith MatXXXXSetPreallocation() paradgm instead of this routine directly. 1920175b88e8SBarry Smith [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation] 1921175b88e8SBarry Smith 1922c464158bSHong Zhang Notes: 19236d6d819aSHong Zhang The number of rows and columns must be divisible by blocksize. 19246d6d819aSHong Zhang This matrix type does not support complex Hermitian operation. 1925c464158bSHong Zhang 1926c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1927c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1928c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1929c464158bSHong Zhang matrices. 1930c464158bSHong Zhang 193149a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 193249a6f317SBarry Smith 1933c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1934c464158bSHong Zhang @*/ 1935be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1936c464158bSHong Zhang { 1937dfbe8321SBarry Smith PetscErrorCode ierr; 1938c464158bSHong Zhang 1939c464158bSHong Zhang PetscFunctionBegin; 1940f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1941f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1942c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1943ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr); 194449b5e25fSSatish Balay PetscFunctionReturn(0); 194549b5e25fSSatish Balay } 194649b5e25fSSatish Balay 19474a2ae208SSatish Balay #undef __FUNCT__ 19484a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1949dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 195049b5e25fSSatish Balay { 195149b5e25fSSatish Balay Mat C; 195249b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 19536849ba73SBarry Smith PetscErrorCode ierr; 1954b40805acSSatish Balay PetscInt i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 195549b5e25fSSatish Balay 195649b5e25fSSatish Balay PetscFunctionBegin; 195729bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 195849b5e25fSSatish Balay 195949b5e25fSSatish Balay *B = 0; 19607adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 1961d0f46423SBarry Smith ierr = MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);CHKERRQ(ierr); 19627adad957SLisandro Dalcin ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 19631d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1964692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1965692f9cbeSHong Zhang 1966273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 196749b5e25fSSatish Balay C->factor = A->factor; 196849b5e25fSSatish Balay c->row = 0; 196949b5e25fSSatish Balay c->icol = 0; 197049b5e25fSSatish Balay c->saved_values = 0; 1971a9817697SBarry Smith c->keepnonzeropattern = a->keepnonzeropattern; 197249b5e25fSSatish Balay C->assembled = PETSC_TRUE; 197349b5e25fSSatish Balay 197426283091SBarry Smith ierr = PetscLayoutCopy(A->rmap,&C->rmap);CHKERRQ(ierr); 197526283091SBarry Smith ierr = PetscLayoutCopy(A->cmap,&C->cmap);CHKERRQ(ierr); 197649b5e25fSSatish Balay c->bs2 = a->bs2; 197749b5e25fSSatish Balay c->mbs = a->mbs; 197849b5e25fSSatish Balay c->nbs = a->nbs; 197949b5e25fSSatish Balay 1980c760cd28SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 1981c760cd28SBarry Smith c->imax = a->imax; 1982c760cd28SBarry Smith c->ilen = a->ilen; 1983c760cd28SBarry Smith c->free_imax_ilen = PETSC_FALSE; 1984c760cd28SBarry Smith } else { 19858777fc3fSSatish Balay ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr); 1986c760cd28SBarry Smith ierr = PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 198749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 198849b5e25fSSatish Balay c->imax[i] = a->imax[i]; 198949b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 199049b5e25fSSatish Balay } 1991c760cd28SBarry Smith c->free_imax_ilen = PETSC_TRUE; 1992c760cd28SBarry Smith } 199349b5e25fSSatish Balay 199449b5e25fSSatish Balay /* allocate the matrix space */ 19954da8f245SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 19964da8f245SBarry Smith ierr = PetscMalloc(bs2*nz*sizeof(MatScalar),&c->a);CHKERRQ(ierr); 19974da8f245SBarry Smith ierr = PetscLogObjectMemory(C,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 19984da8f245SBarry Smith c->singlemalloc = PETSC_FALSE; 19994da8f245SBarry Smith c->free_ij = PETSC_FALSE; 20004da8f245SBarry Smith c->parent = A; 20014da8f245SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 20024da8f245SBarry Smith ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 20034da8f245SBarry Smith ierr = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 20044da8f245SBarry Smith } else { 2005b40805acSSatish Balay ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr); 200613f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 2007b40805acSSatish Balay ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr); 20084da8f245SBarry Smith c->singlemalloc = PETSC_TRUE; 20094da8f245SBarry Smith c->free_ij = PETSC_TRUE; 20104da8f245SBarry Smith } 201149b5e25fSSatish Balay if (mbs > 0) { 20124da8f245SBarry Smith if (cpvalues != MAT_SHARE_NONZERO_PATTERN) { 201313f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 20144da8f245SBarry Smith } 201549b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 201649b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 201749b5e25fSSatish Balay } else { 201849b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 201949b5e25fSSatish Balay } 2020a1c3900fSBarry Smith if (a->jshort) { 20214da8f245SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 20224da8f245SBarry Smith c->jshort = a->jshort; 20234da8f245SBarry Smith c->free_jshort = PETSC_FALSE; 20244da8f245SBarry Smith } else { 2025a1c3900fSBarry Smith ierr = PetscMalloc(nz*sizeof(unsigned short),&c->jshort);CHKERRQ(ierr); 2026c760cd28SBarry Smith ierr = PetscLogObjectMemory(C,nz*sizeof(unsigned short));CHKERRQ(ierr); 2027a1c3900fSBarry Smith ierr = PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));CHKERRQ(ierr); 20284da8f245SBarry Smith c->free_jshort = PETSC_TRUE; 20294da8f245SBarry Smith } 2030a1c3900fSBarry Smith } 203149b5e25fSSatish Balay } 203249b5e25fSSatish Balay 203349b5e25fSSatish Balay c->roworiented = a->roworiented; 203449b5e25fSSatish Balay c->nonew = a->nonew; 203549b5e25fSSatish Balay 203649b5e25fSSatish Balay if (a->diag) { 2037c760cd28SBarry Smith if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { 2038c760cd28SBarry Smith c->diag = a->diag; 2039c760cd28SBarry Smith c->free_diag = PETSC_FALSE; 2040c760cd28SBarry Smith } else { 204113f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 204252e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 204349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 204449b5e25fSSatish Balay c->diag[i] = a->diag[i]; 204549b5e25fSSatish Balay } 2046c760cd28SBarry Smith c->free_diag = PETSC_TRUE; 2047c760cd28SBarry Smith } 204849b5e25fSSatish Balay } else c->diag = 0; 20496c6c5352SBarry Smith c->nz = a->nz; 20506c6c5352SBarry Smith c->maxnz = a->maxnz; 205149b5e25fSSatish Balay c->solve_work = 0; 205249b5e25fSSatish Balay c->mult_work = 0; 2053e6b907acSBarry Smith c->free_a = PETSC_TRUE; 205449b5e25fSSatish Balay *B = C; 20557adad957SLisandro Dalcin ierr = PetscFListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);CHKERRQ(ierr); 205649b5e25fSSatish Balay PetscFunctionReturn(0); 205749b5e25fSSatish Balay } 205849b5e25fSSatish Balay 20594a2ae208SSatish Balay #undef __FUNCT__ 20604a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 2061a313700dSBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, const MatType type,Mat *A) 206249b5e25fSSatish Balay { 206349b5e25fSSatish Balay Mat_SeqSBAIJ *a; 206449b5e25fSSatish Balay Mat B; 20656849ba73SBarry Smith PetscErrorCode ierr; 206613f74950SBarry Smith int fd; 206713f74950SBarry Smith PetscMPIInt size; 206813f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 206913f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 207013f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 207113f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 207287828ca2SBarry Smith PetscScalar *aa; 207349b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 207449b5e25fSSatish Balay 207549b5e25fSSatish Balay PetscFunctionBegin; 2076b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 207749b5e25fSSatish Balay bs2 = bs*bs; 207849b5e25fSSatish Balay 207949b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 208029bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 2081b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 208249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 2083552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 208449b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 208549b5e25fSSatish Balay 208649b5e25fSSatish Balay if (header[3] < 0) { 208729bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 208849b5e25fSSatish Balay } 208949b5e25fSSatish Balay 209029bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 209149b5e25fSSatish Balay 209249b5e25fSSatish Balay /* 209349b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 209449b5e25fSSatish Balay divisible by the blocksize 209549b5e25fSSatish Balay */ 209649b5e25fSSatish Balay mbs = M/bs; 209749b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 209849b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 209949b5e25fSSatish Balay else mbs++; 210049b5e25fSSatish Balay if (extra_rows) { 21011e2582c4SBarry Smith ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr); 210249b5e25fSSatish Balay } 210349b5e25fSSatish Balay 210449b5e25fSSatish Balay /* read in row lengths */ 210513f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 210649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 210749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 210849b5e25fSSatish Balay 210949b5e25fSSatish Balay /* read in column indices */ 211013f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 211149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 211249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 211349b5e25fSSatish Balay 211449b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 211513f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 211613f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 211713f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 211813f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 211949b5e25fSSatish Balay masked = mask + mbs; 212049b5e25fSSatish Balay rowcount = 0; nzcount = 0; 212149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 212249b5e25fSSatish Balay nmask = 0; 212349b5e25fSSatish Balay for (j=0; j<bs; j++) { 212449b5e25fSSatish Balay kmax = rowlengths[rowcount]; 212549b5e25fSSatish Balay for (k=0; k<kmax; k++) { 21262d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 212703630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 212849b5e25fSSatish Balay } 212949b5e25fSSatish Balay rowcount++; 213049b5e25fSSatish Balay } 2131574b2666SHong Zhang s_browlengths[i] += nmask; 2132574b2666SHong Zhang 213349b5e25fSSatish Balay /* zero out the mask elements we set */ 213449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 213549b5e25fSSatish Balay } 213649b5e25fSSatish Balay 213749b5e25fSSatish Balay /* create our matrix */ 2138f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 2139f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr); 21409abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 2141ab93d7beSBarry Smith ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr); 214249b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 214349b5e25fSSatish Balay 214449b5e25fSSatish Balay /* set matrix "i" values */ 214549b5e25fSSatish Balay a->i[0] = 0; 214649b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 2147574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 2148574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 214949b5e25fSSatish Balay } 21506c6c5352SBarry Smith a->nz = a->i[mbs]; 215149b5e25fSSatish Balay 215249b5e25fSSatish Balay /* read in nonzero values */ 215387828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 215449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 215549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 215649b5e25fSSatish Balay 215749b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 215849b5e25fSSatish Balay nzcount = 0; jcount = 0; 215949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 216049b5e25fSSatish Balay nzcountb = nzcount; 216149b5e25fSSatish Balay nmask = 0; 216249b5e25fSSatish Balay for (j=0; j<bs; j++) { 216349b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 216449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 21652d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 216603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 21672d703238SHong Zhang } 21682d703238SHong Zhang } 21692d703238SHong Zhang /* sort the masked values */ 21702d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 21712d703238SHong Zhang 21722d703238SHong Zhang /* set "j" values into matrix */ 21732d703238SHong Zhang maskcount = 1; 21742d703238SHong Zhang for (j=0; j<nmask; j++) { 217549b5e25fSSatish Balay a->j[jcount++] = masked[j]; 217649b5e25fSSatish Balay mask[masked[j]] = maskcount++; 217749b5e25fSSatish Balay } 2178574b2666SHong Zhang 217949b5e25fSSatish Balay /* set "a" values into matrix */ 218049b5e25fSSatish Balay ishift = bs2*a->i[i]; 218149b5e25fSSatish Balay for (j=0; j<bs; j++) { 218249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 218349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 2184574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 2185574b2666SHong Zhang if (tmp >= i){ 218649b5e25fSSatish Balay block = mask[tmp] - 1; 218749b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 218849b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 2189574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 2190574b2666SHong Zhang } 2191574b2666SHong Zhang nzcountb++; 219249b5e25fSSatish Balay } 219349b5e25fSSatish Balay } 219449b5e25fSSatish Balay /* zero out the mask elements we set */ 219549b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 219649b5e25fSSatish Balay } 21976c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 219849b5e25fSSatish Balay 219949b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2200574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 220149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 220249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 220349b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 220449b5e25fSSatish Balay 22059abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 22069abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 220749b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 22089abb65ffSKris Buschelman *A = B; 220949b5e25fSSatish Balay PetscFunctionReturn(0); 221049b5e25fSSatish Balay } 2211574b2666SHong Zhang 2212d06b337dSHong Zhang #undef __FUNCT__ 2213c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays" 2214c75a6043SHong Zhang /*@ 2215c75a6043SHong Zhang MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements 2216c75a6043SHong Zhang (upper triangular entries in CSR format) provided by the user. 2217c75a6043SHong Zhang 2218c75a6043SHong Zhang Collective on MPI_Comm 2219c75a6043SHong Zhang 2220c75a6043SHong Zhang Input Parameters: 2221c75a6043SHong Zhang + comm - must be an MPI communicator of size 1 2222c75a6043SHong Zhang . bs - size of block 2223c75a6043SHong Zhang . m - number of rows 2224c75a6043SHong Zhang . n - number of columns 2225c75a6043SHong Zhang . i - row indices 2226c75a6043SHong Zhang . j - column indices 2227c75a6043SHong Zhang - a - matrix values 2228c75a6043SHong Zhang 2229c75a6043SHong Zhang Output Parameter: 2230c75a6043SHong Zhang . mat - the matrix 2231c75a6043SHong Zhang 2232c75a6043SHong Zhang Level: intermediate 2233c75a6043SHong Zhang 2234c75a6043SHong Zhang Notes: 2235c75a6043SHong Zhang The i, j, and a arrays are not copied by this routine, the user must free these arrays 2236c75a6043SHong Zhang once the matrix is destroyed 2237c75a6043SHong Zhang 2238c75a6043SHong Zhang You cannot set new nonzero locations into this matrix, that will generate an error. 2239c75a6043SHong Zhang 2240c75a6043SHong Zhang The i and j indices are 0 based 2241c75a6043SHong Zhang 2242c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ() 2243c75a6043SHong Zhang 2244c75a6043SHong Zhang @*/ 2245c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) 2246c75a6043SHong Zhang { 2247c75a6043SHong Zhang PetscErrorCode ierr; 2248c75a6043SHong Zhang PetscInt ii; 2249c75a6043SHong Zhang Mat_SeqSBAIJ *sbaij; 2250c75a6043SHong Zhang 2251c75a6043SHong Zhang PetscFunctionBegin; 2252c75a6043SHong Zhang if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs); 2253c75a6043SHong Zhang if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 2254c75a6043SHong Zhang 2255c75a6043SHong Zhang ierr = MatCreate(comm,mat);CHKERRQ(ierr); 2256c75a6043SHong Zhang ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); 2257c75a6043SHong Zhang ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr); 2258c75a6043SHong Zhang ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); 2259c75a6043SHong Zhang sbaij = (Mat_SeqSBAIJ*)(*mat)->data; 2260c75a6043SHong Zhang ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr); 2261c760cd28SBarry Smith ierr = PetscLogObjectMemory(*mat,2*m*sizeof(PetscInt));CHKERRQ(ierr); 2262c75a6043SHong Zhang 2263c75a6043SHong Zhang sbaij->i = i; 2264c75a6043SHong Zhang sbaij->j = j; 2265c75a6043SHong Zhang sbaij->a = a; 2266c75a6043SHong Zhang sbaij->singlemalloc = PETSC_FALSE; 2267c75a6043SHong Zhang sbaij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 2268e6b907acSBarry Smith sbaij->free_a = PETSC_FALSE; 2269e6b907acSBarry Smith sbaij->free_ij = PETSC_FALSE; 2270c75a6043SHong Zhang 2271c75a6043SHong Zhang for (ii=0; ii<m; ii++) { 2272c75a6043SHong Zhang sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii]; 2273c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2274c75a6043SHong 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]); 2275c75a6043SHong Zhang #endif 2276c75a6043SHong Zhang } 2277c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG) 2278c75a6043SHong Zhang for (ii=0; ii<sbaij->i[m]; ii++) { 2279c75a6043SHong Zhang if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); 2280c75a6043SHong Zhang if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); 2281c75a6043SHong Zhang } 2282c75a6043SHong Zhang #endif 2283c75a6043SHong Zhang 2284c75a6043SHong Zhang ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2285c75a6043SHong Zhang ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2286c75a6043SHong Zhang PetscFunctionReturn(0); 2287c75a6043SHong Zhang } 2288d06b337dSHong Zhang 2289d06b337dSHong Zhang 2290d06b337dSHong Zhang 229149b5e25fSSatish Balay 229249b5e25fSSatish Balay 2293