149b5e25fSSatish Balay 249b5e25fSSatish Balay /* 3a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 449b5e25fSSatish Balay matrix storage format. 549b5e25fSSatish Balay */ 69e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 749b5e25fSSatish Balay #include "src/inline/spops.h" 8aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 949b5e25fSSatish Balay 1049b5e25fSSatish Balay #define CHUNKSIZE 10 1149b5e25fSSatish Balay 1249b5e25fSSatish Balay /* 1349b5e25fSSatish Balay Checks for missing diagonals 1449b5e25fSSatish Balay */ 154a2ae208SSatish Balay #undef __FUNCT__ 164a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 17dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A) 1849b5e25fSSatish Balay { 19045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 206849ba73SBarry Smith PetscErrorCode ierr; 2113f74950SBarry Smith PetscInt *diag,*jj = a->j,i; 2249b5e25fSSatish Balay 2349b5e25fSSatish Balay PetscFunctionBegin; 24045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2549b5e25fSSatish Balay diag = a->diag; 2649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2777431f27SBarry Smith if (jj[diag[i]] != i) SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Matrix is missing diagonal number %D",i); 2849b5e25fSSatish Balay } 2949b5e25fSSatish Balay PetscFunctionReturn(0); 3049b5e25fSSatish Balay } 3149b5e25fSSatish Balay 324a2ae208SSatish Balay #undef __FUNCT__ 334a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 34dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A) 3549b5e25fSSatish Balay { 36045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 376849ba73SBarry Smith PetscErrorCode ierr; 3813f74950SBarry Smith PetscInt i,mbs = a->mbs; 3949b5e25fSSatish Balay 4049b5e25fSSatish Balay PetscFunctionBegin; 4149b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4249b5e25fSSatish Balay 4313f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 4452e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 45b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4649b5e25fSSatish Balay PetscFunctionReturn(0); 4749b5e25fSSatish Balay } 4849b5e25fSSatish Balay 494a2ae208SSatish Balay #undef __FUNCT__ 504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 5113f74950SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 5249b5e25fSSatish Balay { 53a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5413f74950SBarry Smith PetscInt n = a->mbs,i; 5549b5e25fSSatish Balay 5649b5e25fSSatish Balay PetscFunctionBegin; 57d3e5a4abSHong Zhang *nn = n; 58a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 59a1373b80SHong Zhang 60a6ece127SHong Zhang if (oshift == 1) { 61a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 6213f74950SBarry Smith PetscInt nz = a->i[n]; 636c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 64a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 65a1373b80SHong Zhang *ia = a->i; *ja = a->j; 66a1373b80SHong Zhang } else { 67a1373b80SHong Zhang *ia = a->i; *ja = a->j; 68a6ece127SHong Zhang } 6949b5e25fSSatish Balay PetscFunctionReturn(0); 7049b5e25fSSatish Balay } 7149b5e25fSSatish Balay 724a2ae208SSatish Balay #undef __FUNCT__ 734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 7413f74950SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 7549b5e25fSSatish Balay { 76b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7713f74950SBarry Smith PetscInt i,n = a->mbs; 78a6ece127SHong Zhang 7949b5e25fSSatish Balay PetscFunctionBegin; 8049b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 81a6ece127SHong Zhang 82a6ece127SHong Zhang if (oshift == 1) { 8313f74950SBarry Smith PetscInt nz = a->i[n]-1; 846c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 85a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 86a6ece127SHong Zhang } 87a6ece127SHong Zhang PetscFunctionReturn(0); 8849b5e25fSSatish Balay } 8949b5e25fSSatish Balay 904a2ae208SSatish Balay #undef __FUNCT__ 914a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 92dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A) 9349b5e25fSSatish Balay { 9449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 95dfbe8321SBarry Smith PetscErrorCode ierr; 9649b5e25fSSatish Balay 9749b5e25fSSatish Balay PetscFunctionBegin; 9877431f27SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->m,a->nz); 9949b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 10049b5e25fSSatish Balay if (!a->singlemalloc) { 10149b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 10263c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 10349b5e25fSSatish Balay } 10449b5e25fSSatish Balay if (a->row) { 10549b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 10649b5e25fSSatish Balay } 10749b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 10849b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 10949b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 11049b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 111d59c15a7SBarry Smith if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 112d59c15a7SBarry Smith if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);} 113d59c15a7SBarry Smith if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 11449b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 115c4319e64SHong Zhang if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 1161a3463dfSHong Zhang 1171a3463dfSHong Zhang if (a->inew){ 1181a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1191a3463dfSHong Zhang a->inew = 0; 1201a3463dfSHong Zhang } 12149b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 122901853e0SKris Buschelman 123901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); 124901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); 125901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); 126901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 127901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); 128901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 12949b5e25fSSatish Balay PetscFunctionReturn(0); 13049b5e25fSSatish Balay } 13149b5e25fSSatish Balay 1324a2ae208SSatish Balay #undef __FUNCT__ 1334a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 134dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op) 13549b5e25fSSatish Balay { 136045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 13749b5e25fSSatish Balay 13849b5e25fSSatish Balay PetscFunctionBegin; 1394d9d31abSKris Buschelman switch (op) { 1404d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1414d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1424d9d31abSKris Buschelman break; 1434d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1444d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1454d9d31abSKris Buschelman break; 1464d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1474d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1484d9d31abSKris Buschelman break; 1494d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1504d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1514d9d31abSKris Buschelman break; 1524d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1534d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1544d9d31abSKris Buschelman break; 1554d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1564d9d31abSKris Buschelman a->nonew = 1; 1574d9d31abSKris Buschelman break; 1584d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1594d9d31abSKris Buschelman a->nonew = -1; 1604d9d31abSKris Buschelman break; 1614d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1624d9d31abSKris Buschelman a->nonew = -2; 1634d9d31abSKris Buschelman break; 1644d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1654d9d31abSKris Buschelman a->nonew = 0; 1664d9d31abSKris Buschelman break; 1674d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1684d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1694d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1704d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1714d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 172b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1734d9d31abSKris Buschelman break; 1744d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 17529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1769a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1779a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1789a4540c5SBarry Smith case MAT_HERMITIAN: 1799a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 18077e54ba9SKris Buschelman case MAT_SYMMETRIC: 18177e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1829a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1839a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1849a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 18577e54ba9SKris Buschelman break; 1864d9d31abSKris Buschelman default: 18729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 18849b5e25fSSatish Balay } 18949b5e25fSSatish Balay PetscFunctionReturn(0); 19049b5e25fSSatish Balay } 19149b5e25fSSatish Balay 1924a2ae208SSatish Balay #undef __FUNCT__ 1934a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 19413f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v) 19549b5e25fSSatish Balay { 19649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1976849ba73SBarry Smith PetscErrorCode ierr; 19813f74950SBarry Smith PetscInt itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 19949b5e25fSSatish Balay MatScalar *aa,*aa_i; 20087828ca2SBarry Smith PetscScalar *v_i; 20149b5e25fSSatish Balay 20249b5e25fSSatish Balay PetscFunctionBegin; 203521d7252SBarry Smith bs = A->bs; 20449b5e25fSSatish Balay ai = a->i; 20549b5e25fSSatish Balay aj = a->j; 20649b5e25fSSatish Balay aa = a->a; 20749b5e25fSSatish Balay bs2 = a->bs2; 20849b5e25fSSatish Balay 20977431f27SBarry Smith if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row); 21049b5e25fSSatish Balay 21149b5e25fSSatish Balay bn = row/bs; /* Block number */ 21249b5e25fSSatish Balay bp = row % bs; /* Block position */ 21349b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 21449b5e25fSSatish Balay *ncols = bs*M; 21549b5e25fSSatish Balay 21649b5e25fSSatish Balay if (v) { 21749b5e25fSSatish Balay *v = 0; 21849b5e25fSSatish Balay if (*ncols) { 21987828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 22049b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22149b5e25fSSatish Balay v_i = *v + i*bs; 22249b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 22349b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 22449b5e25fSSatish Balay } 22549b5e25fSSatish Balay } 22649b5e25fSSatish Balay } 22749b5e25fSSatish Balay 22849b5e25fSSatish Balay if (cols) { 22949b5e25fSSatish Balay *cols = 0; 23049b5e25fSSatish Balay if (*ncols) { 23113f74950SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr); 23249b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23349b5e25fSSatish Balay cols_i = *cols + i*bs; 23449b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 23549b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay } 23849b5e25fSSatish Balay } 23949b5e25fSSatish Balay 24049b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2415ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2425ddb2528SHong Zhang #ifdef column_search 24349b5e25fSSatish Balay v_i = *v + M*bs; 24449b5e25fSSatish Balay cols_i = *cols + M*bs; 24549b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 24649b5e25fSSatish Balay M = ai[i+1] - ai[i]; 24749b5e25fSSatish Balay for (j=0; j<M; j++){ 24849b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 24949b5e25fSSatish Balay if (itmp == bn){ 25049b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 25149b5e25fSSatish Balay for (k=0; k<bs; k++) { 25249b5e25fSSatish Balay *cols_i++ = i*bs+k; 25349b5e25fSSatish Balay *v_i++ = aa_i[k]; 25449b5e25fSSatish Balay } 25549b5e25fSSatish Balay *ncols += bs; 25649b5e25fSSatish Balay break; 25749b5e25fSSatish Balay } 25849b5e25fSSatish Balay } 25949b5e25fSSatish Balay } 2605ddb2528SHong Zhang #endif 26149b5e25fSSatish Balay PetscFunctionReturn(0); 26249b5e25fSSatish Balay } 26349b5e25fSSatish Balay 2644a2ae208SSatish Balay #undef __FUNCT__ 2654a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 26613f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 26749b5e25fSSatish Balay { 268dfbe8321SBarry Smith PetscErrorCode ierr; 26949b5e25fSSatish Balay 27049b5e25fSSatish Balay PetscFunctionBegin; 27149b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 27249b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 27349b5e25fSSatish Balay PetscFunctionReturn(0); 27449b5e25fSSatish Balay } 27549b5e25fSSatish Balay 2764a2ae208SSatish Balay #undef __FUNCT__ 2774a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 278dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 27949b5e25fSSatish Balay { 280dfbe8321SBarry Smith PetscErrorCode ierr; 28149b5e25fSSatish Balay PetscFunctionBegin; 282999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2838115998fSBarry Smith PetscFunctionReturn(0); 28449b5e25fSSatish Balay } 28549b5e25fSSatish Balay 2864a2ae208SSatish Balay #undef __FUNCT__ 2874a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 2886849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 28949b5e25fSSatish Balay { 29049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 291dfbe8321SBarry Smith PetscErrorCode ierr; 292521d7252SBarry Smith PetscInt i,j,bs = A->bs,k,l,bs2=a->bs2; 293fb9695e5SSatish Balay char *name; 294f3ef73ceSBarry Smith PetscViewerFormat format; 29549b5e25fSSatish Balay 29649b5e25fSSatish Balay PetscFunctionBegin; 29780fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 298b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 299456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 30077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 301fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 30229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 303fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 304b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 30549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30649b5e25fSSatish Balay for (j=0; j<bs; j++) { 30777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 30849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 30949b5e25fSSatish Balay for (l=0; l<bs; l++) { 31049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 31149b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 31349b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 31449b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 31649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 31749b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 31949b5e25fSSatish Balay } 32049b5e25fSSatish Balay #else 32149b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 32277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 32349b5e25fSSatish Balay } 32449b5e25fSSatish Balay #endif 32549b5e25fSSatish Balay } 32649b5e25fSSatish Balay } 327b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 32849b5e25fSSatish Balay } 32949b5e25fSSatish Balay } 330b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 33149b5e25fSSatish Balay } else { 332b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 33349b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33449b5e25fSSatish Balay for (j=0; j<bs; j++) { 33577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 33649b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33749b5e25fSSatish Balay for (l=0; l<bs; l++) { 33849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 33949b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 34077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 34149b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34249b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 34377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 34449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34549b5e25fSSatish Balay } else { 34677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34749b5e25fSSatish Balay } 34849b5e25fSSatish Balay #else 34977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35049b5e25fSSatish Balay #endif 35149b5e25fSSatish Balay } 35249b5e25fSSatish Balay } 353b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 35449b5e25fSSatish Balay } 35549b5e25fSSatish Balay } 356b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 35749b5e25fSSatish Balay } 358b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 35949b5e25fSSatish Balay PetscFunctionReturn(0); 36049b5e25fSSatish Balay } 36149b5e25fSSatish Balay 3624a2ae208SSatish Balay #undef __FUNCT__ 3634a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 3646849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 36549b5e25fSSatish Balay { 36649b5e25fSSatish Balay Mat A = (Mat) Aa; 36749b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 3686849ba73SBarry Smith PetscErrorCode ierr; 36913f74950SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2; 37013f74950SBarry Smith PetscMPIInt rank; 37149b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 37249b5e25fSSatish Balay MatScalar *aa; 37349b5e25fSSatish Balay MPI_Comm comm; 374b0a32e0cSBarry Smith PetscViewer viewer; 37549b5e25fSSatish Balay 37649b5e25fSSatish Balay PetscFunctionBegin; 37749b5e25fSSatish Balay /* 37849b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 37949b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 38049b5e25fSSatish Balay rest should return immediately. 38149b5e25fSSatish Balay */ 38249b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 38349b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 38449b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 38549b5e25fSSatish Balay 38649b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 38749b5e25fSSatish Balay 388b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 389b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 39049b5e25fSSatish Balay 39149b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 392b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 39349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 39449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 395c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 39649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 39749b5e25fSSatish Balay aa = a->a + j*bs2; 39849b5e25fSSatish Balay for (k=0; k<bs; k++) { 39949b5e25fSSatish Balay for (l=0; l<bs; l++) { 40049b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 401b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 40249b5e25fSSatish Balay } 40349b5e25fSSatish Balay } 40449b5e25fSSatish Balay } 40549b5e25fSSatish Balay } 406b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 40749b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 40849b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 409c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 41049b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 41149b5e25fSSatish Balay aa = a->a + j*bs2; 41249b5e25fSSatish Balay for (k=0; k<bs; k++) { 41349b5e25fSSatish Balay for (l=0; l<bs; l++) { 41449b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 415b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 41649b5e25fSSatish Balay } 41749b5e25fSSatish Balay } 41849b5e25fSSatish Balay } 41949b5e25fSSatish Balay } 42049b5e25fSSatish Balay 421b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 42249b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 42349b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 424c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 42549b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42649b5e25fSSatish Balay aa = a->a + j*bs2; 42749b5e25fSSatish Balay for (k=0; k<bs; k++) { 42849b5e25fSSatish Balay for (l=0; l<bs; l++) { 42949b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 430b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 43149b5e25fSSatish Balay } 43249b5e25fSSatish Balay } 43349b5e25fSSatish Balay } 43449b5e25fSSatish Balay } 43549b5e25fSSatish Balay PetscFunctionReturn(0); 43649b5e25fSSatish Balay } 43749b5e25fSSatish Balay 4384a2ae208SSatish Balay #undef __FUNCT__ 4394a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4406849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 44149b5e25fSSatish Balay { 442dfbe8321SBarry Smith PetscErrorCode ierr; 44349b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 444b0a32e0cSBarry Smith PetscDraw draw; 44549b5e25fSSatish Balay PetscTruth isnull; 44649b5e25fSSatish Balay 44749b5e25fSSatish Balay PetscFunctionBegin; 44849b5e25fSSatish Balay 449b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 450b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 45149b5e25fSSatish Balay 45249b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 453c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 45449b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 455b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 456b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 45749b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 45849b5e25fSSatish Balay PetscFunctionReturn(0); 45949b5e25fSSatish Balay } 46049b5e25fSSatish Balay 4614a2ae208SSatish Balay #undef __FUNCT__ 4624a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 463dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 46449b5e25fSSatish Balay { 465dfbe8321SBarry Smith PetscErrorCode ierr; 46632077d6dSBarry Smith PetscTruth iascii,isdraw; 46749b5e25fSSatish Balay 46849b5e25fSSatish Balay PetscFunctionBegin; 46932077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 470fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 47132077d6dSBarry Smith if (iascii){ 47249b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 47349b5e25fSSatish Balay } else if (isdraw) { 47449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 47549b5e25fSSatish Balay } else { 476a5e6ed63SBarry Smith Mat B; 477*ceb03754SKris Buschelman ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 478a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 479a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 48049b5e25fSSatish Balay } 48149b5e25fSSatish Balay PetscFunctionReturn(0); 48249b5e25fSSatish Balay } 48349b5e25fSSatish Balay 48449b5e25fSSatish Balay 4854a2ae208SSatish Balay #undef __FUNCT__ 4864a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 48713f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 48849b5e25fSSatish Balay { 489045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 49013f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 49113f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 49213f74950SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2; 49349b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 49449b5e25fSSatish Balay 49549b5e25fSSatish Balay PetscFunctionBegin; 49649b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 49749b5e25fSSatish Balay row = im[k]; brow = row/bs; 49877431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 49977431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 50049b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 50149b5e25fSSatish Balay nrow = ailen[brow]; 50249b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 50377431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 50477431f27SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1); 50549b5e25fSSatish Balay col = in[l] ; 50649b5e25fSSatish Balay bcol = col/bs; 50749b5e25fSSatish Balay cidx = col%bs; 50849b5e25fSSatish Balay ridx = row%bs; 50949b5e25fSSatish Balay high = nrow; 51049b5e25fSSatish Balay low = 0; /* assume unsorted */ 51149b5e25fSSatish Balay while (high-low > 5) { 51249b5e25fSSatish Balay t = (low+high)/2; 51349b5e25fSSatish Balay if (rp[t] > bcol) high = t; 51449b5e25fSSatish Balay else low = t; 51549b5e25fSSatish Balay } 51649b5e25fSSatish Balay for (i=low; i<high; i++) { 51749b5e25fSSatish Balay if (rp[i] > bcol) break; 51849b5e25fSSatish Balay if (rp[i] == bcol) { 51949b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 52049b5e25fSSatish Balay goto finished; 52149b5e25fSSatish Balay } 52249b5e25fSSatish Balay } 52349b5e25fSSatish Balay *v++ = zero; 52449b5e25fSSatish Balay finished:; 52549b5e25fSSatish Balay } 52649b5e25fSSatish Balay } 52749b5e25fSSatish Balay PetscFunctionReturn(0); 52849b5e25fSSatish Balay } 52949b5e25fSSatish Balay 53049b5e25fSSatish Balay 5314a2ae208SSatish Balay #undef __FUNCT__ 5324a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 53313f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 53449b5e25fSSatish Balay { 5350880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5366849ba73SBarry Smith PetscErrorCode ierr; 537e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1; 53813f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 53913f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval; 5400880e062SHong Zhang PetscTruth roworiented=a->roworiented; 541f15d580aSBarry Smith const MatScalar *value = v; 542f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5430880e062SHong Zhang 54449b5e25fSSatish Balay PetscFunctionBegin; 5450880e062SHong Zhang if (roworiented) { 5460880e062SHong Zhang stepval = (n-1)*bs; 5470880e062SHong Zhang } else { 5480880e062SHong Zhang stepval = (m-1)*bs; 5490880e062SHong Zhang } 5500880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5510880e062SHong Zhang row = im[k]; 5520880e062SHong Zhang if (row < 0) continue; 5532515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 55477431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5550880e062SHong Zhang #endif 5560880e062SHong Zhang rp = aj + ai[row]; 5570880e062SHong Zhang ap = aa + bs2*ai[row]; 5580880e062SHong Zhang rmax = imax[row]; 5590880e062SHong Zhang nrow = ailen[row]; 5600880e062SHong Zhang low = 0; 561818f2c47SBarry Smith high = nrow; 5620880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5630880e062SHong Zhang if (in[l] < 0) continue; 5640880e062SHong Zhang col = in[l]; 5652515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 56677431f27SBarry Smith if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1); 567b1823623SSatish Balay #endif 5680880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 5690880e062SHong Zhang if (roworiented) { 5700880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 5710880e062SHong Zhang } else { 5720880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 5730880e062SHong Zhang } 574818f2c47SBarry Smith if (col < lastcol) low = 0; else high = nrow; 575e2ee6c50SBarry Smith lastcol = col; 5760880e062SHong Zhang while (high-low > 7) { 5770880e062SHong Zhang t = (low+high)/2; 5780880e062SHong Zhang if (rp[t] > col) high = t; 5790880e062SHong Zhang else low = t; 5800880e062SHong Zhang } 5810880e062SHong Zhang for (i=low; i<high; i++) { 5820880e062SHong Zhang if (rp[i] > col) break; 5830880e062SHong Zhang if (rp[i] == col) { 5840880e062SHong Zhang bap = ap + bs2*i; 5850880e062SHong Zhang if (roworiented) { 5860880e062SHong Zhang if (is == ADD_VALUES) { 5870880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5880880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5890880e062SHong Zhang bap[jj] += *value++; 5900880e062SHong Zhang } 5910880e062SHong Zhang } 5920880e062SHong Zhang } else { 5930880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5940880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5950880e062SHong Zhang bap[jj] = *value++; 5960880e062SHong Zhang } 5970880e062SHong Zhang } 5980880e062SHong Zhang } 5990880e062SHong Zhang } else { 6000880e062SHong Zhang if (is == ADD_VALUES) { 6010880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6020880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6030880e062SHong Zhang *bap++ += *value++; 6040880e062SHong Zhang } 6050880e062SHong Zhang } 6060880e062SHong Zhang } else { 6070880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6080880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6090880e062SHong Zhang *bap++ = *value++; 6100880e062SHong Zhang } 6110880e062SHong Zhang } 6120880e062SHong Zhang } 6130880e062SHong Zhang } 6140880e062SHong Zhang goto noinsert2; 6150880e062SHong Zhang } 6160880e062SHong Zhang } 6170880e062SHong Zhang if (nonew == 1) goto noinsert2; 61877431f27SBarry Smith else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 6190880e062SHong Zhang if (nrow >= rmax) { 6200880e062SHong Zhang /* there is no extra room in row, therefore enlarge */ 62113f74950SBarry Smith PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6220880e062SHong Zhang MatScalar *new_a; 6230880e062SHong Zhang 62477431f27SBarry Smith if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 6250880e062SHong Zhang 6260880e062SHong Zhang /* malloc new storage space */ 62713f74950SBarry Smith len = new_nz*(sizeof(PetscInt)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(PetscInt); 6280880e062SHong Zhang ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 62913f74950SBarry Smith new_j = (PetscInt*)(new_a + bs2*new_nz); 6300880e062SHong Zhang new_i = new_j + new_nz; 6310880e062SHong Zhang 6320880e062SHong Zhang /* copy over old data into new slots */ 6330880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6340880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 63513f74950SBarry Smith ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 6360880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 63713f74950SBarry Smith ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 6380880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 6390880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 6400880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 6410880e062SHong Zhang /* free up old matrix storage */ 6420880e062SHong Zhang ierr = PetscFree(a->a);CHKERRQ(ierr); 6430880e062SHong Zhang if (!a->singlemalloc) { 6440880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 6450880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 6460880e062SHong Zhang } 6470880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 6480880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 6490880e062SHong Zhang 6500880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 6510880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 65252e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar)));CHKERRQ(ierr); 6536c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 6540880e062SHong Zhang a->reallocs++; 6556c6c5352SBarry Smith a->nz++; 6560880e062SHong Zhang } 6570880e062SHong Zhang N = nrow++ - 1; 6580880e062SHong Zhang /* shift up all the later entries in this row */ 6590880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6600880e062SHong Zhang rp[ii+1] = rp[ii]; 6610880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6620880e062SHong Zhang } 6630880e062SHong Zhang if (N >= i) { 6640880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6650880e062SHong Zhang } 6660880e062SHong Zhang rp[i] = col; 6670880e062SHong Zhang bap = ap + bs2*i; 6680880e062SHong Zhang if (roworiented) { 6690880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6700880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6710880e062SHong Zhang bap[jj] = *value++; 6720880e062SHong Zhang } 6730880e062SHong Zhang } 6740880e062SHong Zhang } else { 6750880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6760880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6770880e062SHong Zhang *bap++ = *value++; 6780880e062SHong Zhang } 6790880e062SHong Zhang } 6800880e062SHong Zhang } 6810880e062SHong Zhang noinsert2:; 6820880e062SHong Zhang low = i; 6830880e062SHong Zhang } 6840880e062SHong Zhang ailen[row] = nrow; 6850880e062SHong Zhang } 6860880e062SHong Zhang PetscFunctionReturn(0); 68749b5e25fSSatish Balay } 68849b5e25fSSatish Balay 6894a2ae208SSatish Balay #undef __FUNCT__ 6904a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 691dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 69249b5e25fSSatish Balay { 69349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6946849ba73SBarry Smith PetscErrorCode ierr; 69513f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 69613f74950SBarry Smith PetscInt m = A->m,*ip,N,*ailen = a->ilen; 69713f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 69849b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 69949b5e25fSSatish Balay 70049b5e25fSSatish Balay PetscFunctionBegin; 70149b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 70249b5e25fSSatish Balay 70349b5e25fSSatish Balay if (m) rmax = ailen[0]; 70449b5e25fSSatish Balay for (i=1; i<mbs; i++) { 70549b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 70649b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 70749b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 70849b5e25fSSatish Balay if (fshift) { 70949b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 71049b5e25fSSatish Balay N = ailen[i]; 71149b5e25fSSatish Balay for (j=0; j<N; j++) { 71249b5e25fSSatish Balay ip[j-fshift] = ip[j]; 71349b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 71449b5e25fSSatish Balay } 71549b5e25fSSatish Balay } 71649b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 71749b5e25fSSatish Balay } 71849b5e25fSSatish Balay if (mbs) { 71949b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 72049b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 72149b5e25fSSatish Balay } 72249b5e25fSSatish Balay /* reset ilen and imax for each row */ 72349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 72449b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 72549b5e25fSSatish Balay } 7266c6c5352SBarry Smith a->nz = ai[mbs]; 72749b5e25fSSatish Balay 728b424e231SHong Zhang /* diagonals may have moved, reset it */ 729b424e231SHong Zhang if (a->diag) { 73013f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 73149b5e25fSSatish Balay } 73277431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n", 733521d7252SBarry Smith m,A->m,A->bs,fshift*bs2,a->nz*bs2); 73477431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %D\n", 73549b5e25fSSatish Balay a->reallocs); 73677431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %D\n",rmax); 73749b5e25fSSatish Balay a->reallocs = 0; 73849b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 73949b5e25fSSatish Balay PetscFunctionReturn(0); 74049b5e25fSSatish Balay } 74149b5e25fSSatish Balay 74249b5e25fSSatish Balay /* 74349b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 74449b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 74549b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 74649b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 74749b5e25fSSatish Balay */ 7484a2ae208SSatish Balay #undef __FUNCT__ 7494a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 75013f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max) 75149b5e25fSSatish Balay { 75213f74950SBarry Smith PetscInt i,j,k,row; 75349b5e25fSSatish Balay PetscTruth flg; 75449b5e25fSSatish Balay 75549b5e25fSSatish Balay PetscFunctionBegin; 75649b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 75749b5e25fSSatish Balay row = idx[i]; 75849b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 75949b5e25fSSatish Balay sizes[j] = 1; 76049b5e25fSSatish Balay i++; 76149b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 76249b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 76349b5e25fSSatish Balay i++; 76449b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 76549b5e25fSSatish Balay flg = PETSC_TRUE; 76649b5e25fSSatish Balay for (k=1; k<bs; k++) { 76749b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 76849b5e25fSSatish Balay flg = PETSC_FALSE; 76949b5e25fSSatish Balay break; 77049b5e25fSSatish Balay } 77149b5e25fSSatish Balay } 772abc0a331SBarry Smith if (flg) { /* No break in the bs */ 77349b5e25fSSatish Balay sizes[j] = bs; 77449b5e25fSSatish Balay i+= bs; 77549b5e25fSSatish Balay } else { 77649b5e25fSSatish Balay sizes[j] = 1; 77749b5e25fSSatish Balay i++; 77849b5e25fSSatish Balay } 77949b5e25fSSatish Balay } 78049b5e25fSSatish Balay } 78149b5e25fSSatish Balay *bs_max = j; 78249b5e25fSSatish Balay PetscFunctionReturn(0); 78349b5e25fSSatish Balay } 78449b5e25fSSatish Balay 78549b5e25fSSatish Balay 78649b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 78749b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 78849b5e25fSSatish Balay */ 78949b5e25fSSatish Balay 7904a2ae208SSatish Balay #undef __FUNCT__ 7914a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 79213f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 79349b5e25fSSatish Balay { 79449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 7956849ba73SBarry Smith PetscErrorCode ierr; 796e2ee6c50SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1; 79713f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 79813f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs=A->bs,brow,bcol; 79913f74950SBarry Smith PetscInt ridx,cidx,bs2=a->bs2; 80049b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 80149b5e25fSSatish Balay 80249b5e25fSSatish Balay PetscFunctionBegin; 80349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 80449b5e25fSSatish Balay row = im[k]; /* row number */ 80549b5e25fSSatish Balay brow = row/bs; /* block row number */ 80649b5e25fSSatish Balay if (row < 0) continue; 8072515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 80877431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 80949b5e25fSSatish Balay #endif 81049b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 81149b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 81249b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 81349b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 81449b5e25fSSatish Balay low = 0; 81549b5e25fSSatish Balay 81649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 81749b5e25fSSatish Balay if (in[l] < 0) continue; 8182515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 81977431f27SBarry Smith if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1); 82049b5e25fSSatish Balay #endif 82149b5e25fSSatish Balay col = in[l]; 82249b5e25fSSatish Balay bcol = col/bs; /* block col number */ 82349b5e25fSSatish Balay 824f4989cb3SHong Zhang if (brow > bcol) continue; /* ignore lower triangular values */ 825f4989cb3SHong Zhang 82649b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8278549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 82849b5e25fSSatish Balay /* element value a(k,l) */ 82949b5e25fSSatish Balay if (roworiented) { 83049b5e25fSSatish Balay value = v[l + k*n]; 83149b5e25fSSatish Balay } else { 83249b5e25fSSatish Balay value = v[k + l*m]; 83349b5e25fSSatish Balay } 83449b5e25fSSatish Balay 83549b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 836e2ee6c50SBarry Smith if (col < lastcol) low = 0; high = nrow; 837e2ee6c50SBarry Smith lastcol = col; 83849b5e25fSSatish Balay while (high-low > 7) { 83949b5e25fSSatish Balay t = (low+high)/2; 84049b5e25fSSatish Balay if (rp[t] > bcol) high = t; 84149b5e25fSSatish Balay else low = t; 84249b5e25fSSatish Balay } 84349b5e25fSSatish Balay for (i=low; i<high; i++) { 84449b5e25fSSatish Balay if (rp[i] > bcol) break; 84549b5e25fSSatish Balay if (rp[i] == bcol) { 84649b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 84749b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 84849b5e25fSSatish Balay else *bap = value; 8498549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8508549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8518549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8528549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8538549e402SHong Zhang else *bap = value; 8548549e402SHong Zhang } 85549b5e25fSSatish Balay goto noinsert1; 85649b5e25fSSatish Balay } 85749b5e25fSSatish Balay } 85849b5e25fSSatish Balay 85949b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 86077431f27SBarry Smith else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 86149b5e25fSSatish Balay if (nrow >= rmax) { 86249b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 86313f74950SBarry Smith PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 86449b5e25fSSatish Balay MatScalar *new_a; 86549b5e25fSSatish Balay 86677431f27SBarry Smith if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 86749b5e25fSSatish Balay 86849b5e25fSSatish Balay /* Malloc new storage space */ 86913f74950SBarry Smith len = new_nz*(sizeof(PetscInt)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(PetscInt); 87082502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 87113f74950SBarry Smith new_j = (PetscInt*)(new_a + bs2*new_nz); 87249b5e25fSSatish Balay new_i = new_j + new_nz; 87349b5e25fSSatish Balay 87449b5e25fSSatish Balay /* copy over old data into new slots */ 87549b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 87649b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 87713f74950SBarry Smith ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 87849b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 87913f74950SBarry Smith ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 88049b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 88149b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 88249b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 88349b5e25fSSatish Balay /* free up old matrix storage */ 88449b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 88549b5e25fSSatish Balay if (!a->singlemalloc) { 88649b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 88749b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 88849b5e25fSSatish Balay } 88949b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 89049b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 89149b5e25fSSatish Balay 89249b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 89349b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 89452e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar)));CHKERRQ(ierr); 8956c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 89649b5e25fSSatish Balay a->reallocs++; 8976c6c5352SBarry Smith a->nz++; 89849b5e25fSSatish Balay } 89949b5e25fSSatish Balay 90049b5e25fSSatish Balay N = nrow++ - 1; 90149b5e25fSSatish Balay /* shift up all the later entries in this row */ 90249b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 90349b5e25fSSatish Balay rp[ii+1] = rp[ii]; 90449b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 90549b5e25fSSatish Balay } 90649b5e25fSSatish Balay if (N>=i) { 90749b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 90849b5e25fSSatish Balay } 90949b5e25fSSatish Balay rp[i] = bcol; 91049b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 91149b5e25fSSatish Balay noinsert1:; 91249b5e25fSSatish Balay low = i; 9138549e402SHong Zhang } 91449b5e25fSSatish Balay } /* end of loop over added columns */ 91549b5e25fSSatish Balay ailen[brow] = nrow; 91649b5e25fSSatish Balay } /* end of loop over added rows */ 91749b5e25fSSatish Balay PetscFunctionReturn(0); 91849b5e25fSSatish Balay } 91949b5e25fSSatish Balay 9206849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 9216849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 92213f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt); 92313f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 92413f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]); 9256849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat); 9266849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 9276849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 9286849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec); 9296849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 9306849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 9316849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat); 9326849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec); 93349b5e25fSSatish Balay 9346849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 9356849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 9366849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 9376849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 9386849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 9396849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 9406849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 9416849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 94249b5e25fSSatish Balay 9436849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 944d59c15a7SBarry Smith 9456849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 9466849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 9476849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 9486849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 9496849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 9506849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 9516849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 9526849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 95307f98182SSatish Balay 954af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,MatFactorInfo*,Mat*); 955af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,MatFactorInfo*,Mat*); 956af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,MatFactorInfo*,Mat*); 957af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,MatFactorInfo*,Mat*); 958af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,MatFactorInfo*,Mat*); 959af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,MatFactorInfo*,Mat*); 960af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,MatFactorInfo*,Mat*); 961af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,MatFactorInfo*,Mat*); 96213f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*); 96349b5e25fSSatish Balay 9646849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 9656849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 9666849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 9676849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 9686849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 9696849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 9706849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 9716849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 97249b5e25fSSatish Balay 9736849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 9746849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 9756849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 9766849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 9776849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 9786849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 9796849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 9806849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 98149b5e25fSSatish Balay 9824d101231SSatish Balay #ifdef HAVE_MatICCFactor 9836849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work! */ 9844a2ae208SSatish Balay #undef __FUNCT__ 9854d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 986dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 98749b5e25fSSatish Balay { 9884ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 98949b5e25fSSatish Balay Mat outA; 990dfbe8321SBarry Smith PetscErrorCode ierr; 99149b5e25fSSatish Balay PetscTruth row_identity,col_identity; 99249b5e25fSSatish Balay 99349b5e25fSSatish Balay PetscFunctionBegin; 99449b5e25fSSatish Balay outA = inA; 9951260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 99649b5e25fSSatish Balay 99749b5e25fSSatish Balay if (!a->diag) { 9981a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 99949b5e25fSSatish Balay } 100049b5e25fSSatish Balay /* 100149b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 100249b5e25fSSatish Balay for ILU(0) factorization with natural ordering 100349b5e25fSSatish Balay */ 100449b5e25fSSatish Balay switch (a->bs) { 100549b5e25fSSatish Balay case 1: 10060fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 100707f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1008d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 10094d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 101049b5e25fSSatish Balay case 2: 10111a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10121a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10131a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10144d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 101549b5e25fSSatish Balay break; 101649b5e25fSSatish Balay case 3: 10171a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10181a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10191a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10204d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 102149b5e25fSSatish Balay break; 102249b5e25fSSatish Balay case 4: 10231a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10241a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10251a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10264d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 102749b5e25fSSatish Balay break; 102849b5e25fSSatish Balay case 5: 10291a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10301a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10311a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10324d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 103349b5e25fSSatish Balay break; 103449b5e25fSSatish Balay case 6: 10351a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10361a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10371a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10384d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 103949b5e25fSSatish Balay break; 104049b5e25fSSatish Balay case 7: 10411a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10421a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10431a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10444d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 104549b5e25fSSatish Balay break; 104649b5e25fSSatish Balay default: 104749b5e25fSSatish Balay a->row = row; 10481a3463dfSHong Zhang a->icol = col; 104949b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 105049b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 105149b5e25fSSatish Balay 105249b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 105349b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 105452e6d16bSBarry Smith ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); 105549b5e25fSSatish Balay 105649b5e25fSSatish Balay if (!a->solve_work) { 105787828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 105852e6d16bSBarry Smith ierr = PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));CHKERRQ(ierr); 105949b5e25fSSatish Balay } 106049b5e25fSSatish Balay } 106149b5e25fSSatish Balay 1062af281ebdSHong Zhang ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr); 106349b5e25fSSatish Balay PetscFunctionReturn(0); 106449b5e25fSSatish Balay } 1065950f1e5bSHong Zhang #endif 1066950f1e5bSHong Zhang 10674a2ae208SSatish Balay #undef __FUNCT__ 10684a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 1069dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A) 107049b5e25fSSatish Balay { 107149b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 107249b5e25fSSatish Balay MPI_Comm comm = A->comm; 1073dfbe8321SBarry Smith PetscErrorCode ierr; 107449b5e25fSSatish Balay 107549b5e25fSSatish Balay PetscFunctionBegin; 107649b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 107749b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 107849b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 107949b5e25fSSatish Balay PetscFunctionReturn(0); 108049b5e25fSSatish Balay } 108149b5e25fSSatish Balay 108249b5e25fSSatish Balay EXTERN_C_BEGIN 10834a2ae208SSatish Balay #undef __FUNCT__ 10844a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 108513f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 108649b5e25fSSatish Balay { 1087045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 108813f74950SBarry Smith PetscInt i,nz,n; 108949b5e25fSSatish Balay 109049b5e25fSSatish Balay PetscFunctionBegin; 10916c6c5352SBarry Smith nz = baij->maxnz; 1092c464158bSHong Zhang n = mat->n; 109349b5e25fSSatish Balay for (i=0; i<nz; i++) { 109449b5e25fSSatish Balay baij->j[i] = indices[i]; 109549b5e25fSSatish Balay } 10966c6c5352SBarry Smith baij->nz = nz; 109749b5e25fSSatish Balay for (i=0; i<n; i++) { 109849b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 109949b5e25fSSatish Balay } 110049b5e25fSSatish Balay PetscFunctionReturn(0); 110149b5e25fSSatish Balay } 110249b5e25fSSatish Balay EXTERN_C_END 110349b5e25fSSatish Balay 11044a2ae208SSatish Balay #undef __FUNCT__ 11054a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 110649b5e25fSSatish Balay /*@ 110719585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 110849b5e25fSSatish Balay in the matrix. 110949b5e25fSSatish Balay 111049b5e25fSSatish Balay Input Parameters: 111119585528SSatish Balay + mat - the SeqSBAIJ matrix 111249b5e25fSSatish Balay - indices - the column indices 111349b5e25fSSatish Balay 111449b5e25fSSatish Balay Level: advanced 111549b5e25fSSatish Balay 111649b5e25fSSatish Balay Notes: 111749b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 111849b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 111949b5e25fSSatish Balay of the MatSetValues() operation. 112049b5e25fSSatish Balay 112149b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 1122d1be2dadSMatthew Knepley MatCreateSeqSBAIJ(), and the columns indices MUST be sorted. 112349b5e25fSSatish Balay 1124ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 112549b5e25fSSatish Balay 1126ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 112749b5e25fSSatish Balay @*/ 112813f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 112949b5e25fSSatish Balay { 113013f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 113149b5e25fSSatish Balay 113249b5e25fSSatish Balay PetscFunctionBegin; 11334482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 11344482741eSBarry Smith PetscValidPointer(indices,2); 1135c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 113649b5e25fSSatish Balay if (f) { 113749b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 113849b5e25fSSatish Balay } else { 1139e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 114049b5e25fSSatish Balay } 114149b5e25fSSatish Balay PetscFunctionReturn(0); 114249b5e25fSSatish Balay } 114349b5e25fSSatish Balay 11444a2ae208SSatish Balay #undef __FUNCT__ 11454a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1146dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1147273d9f13SBarry Smith { 1148dfbe8321SBarry Smith PetscErrorCode ierr; 1149273d9f13SBarry Smith 1150273d9f13SBarry Smith PetscFunctionBegin; 1151273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1152273d9f13SBarry Smith PetscFunctionReturn(0); 1153273d9f13SBarry Smith } 1154273d9f13SBarry Smith 1155a6ece127SHong Zhang #undef __FUNCT__ 1156a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1157dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1158a6ece127SHong Zhang { 1159a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1160a6ece127SHong Zhang PetscFunctionBegin; 1161a6ece127SHong Zhang *array = a->a; 1162a6ece127SHong Zhang PetscFunctionReturn(0); 1163a6ece127SHong Zhang } 1164a6ece127SHong Zhang 1165a6ece127SHong Zhang #undef __FUNCT__ 1166a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1167dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1168a6ece127SHong Zhang { 1169a6ece127SHong Zhang PetscFunctionBegin; 1170a6ece127SHong Zhang PetscFunctionReturn(0); 1171a6ece127SHong Zhang } 1172a6ece127SHong Zhang 117342ee4b1aSHong Zhang #include "petscblaslapack.h" 117442ee4b1aSHong Zhang #undef __FUNCT__ 117542ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1176dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 117742ee4b1aSHong Zhang { 117842ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1179dfbe8321SBarry Smith PetscErrorCode ierr; 1180521d7252SBarry Smith PetscInt i,bs=Y->bs,bs2,j; 11814ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 118242ee4b1aSHong Zhang 118342ee4b1aSHong Zhang PetscFunctionBegin; 118442ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 118571044d3cSBarry Smith BLASaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one); 1186c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1187c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1188c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1189c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1190c537a176SHong Zhang } 1191c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1192c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1193c4319e64SHong Zhang y->XtoY = X; 1194c537a176SHong Zhang } 1195c4319e64SHong Zhang bs2 = bs*bs; 11966c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1197c4319e64SHong Zhang j = 0; 1198c4319e64SHong Zhang while (j < bs2){ 11996550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1200c4319e64SHong Zhang j++; 1201c537a176SHong Zhang } 1202c4319e64SHong Zhang } 120377431f27SBarry Smith PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %D/%D = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz)); 120442ee4b1aSHong Zhang } else { 120542ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 120642ee4b1aSHong Zhang } 120742ee4b1aSHong Zhang PetscFunctionReturn(0); 120842ee4b1aSHong Zhang } 120942ee4b1aSHong Zhang 1210efcf0fc3SBarry Smith #undef __FUNCT__ 1211efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1212dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1213efcf0fc3SBarry Smith { 1214efcf0fc3SBarry Smith PetscFunctionBegin; 1215efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1216efcf0fc3SBarry Smith PetscFunctionReturn(0); 1217efcf0fc3SBarry Smith } 1218efcf0fc3SBarry Smith 1219efcf0fc3SBarry Smith #undef __FUNCT__ 1220efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1221dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1222efcf0fc3SBarry Smith { 1223efcf0fc3SBarry Smith PetscFunctionBegin; 1224efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1225efcf0fc3SBarry Smith PetscFunctionReturn(0); 1226efcf0fc3SBarry Smith } 1227efcf0fc3SBarry Smith 1228efcf0fc3SBarry Smith #undef __FUNCT__ 1229efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1230dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1231efcf0fc3SBarry Smith { 1232efcf0fc3SBarry Smith PetscFunctionBegin; 1233efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1234efcf0fc3SBarry Smith PetscFunctionReturn(0); 1235efcf0fc3SBarry Smith } 1236efcf0fc3SBarry Smith 123749b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 123849b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 123949b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 124049b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 124149b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 124297304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1243e005ede5SBarry Smith MatMult_SeqSBAIJ_N, 1244e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 124549b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 124649b5e25fSSatish Balay 0, 124749b5e25fSSatish Balay 0, 124897304618SKris Buschelman /*10*/ 0, 124949b5e25fSSatish Balay 0, 12505f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1251d06b337dSHong Zhang MatRelax_SeqSBAIJ, 125249b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 125397304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 125449b5e25fSSatish Balay MatEqual_SeqSBAIJ, 125549b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 125649b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 125749b5e25fSSatish Balay MatNorm_SeqSBAIJ, 125897304618SKris Buschelman /*20*/ 0, 125949b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 126049b5e25fSSatish Balay 0, 126149b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 126249b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1263dcf5cc72SBarry Smith /*25*/ 0, 126449b5e25fSSatish Balay 0, 126549b5e25fSSatish Balay 0, 12665f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12675f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 126897304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1269c464158bSHong Zhang 0, 12704d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1271a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1272a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 127397304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 127449b5e25fSSatish Balay 0, 127549b5e25fSSatish Balay 0, 127649b5e25fSSatish Balay 0, 1277950f1e5bSHong Zhang 0, 127897304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 127949b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 128049b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 128149b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 128249b5e25fSSatish Balay 0, 128397304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 128449b5e25fSSatish Balay MatScale_SeqSBAIJ, 128549b5e25fSSatish Balay 0, 128649b5e25fSSatish Balay 0, 128749b5e25fSSatish Balay 0, 1288521d7252SBarry Smith /*50*/ 0, 128949b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 129049b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 129149b5e25fSSatish Balay 0, 129249b5e25fSSatish Balay 0, 129397304618SKris Buschelman /*55*/ 0, 129449b5e25fSSatish Balay 0, 129549b5e25fSSatish Balay 0, 129649b5e25fSSatish Balay 0, 129749b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 129897304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 129949b5e25fSSatish Balay 0, 130049b5e25fSSatish Balay 0, 13018a124369SBarry Smith MatGetPetscMaps_Petsc, 1302d959ec07SHong Zhang 0, 130397304618SKris Buschelman /*65*/ 0, 1304d959ec07SHong Zhang 0, 1305d959ec07SHong Zhang 0, 1306d959ec07SHong Zhang 0, 1307d959ec07SHong Zhang 0, 130897304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13093e0d88b5SBarry Smith 0, 13103e0d88b5SBarry Smith 0, 13113e0d88b5SBarry Smith 0, 13123e0d88b5SBarry Smith 0, 131397304618SKris Buschelman /*75*/ 0, 13143e0d88b5SBarry Smith 0, 13153e0d88b5SBarry Smith 0, 13163e0d88b5SBarry Smith 0, 13173e0d88b5SBarry Smith 0, 131897304618SKris Buschelman /*80*/ 0, 13193e0d88b5SBarry Smith 0, 13203e0d88b5SBarry Smith 0, 13213e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 132297304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13233e0d88b5SBarry Smith #else 132497304618SKris Buschelman 0, 13253e0d88b5SBarry Smith #endif 1326865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1327865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1328865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1329efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1330865e5f61SKris Buschelman 0, 1331865e5f61SKris Buschelman 0, 1332865e5f61SKris Buschelman /*90*/ 0, 1333865e5f61SKris Buschelman 0, 1334865e5f61SKris Buschelman 0, 1335865e5f61SKris Buschelman 0, 1336865e5f61SKris Buschelman 0, 1337865e5f61SKris Buschelman /*95*/ 0, 1338865e5f61SKris Buschelman 0, 1339865e5f61SKris Buschelman 0, 1340865e5f61SKris Buschelman 0}; 134149b5e25fSSatish Balay 134249b5e25fSSatish Balay EXTERN_C_BEGIN 13434a2ae208SSatish Balay #undef __FUNCT__ 13444a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1345dfbe8321SBarry Smith PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat) 134649b5e25fSSatish Balay { 13474afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1348521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 1349dfbe8321SBarry Smith PetscErrorCode ierr; 135049b5e25fSSatish Balay 135149b5e25fSSatish Balay PetscFunctionBegin; 135249b5e25fSSatish Balay if (aij->nonew != 1) { 1353e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 135449b5e25fSSatish Balay } 135549b5e25fSSatish Balay 135649b5e25fSSatish Balay /* allocate space for values if not already there */ 135749b5e25fSSatish Balay if (!aij->saved_values) { 135887828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 135949b5e25fSSatish Balay } 136049b5e25fSSatish Balay 136149b5e25fSSatish Balay /* copy values over */ 136287828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 136349b5e25fSSatish Balay PetscFunctionReturn(0); 136449b5e25fSSatish Balay } 136549b5e25fSSatish Balay EXTERN_C_END 136649b5e25fSSatish Balay 136749b5e25fSSatish Balay EXTERN_C_BEGIN 13684a2ae208SSatish Balay #undef __FUNCT__ 13694a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1370dfbe8321SBarry Smith PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat) 137149b5e25fSSatish Balay { 13724afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13736849ba73SBarry Smith PetscErrorCode ierr; 1374521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 137549b5e25fSSatish Balay 137649b5e25fSSatish Balay PetscFunctionBegin; 137749b5e25fSSatish Balay if (aij->nonew != 1) { 1378e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 137949b5e25fSSatish Balay } 138049b5e25fSSatish Balay if (!aij->saved_values) { 1381e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 138249b5e25fSSatish Balay } 138349b5e25fSSatish Balay 138449b5e25fSSatish Balay /* copy values over */ 138587828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 138649b5e25fSSatish Balay PetscFunctionReturn(0); 138749b5e25fSSatish Balay } 138849b5e25fSSatish Balay EXTERN_C_END 138949b5e25fSSatish Balay 13908549e402SHong Zhang EXTERN_C_BEGIN 13914a2ae208SSatish Balay #undef __FUNCT__ 1392a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 139313f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 139449b5e25fSSatish Balay { 1395c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13966849ba73SBarry Smith PetscErrorCode ierr; 139713f74950SBarry Smith PetscInt i,len,mbs,bs2; 139849b5e25fSSatish Balay PetscTruth flg; 139949b5e25fSSatish Balay 140049b5e25fSSatish Balay PetscFunctionBegin; 1401273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1402e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1403c464158bSHong Zhang mbs = B->m/bs; 140449b5e25fSSatish Balay bs2 = bs*bs; 140549b5e25fSSatish Balay 1406c464158bSHong Zhang if (mbs*bs != B->m) { 140729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 140849b5e25fSSatish Balay } 140949b5e25fSSatish Balay 1410435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 141177431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 141249b5e25fSSatish Balay if (nnz) { 141349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 141477431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 141577431f27SBarry 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); 141649b5e25fSSatish Balay } 141749b5e25fSSatish Balay } 141849b5e25fSSatish Balay 1419e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 142049b5e25fSSatish Balay if (!flg) { 142149b5e25fSSatish Balay switch (bs) { 142249b5e25fSSatish Balay case 1: 14234ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 142449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1425d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1426e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 142749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 142849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 142949b5e25fSSatish Balay break; 143049b5e25fSSatish Balay case 2: 14314ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 143249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1433e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 143449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 143549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 143649b5e25fSSatish Balay break; 143749b5e25fSSatish Balay case 3: 14385f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 143949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1440e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 144149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 144249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 144349b5e25fSSatish Balay break; 144449b5e25fSSatish Balay case 4: 14455f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 144649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1447e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 144849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 144949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 145049b5e25fSSatish Balay break; 145149b5e25fSSatish Balay case 5: 14525f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 145349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1454e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 145549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 145649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 145749b5e25fSSatish Balay break; 145849b5e25fSSatish Balay case 6: 14595f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 146049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1461e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 146249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 146349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 146449b5e25fSSatish Balay break; 146549b5e25fSSatish Balay case 7: 1466de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 146749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1468e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1469de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 147049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 147149b5e25fSSatish Balay break; 147249b5e25fSSatish Balay } 147349b5e25fSSatish Balay } 147449b5e25fSSatish Balay 147549b5e25fSSatish Balay b->mbs = mbs; 14764afc71dfSHong Zhang b->nbs = mbs; 147713f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->imax);CHKERRQ(ierr); 147849b5e25fSSatish Balay if (!nnz) { 1479435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 148049b5e25fSSatish Balay else if (nz <= 0) nz = 1; 148149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14828cef66ccSHong Zhang b->imax[i] = nz; 148349b5e25fSSatish Balay } 1484153ea458SHong Zhang nz = nz*mbs; /* total nz */ 148549b5e25fSSatish Balay } else { 148649b5e25fSSatish Balay nz = 0; 14878cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 148849b5e25fSSatish Balay } 14896c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 149049b5e25fSSatish Balay 149149b5e25fSSatish Balay /* allocate the matrix space */ 149213f74950SBarry Smith len = nz*sizeof(PetscInt) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(PetscInt); 149382502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 14946c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 149513f74950SBarry Smith b->j = (PetscInt*)(b->a + nz*bs2); 149613f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 14976c6c5352SBarry Smith b->i = b->j + nz; 149849b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 149949b5e25fSSatish Balay 150049b5e25fSSatish Balay /* pointer to beginning of each row */ 150149b5e25fSSatish Balay b->i[0] = 0; 150249b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 150349b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 150449b5e25fSSatish Balay } 150549b5e25fSSatish Balay 150649b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 150713f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->ilen);CHKERRQ(ierr); 150852e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 150949b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 151049b5e25fSSatish Balay 1511521d7252SBarry Smith B->bs = bs; 151249b5e25fSSatish Balay b->bs2 = bs2; 15136c6c5352SBarry Smith b->nz = 0; 15146c6c5352SBarry Smith b->maxnz = nz*bs2; 1515153ea458SHong Zhang 151616cdd363SHong Zhang b->inew = 0; 151716cdd363SHong Zhang b->jnew = 0; 151816cdd363SHong Zhang b->anew = 0; 151916cdd363SHong Zhang b->a2anew = 0; 15201a3463dfSHong Zhang b->permute = PETSC_FALSE; 1521c464158bSHong Zhang PetscFunctionReturn(0); 1522c464158bSHong Zhang } 1523a23d5eceSKris Buschelman EXTERN_C_END 1524153ea458SHong Zhang 1525d769727bSBarry Smith EXTERN_C_BEGIN 1526*ceb03754SKris Buschelman EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,MatReuse,Mat*); 1527*ceb03754SKris Buschelman EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*); 1528d769727bSBarry Smith EXTERN_C_END 1529d769727bSBarry Smith 15300bad9183SKris Buschelman /*MC 1531fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15320bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15330bad9183SKris Buschelman 15340bad9183SKris Buschelman Options Database Keys: 15350bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15360bad9183SKris Buschelman 15370bad9183SKris Buschelman Level: beginner 15380bad9183SKris Buschelman 15390bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15400bad9183SKris Buschelman M*/ 15410bad9183SKris Buschelman 1542a23d5eceSKris Buschelman EXTERN_C_BEGIN 1543a23d5eceSKris Buschelman #undef __FUNCT__ 1544a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1545dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqSBAIJ(Mat B) 1546a23d5eceSKris Buschelman { 1547a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1548dfbe8321SBarry Smith PetscErrorCode ierr; 154913f74950SBarry Smith PetscMPIInt size; 1550a23d5eceSKris Buschelman 1551a23d5eceSKris Buschelman PetscFunctionBegin; 1552a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1553a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1554a23d5eceSKris Buschelman B->m = B->M = PetscMax(B->m,B->M); 1555a23d5eceSKris Buschelman B->n = B->N = PetscMax(B->n,B->N); 1556a23d5eceSKris Buschelman 1557a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1558a23d5eceSKris Buschelman B->data = (void*)b; 1559a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1560a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1561a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1562a23d5eceSKris Buschelman B->factor = 0; 1563a23d5eceSKris Buschelman B->lupivotthreshold = 1.0; 1564a23d5eceSKris Buschelman B->mapping = 0; 1565a23d5eceSKris Buschelman b->row = 0; 1566a23d5eceSKris Buschelman b->icol = 0; 1567a23d5eceSKris Buschelman b->reallocs = 0; 1568a23d5eceSKris Buschelman b->saved_values = 0; 1569a23d5eceSKris Buschelman 1570a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1571a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1572a23d5eceSKris Buschelman 1573a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1574a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1575a23d5eceSKris Buschelman b->nonew = 0; 1576a23d5eceSKris Buschelman b->diag = 0; 1577a23d5eceSKris Buschelman b->solve_work = 0; 1578a23d5eceSKris Buschelman b->mult_work = 0; 1579a23d5eceSKris Buschelman B->spptr = 0; 1580a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1581a23d5eceSKris Buschelman b->xtoy = 0; 1582a23d5eceSKris Buschelman b->XtoY = 0; 1583a23d5eceSKris Buschelman 1584a23d5eceSKris Buschelman b->inew = 0; 1585a23d5eceSKris Buschelman b->jnew = 0; 1586a23d5eceSKris Buschelman b->anew = 0; 1587a23d5eceSKris Buschelman b->a2anew = 0; 1588a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1589a23d5eceSKris Buschelman 1590a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1591a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1592a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1593a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1594a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1595a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1596a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1597a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1598a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 15994e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16004e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16014e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1602a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1603a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1604a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1605a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1606a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1607a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 160823ce1328SBarry Smith 160923ce1328SBarry Smith B->symmetric = PETSC_TRUE; 161023ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 161123ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 161223ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 1613a23d5eceSKris Buschelman PetscFunctionReturn(0); 1614a23d5eceSKris Buschelman } 1615a23d5eceSKris Buschelman EXTERN_C_END 1616a23d5eceSKris Buschelman 1617a23d5eceSKris Buschelman #undef __FUNCT__ 1618a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1619a23d5eceSKris Buschelman /*@C 1620a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1621a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1622a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1623a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1624a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1625a23d5eceSKris Buschelman 1626a23d5eceSKris Buschelman Collective on Mat 1627a23d5eceSKris Buschelman 1628a23d5eceSKris Buschelman Input Parameters: 1629a23d5eceSKris Buschelman + A - the symmetric matrix 1630a23d5eceSKris Buschelman . bs - size of block 1631a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1632a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1633a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1634a23d5eceSKris Buschelman 1635a23d5eceSKris Buschelman Options Database Keys: 1636a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1637a23d5eceSKris Buschelman block calculations (much slower) 1638a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1639a23d5eceSKris Buschelman 1640a23d5eceSKris Buschelman Level: intermediate 1641a23d5eceSKris Buschelman 1642a23d5eceSKris Buschelman Notes: 1643a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1644a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1645a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1646a23d5eceSKris Buschelman matrices. 1647a23d5eceSKris Buschelman 164849a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 164949a6f317SBarry Smith 165049a6f317SBarry Smith 1651a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1652a23d5eceSKris Buschelman @*/ 165313f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 165413f74950SBarry Smith { 165513f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1656a23d5eceSKris Buschelman 1657a23d5eceSKris Buschelman PetscFunctionBegin; 1658a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1659a23d5eceSKris Buschelman if (f) { 1660a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1661a23d5eceSKris Buschelman } 1662a23d5eceSKris Buschelman PetscFunctionReturn(0); 1663a23d5eceSKris Buschelman } 166449b5e25fSSatish Balay 16654a2ae208SSatish Balay #undef __FUNCT__ 16664a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1667c464158bSHong Zhang /*@C 1668c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1669c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1670c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1671c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1672c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 167349b5e25fSSatish Balay 1674c464158bSHong Zhang Collective on MPI_Comm 1675c464158bSHong Zhang 1676c464158bSHong Zhang Input Parameters: 1677c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1678c464158bSHong Zhang . bs - size of block 1679c464158bSHong Zhang . m - number of rows, or number of columns 1680c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1681744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1682744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1683c464158bSHong Zhang 1684c464158bSHong Zhang Output Parameter: 1685c464158bSHong Zhang . A - the symmetric matrix 1686c464158bSHong Zhang 1687c464158bSHong Zhang Options Database Keys: 1688c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1689c464158bSHong Zhang block calculations (much slower) 1690c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1691c464158bSHong Zhang 1692c464158bSHong Zhang Level: intermediate 1693c464158bSHong Zhang 1694c464158bSHong Zhang Notes: 1695c464158bSHong Zhang 1696c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1697c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1698c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1699c464158bSHong Zhang matrices. 1700c464158bSHong Zhang 170149a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 170249a6f317SBarry Smith 1703c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1704c464158bSHong Zhang @*/ 170513f74950SBarry Smith PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1706c464158bSHong Zhang { 1707dfbe8321SBarry Smith PetscErrorCode ierr; 1708c464158bSHong Zhang 1709c464158bSHong Zhang PetscFunctionBegin; 1710c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1711c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1712273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 171349b5e25fSSatish Balay PetscFunctionReturn(0); 171449b5e25fSSatish Balay } 171549b5e25fSSatish Balay 17164a2ae208SSatish Balay #undef __FUNCT__ 17174a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1718dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 171949b5e25fSSatish Balay { 172049b5e25fSSatish Balay Mat C; 172149b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17226849ba73SBarry Smith PetscErrorCode ierr; 172313f74950SBarry Smith PetscInt i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 172449b5e25fSSatish Balay 172549b5e25fSSatish Balay PetscFunctionBegin; 172629bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 172749b5e25fSSatish Balay 172849b5e25fSSatish Balay *B = 0; 1729692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1730be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17311d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1732692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1733692f9cbeSHong Zhang 1734273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 173549b5e25fSSatish Balay C->factor = A->factor; 173649b5e25fSSatish Balay c->row = 0; 173749b5e25fSSatish Balay c->icol = 0; 173849b5e25fSSatish Balay c->saved_values = 0; 173949b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 174049b5e25fSSatish Balay C->assembled = PETSC_TRUE; 174149b5e25fSSatish Balay 174282327fa8SHong Zhang C->M = A->M; 174382327fa8SHong Zhang C->N = A->N; 1744521d7252SBarry Smith C->bs = A->bs; 174549b5e25fSSatish Balay c->bs2 = a->bs2; 174649b5e25fSSatish Balay c->mbs = a->mbs; 174749b5e25fSSatish Balay c->nbs = a->nbs; 174849b5e25fSSatish Balay 174913f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr); 175013f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr); 175149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 175249b5e25fSSatish Balay c->imax[i] = a->imax[i]; 175349b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 175449b5e25fSSatish Balay } 175549b5e25fSSatish Balay 175649b5e25fSSatish Balay /* allocate the matrix space */ 175749b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 175813f74950SBarry Smith len = (mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)); 175982502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 176013f74950SBarry Smith c->j = (PetscInt*)(c->a + nz*bs2); 176149b5e25fSSatish Balay c->i = c->j + nz; 176213f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 176349b5e25fSSatish Balay if (mbs > 0) { 176413f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 176549b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 176649b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 176749b5e25fSSatish Balay } else { 176849b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 176949b5e25fSSatish Balay } 177049b5e25fSSatish Balay } 177149b5e25fSSatish Balay 177252e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 177349b5e25fSSatish Balay c->sorted = a->sorted; 177449b5e25fSSatish Balay c->roworiented = a->roworiented; 177549b5e25fSSatish Balay c->nonew = a->nonew; 177649b5e25fSSatish Balay 177749b5e25fSSatish Balay if (a->diag) { 177813f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 177952e6d16bSBarry Smith ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 178049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 178149b5e25fSSatish Balay c->diag[i] = a->diag[i]; 178249b5e25fSSatish Balay } 178349b5e25fSSatish Balay } else c->diag = 0; 17846c6c5352SBarry Smith c->nz = a->nz; 17856c6c5352SBarry Smith c->maxnz = a->maxnz; 178649b5e25fSSatish Balay c->solve_work = 0; 178749b5e25fSSatish Balay c->mult_work = 0; 178849b5e25fSSatish Balay *B = C; 1789b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 179049b5e25fSSatish Balay PetscFunctionReturn(0); 179149b5e25fSSatish Balay } 179249b5e25fSSatish Balay 17934a2ae208SSatish Balay #undef __FUNCT__ 17944a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1795dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A) 179649b5e25fSSatish Balay { 179749b5e25fSSatish Balay Mat_SeqSBAIJ *a; 179849b5e25fSSatish Balay Mat B; 17996849ba73SBarry Smith PetscErrorCode ierr; 180013f74950SBarry Smith int fd; 180113f74950SBarry Smith PetscMPIInt size; 180213f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 180313f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 180413f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 180513f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 180687828ca2SBarry Smith PetscScalar *aa; 180749b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 180849b5e25fSSatish Balay 180949b5e25fSSatish Balay PetscFunctionBegin; 1810b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 181149b5e25fSSatish Balay bs2 = bs*bs; 181249b5e25fSSatish Balay 181349b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 181429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1815b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 181649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1817552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 181849b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 181949b5e25fSSatish Balay 182049b5e25fSSatish Balay if (header[3] < 0) { 182129bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 182249b5e25fSSatish Balay } 182349b5e25fSSatish Balay 182429bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 182549b5e25fSSatish Balay 182649b5e25fSSatish Balay /* 182749b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 182849b5e25fSSatish Balay divisible by the blocksize 182949b5e25fSSatish Balay */ 183049b5e25fSSatish Balay mbs = M/bs; 183149b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 183249b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 183349b5e25fSSatish Balay else mbs++; 183449b5e25fSSatish Balay if (extra_rows) { 1835b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 183649b5e25fSSatish Balay } 183749b5e25fSSatish Balay 183849b5e25fSSatish Balay /* read in row lengths */ 183913f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 184049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 184149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 184249b5e25fSSatish Balay 184349b5e25fSSatish Balay /* read in column indices */ 184413f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 184549b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 184649b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 184749b5e25fSSatish Balay 184849b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 184913f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 185013f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 185113f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 185213f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 185349b5e25fSSatish Balay masked = mask + mbs; 185449b5e25fSSatish Balay rowcount = 0; nzcount = 0; 185549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 185649b5e25fSSatish Balay nmask = 0; 185749b5e25fSSatish Balay for (j=0; j<bs; j++) { 185849b5e25fSSatish Balay kmax = rowlengths[rowcount]; 185949b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18602d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 186103630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 186249b5e25fSSatish Balay } 186349b5e25fSSatish Balay rowcount++; 186449b5e25fSSatish Balay } 1865574b2666SHong Zhang s_browlengths[i] += nmask; 1866574b2666SHong Zhang 186749b5e25fSSatish Balay /* zero out the mask elements we set */ 186849b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 186949b5e25fSSatish Balay } 187049b5e25fSSatish Balay 187149b5e25fSSatish Balay /* create our matrix */ 18729abb65ffSKris Buschelman ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr); 18739abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 187478473fd7SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr); 187549b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 187649b5e25fSSatish Balay 187749b5e25fSSatish Balay /* set matrix "i" values */ 187849b5e25fSSatish Balay a->i[0] = 0; 187949b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1880574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1881574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 188249b5e25fSSatish Balay } 18836c6c5352SBarry Smith a->nz = a->i[mbs]; 188449b5e25fSSatish Balay 188549b5e25fSSatish Balay /* read in nonzero values */ 188687828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 188749b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 188849b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 188949b5e25fSSatish Balay 189049b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 189149b5e25fSSatish Balay nzcount = 0; jcount = 0; 189249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 189349b5e25fSSatish Balay nzcountb = nzcount; 189449b5e25fSSatish Balay nmask = 0; 189549b5e25fSSatish Balay for (j=0; j<bs; j++) { 189649b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 189749b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18982d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 189903630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19002d703238SHong Zhang } 19012d703238SHong Zhang } 19022d703238SHong Zhang /* sort the masked values */ 19032d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19042d703238SHong Zhang 19052d703238SHong Zhang /* set "j" values into matrix */ 19062d703238SHong Zhang maskcount = 1; 19072d703238SHong Zhang for (j=0; j<nmask; j++) { 190849b5e25fSSatish Balay a->j[jcount++] = masked[j]; 190949b5e25fSSatish Balay mask[masked[j]] = maskcount++; 191049b5e25fSSatish Balay } 1911574b2666SHong Zhang 191249b5e25fSSatish Balay /* set "a" values into matrix */ 191349b5e25fSSatish Balay ishift = bs2*a->i[i]; 191449b5e25fSSatish Balay for (j=0; j<bs; j++) { 191549b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 191649b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1917574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1918574b2666SHong Zhang if (tmp >= i){ 191949b5e25fSSatish Balay block = mask[tmp] - 1; 192049b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 192149b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1922574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1923574b2666SHong Zhang } 1924574b2666SHong Zhang nzcountb++; 192549b5e25fSSatish Balay } 192649b5e25fSSatish Balay } 192749b5e25fSSatish Balay /* zero out the mask elements we set */ 192849b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 192949b5e25fSSatish Balay } 19306c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 193149b5e25fSSatish Balay 193249b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1933574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 193449b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 193549b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 193649b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 193749b5e25fSSatish Balay 19389abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19399abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 194049b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19419abb65ffSKris Buschelman *A = B; 194249b5e25fSSatish Balay PetscFunctionReturn(0); 194349b5e25fSSatish Balay } 1944574b2666SHong Zhang 1945d06b337dSHong Zhang #undef __FUNCT__ 1946d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 194713f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1948d06b337dSHong Zhang { 1949d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1950d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1951d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19526849ba73SBarry Smith PetscErrorCode ierr; 1953521d7252SBarry Smith PetscInt m=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j; 1954521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 1955d06b337dSHong Zhang 1956d06b337dSHong Zhang PetscFunctionBegin; 1957b965ef7fSBarry Smith its = its*lits; 195877431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1959d06b337dSHong Zhang 1960d06b337dSHong Zhang if (bs > 1) 1961d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1962d06b337dSHong Zhang 19631ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1964d06b337dSHong Zhang if (xx != bb) { 19651ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1966d06b337dSHong Zhang } else { 1967d06b337dSHong Zhang b = x; 1968d06b337dSHong Zhang } 1969d06b337dSHong Zhang 1970d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1971d06b337dSHong Zhang 1972d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19732798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1974d06b337dSHong Zhang for (i=0; i<m; i++) 1975d06b337dSHong Zhang t[i] = b[i]; 1976d06b337dSHong Zhang 1977d06b337dSHong Zhang for (i=0; i<m; i++){ 197844706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1979d06b337dSHong Zhang v = aa + ai[i] + 1; 1980d06b337dSHong Zhang vj = aj + ai[i] + 1; 1981d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1982d06b337dSHong Zhang x[i] = omega*t[i]/d; 1983d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 198444706875SHong Zhang PetscLogFlops(2*nz-1); 1985d06b337dSHong Zhang } 1986d06b337dSHong Zhang } 1987d06b337dSHong Zhang 19882798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1989d06b337dSHong Zhang for (i=0; i<m; i++) 1990d06b337dSHong Zhang t[i] = b[i]; 1991d06b337dSHong Zhang 1992d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1993d06b337dSHong Zhang v = aa + ai[i] + 1; 1994d06b337dSHong Zhang vj = aj + ai[i] + 1; 1995d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1996d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 199744706875SHong Zhang PetscLogFlops(2*nz-1); 1998d06b337dSHong Zhang } 1999d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2000d06b337dSHong Zhang d = *(aa + ai[i]); 2001d06b337dSHong Zhang v = aa + ai[i] + 1; 2002d06b337dSHong Zhang vj = aj + ai[i] + 1; 2003d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2004d06b337dSHong Zhang sum = t[i]; 2005d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 200644706875SHong Zhang PetscLogFlops(2*nz-1); 2007d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2008d06b337dSHong Zhang } 2009d06b337dSHong Zhang } 2010d06b337dSHong Zhang its--; 2011d06b337dSHong Zhang } 2012d06b337dSHong Zhang 2013d06b337dSHong Zhang while (its--) { 201444706875SHong Zhang /* 201544706875SHong Zhang forward sweep: 201644706875SHong Zhang for i=0,...,m-1: 201744706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 201844706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 201944706875SHong Zhang b = b - x[i]*U^T(i,:); 2020d06b337dSHong Zhang 202144706875SHong Zhang */ 20222798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2023d06b337dSHong Zhang for (i=0; i<m; i++) 2024d06b337dSHong Zhang t[i] = b[i]; 2025d06b337dSHong Zhang 2026d06b337dSHong Zhang for (i=0; i<m; i++){ 202744706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2028d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2029d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2030d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2031d06b337dSHong Zhang sum = t[i]; 2032d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2033d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2034d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 203544706875SHong Zhang PetscLogFlops(4*nz-2); 2036d06b337dSHong Zhang } 2037d06b337dSHong Zhang } 2038d06b337dSHong Zhang 20392798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 204044706875SHong Zhang /* 204144706875SHong Zhang backward sweep: 204244706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 204344706875SHong Zhang for i=m-1,...,0: 204444706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 204544706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 204644706875SHong Zhang */ 2047d06b337dSHong Zhang for (i=0; i<m; i++) 2048d06b337dSHong Zhang t[i] = b[i]; 2049d06b337dSHong Zhang 2050d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2051d06b337dSHong Zhang v = aa + ai[i] + 1; 2052d06b337dSHong Zhang vj = aj + ai[i] + 1; 2053d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2054d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 205544706875SHong Zhang PetscLogFlops(2*nz-1); 2056d06b337dSHong Zhang } 2057d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2058d06b337dSHong Zhang d = *(aa + ai[i]); 2059d06b337dSHong Zhang v = aa + ai[i] + 1; 2060d06b337dSHong Zhang vj = aj + ai[i] + 1; 2061d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2062d06b337dSHong Zhang sum = t[i]; 2063d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 206444706875SHong Zhang PetscLogFlops(2*nz-1); 2065d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2066d06b337dSHong Zhang } 2067d06b337dSHong Zhang } 2068d06b337dSHong Zhang } 2069d06b337dSHong Zhang 2070d06b337dSHong Zhang ierr = PetscFree(t);CHKERRQ(ierr); 20711ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2072d06b337dSHong Zhang if (bb != xx) { 20731ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2074d06b337dSHong Zhang } 2075d06b337dSHong Zhang PetscFunctionReturn(0); 2076d06b337dSHong Zhang } 2077d06b337dSHong Zhang 2078d06b337dSHong Zhang 2079d06b337dSHong Zhang 2080d06b337dSHong Zhang 208149b5e25fSSatish Balay 208249b5e25fSSatish Balay 2083