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); 4413f74950SBarry Smith PetscLogObjectMemory(A,(mbs+1)*sizeof(PetscInt)); 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 26249b5e25fSSatish Balay PetscFunctionReturn(0); 26349b5e25fSSatish Balay } 26449b5e25fSSatish Balay 2654a2ae208SSatish Balay #undef __FUNCT__ 2664a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 26713f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 26849b5e25fSSatish Balay { 269dfbe8321SBarry Smith PetscErrorCode ierr; 27049b5e25fSSatish Balay 27149b5e25fSSatish Balay PetscFunctionBegin; 27249b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 27349b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 27449b5e25fSSatish Balay PetscFunctionReturn(0); 27549b5e25fSSatish Balay } 27649b5e25fSSatish Balay 2774a2ae208SSatish Balay #undef __FUNCT__ 2784a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 279dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B) 28049b5e25fSSatish Balay { 281dfbe8321SBarry Smith PetscErrorCode ierr; 28249b5e25fSSatish Balay PetscFunctionBegin; 283999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2848115998fSBarry Smith PetscFunctionReturn(0); 28549b5e25fSSatish Balay } 28649b5e25fSSatish Balay 2874a2ae208SSatish Balay #undef __FUNCT__ 2884a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 2896849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 29049b5e25fSSatish Balay { 29149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 292dfbe8321SBarry Smith PetscErrorCode ierr; 293521d7252SBarry Smith PetscInt i,j,bs = A->bs,k,l,bs2=a->bs2; 294fb9695e5SSatish Balay char *name; 295f3ef73ceSBarry Smith PetscViewerFormat format; 29649b5e25fSSatish Balay 29749b5e25fSSatish Balay PetscFunctionBegin; 29880fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 299b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 300456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 30177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %D\n",bs);CHKERRQ(ierr); 302fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 30329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 304fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 305b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 30649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30749b5e25fSSatish Balay for (j=0; j<bs; j++) { 30877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 30949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31049b5e25fSSatish Balay for (l=0; l<bs; l++) { 31149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 31249b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 31449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 31549b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 31749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 31849b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32049b5e25fSSatish Balay } 32149b5e25fSSatish Balay #else 32249b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 32377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 32449b5e25fSSatish Balay } 32549b5e25fSSatish Balay #endif 32649b5e25fSSatish Balay } 32749b5e25fSSatish Balay } 328b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 32949b5e25fSSatish Balay } 33049b5e25fSSatish Balay } 331b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 33249b5e25fSSatish Balay } else { 333b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 33449b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33549b5e25fSSatish Balay for (j=0; j<bs; j++) { 33677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr); 33749b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33849b5e25fSSatish Balay for (l=0; l<bs; l++) { 33949b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 34049b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 34177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l, 34249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34349b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 34477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l, 34549b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34649b5e25fSSatish Balay } else { 34777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34849b5e25fSSatish Balay } 34949b5e25fSSatish Balay #else 35077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35149b5e25fSSatish Balay #endif 35249b5e25fSSatish Balay } 35349b5e25fSSatish Balay } 354b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 35549b5e25fSSatish Balay } 35649b5e25fSSatish Balay } 357b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 35849b5e25fSSatish Balay } 359b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 36049b5e25fSSatish Balay PetscFunctionReturn(0); 36149b5e25fSSatish Balay } 36249b5e25fSSatish Balay 3634a2ae208SSatish Balay #undef __FUNCT__ 3644a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 3656849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 36649b5e25fSSatish Balay { 36749b5e25fSSatish Balay Mat A = (Mat) Aa; 36849b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 3696849ba73SBarry Smith PetscErrorCode ierr; 37013f74950SBarry Smith PetscInt row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2; 37113f74950SBarry Smith PetscMPIInt rank; 37249b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 37349b5e25fSSatish Balay MatScalar *aa; 37449b5e25fSSatish Balay MPI_Comm comm; 375b0a32e0cSBarry Smith PetscViewer viewer; 37649b5e25fSSatish Balay 37749b5e25fSSatish Balay PetscFunctionBegin; 37849b5e25fSSatish Balay /* 37949b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 38049b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 38149b5e25fSSatish Balay rest should return immediately. 38249b5e25fSSatish Balay */ 38349b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 38449b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 38549b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 38649b5e25fSSatish Balay 38749b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 38849b5e25fSSatish Balay 389b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 390b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 39149b5e25fSSatish Balay 39249b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 393b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 39449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 39549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 396c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 39749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 39849b5e25fSSatish Balay aa = a->a + j*bs2; 39949b5e25fSSatish Balay for (k=0; k<bs; k++) { 40049b5e25fSSatish Balay for (l=0; l<bs; l++) { 40149b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 402b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 40349b5e25fSSatish Balay } 40449b5e25fSSatish Balay } 40549b5e25fSSatish Balay } 40649b5e25fSSatish Balay } 407b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 40849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 40949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 410c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 41149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 41249b5e25fSSatish Balay aa = a->a + j*bs2; 41349b5e25fSSatish Balay for (k=0; k<bs; k++) { 41449b5e25fSSatish Balay for (l=0; l<bs; l++) { 41549b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 416b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 41749b5e25fSSatish Balay } 41849b5e25fSSatish Balay } 41949b5e25fSSatish Balay } 42049b5e25fSSatish Balay } 42149b5e25fSSatish Balay 422b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 42349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 42449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 425c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 42649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42749b5e25fSSatish Balay aa = a->a + j*bs2; 42849b5e25fSSatish Balay for (k=0; k<bs; k++) { 42949b5e25fSSatish Balay for (l=0; l<bs; l++) { 43049b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 431b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 43249b5e25fSSatish Balay } 43349b5e25fSSatish Balay } 43449b5e25fSSatish Balay } 43549b5e25fSSatish Balay } 43649b5e25fSSatish Balay PetscFunctionReturn(0); 43749b5e25fSSatish Balay } 43849b5e25fSSatish Balay 4394a2ae208SSatish Balay #undef __FUNCT__ 4404a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 4416849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 44249b5e25fSSatish Balay { 443dfbe8321SBarry Smith PetscErrorCode ierr; 44449b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 445b0a32e0cSBarry Smith PetscDraw draw; 44649b5e25fSSatish Balay PetscTruth isnull; 44749b5e25fSSatish Balay 44849b5e25fSSatish Balay PetscFunctionBegin; 44949b5e25fSSatish Balay 450b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 451b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 45249b5e25fSSatish Balay 45349b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 454c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 45549b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 456b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 457b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 45849b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 45949b5e25fSSatish Balay PetscFunctionReturn(0); 46049b5e25fSSatish Balay } 46149b5e25fSSatish Balay 4624a2ae208SSatish Balay #undef __FUNCT__ 4634a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 464dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 46549b5e25fSSatish Balay { 466dfbe8321SBarry Smith PetscErrorCode ierr; 46732077d6dSBarry Smith PetscTruth iascii,isdraw; 46849b5e25fSSatish Balay 46949b5e25fSSatish Balay PetscFunctionBegin; 47032077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 471fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 47232077d6dSBarry Smith if (iascii){ 47349b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 47449b5e25fSSatish Balay } else if (isdraw) { 47549b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 47649b5e25fSSatish Balay } else { 477a5e6ed63SBarry Smith Mat B; 478a5e6ed63SBarry Smith ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr); 479a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 480a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 48149b5e25fSSatish Balay } 48249b5e25fSSatish Balay PetscFunctionReturn(0); 48349b5e25fSSatish Balay } 48449b5e25fSSatish Balay 48549b5e25fSSatish Balay 4864a2ae208SSatish Balay #undef __FUNCT__ 4874a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 48813f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) 48949b5e25fSSatish Balay { 490045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 49113f74950SBarry Smith PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 49213f74950SBarry Smith PetscInt *ai = a->i,*ailen = a->ilen; 49313f74950SBarry Smith PetscInt brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2; 49449b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 49549b5e25fSSatish Balay 49649b5e25fSSatish Balay PetscFunctionBegin; 49749b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 49849b5e25fSSatish Balay row = im[k]; brow = row/bs; 49977431f27SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); 50077431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 50149b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 50249b5e25fSSatish Balay nrow = ailen[brow]; 50349b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 50477431f27SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); 50577431f27SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1); 50649b5e25fSSatish Balay col = in[l] ; 50749b5e25fSSatish Balay bcol = col/bs; 50849b5e25fSSatish Balay cidx = col%bs; 50949b5e25fSSatish Balay ridx = row%bs; 51049b5e25fSSatish Balay high = nrow; 51149b5e25fSSatish Balay low = 0; /* assume unsorted */ 51249b5e25fSSatish Balay while (high-low > 5) { 51349b5e25fSSatish Balay t = (low+high)/2; 51449b5e25fSSatish Balay if (rp[t] > bcol) high = t; 51549b5e25fSSatish Balay else low = t; 51649b5e25fSSatish Balay } 51749b5e25fSSatish Balay for (i=low; i<high; i++) { 51849b5e25fSSatish Balay if (rp[i] > bcol) break; 51949b5e25fSSatish Balay if (rp[i] == bcol) { 52049b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 52149b5e25fSSatish Balay goto finished; 52249b5e25fSSatish Balay } 52349b5e25fSSatish Balay } 52449b5e25fSSatish Balay *v++ = zero; 52549b5e25fSSatish Balay finished:; 52649b5e25fSSatish Balay } 52749b5e25fSSatish Balay } 52849b5e25fSSatish Balay PetscFunctionReturn(0); 52949b5e25fSSatish Balay } 53049b5e25fSSatish Balay 53149b5e25fSSatish Balay 5324a2ae208SSatish Balay #undef __FUNCT__ 5334a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 53413f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) 53549b5e25fSSatish Balay { 5360880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5376849ba73SBarry Smith PetscErrorCode ierr; 53813f74950SBarry Smith PetscInt *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 53913f74950SBarry Smith PetscInt *imax=a->imax,*ai=a->i,*ailen=a->ilen; 54013f74950SBarry Smith PetscInt *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval; 5410880e062SHong Zhang PetscTruth roworiented=a->roworiented; 542f15d580aSBarry Smith const MatScalar *value = v; 543f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5440880e062SHong Zhang 54549b5e25fSSatish Balay PetscFunctionBegin; 5460880e062SHong Zhang if (roworiented) { 5470880e062SHong Zhang stepval = (n-1)*bs; 5480880e062SHong Zhang } else { 5490880e062SHong Zhang stepval = (m-1)*bs; 5500880e062SHong Zhang } 5510880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5520880e062SHong Zhang row = im[k]; 5530880e062SHong Zhang if (row < 0) continue; 5540880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 55577431f27SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1); 5560880e062SHong Zhang #endif 5570880e062SHong Zhang rp = aj + ai[row]; 5580880e062SHong Zhang ap = aa + bs2*ai[row]; 5590880e062SHong Zhang rmax = imax[row]; 5600880e062SHong Zhang nrow = ailen[row]; 5610880e062SHong Zhang low = 0; 5620880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5630880e062SHong Zhang if (in[l] < 0) continue; 5640880e062SHong Zhang col = in[l]; 565b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g) 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 } 5740880e062SHong Zhang if (!sorted) low = 0; high = nrow; 5750880e062SHong Zhang while (high-low > 7) { 5760880e062SHong Zhang t = (low+high)/2; 5770880e062SHong Zhang if (rp[t] > col) high = t; 5780880e062SHong Zhang else low = t; 5790880e062SHong Zhang } 5800880e062SHong Zhang for (i=low; i<high; i++) { 5810880e062SHong Zhang if (rp[i] > col) break; 5820880e062SHong Zhang if (rp[i] == col) { 5830880e062SHong Zhang bap = ap + bs2*i; 5840880e062SHong Zhang if (roworiented) { 5850880e062SHong Zhang if (is == ADD_VALUES) { 5860880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5870880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5880880e062SHong Zhang bap[jj] += *value++; 5890880e062SHong Zhang } 5900880e062SHong Zhang } 5910880e062SHong Zhang } else { 5920880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5930880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5940880e062SHong Zhang bap[jj] = *value++; 5950880e062SHong Zhang } 5960880e062SHong Zhang } 5970880e062SHong Zhang } 5980880e062SHong Zhang } else { 5990880e062SHong Zhang if (is == ADD_VALUES) { 6000880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6010880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6020880e062SHong Zhang *bap++ += *value++; 6030880e062SHong Zhang } 6040880e062SHong Zhang } 6050880e062SHong Zhang } else { 6060880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6070880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6080880e062SHong Zhang *bap++ = *value++; 6090880e062SHong Zhang } 6100880e062SHong Zhang } 6110880e062SHong Zhang } 6120880e062SHong Zhang } 6130880e062SHong Zhang goto noinsert2; 6140880e062SHong Zhang } 6150880e062SHong Zhang } 6160880e062SHong Zhang if (nonew == 1) goto noinsert2; 61777431f27SBarry Smith else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 6180880e062SHong Zhang if (nrow >= rmax) { 6190880e062SHong Zhang /* there is no extra room in row, therefore enlarge */ 62013f74950SBarry Smith PetscInt new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6210880e062SHong Zhang MatScalar *new_a; 6220880e062SHong Zhang 62377431f27SBarry Smith if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col); 6240880e062SHong Zhang 6250880e062SHong Zhang /* malloc new storage space */ 62613f74950SBarry Smith len = new_nz*(sizeof(PetscInt)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(PetscInt); 6270880e062SHong Zhang ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 62813f74950SBarry Smith new_j = (PetscInt*)(new_a + bs2*new_nz); 6290880e062SHong Zhang new_i = new_j + new_nz; 6300880e062SHong Zhang 6310880e062SHong Zhang /* copy over old data into new slots */ 6320880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6330880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 63413f74950SBarry Smith ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(PetscInt));CHKERRQ(ierr); 6350880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 63613f74950SBarry Smith ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(PetscInt));CHKERRQ(ierr); 6370880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 6380880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 6390880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 6400880e062SHong Zhang /* free up old matrix storage */ 6410880e062SHong Zhang ierr = PetscFree(a->a);CHKERRQ(ierr); 6420880e062SHong Zhang if (!a->singlemalloc) { 6430880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 6440880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 6450880e062SHong Zhang } 6460880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 6470880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 6480880e062SHong Zhang 6490880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 6500880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 65113f74950SBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar))); 6526c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 6530880e062SHong Zhang a->reallocs++; 6546c6c5352SBarry Smith a->nz++; 6550880e062SHong Zhang } 6560880e062SHong Zhang N = nrow++ - 1; 6570880e062SHong Zhang /* shift up all the later entries in this row */ 6580880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6590880e062SHong Zhang rp[ii+1] = rp[ii]; 6600880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6610880e062SHong Zhang } 6620880e062SHong Zhang if (N >= i) { 6630880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6640880e062SHong Zhang } 6650880e062SHong Zhang rp[i] = col; 6660880e062SHong Zhang bap = ap + bs2*i; 6670880e062SHong Zhang if (roworiented) { 6680880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6690880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6700880e062SHong Zhang bap[jj] = *value++; 6710880e062SHong Zhang } 6720880e062SHong Zhang } 6730880e062SHong Zhang } else { 6740880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6750880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6760880e062SHong Zhang *bap++ = *value++; 6770880e062SHong Zhang } 6780880e062SHong Zhang } 6790880e062SHong Zhang } 6800880e062SHong Zhang noinsert2:; 6810880e062SHong Zhang low = i; 6820880e062SHong Zhang } 6830880e062SHong Zhang ailen[row] = nrow; 6840880e062SHong Zhang } 6850880e062SHong Zhang PetscFunctionReturn(0); 68649b5e25fSSatish Balay } 68749b5e25fSSatish Balay 6884a2ae208SSatish Balay #undef __FUNCT__ 6894a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 690dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 69149b5e25fSSatish Balay { 69249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6936849ba73SBarry Smith PetscErrorCode ierr; 69413f74950SBarry Smith PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 69513f74950SBarry Smith PetscInt m = A->m,*ip,N,*ailen = a->ilen; 69613f74950SBarry Smith PetscInt mbs = a->mbs,bs2 = a->bs2,rmax = 0; 69749b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 69849b5e25fSSatish Balay 69949b5e25fSSatish Balay PetscFunctionBegin; 70049b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 70149b5e25fSSatish Balay 70249b5e25fSSatish Balay if (m) rmax = ailen[0]; 70349b5e25fSSatish Balay for (i=1; i<mbs; i++) { 70449b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 70549b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 70649b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 70749b5e25fSSatish Balay if (fshift) { 70849b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 70949b5e25fSSatish Balay N = ailen[i]; 71049b5e25fSSatish Balay for (j=0; j<N; j++) { 71149b5e25fSSatish Balay ip[j-fshift] = ip[j]; 71249b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 71349b5e25fSSatish Balay } 71449b5e25fSSatish Balay } 71549b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 71649b5e25fSSatish Balay } 71749b5e25fSSatish Balay if (mbs) { 71849b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 71949b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 72049b5e25fSSatish Balay } 72149b5e25fSSatish Balay /* reset ilen and imax for each row */ 72249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 72349b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 72449b5e25fSSatish Balay } 7256c6c5352SBarry Smith a->nz = ai[mbs]; 72649b5e25fSSatish Balay 727b424e231SHong Zhang /* diagonals may have moved, reset it */ 728b424e231SHong Zhang if (a->diag) { 72913f74950SBarry Smith ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 73049b5e25fSSatish Balay } 73177431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n", 732521d7252SBarry Smith m,A->m,A->bs,fshift*bs2,a->nz*bs2); 73377431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %D\n", 73449b5e25fSSatish Balay a->reallocs); 73577431f27SBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %D\n",rmax); 73649b5e25fSSatish Balay a->reallocs = 0; 73749b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 73849b5e25fSSatish Balay 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 } 77249b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* 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; 79613f74950SBarry Smith PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 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 80449b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 80549b5e25fSSatish Balay row = im[k]; /* row number */ 80649b5e25fSSatish Balay brow = row/bs; /* block row number */ 80749b5e25fSSatish Balay if (row < 0) continue; 80849b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 80977431f27SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1); 81049b5e25fSSatish Balay #endif 81149b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 81249b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 81349b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 81449b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 81549b5e25fSSatish Balay low = 0; 81649b5e25fSSatish Balay 81749b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 81849b5e25fSSatish Balay if (in[l] < 0) continue; 81949b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 82077431f27SBarry Smith if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1); 82149b5e25fSSatish Balay #endif 82249b5e25fSSatish Balay col = in[l]; 82349b5e25fSSatish Balay bcol = col/bs; /* block col number */ 82449b5e25fSSatish Balay 82549b5e25fSSatish Balay if (brow <= bcol){ 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 */ 83649b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 83749b5e25fSSatish Balay while (high-low > 7) { 83849b5e25fSSatish Balay t = (low+high)/2; 83949b5e25fSSatish Balay if (rp[t] > bcol) high = t; 84049b5e25fSSatish Balay else low = t; 84149b5e25fSSatish Balay } 84249b5e25fSSatish Balay for (i=low; i<high; i++) { 84377431f27SBarry Smith /* printf("The loop of i=low.., rp[%D]=%D\n",i,rp[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; 89413f74950SBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(PetscInt) + bs2*sizeof(MatScalar))); 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; 91349b5e25fSSatish Balay /* } */ 9148549e402SHong Zhang } 91549b5e25fSSatish Balay } /* end of if .. if.. */ 91649b5e25fSSatish Balay } /* end of loop over added columns */ 91749b5e25fSSatish Balay ailen[brow] = nrow; 91849b5e25fSSatish Balay } /* end of loop over added rows */ 91949b5e25fSSatish Balay 92049b5e25fSSatish Balay PetscFunctionReturn(0); 92149b5e25fSSatish Balay } 92249b5e25fSSatish Balay 9236849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 9246849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 92513f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt); 92613f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*); 92713f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]); 9286849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat); 9296849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 9306849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 9316849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec); 9326849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 9336849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 9346849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat); 9356849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec); 93649b5e25fSSatish Balay 9376849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 9386849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 9396849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 9406849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 9416849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 9426849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 9436849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 9446849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 94549b5e25fSSatish Balay 9466849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 947d59c15a7SBarry Smith 9486849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 9496849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 9506849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 9516849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 9526849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 9536849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 9546849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 9556849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 95607f98182SSatish Balay 9576849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9586849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9596849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9606849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9616849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9626849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9636849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 9646849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 96513f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*); 96649b5e25fSSatish Balay 9676849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 9686849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 9696849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 9706849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 9716849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 9726849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 9736849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 9746849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 97549b5e25fSSatish Balay 9766849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 9776849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 9786849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 9796849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 9806849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 9816849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 9826849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 9836849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 98449b5e25fSSatish Balay 9854d101231SSatish Balay #ifdef HAVE_MatICCFactor 9866849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work! */ 9874a2ae208SSatish Balay #undef __FUNCT__ 9884d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 989dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 99049b5e25fSSatish Balay { 9914ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 99249b5e25fSSatish Balay Mat outA; 993dfbe8321SBarry Smith PetscErrorCode ierr; 99449b5e25fSSatish Balay PetscTruth row_identity,col_identity; 99549b5e25fSSatish Balay 99649b5e25fSSatish Balay PetscFunctionBegin; 99749b5e25fSSatish Balay outA = inA; 9981260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 99949b5e25fSSatish Balay 100049b5e25fSSatish Balay if (!a->diag) { 10011a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 100249b5e25fSSatish Balay } 100349b5e25fSSatish Balay /* 100449b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 100549b5e25fSSatish Balay for ILU(0) factorization with natural ordering 100649b5e25fSSatish Balay */ 100749b5e25fSSatish Balay switch (a->bs) { 100849b5e25fSSatish Balay case 1: 10090fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 101007f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1011d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 10124d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 101349b5e25fSSatish Balay case 2: 10141a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10151a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10161a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10174d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 101849b5e25fSSatish Balay break; 101949b5e25fSSatish Balay case 3: 10201a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10211a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10221a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10234d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 102449b5e25fSSatish Balay break; 102549b5e25fSSatish Balay case 4: 10261a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10271a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10281a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10294d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 103049b5e25fSSatish Balay break; 103149b5e25fSSatish Balay case 5: 10321a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10331a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10341a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10354d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 103649b5e25fSSatish Balay break; 103749b5e25fSSatish Balay case 6: 10381a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10391a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10401a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10414d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 104249b5e25fSSatish Balay break; 104349b5e25fSSatish Balay case 7: 10441a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10451a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10461a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10474d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 104849b5e25fSSatish Balay break; 104949b5e25fSSatish Balay default: 105049b5e25fSSatish Balay a->row = row; 10511a3463dfSHong Zhang a->icol = col; 105249b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 105349b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 105449b5e25fSSatish Balay 105549b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 105649b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1057b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 105849b5e25fSSatish Balay 105949b5e25fSSatish Balay if (!a->solve_work) { 106087828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 106187828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 106249b5e25fSSatish Balay } 106349b5e25fSSatish Balay } 106449b5e25fSSatish Balay 10651a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 106649b5e25fSSatish Balay 106749b5e25fSSatish Balay PetscFunctionReturn(0); 106849b5e25fSSatish Balay } 1069950f1e5bSHong Zhang #endif 1070950f1e5bSHong Zhang 10714a2ae208SSatish Balay #undef __FUNCT__ 10724a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 1073dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A) 107449b5e25fSSatish Balay { 107549b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 107649b5e25fSSatish Balay MPI_Comm comm = A->comm; 1077dfbe8321SBarry Smith PetscErrorCode ierr; 107849b5e25fSSatish Balay 107949b5e25fSSatish Balay PetscFunctionBegin; 108049b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 108149b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 108249b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 108349b5e25fSSatish Balay PetscFunctionReturn(0); 108449b5e25fSSatish Balay } 108549b5e25fSSatish Balay 108649b5e25fSSatish Balay EXTERN_C_BEGIN 10874a2ae208SSatish Balay #undef __FUNCT__ 10884a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 108913f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices) 109049b5e25fSSatish Balay { 1091045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 109213f74950SBarry Smith PetscInt i,nz,n; 109349b5e25fSSatish Balay 109449b5e25fSSatish Balay PetscFunctionBegin; 10956c6c5352SBarry Smith nz = baij->maxnz; 1096c464158bSHong Zhang n = mat->n; 109749b5e25fSSatish Balay for (i=0; i<nz; i++) { 109849b5e25fSSatish Balay baij->j[i] = indices[i]; 109949b5e25fSSatish Balay } 11006c6c5352SBarry Smith baij->nz = nz; 110149b5e25fSSatish Balay for (i=0; i<n; i++) { 110249b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 110349b5e25fSSatish Balay } 110449b5e25fSSatish Balay 110549b5e25fSSatish Balay PetscFunctionReturn(0); 110649b5e25fSSatish Balay } 110749b5e25fSSatish Balay EXTERN_C_END 110849b5e25fSSatish Balay 11094a2ae208SSatish Balay #undef __FUNCT__ 11104a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 111149b5e25fSSatish Balay /*@ 111219585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 111349b5e25fSSatish Balay in the matrix. 111449b5e25fSSatish Balay 111549b5e25fSSatish Balay Input Parameters: 111619585528SSatish Balay + mat - the SeqSBAIJ matrix 111749b5e25fSSatish Balay - indices - the column indices 111849b5e25fSSatish Balay 111949b5e25fSSatish Balay Level: advanced 112049b5e25fSSatish Balay 112149b5e25fSSatish Balay Notes: 112249b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 112349b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 112449b5e25fSSatish Balay of the MatSetValues() operation. 112549b5e25fSSatish Balay 112649b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 112719585528SSatish Balay MatCreateSeqSBAIJ(). 112849b5e25fSSatish Balay 1129ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 113049b5e25fSSatish Balay 1131ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 113249b5e25fSSatish Balay @*/ 113313f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices) 113449b5e25fSSatish Balay { 113513f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt *); 113649b5e25fSSatish Balay 113749b5e25fSSatish Balay PetscFunctionBegin; 11384482741eSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 11394482741eSBarry Smith PetscValidPointer(indices,2); 1140c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 114149b5e25fSSatish Balay if (f) { 114249b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 114349b5e25fSSatish Balay } else { 1144e005ede5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); 114549b5e25fSSatish Balay } 114649b5e25fSSatish Balay PetscFunctionReturn(0); 114749b5e25fSSatish Balay } 114849b5e25fSSatish Balay 11494a2ae208SSatish Balay #undef __FUNCT__ 11504a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1151dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A) 1152273d9f13SBarry Smith { 1153dfbe8321SBarry Smith PetscErrorCode ierr; 1154273d9f13SBarry Smith 1155273d9f13SBarry Smith PetscFunctionBegin; 1156273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1157273d9f13SBarry Smith PetscFunctionReturn(0); 1158273d9f13SBarry Smith } 1159273d9f13SBarry Smith 1160a6ece127SHong Zhang #undef __FUNCT__ 1161a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1162dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1163a6ece127SHong Zhang { 1164a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1165a6ece127SHong Zhang PetscFunctionBegin; 1166a6ece127SHong Zhang *array = a->a; 1167a6ece127SHong Zhang PetscFunctionReturn(0); 1168a6ece127SHong Zhang } 1169a6ece127SHong Zhang 1170a6ece127SHong Zhang #undef __FUNCT__ 1171a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1172dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1173a6ece127SHong Zhang { 1174a6ece127SHong Zhang PetscFunctionBegin; 1175a6ece127SHong Zhang PetscFunctionReturn(0); 1176a6ece127SHong Zhang } 1177a6ece127SHong Zhang 117842ee4b1aSHong Zhang #include "petscblaslapack.h" 117942ee4b1aSHong Zhang #undef __FUNCT__ 118042ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1181dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 118242ee4b1aSHong Zhang { 118342ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1184dfbe8321SBarry Smith PetscErrorCode ierr; 1185521d7252SBarry Smith PetscInt i,bs=Y->bs,bs2,j; 11864ce68768SBarry Smith PetscBLASInt bnz = (PetscBLASInt)x->nz,one = 1; 118742ee4b1aSHong Zhang 118842ee4b1aSHong Zhang PetscFunctionBegin; 118942ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 11904ce68768SBarry Smith BLaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one); 1191c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1192c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1193c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1194c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1195c537a176SHong Zhang } 1196c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1197c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1198c4319e64SHong Zhang y->XtoY = X; 1199c537a176SHong Zhang } 1200c4319e64SHong Zhang bs2 = bs*bs; 12016c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1202c4319e64SHong Zhang j = 0; 1203c4319e64SHong Zhang while (j < bs2){ 12046550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1205c4319e64SHong Zhang j++; 1206c537a176SHong Zhang } 1207c4319e64SHong Zhang } 120877431f27SBarry 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)); 120942ee4b1aSHong Zhang } else { 121042ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 121142ee4b1aSHong Zhang } 121242ee4b1aSHong Zhang PetscFunctionReturn(0); 121342ee4b1aSHong Zhang } 121442ee4b1aSHong Zhang 1215efcf0fc3SBarry Smith #undef __FUNCT__ 1216efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1217dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg) 1218efcf0fc3SBarry Smith { 1219efcf0fc3SBarry Smith PetscFunctionBegin; 1220efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1221efcf0fc3SBarry Smith PetscFunctionReturn(0); 1222efcf0fc3SBarry Smith } 1223efcf0fc3SBarry Smith 1224efcf0fc3SBarry Smith #undef __FUNCT__ 1225efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1226dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1227efcf0fc3SBarry Smith { 1228efcf0fc3SBarry Smith PetscFunctionBegin; 1229efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1230efcf0fc3SBarry Smith PetscFunctionReturn(0); 1231efcf0fc3SBarry Smith } 1232efcf0fc3SBarry Smith 1233efcf0fc3SBarry Smith #undef __FUNCT__ 1234efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1235dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1236efcf0fc3SBarry Smith { 1237efcf0fc3SBarry Smith PetscFunctionBegin; 1238efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1239efcf0fc3SBarry Smith PetscFunctionReturn(0); 1240efcf0fc3SBarry Smith } 1241efcf0fc3SBarry Smith 124249b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 124349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 124449b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 124549b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 124649b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 124797304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 1248e005ede5SBarry Smith MatMult_SeqSBAIJ_N, 1249e005ede5SBarry Smith MatMultAdd_SeqSBAIJ_N, 125049b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 125149b5e25fSSatish Balay 0, 125249b5e25fSSatish Balay 0, 125397304618SKris Buschelman /*10*/ 0, 125449b5e25fSSatish Balay 0, 12555f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1256d06b337dSHong Zhang MatRelax_SeqSBAIJ, 125749b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 125897304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 125949b5e25fSSatish Balay MatEqual_SeqSBAIJ, 126049b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 126149b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 126249b5e25fSSatish Balay MatNorm_SeqSBAIJ, 126397304618SKris Buschelman /*20*/ 0, 126449b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 126549b5e25fSSatish Balay 0, 126649b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 126749b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 1268*dcf5cc72SBarry Smith /*25*/ 0, 126949b5e25fSSatish Balay 0, 127049b5e25fSSatish Balay 0, 12715f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12725f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 127397304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1274c464158bSHong Zhang 0, 12754d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1276a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1277a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 127897304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 127949b5e25fSSatish Balay 0, 128049b5e25fSSatish Balay 0, 128149b5e25fSSatish Balay 0, 1282950f1e5bSHong Zhang 0, 128397304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 128449b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 128549b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 128649b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 128749b5e25fSSatish Balay 0, 128897304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 128949b5e25fSSatish Balay MatScale_SeqSBAIJ, 129049b5e25fSSatish Balay 0, 129149b5e25fSSatish Balay 0, 129249b5e25fSSatish Balay 0, 1293521d7252SBarry Smith /*50*/ 0, 129449b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 129549b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 129649b5e25fSSatish Balay 0, 129749b5e25fSSatish Balay 0, 129897304618SKris Buschelman /*55*/ 0, 129949b5e25fSSatish Balay 0, 130049b5e25fSSatish Balay 0, 130149b5e25fSSatish Balay 0, 130249b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 130397304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 130449b5e25fSSatish Balay 0, 130549b5e25fSSatish Balay 0, 13068a124369SBarry Smith MatGetPetscMaps_Petsc, 1307d959ec07SHong Zhang 0, 130897304618SKris Buschelman /*65*/ 0, 1309d959ec07SHong Zhang 0, 1310d959ec07SHong Zhang 0, 1311d959ec07SHong Zhang 0, 1312d959ec07SHong Zhang 0, 131397304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13143e0d88b5SBarry Smith 0, 13153e0d88b5SBarry Smith 0, 13163e0d88b5SBarry Smith 0, 13173e0d88b5SBarry Smith 0, 131897304618SKris Buschelman /*75*/ 0, 13193e0d88b5SBarry Smith 0, 13203e0d88b5SBarry Smith 0, 13213e0d88b5SBarry Smith 0, 13223e0d88b5SBarry Smith 0, 132397304618SKris Buschelman /*80*/ 0, 13243e0d88b5SBarry Smith 0, 13253e0d88b5SBarry Smith 0, 13263e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 132797304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13283e0d88b5SBarry Smith #else 132997304618SKris Buschelman 0, 13303e0d88b5SBarry Smith #endif 1331865e5f61SKris Buschelman MatLoad_SeqSBAIJ, 1332865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ, 1333865e5f61SKris Buschelman MatIsHermitian_SeqSBAIJ, 1334efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1335865e5f61SKris Buschelman 0, 1336865e5f61SKris Buschelman 0, 1337865e5f61SKris Buschelman /*90*/ 0, 1338865e5f61SKris Buschelman 0, 1339865e5f61SKris Buschelman 0, 1340865e5f61SKris Buschelman 0, 1341865e5f61SKris Buschelman 0, 1342865e5f61SKris Buschelman /*95*/ 0, 1343865e5f61SKris Buschelman 0, 1344865e5f61SKris Buschelman 0, 1345865e5f61SKris Buschelman 0}; 134649b5e25fSSatish Balay 134749b5e25fSSatish Balay EXTERN_C_BEGIN 13484a2ae208SSatish Balay #undef __FUNCT__ 13494a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 1350dfbe8321SBarry Smith PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat) 135149b5e25fSSatish Balay { 13524afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1353521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 1354dfbe8321SBarry Smith PetscErrorCode ierr; 135549b5e25fSSatish Balay 135649b5e25fSSatish Balay PetscFunctionBegin; 135749b5e25fSSatish Balay if (aij->nonew != 1) { 1358e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 135949b5e25fSSatish Balay } 136049b5e25fSSatish Balay 136149b5e25fSSatish Balay /* allocate space for values if not already there */ 136249b5e25fSSatish Balay if (!aij->saved_values) { 136387828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 136449b5e25fSSatish Balay } 136549b5e25fSSatish Balay 136649b5e25fSSatish Balay /* copy values over */ 136787828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 136849b5e25fSSatish Balay PetscFunctionReturn(0); 136949b5e25fSSatish Balay } 137049b5e25fSSatish Balay EXTERN_C_END 137149b5e25fSSatish Balay 137249b5e25fSSatish Balay EXTERN_C_BEGIN 13734a2ae208SSatish Balay #undef __FUNCT__ 13744a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 1375dfbe8321SBarry Smith PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat) 137649b5e25fSSatish Balay { 13774afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 13786849ba73SBarry Smith PetscErrorCode ierr; 1379521d7252SBarry Smith PetscInt nz = aij->i[mat->m]*mat->bs*aij->bs2; 138049b5e25fSSatish Balay 138149b5e25fSSatish Balay PetscFunctionBegin; 138249b5e25fSSatish Balay if (aij->nonew != 1) { 1383e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 138449b5e25fSSatish Balay } 138549b5e25fSSatish Balay if (!aij->saved_values) { 1386e005ede5SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); 138749b5e25fSSatish Balay } 138849b5e25fSSatish Balay 138949b5e25fSSatish Balay /* copy values over */ 139087828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 139149b5e25fSSatish Balay PetscFunctionReturn(0); 139249b5e25fSSatish Balay } 139349b5e25fSSatish Balay EXTERN_C_END 139449b5e25fSSatish Balay 13958549e402SHong Zhang EXTERN_C_BEGIN 13964a2ae208SSatish Balay #undef __FUNCT__ 1397a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 139813f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz) 139949b5e25fSSatish Balay { 1400c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14016849ba73SBarry Smith PetscErrorCode ierr; 140213f74950SBarry Smith PetscInt i,len,mbs,bs2; 140349b5e25fSSatish Balay PetscTruth flg; 140449b5e25fSSatish Balay 140549b5e25fSSatish Balay PetscFunctionBegin; 1406273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1407e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1408c464158bSHong Zhang mbs = B->m/bs; 140949b5e25fSSatish Balay bs2 = bs*bs; 141049b5e25fSSatish Balay 1411c464158bSHong Zhang if (mbs*bs != B->m) { 141229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 141349b5e25fSSatish Balay } 141449b5e25fSSatish Balay 1415435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 141677431f27SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz); 141749b5e25fSSatish Balay if (nnz) { 141849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 141977431f27SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]); 142077431f27SBarry 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); 142149b5e25fSSatish Balay } 142249b5e25fSSatish Balay } 142349b5e25fSSatish Balay 1424e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 142549b5e25fSSatish Balay if (!flg) { 142649b5e25fSSatish Balay switch (bs) { 142749b5e25fSSatish Balay case 1: 14284ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 142949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1430d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 1431e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 143249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 143349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 143449b5e25fSSatish Balay break; 143549b5e25fSSatish Balay case 2: 14364ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 143749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 1438e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_2; 143949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 144049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 144149b5e25fSSatish Balay break; 144249b5e25fSSatish Balay case 3: 14435f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 144449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 1445e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_3; 144649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 144749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 144849b5e25fSSatish Balay break; 144949b5e25fSSatish Balay case 4: 14505f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 145149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 1452e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_4; 145349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 145449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 145549b5e25fSSatish Balay break; 145649b5e25fSSatish Balay case 5: 14575f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 145849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 1459e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_5; 146049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 146149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 146249b5e25fSSatish Balay break; 146349b5e25fSSatish Balay case 6: 14645f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 146549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 1466e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_6; 146749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 146849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 146949b5e25fSSatish Balay break; 147049b5e25fSSatish Balay case 7: 1471de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 147249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 1473e005ede5SBarry Smith B->ops->solvetranspose = MatSolve_SeqSBAIJ_7; 1474de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 147549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 147649b5e25fSSatish Balay break; 147749b5e25fSSatish Balay } 147849b5e25fSSatish Balay } 147949b5e25fSSatish Balay 148049b5e25fSSatish Balay b->mbs = mbs; 14814afc71dfSHong Zhang b->nbs = mbs; 148213f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->imax);CHKERRQ(ierr); 148349b5e25fSSatish Balay if (!nnz) { 1484435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 148549b5e25fSSatish Balay else if (nz <= 0) nz = 1; 148649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14878cef66ccSHong Zhang b->imax[i] = nz; 148849b5e25fSSatish Balay } 1489153ea458SHong Zhang nz = nz*mbs; /* total nz */ 149049b5e25fSSatish Balay } else { 149149b5e25fSSatish Balay nz = 0; 14928cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 149349b5e25fSSatish Balay } 14946c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 149549b5e25fSSatish Balay 149649b5e25fSSatish Balay /* allocate the matrix space */ 149713f74950SBarry Smith len = nz*sizeof(PetscInt) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(PetscInt); 149882502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 14996c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 150013f74950SBarry Smith b->j = (PetscInt*)(b->a + nz*bs2); 150113f74950SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 15026c6c5352SBarry Smith b->i = b->j + nz; 150349b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 150449b5e25fSSatish Balay 150549b5e25fSSatish Balay /* pointer to beginning of each row */ 150649b5e25fSSatish Balay b->i[0] = 0; 150749b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 150849b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 150949b5e25fSSatish Balay } 151049b5e25fSSatish Balay 151149b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 151213f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&b->ilen);CHKERRQ(ierr); 151313f74950SBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 151449b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 151549b5e25fSSatish Balay 1516521d7252SBarry Smith B->bs = bs; 151749b5e25fSSatish Balay b->bs2 = bs2; 15186c6c5352SBarry Smith b->nz = 0; 15196c6c5352SBarry Smith b->maxnz = nz*bs2; 1520153ea458SHong Zhang 152116cdd363SHong Zhang b->inew = 0; 152216cdd363SHong Zhang b->jnew = 0; 152316cdd363SHong Zhang b->anew = 0; 152416cdd363SHong Zhang b->a2anew = 0; 15251a3463dfSHong Zhang b->permute = PETSC_FALSE; 1526c464158bSHong Zhang PetscFunctionReturn(0); 1527c464158bSHong Zhang } 1528a23d5eceSKris Buschelman EXTERN_C_END 1529153ea458SHong Zhang 1530d769727bSBarry Smith EXTERN_C_BEGIN 1531dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,Mat*); 1532dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,Mat*); 1533d769727bSBarry Smith EXTERN_C_END 1534d769727bSBarry Smith 15350bad9183SKris Buschelman /*MC 1536fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15370bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15380bad9183SKris Buschelman 15390bad9183SKris Buschelman Options Database Keys: 15400bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15410bad9183SKris Buschelman 15420bad9183SKris Buschelman Level: beginner 15430bad9183SKris Buschelman 15440bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15450bad9183SKris Buschelman M*/ 15460bad9183SKris Buschelman 1547a23d5eceSKris Buschelman EXTERN_C_BEGIN 1548a23d5eceSKris Buschelman #undef __FUNCT__ 1549a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1550dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqSBAIJ(Mat B) 1551a23d5eceSKris Buschelman { 1552a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1553dfbe8321SBarry Smith PetscErrorCode ierr; 155413f74950SBarry Smith PetscMPIInt size; 1555a23d5eceSKris Buschelman 1556a23d5eceSKris Buschelman PetscFunctionBegin; 1557a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1558a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1559a23d5eceSKris Buschelman B->m = B->M = PetscMax(B->m,B->M); 1560a23d5eceSKris Buschelman B->n = B->N = PetscMax(B->n,B->N); 1561a23d5eceSKris Buschelman 1562a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1563a23d5eceSKris Buschelman B->data = (void*)b; 1564a23d5eceSKris Buschelman ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1565a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1566a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1567a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1568a23d5eceSKris Buschelman B->factor = 0; 1569a23d5eceSKris Buschelman B->lupivotthreshold = 1.0; 1570a23d5eceSKris Buschelman B->mapping = 0; 1571a23d5eceSKris Buschelman b->row = 0; 1572a23d5eceSKris Buschelman b->icol = 0; 1573a23d5eceSKris Buschelman b->reallocs = 0; 1574a23d5eceSKris Buschelman b->saved_values = 0; 1575a23d5eceSKris Buschelman 1576a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1577a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1578a23d5eceSKris Buschelman 1579a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1580a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1581a23d5eceSKris Buschelman b->nonew = 0; 1582a23d5eceSKris Buschelman b->diag = 0; 1583a23d5eceSKris Buschelman b->solve_work = 0; 1584a23d5eceSKris Buschelman b->mult_work = 0; 1585a23d5eceSKris Buschelman B->spptr = 0; 1586a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1587a23d5eceSKris Buschelman b->xtoy = 0; 1588a23d5eceSKris Buschelman b->XtoY = 0; 1589a23d5eceSKris Buschelman 1590a23d5eceSKris Buschelman b->inew = 0; 1591a23d5eceSKris Buschelman b->jnew = 0; 1592a23d5eceSKris Buschelman b->anew = 0; 1593a23d5eceSKris Buschelman b->a2anew = 0; 1594a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1595a23d5eceSKris Buschelman 1596a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1597a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1598a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1599a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1600a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1601a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1602a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1603a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1604a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16054e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16064e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16074e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1608a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1609a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1610a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1611a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1612a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1613a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 161423ce1328SBarry Smith 161523ce1328SBarry Smith B->symmetric = PETSC_TRUE; 161623ce1328SBarry Smith B->structurally_symmetric = PETSC_TRUE; 161723ce1328SBarry Smith B->symmetric_set = PETSC_TRUE; 161823ce1328SBarry Smith B->structurally_symmetric_set = PETSC_TRUE; 1619a23d5eceSKris Buschelman PetscFunctionReturn(0); 1620a23d5eceSKris Buschelman } 1621a23d5eceSKris Buschelman EXTERN_C_END 1622a23d5eceSKris Buschelman 1623a23d5eceSKris Buschelman #undef __FUNCT__ 1624a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1625a23d5eceSKris Buschelman /*@C 1626a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1627a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1628a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1629a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1630a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1631a23d5eceSKris Buschelman 1632a23d5eceSKris Buschelman Collective on Mat 1633a23d5eceSKris Buschelman 1634a23d5eceSKris Buschelman Input Parameters: 1635a23d5eceSKris Buschelman + A - the symmetric matrix 1636a23d5eceSKris Buschelman . bs - size of block 1637a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1638a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1639a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1640a23d5eceSKris Buschelman 1641a23d5eceSKris Buschelman Options Database Keys: 1642a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1643a23d5eceSKris Buschelman block calculations (much slower) 1644a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1645a23d5eceSKris Buschelman 1646a23d5eceSKris Buschelman Level: intermediate 1647a23d5eceSKris Buschelman 1648a23d5eceSKris Buschelman Notes: 1649a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1650a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1651a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1652a23d5eceSKris Buschelman matrices. 1653a23d5eceSKris Buschelman 165449a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 165549a6f317SBarry Smith 165649a6f317SBarry Smith 1657a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1658a23d5eceSKris Buschelman @*/ 165913f74950SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[]) 166013f74950SBarry Smith { 166113f74950SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]); 1662a23d5eceSKris Buschelman 1663a23d5eceSKris Buschelman PetscFunctionBegin; 1664a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1665a23d5eceSKris Buschelman if (f) { 1666a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1667a23d5eceSKris Buschelman } 1668a23d5eceSKris Buschelman PetscFunctionReturn(0); 1669a23d5eceSKris Buschelman } 167049b5e25fSSatish Balay 16714a2ae208SSatish Balay #undef __FUNCT__ 16724a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1673c464158bSHong Zhang /*@C 1674c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1675c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1676c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1677c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1678c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 167949b5e25fSSatish Balay 1680c464158bSHong Zhang Collective on MPI_Comm 1681c464158bSHong Zhang 1682c464158bSHong Zhang Input Parameters: 1683c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1684c464158bSHong Zhang . bs - size of block 1685c464158bSHong Zhang . m - number of rows, or number of columns 1686c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1687744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1688744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1689c464158bSHong Zhang 1690c464158bSHong Zhang Output Parameter: 1691c464158bSHong Zhang . A - the symmetric matrix 1692c464158bSHong Zhang 1693c464158bSHong Zhang Options Database Keys: 1694c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1695c464158bSHong Zhang block calculations (much slower) 1696c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1697c464158bSHong Zhang 1698c464158bSHong Zhang Level: intermediate 1699c464158bSHong Zhang 1700c464158bSHong Zhang Notes: 1701c464158bSHong Zhang 1702c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1703c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1704c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1705c464158bSHong Zhang matrices. 1706c464158bSHong Zhang 170749a6f317SBarry Smith If the nnz parameter is given then the nz parameter is ignored 170849a6f317SBarry Smith 1709c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1710c464158bSHong Zhang @*/ 171113f74950SBarry Smith PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) 1712c464158bSHong Zhang { 1713dfbe8321SBarry Smith PetscErrorCode ierr; 1714c464158bSHong Zhang 1715c464158bSHong Zhang PetscFunctionBegin; 1716c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1717c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1718273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 171949b5e25fSSatish Balay PetscFunctionReturn(0); 172049b5e25fSSatish Balay } 172149b5e25fSSatish Balay 17224a2ae208SSatish Balay #undef __FUNCT__ 17234a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 1724dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 172549b5e25fSSatish Balay { 172649b5e25fSSatish Balay Mat C; 172749b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17286849ba73SBarry Smith PetscErrorCode ierr; 172913f74950SBarry Smith PetscInt i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 173049b5e25fSSatish Balay 173149b5e25fSSatish Balay PetscFunctionBegin; 173229bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 173349b5e25fSSatish Balay 173449b5e25fSSatish Balay *B = 0; 1735692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1736be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 17371d5dac46SHong Zhang ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1738692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1739692f9cbeSHong Zhang 1740273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 174149b5e25fSSatish Balay C->factor = A->factor; 174249b5e25fSSatish Balay c->row = 0; 174349b5e25fSSatish Balay c->icol = 0; 174449b5e25fSSatish Balay c->saved_values = 0; 174549b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 174649b5e25fSSatish Balay C->assembled = PETSC_TRUE; 174749b5e25fSSatish Balay 174882327fa8SHong Zhang C->M = A->M; 174982327fa8SHong Zhang C->N = A->N; 1750521d7252SBarry Smith C->bs = A->bs; 175149b5e25fSSatish Balay c->bs2 = a->bs2; 175249b5e25fSSatish Balay c->mbs = a->mbs; 175349b5e25fSSatish Balay c->nbs = a->nbs; 175449b5e25fSSatish Balay 175513f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr); 175613f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr); 175749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 175849b5e25fSSatish Balay c->imax[i] = a->imax[i]; 175949b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 176049b5e25fSSatish Balay } 176149b5e25fSSatish Balay 176249b5e25fSSatish Balay /* allocate the matrix space */ 176349b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 176413f74950SBarry Smith len = (mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)); 176582502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 176613f74950SBarry Smith c->j = (PetscInt*)(c->a + nz*bs2); 176749b5e25fSSatish Balay c->i = c->j + nz; 176813f74950SBarry Smith ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr); 176949b5e25fSSatish Balay if (mbs > 0) { 177013f74950SBarry Smith ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr); 177149b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 177249b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 177349b5e25fSSatish Balay } else { 177449b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 177549b5e25fSSatish Balay } 177649b5e25fSSatish Balay } 177749b5e25fSSatish Balay 177813f74950SBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 177949b5e25fSSatish Balay c->sorted = a->sorted; 178049b5e25fSSatish Balay c->roworiented = a->roworiented; 178149b5e25fSSatish Balay c->nonew = a->nonew; 178249b5e25fSSatish Balay 178349b5e25fSSatish Balay if (a->diag) { 178413f74950SBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); 178513f74950SBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt)); 178649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 178749b5e25fSSatish Balay c->diag[i] = a->diag[i]; 178849b5e25fSSatish Balay } 178949b5e25fSSatish Balay } else c->diag = 0; 17906c6c5352SBarry Smith c->nz = a->nz; 17916c6c5352SBarry Smith c->maxnz = a->maxnz; 179249b5e25fSSatish Balay c->solve_work = 0; 179349b5e25fSSatish Balay c->mult_work = 0; 179449b5e25fSSatish Balay *B = C; 1795b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 179649b5e25fSSatish Balay PetscFunctionReturn(0); 179749b5e25fSSatish Balay } 179849b5e25fSSatish Balay 17994a2ae208SSatish Balay #undef __FUNCT__ 18004a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1801dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A) 180249b5e25fSSatish Balay { 180349b5e25fSSatish Balay Mat_SeqSBAIJ *a; 180449b5e25fSSatish Balay Mat B; 18056849ba73SBarry Smith PetscErrorCode ierr; 180613f74950SBarry Smith int fd; 180713f74950SBarry Smith PetscMPIInt size; 180813f74950SBarry Smith PetscInt i,nz,header[4],*rowlengths=0,M,N,bs=1; 180913f74950SBarry Smith PetscInt *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 181013f74950SBarry Smith PetscInt kmax,jcount,block,idx,point,nzcountb,extra_rows; 181113f74950SBarry Smith PetscInt *masked,nmask,tmp,bs2,ishift; 181287828ca2SBarry Smith PetscScalar *aa; 181349b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 181449b5e25fSSatish Balay 181549b5e25fSSatish Balay PetscFunctionBegin; 1816b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 181749b5e25fSSatish Balay bs2 = bs*bs; 181849b5e25fSSatish Balay 181949b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 182029bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1821b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 182249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1823552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 182449b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 182549b5e25fSSatish Balay 182649b5e25fSSatish Balay if (header[3] < 0) { 182729bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 182849b5e25fSSatish Balay } 182949b5e25fSSatish Balay 183029bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 183149b5e25fSSatish Balay 183249b5e25fSSatish Balay /* 183349b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 183449b5e25fSSatish Balay divisible by the blocksize 183549b5e25fSSatish Balay */ 183649b5e25fSSatish Balay mbs = M/bs; 183749b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 183849b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 183949b5e25fSSatish Balay else mbs++; 184049b5e25fSSatish Balay if (extra_rows) { 1841b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 184249b5e25fSSatish Balay } 184349b5e25fSSatish Balay 184449b5e25fSSatish Balay /* read in row lengths */ 184513f74950SBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 184649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 184749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 184849b5e25fSSatish Balay 184949b5e25fSSatish Balay /* read in column indices */ 185013f74950SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr); 185149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 185249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 185349b5e25fSSatish Balay 185449b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 185513f74950SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr); 185613f74950SBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr); 185713f74950SBarry Smith ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr); 185813f74950SBarry Smith ierr = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr); 185949b5e25fSSatish Balay masked = mask + mbs; 186049b5e25fSSatish Balay rowcount = 0; nzcount = 0; 186149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 186249b5e25fSSatish Balay nmask = 0; 186349b5e25fSSatish Balay for (j=0; j<bs; j++) { 186449b5e25fSSatish Balay kmax = rowlengths[rowcount]; 186549b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18662d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 186703630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 186849b5e25fSSatish Balay } 186949b5e25fSSatish Balay rowcount++; 187049b5e25fSSatish Balay } 1871574b2666SHong Zhang s_browlengths[i] += nmask; 1872574b2666SHong Zhang 187349b5e25fSSatish Balay /* zero out the mask elements we set */ 187449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 187549b5e25fSSatish Balay } 187649b5e25fSSatish Balay 187749b5e25fSSatish Balay /* create our matrix */ 18789abb65ffSKris Buschelman ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr); 18799abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 188078473fd7SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr); 188149b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 188249b5e25fSSatish Balay 188349b5e25fSSatish Balay /* set matrix "i" values */ 188449b5e25fSSatish Balay a->i[0] = 0; 188549b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1886574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1887574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 188849b5e25fSSatish Balay } 18896c6c5352SBarry Smith a->nz = a->i[mbs]; 189049b5e25fSSatish Balay 189149b5e25fSSatish Balay /* read in nonzero values */ 189287828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 189349b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 189449b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 189549b5e25fSSatish Balay 189649b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 189749b5e25fSSatish Balay nzcount = 0; jcount = 0; 189849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 189949b5e25fSSatish Balay nzcountb = nzcount; 190049b5e25fSSatish Balay nmask = 0; 190149b5e25fSSatish Balay for (j=0; j<bs; j++) { 190249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 190349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 19042d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 190503630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 19062d703238SHong Zhang } 19072d703238SHong Zhang } 19082d703238SHong Zhang /* sort the masked values */ 19092d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19102d703238SHong Zhang 19112d703238SHong Zhang /* set "j" values into matrix */ 19122d703238SHong Zhang maskcount = 1; 19132d703238SHong Zhang for (j=0; j<nmask; j++) { 191449b5e25fSSatish Balay a->j[jcount++] = masked[j]; 191549b5e25fSSatish Balay mask[masked[j]] = maskcount++; 191649b5e25fSSatish Balay } 1917574b2666SHong Zhang 191849b5e25fSSatish Balay /* set "a" values into matrix */ 191949b5e25fSSatish Balay ishift = bs2*a->i[i]; 192049b5e25fSSatish Balay for (j=0; j<bs; j++) { 192149b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 192249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1923574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1924574b2666SHong Zhang if (tmp >= i){ 192549b5e25fSSatish Balay block = mask[tmp] - 1; 192649b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 192749b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1928574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1929574b2666SHong Zhang } 1930574b2666SHong Zhang nzcountb++; 193149b5e25fSSatish Balay } 193249b5e25fSSatish Balay } 193349b5e25fSSatish Balay /* zero out the mask elements we set */ 193449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 193549b5e25fSSatish Balay } 19366c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 193749b5e25fSSatish Balay 193849b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1939574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 194049b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 194149b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 194249b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 194349b5e25fSSatish Balay 19449abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19459abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 194649b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19479abb65ffSKris Buschelman *A = B; 194849b5e25fSSatish Balay PetscFunctionReturn(0); 194949b5e25fSSatish Balay } 1950574b2666SHong Zhang 1951d06b337dSHong Zhang #undef __FUNCT__ 1952d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 195313f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 1954d06b337dSHong Zhang { 1955d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1956d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1957d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 19586849ba73SBarry Smith PetscErrorCode ierr; 1959521d7252SBarry Smith PetscInt m=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j; 1960521d7252SBarry Smith PetscInt nz,nz1,*vj,*vj1,i; 1961d06b337dSHong Zhang 1962d06b337dSHong Zhang PetscFunctionBegin; 1963b965ef7fSBarry Smith its = its*lits; 196477431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 1965d06b337dSHong Zhang 1966d06b337dSHong Zhang if (bs > 1) 1967d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1968d06b337dSHong Zhang 19691ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1970d06b337dSHong Zhang if (xx != bb) { 19711ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1972d06b337dSHong Zhang } else { 1973d06b337dSHong Zhang b = x; 1974d06b337dSHong Zhang } 1975d06b337dSHong Zhang 1976d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1977d06b337dSHong Zhang 1978d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19792798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1980d06b337dSHong Zhang for (i=0; i<m; i++) 1981d06b337dSHong Zhang t[i] = b[i]; 1982d06b337dSHong Zhang 1983d06b337dSHong Zhang for (i=0; i<m; i++){ 198444706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1985d06b337dSHong Zhang v = aa + ai[i] + 1; 1986d06b337dSHong Zhang vj = aj + ai[i] + 1; 1987d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1988d06b337dSHong Zhang x[i] = omega*t[i]/d; 1989d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 199044706875SHong Zhang PetscLogFlops(2*nz-1); 1991d06b337dSHong Zhang } 1992d06b337dSHong Zhang } 1993d06b337dSHong Zhang 19942798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1995d06b337dSHong Zhang for (i=0; i<m; i++) 1996d06b337dSHong Zhang t[i] = b[i]; 1997d06b337dSHong Zhang 1998d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1999d06b337dSHong Zhang v = aa + ai[i] + 1; 2000d06b337dSHong Zhang vj = aj + ai[i] + 1; 2001d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2002d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 200344706875SHong Zhang PetscLogFlops(2*nz-1); 2004d06b337dSHong Zhang } 2005d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2006d06b337dSHong Zhang d = *(aa + ai[i]); 2007d06b337dSHong Zhang v = aa + ai[i] + 1; 2008d06b337dSHong Zhang vj = aj + ai[i] + 1; 2009d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2010d06b337dSHong Zhang sum = t[i]; 2011d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 201244706875SHong Zhang PetscLogFlops(2*nz-1); 2013d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2014d06b337dSHong Zhang } 2015d06b337dSHong Zhang } 2016d06b337dSHong Zhang its--; 2017d06b337dSHong Zhang } 2018d06b337dSHong Zhang 2019d06b337dSHong Zhang while (its--) { 202044706875SHong Zhang /* 202144706875SHong Zhang forward sweep: 202244706875SHong Zhang for i=0,...,m-1: 202344706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 202444706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 202544706875SHong Zhang b = b - x[i]*U^T(i,:); 2026d06b337dSHong Zhang 202744706875SHong Zhang */ 20282798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2029d06b337dSHong Zhang for (i=0; i<m; i++) 2030d06b337dSHong Zhang t[i] = b[i]; 2031d06b337dSHong Zhang 2032d06b337dSHong Zhang for (i=0; i<m; i++){ 203344706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2034d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2035d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2036d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2037d06b337dSHong Zhang sum = t[i]; 2038d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2039d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2040d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 204144706875SHong Zhang PetscLogFlops(4*nz-2); 2042d06b337dSHong Zhang } 2043d06b337dSHong Zhang } 2044d06b337dSHong Zhang 20452798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 204644706875SHong Zhang /* 204744706875SHong Zhang backward sweep: 204844706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 204944706875SHong Zhang for i=m-1,...,0: 205044706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 205144706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 205244706875SHong Zhang */ 2053d06b337dSHong Zhang for (i=0; i<m; i++) 2054d06b337dSHong Zhang t[i] = b[i]; 2055d06b337dSHong Zhang 2056d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2057d06b337dSHong Zhang v = aa + ai[i] + 1; 2058d06b337dSHong Zhang vj = aj + ai[i] + 1; 2059d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2060d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 206144706875SHong Zhang PetscLogFlops(2*nz-1); 2062d06b337dSHong Zhang } 2063d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2064d06b337dSHong Zhang d = *(aa + ai[i]); 2065d06b337dSHong Zhang v = aa + ai[i] + 1; 2066d06b337dSHong Zhang vj = aj + ai[i] + 1; 2067d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2068d06b337dSHong Zhang sum = t[i]; 2069d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 207044706875SHong Zhang PetscLogFlops(2*nz-1); 2071d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2072d06b337dSHong Zhang } 2073d06b337dSHong Zhang } 2074d06b337dSHong Zhang } 2075d06b337dSHong Zhang 2076d06b337dSHong Zhang ierr = PetscFree(t);CHKERRQ(ierr); 20771ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2078d06b337dSHong Zhang if (bb != xx) { 20791ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2080d06b337dSHong Zhang } 20812798e883SHong Zhang 2082d06b337dSHong Zhang PetscFunctionReturn(0); 2083d06b337dSHong Zhang } 2084d06b337dSHong Zhang 2085d06b337dSHong Zhang 2086d06b337dSHong Zhang 2087d06b337dSHong Zhang 208849b5e25fSSatish Balay 208949b5e25fSSatish Balay 2090