173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 4a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/vec/vecimpl.h" 949b5e25fSSatish Balay #include "src/inline/spops.h" 10aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1149b5e25fSSatish Balay 1249b5e25fSSatish Balay #define CHUNKSIZE 10 1349b5e25fSSatish Balay 1449b5e25fSSatish Balay /* 1549b5e25fSSatish Balay Checks for missing diagonals 1649b5e25fSSatish Balay */ 174a2ae208SSatish Balay #undef __FUNCT__ 184a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 1949b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2049b5e25fSSatish Balay { 21045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2249b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2349b5e25fSSatish Balay 2449b5e25fSSatish Balay PetscFunctionBegin; 25045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2649b5e25fSSatish Balay diag = a->diag; 2749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 28e11b0e14SHong Zhang if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i); 2949b5e25fSSatish Balay } 3049b5e25fSSatish Balay PetscFunctionReturn(0); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay 334a2ae208SSatish Balay #undef __FUNCT__ 344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 3549b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3649b5e25fSSatish Balay { 37045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 38b424e231SHong Zhang int i,mbs = a->mbs,ierr; 3949b5e25fSSatish Balay 4049b5e25fSSatish Balay PetscFunctionBegin; 4149b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4249b5e25fSSatish Balay 43b424e231SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr); 44b424e231SHong Zhang PetscLogObjectMemory(A,(mbs+1)*sizeof(int)); 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" 5149b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 5249b5e25fSSatish Balay { 53a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 549e37e433SSatish Balay int 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 */ 62a6ece127SHong Zhang int s_nz = a->i[n]; 63a6ece127SHong Zhang for (i=0; i<s_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" 74045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 7549b5e25fSSatish Balay { 76b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 77a15ac5e4SHong Zhang int i,n = a->mbs; 78a6ece127SHong Zhang 7949b5e25fSSatish Balay PetscFunctionBegin; 8049b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 81a6ece127SHong Zhang 82a6ece127SHong Zhang if (oshift == 1) { 83a6ece127SHong Zhang int s_nz = a->i[n]-1; 84a6ece127SHong Zhang for (i=0; i<s_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__ "MatGetBlockSize_SeqSBAIJ" 9249b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9349b5e25fSSatish Balay { 9449b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 9549b5e25fSSatish Balay 9649b5e25fSSatish Balay PetscFunctionBegin; 9749b5e25fSSatish Balay *bs = sbaij->bs; 9849b5e25fSSatish Balay PetscFunctionReturn(0); 9949b5e25fSSatish Balay } 10049b5e25fSSatish Balay 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 10349b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 10449b5e25fSSatish Balay { 10549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10649b5e25fSSatish Balay int ierr; 10749b5e25fSSatish Balay 10849b5e25fSSatish Balay PetscFunctionBegin; 10949b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 110b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 11149b5e25fSSatish Balay #endif 11249b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 11349b5e25fSSatish Balay if (!a->singlemalloc) { 11449b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 11563c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 11649b5e25fSSatish Balay } 11749b5e25fSSatish Balay if (a->row) { 11849b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 11949b5e25fSSatish Balay } 12049b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 12149b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 12249b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 12349b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 124d59c15a7SBarry Smith if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 125d59c15a7SBarry Smith if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);} 126d59c15a7SBarry Smith if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 12749b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1281a3463dfSHong Zhang 1291a3463dfSHong Zhang if (a->inew){ 1301a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1311a3463dfSHong Zhang a->inew = 0; 1321a3463dfSHong Zhang } 13349b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 13449b5e25fSSatish Balay PetscFunctionReturn(0); 13549b5e25fSSatish Balay } 13649b5e25fSSatish Balay 1374a2ae208SSatish Balay #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 13949b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 14049b5e25fSSatish Balay { 141045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 14249b5e25fSSatish Balay 14349b5e25fSSatish Balay PetscFunctionBegin; 1444d9d31abSKris Buschelman switch (op) { 1454d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1464d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1474d9d31abSKris Buschelman break; 1484d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1494d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1504d9d31abSKris Buschelman break; 1514d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1524d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1554d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1584d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1614d9d31abSKris Buschelman a->nonew = 1; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1644d9d31abSKris Buschelman a->nonew = -1; 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1674d9d31abSKris Buschelman a->nonew = -2; 1684d9d31abSKris Buschelman break; 1694d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1704d9d31abSKris Buschelman a->nonew = 0; 1714d9d31abSKris Buschelman break; 1724d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1734d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1744d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1754d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1764d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 177b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1784d9d31abSKris Buschelman break; 1794d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 18029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1814d9d31abSKris Buschelman default: 18229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 18349b5e25fSSatish Balay } 18449b5e25fSSatish Balay PetscFunctionReturn(0); 18549b5e25fSSatish Balay } 18649b5e25fSSatish Balay 1874a2ae208SSatish Balay #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 18987828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v) 19049b5e25fSSatish Balay { 19149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 19282502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr; 19349b5e25fSSatish Balay MatScalar *aa,*aa_i; 19487828ca2SBarry Smith PetscScalar *v_i; 19549b5e25fSSatish Balay 19649b5e25fSSatish Balay PetscFunctionBegin; 19749b5e25fSSatish Balay bs = a->bs; 19849b5e25fSSatish Balay ai = a->i; 19949b5e25fSSatish Balay aj = a->j; 20049b5e25fSSatish Balay aa = a->a; 20149b5e25fSSatish Balay bs2 = a->bs2; 20249b5e25fSSatish Balay 203c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 20449b5e25fSSatish Balay 20549b5e25fSSatish Balay bn = row/bs; /* Block number */ 20649b5e25fSSatish Balay bp = row % bs; /* Block position */ 20749b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 20849b5e25fSSatish Balay *ncols = bs*M; 20949b5e25fSSatish Balay 21049b5e25fSSatish Balay if (v) { 21149b5e25fSSatish Balay *v = 0; 21249b5e25fSSatish Balay if (*ncols) { 21387828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 21449b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 21549b5e25fSSatish Balay v_i = *v + i*bs; 21649b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 21749b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 21849b5e25fSSatish Balay } 21949b5e25fSSatish Balay } 22049b5e25fSSatish Balay } 22149b5e25fSSatish Balay 22249b5e25fSSatish Balay if (cols) { 22349b5e25fSSatish Balay *cols = 0; 22449b5e25fSSatish Balay if (*ncols) { 22582502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 22649b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22749b5e25fSSatish Balay cols_i = *cols + i*bs; 22849b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 22949b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 23049b5e25fSSatish Balay } 23149b5e25fSSatish Balay } 23249b5e25fSSatish Balay } 23349b5e25fSSatish Balay 23449b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2355ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2365ddb2528SHong Zhang #ifdef column_search 23749b5e25fSSatish Balay v_i = *v + M*bs; 23849b5e25fSSatish Balay cols_i = *cols + M*bs; 23949b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 24049b5e25fSSatish Balay M = ai[i+1] - ai[i]; 24149b5e25fSSatish Balay for (j=0; j<M; j++){ 24249b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 24349b5e25fSSatish Balay if (itmp == bn){ 24449b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 24549b5e25fSSatish Balay for (k=0; k<bs; k++) { 24649b5e25fSSatish Balay *cols_i++ = i*bs+k; 24749b5e25fSSatish Balay *v_i++ = aa_i[k]; 24849b5e25fSSatish Balay } 24949b5e25fSSatish Balay *ncols += bs; 25049b5e25fSSatish Balay break; 25149b5e25fSSatish Balay } 25249b5e25fSSatish Balay } 25349b5e25fSSatish Balay } 2545ddb2528SHong Zhang #endif 25549b5e25fSSatish Balay 25649b5e25fSSatish Balay PetscFunctionReturn(0); 25749b5e25fSSatish Balay } 25849b5e25fSSatish Balay 2594a2ae208SSatish Balay #undef __FUNCT__ 2604a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 26187828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 26249b5e25fSSatish Balay { 26349b5e25fSSatish Balay int ierr; 26449b5e25fSSatish Balay 26549b5e25fSSatish Balay PetscFunctionBegin; 26649b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 26749b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 26849b5e25fSSatish Balay PetscFunctionReturn(0); 26949b5e25fSSatish Balay } 27049b5e25fSSatish Balay 2714a2ae208SSatish Balay #undef __FUNCT__ 2724a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 27349b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 27449b5e25fSSatish Balay { 2758115998fSBarry Smith int ierr; 27649b5e25fSSatish Balay PetscFunctionBegin; 277999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2788115998fSBarry Smith PetscFunctionReturn(0); 27949b5e25fSSatish Balay } 28049b5e25fSSatish Balay 2814a2ae208SSatish Balay #undef __FUNCT__ 2824a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary" 283b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer) 28449b5e25fSSatish Balay { 28549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 28649b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 28787828ca2SBarry Smith PetscScalar *aa; 28849b5e25fSSatish Balay FILE *file; 28949b5e25fSSatish Balay 29049b5e25fSSatish Balay PetscFunctionBegin; 291b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 29282502324SSatish Balay ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 293552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 29449b5e25fSSatish Balay 295c464158bSHong Zhang col_lens[1] = A->m; 296c464158bSHong Zhang col_lens[2] = A->m; 29749b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 29849b5e25fSSatish Balay 29949b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 30049b5e25fSSatish Balay count = 0; 30149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30249b5e25fSSatish Balay for (j=0; j<bs; j++) { 30349b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 30449b5e25fSSatish Balay } 30549b5e25fSSatish Balay } 306c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 30749b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 30849b5e25fSSatish Balay 30949b5e25fSSatish Balay /* store column indices (zero start index) */ 31082502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 31149b5e25fSSatish Balay count = 0; 31249b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31349b5e25fSSatish Balay for (j=0; j<bs; j++) { 31449b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31549b5e25fSSatish Balay for (l=0; l<bs; l++) { 31649b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 31749b5e25fSSatish Balay } 31849b5e25fSSatish Balay } 31949b5e25fSSatish Balay } 32049b5e25fSSatish Balay } 32149b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 32249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 32349b5e25fSSatish Balay 32449b5e25fSSatish Balay /* store nonzero values */ 32587828ca2SBarry Smith ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 32649b5e25fSSatish Balay count = 0; 32749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 32849b5e25fSSatish Balay for (j=0; j<bs; j++) { 32949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33049b5e25fSSatish Balay for (l=0; l<bs; l++) { 33149b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 33249b5e25fSSatish Balay } 33349b5e25fSSatish Balay } 33449b5e25fSSatish Balay } 33549b5e25fSSatish Balay } 33649b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 33749b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 33849b5e25fSSatish Balay 339b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 34049b5e25fSSatish Balay if (file) { 34149b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 34249b5e25fSSatish Balay } 34349b5e25fSSatish Balay PetscFunctionReturn(0); 34449b5e25fSSatish Balay } 34549b5e25fSSatish Balay 3464a2ae208SSatish Balay #undef __FUNCT__ 3474a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 348b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 34949b5e25fSSatish Balay { 35049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 351fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 352fb9695e5SSatish Balay char *name; 353f3ef73ceSBarry Smith PetscViewerFormat format; 35449b5e25fSSatish Balay 35549b5e25fSSatish Balay PetscFunctionBegin; 35680fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 357b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 358456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 359b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 360fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 36129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 362fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 363b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36449b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 36549b5e25fSSatish Balay for (j=0; j<bs; j++) { 366b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 36749b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 36849b5e25fSSatish Balay for (l=0; l<bs; l++) { 36949b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37049b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37162b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 37249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37349b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 37549b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37649b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37849b5e25fSSatish Balay } 37949b5e25fSSatish Balay #else 38049b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 38162b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38249b5e25fSSatish Balay } 38349b5e25fSSatish Balay #endif 38449b5e25fSSatish Balay } 38549b5e25fSSatish Balay } 386b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 38749b5e25fSSatish Balay } 38849b5e25fSSatish Balay } 389b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 39049b5e25fSSatish Balay } else { 391b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 39249b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 39349b5e25fSSatish Balay for (j=0; j<bs; j++) { 394b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 39549b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 39649b5e25fSSatish Balay for (l=0; l<bs; l++) { 39749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 39849b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 39962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 40049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40149b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 40262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 40349b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40449b5e25fSSatish Balay } else { 40562b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40649b5e25fSSatish Balay } 40749b5e25fSSatish Balay #else 408b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 40949b5e25fSSatish Balay #endif 41049b5e25fSSatish Balay } 41149b5e25fSSatish Balay } 412b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 41349b5e25fSSatish Balay } 41449b5e25fSSatish Balay } 415b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 41649b5e25fSSatish Balay } 417b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 41849b5e25fSSatish Balay PetscFunctionReturn(0); 41949b5e25fSSatish Balay } 42049b5e25fSSatish Balay 4214a2ae208SSatish Balay #undef __FUNCT__ 4224a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 423b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 42449b5e25fSSatish Balay { 42549b5e25fSSatish Balay Mat A = (Mat) Aa; 42649b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 42749b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 42849b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 42949b5e25fSSatish Balay MatScalar *aa; 43049b5e25fSSatish Balay MPI_Comm comm; 431b0a32e0cSBarry Smith PetscViewer viewer; 43249b5e25fSSatish Balay 43349b5e25fSSatish Balay PetscFunctionBegin; 43449b5e25fSSatish Balay /* 43549b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 43649b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 43749b5e25fSSatish Balay rest should return immediately. 43849b5e25fSSatish Balay */ 43949b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 44049b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 44149b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 44249b5e25fSSatish Balay 44349b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 44449b5e25fSSatish Balay 445b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 446b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 44749b5e25fSSatish Balay 44849b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 449b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 45049b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45149b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 452c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45349b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45449b5e25fSSatish Balay aa = a->a + j*bs2; 45549b5e25fSSatish Balay for (k=0; k<bs; k++) { 45649b5e25fSSatish Balay for (l=0; l<bs; l++) { 45749b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 458b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 45949b5e25fSSatish Balay } 46049b5e25fSSatish Balay } 46149b5e25fSSatish Balay } 46249b5e25fSSatish Balay } 463b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 46449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 466c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 46749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46849b5e25fSSatish Balay aa = a->a + j*bs2; 46949b5e25fSSatish Balay for (k=0; k<bs; k++) { 47049b5e25fSSatish Balay for (l=0; l<bs; l++) { 47149b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 472b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 47349b5e25fSSatish Balay } 47449b5e25fSSatish Balay } 47549b5e25fSSatish Balay } 47649b5e25fSSatish Balay } 47749b5e25fSSatish Balay 478b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 47949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 48049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 481c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 48249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 48349b5e25fSSatish Balay aa = a->a + j*bs2; 48449b5e25fSSatish Balay for (k=0; k<bs; k++) { 48549b5e25fSSatish Balay for (l=0; l<bs; l++) { 48649b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 487b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 48849b5e25fSSatish Balay } 48949b5e25fSSatish Balay } 49049b5e25fSSatish Balay } 49149b5e25fSSatish Balay } 49249b5e25fSSatish Balay PetscFunctionReturn(0); 49349b5e25fSSatish Balay } 49449b5e25fSSatish Balay 4954a2ae208SSatish Balay #undef __FUNCT__ 4964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 497b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 49849b5e25fSSatish Balay { 49949b5e25fSSatish Balay int ierr; 50049b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 501b0a32e0cSBarry Smith PetscDraw draw; 50249b5e25fSSatish Balay PetscTruth isnull; 50349b5e25fSSatish Balay 50449b5e25fSSatish Balay PetscFunctionBegin; 50549b5e25fSSatish Balay 506b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 507b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 50849b5e25fSSatish Balay 50949b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 510c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 51149b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 512b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 513b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 51449b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 51549b5e25fSSatish Balay PetscFunctionReturn(0); 51649b5e25fSSatish Balay } 51749b5e25fSSatish Balay 5184a2ae208SSatish Balay #undef __FUNCT__ 5194a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 520b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 52149b5e25fSSatish Balay { 52249b5e25fSSatish Balay int ierr; 52349b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 52449b5e25fSSatish Balay 52549b5e25fSSatish Balay PetscFunctionBegin; 526b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 527b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 528fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 529fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 53049b5e25fSSatish Balay if (issocket) { 53129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 53249b5e25fSSatish Balay } else if (isascii){ 53349b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 53449b5e25fSSatish Balay } else if (isbinary) { 53549b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 53649b5e25fSSatish Balay } else if (isdraw) { 53749b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 53849b5e25fSSatish Balay } else { 539b7aaefc3SHong Zhang SETERRQ1(1,"Viewer type %s not supported by SeqSBAIJ matrices",((PetscObject)viewer)->type_name); 54049b5e25fSSatish Balay } 54149b5e25fSSatish Balay PetscFunctionReturn(0); 54249b5e25fSSatish Balay } 54349b5e25fSSatish Balay 54449b5e25fSSatish Balay 5454a2ae208SSatish Balay #undef __FUNCT__ 5464a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 54787828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v) 54849b5e25fSSatish Balay { 549045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 55049b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 55149b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 55249b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 55349b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 55449b5e25fSSatish Balay 55549b5e25fSSatish Balay PetscFunctionBegin; 55649b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 55749b5e25fSSatish Balay row = im[k]; brow = row/bs; 55829bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 559c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56049b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 56149b5e25fSSatish Balay nrow = ailen[brow]; 56249b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 56329bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 564c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 56549b5e25fSSatish Balay col = in[l] ; 56649b5e25fSSatish Balay bcol = col/bs; 56749b5e25fSSatish Balay cidx = col%bs; 56849b5e25fSSatish Balay ridx = row%bs; 56949b5e25fSSatish Balay high = nrow; 57049b5e25fSSatish Balay low = 0; /* assume unsorted */ 57149b5e25fSSatish Balay while (high-low > 5) { 57249b5e25fSSatish Balay t = (low+high)/2; 57349b5e25fSSatish Balay if (rp[t] > bcol) high = t; 57449b5e25fSSatish Balay else low = t; 57549b5e25fSSatish Balay } 57649b5e25fSSatish Balay for (i=low; i<high; i++) { 57749b5e25fSSatish Balay if (rp[i] > bcol) break; 57849b5e25fSSatish Balay if (rp[i] == bcol) { 57949b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 58049b5e25fSSatish Balay goto finished; 58149b5e25fSSatish Balay } 58249b5e25fSSatish Balay } 58349b5e25fSSatish Balay *v++ = zero; 58449b5e25fSSatish Balay finished:; 58549b5e25fSSatish Balay } 58649b5e25fSSatish Balay } 58749b5e25fSSatish Balay PetscFunctionReturn(0); 58849b5e25fSSatish Balay } 58949b5e25fSSatish Balay 59049b5e25fSSatish Balay 5914a2ae208SSatish Balay #undef __FUNCT__ 5924a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 59387828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 59449b5e25fSSatish Balay { 5950880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5960880e062SHong Zhang int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 5970880e062SHong Zhang int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 5980880e062SHong Zhang int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 5990880e062SHong Zhang PetscTruth roworiented=a->roworiented; 6000880e062SHong Zhang MatScalar *value = v,*ap,*aa = a->a,*bap; 6010880e062SHong Zhang 60249b5e25fSSatish Balay PetscFunctionBegin; 6030880e062SHong Zhang if (roworiented) { 6040880e062SHong Zhang stepval = (n-1)*bs; 6050880e062SHong Zhang } else { 6060880e062SHong Zhang stepval = (m-1)*bs; 6070880e062SHong Zhang } 6080880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 6090880e062SHong Zhang row = im[k]; 6100880e062SHong Zhang if (row < 0) continue; 6110880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 6120880e062SHong Zhang if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 6130880e062SHong Zhang #endif 6140880e062SHong Zhang rp = aj + ai[row]; 6150880e062SHong Zhang ap = aa + bs2*ai[row]; 6160880e062SHong Zhang rmax = imax[row]; 6170880e062SHong Zhang nrow = ailen[row]; 6180880e062SHong Zhang low = 0; 6190880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6200880e062SHong Zhang if (in[l] < 0) continue; 6210880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 6220880e062SHong Zhang if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 6230880e062SHong Zhang #endif 6240880e062SHong Zhang col = in[l]; 6250880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6260880e062SHong Zhang if (roworiented) { 6270880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6280880e062SHong Zhang } else { 6290880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6300880e062SHong Zhang } 6310880e062SHong Zhang if (!sorted) low = 0; high = nrow; 6320880e062SHong Zhang while (high-low > 7) { 6330880e062SHong Zhang t = (low+high)/2; 6340880e062SHong Zhang if (rp[t] > col) high = t; 6350880e062SHong Zhang else low = t; 6360880e062SHong Zhang } 6370880e062SHong Zhang for (i=low; i<high; i++) { 6380880e062SHong Zhang if (rp[i] > col) break; 6390880e062SHong Zhang if (rp[i] == col) { 6400880e062SHong Zhang bap = ap + bs2*i; 6410880e062SHong Zhang if (roworiented) { 6420880e062SHong Zhang if (is == ADD_VALUES) { 6430880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6440880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6450880e062SHong Zhang bap[jj] += *value++; 6460880e062SHong Zhang } 6470880e062SHong Zhang } 6480880e062SHong Zhang } else { 6490880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6500880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6510880e062SHong Zhang bap[jj] = *value++; 6520880e062SHong Zhang } 6530880e062SHong Zhang } 6540880e062SHong Zhang } 6550880e062SHong Zhang } else { 6560880e062SHong Zhang if (is == ADD_VALUES) { 6570880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6580880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6590880e062SHong Zhang *bap++ += *value++; 6600880e062SHong Zhang } 6610880e062SHong Zhang } 6620880e062SHong Zhang } else { 6630880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6640880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6650880e062SHong Zhang *bap++ = *value++; 6660880e062SHong Zhang } 6670880e062SHong Zhang } 6680880e062SHong Zhang } 6690880e062SHong Zhang } 6700880e062SHong Zhang goto noinsert2; 6710880e062SHong Zhang } 6720880e062SHong Zhang } 6730880e062SHong Zhang if (nonew == 1) goto noinsert2; 6740880e062SHong Zhang else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 6750880e062SHong Zhang if (nrow >= rmax) { 6760880e062SHong Zhang /* there is no extra room in row, therefore enlarge */ 6770880e062SHong Zhang int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6780880e062SHong Zhang MatScalar *new_a; 6790880e062SHong Zhang 6800880e062SHong Zhang if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 6810880e062SHong Zhang 6820880e062SHong Zhang /* malloc new storage space */ 6830880e062SHong Zhang len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 6840880e062SHong Zhang ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 6850880e062SHong Zhang new_j = (int*)(new_a + bs2*new_nz); 6860880e062SHong Zhang new_i = new_j + new_nz; 6870880e062SHong Zhang 6880880e062SHong Zhang /* copy over old data into new slots */ 6890880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6900880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 6910880e062SHong Zhang ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 6920880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 6930880e062SHong Zhang ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 6940880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 6950880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 6960880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 6970880e062SHong Zhang /* free up old matrix storage */ 6980880e062SHong Zhang ierr = PetscFree(a->a);CHKERRQ(ierr); 6990880e062SHong Zhang if (!a->singlemalloc) { 7000880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 7010880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 7020880e062SHong Zhang } 7030880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 7040880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 7050880e062SHong Zhang 7060880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 7070880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 7080880e062SHong Zhang PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 7090880e062SHong Zhang a->s_maxnz += bs2*CHUNKSIZE; 7100880e062SHong Zhang a->reallocs++; 7110880e062SHong Zhang a->s_nz++; 7120880e062SHong Zhang } 7130880e062SHong Zhang N = nrow++ - 1; 7140880e062SHong Zhang /* shift up all the later entries in this row */ 7150880e062SHong Zhang for (ii=N; ii>=i; ii--) { 7160880e062SHong Zhang rp[ii+1] = rp[ii]; 7170880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 7180880e062SHong Zhang } 7190880e062SHong Zhang if (N >= i) { 7200880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 7210880e062SHong Zhang } 7220880e062SHong Zhang rp[i] = col; 7230880e062SHong Zhang bap = ap + bs2*i; 7240880e062SHong Zhang if (roworiented) { 7250880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7260880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 7270880e062SHong Zhang bap[jj] = *value++; 7280880e062SHong Zhang } 7290880e062SHong Zhang } 7300880e062SHong Zhang } else { 7310880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7320880e062SHong Zhang for (jj=0; jj<bs; jj++) { 7330880e062SHong Zhang *bap++ = *value++; 7340880e062SHong Zhang } 7350880e062SHong Zhang } 7360880e062SHong Zhang } 7370880e062SHong Zhang noinsert2:; 7380880e062SHong Zhang low = i; 7390880e062SHong Zhang } 7400880e062SHong Zhang ailen[row] = nrow; 7410880e062SHong Zhang } 7420880e062SHong Zhang PetscFunctionReturn(0); 74349b5e25fSSatish Balay } 74449b5e25fSSatish Balay 7454a2ae208SSatish Balay #undef __FUNCT__ 7464a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 74749b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 74849b5e25fSSatish Balay { 74949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 75049b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 751c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 75249b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 75349b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 754cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 755cf4441caSHong Zhang PetscTruth flag; 756cf4441caSHong Zhang #endif 75749b5e25fSSatish Balay 75849b5e25fSSatish Balay PetscFunctionBegin; 75949b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 76049b5e25fSSatish Balay 76149b5e25fSSatish Balay if (m) rmax = ailen[0]; 76249b5e25fSSatish Balay for (i=1; i<mbs; i++) { 76349b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 76449b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 76549b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 76649b5e25fSSatish Balay if (fshift) { 76749b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 76849b5e25fSSatish Balay N = ailen[i]; 76949b5e25fSSatish Balay for (j=0; j<N; j++) { 77049b5e25fSSatish Balay ip[j-fshift] = ip[j]; 77149b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 77249b5e25fSSatish Balay } 77349b5e25fSSatish Balay } 77449b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 77549b5e25fSSatish Balay } 77649b5e25fSSatish Balay if (mbs) { 77749b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 77849b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 77949b5e25fSSatish Balay } 78049b5e25fSSatish Balay /* reset ilen and imax for each row */ 78149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 78249b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 78349b5e25fSSatish Balay } 78449b5e25fSSatish Balay a->s_nz = ai[mbs]; 78549b5e25fSSatish Balay 786b424e231SHong Zhang /* diagonals may have moved, reset it */ 787b424e231SHong Zhang if (a->diag) { 788b424e231SHong Zhang ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr); 78949b5e25fSSatish Balay } 790b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 791c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 792b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 79349b5e25fSSatish Balay a->reallocs); 794b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 79549b5e25fSSatish Balay a->reallocs = 0; 79649b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 79749b5e25fSSatish Balay 798cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 7992c535e4dSHong Zhang ierr = PetscOptionsHasName(A->prefix,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr); 800cf4441caSHong Zhang if (flag) { ierr = MatUseSpooles_SeqSBAIJ(A);CHKERRQ(ierr); } 801cf4441caSHong Zhang #endif 802cf4441caSHong Zhang 80349b5e25fSSatish Balay PetscFunctionReturn(0); 80449b5e25fSSatish Balay } 80549b5e25fSSatish Balay 80649b5e25fSSatish Balay /* 80749b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 80849b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 80949b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 81049b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 81149b5e25fSSatish Balay */ 8124a2ae208SSatish Balay #undef __FUNCT__ 8134a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 814db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 81549b5e25fSSatish Balay { 81649b5e25fSSatish Balay int i,j,k,row; 81749b5e25fSSatish Balay PetscTruth flg; 81849b5e25fSSatish Balay 81949b5e25fSSatish Balay PetscFunctionBegin; 82049b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 82149b5e25fSSatish Balay row = idx[i]; 82249b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 82349b5e25fSSatish Balay sizes[j] = 1; 82449b5e25fSSatish Balay i++; 82549b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 82649b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 82749b5e25fSSatish Balay i++; 82849b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 82949b5e25fSSatish Balay flg = PETSC_TRUE; 83049b5e25fSSatish Balay for (k=1; k<bs; k++) { 83149b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 83249b5e25fSSatish Balay flg = PETSC_FALSE; 83349b5e25fSSatish Balay break; 83449b5e25fSSatish Balay } 83549b5e25fSSatish Balay } 83649b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 83749b5e25fSSatish Balay sizes[j] = bs; 83849b5e25fSSatish Balay i+= bs; 83949b5e25fSSatish Balay } else { 84049b5e25fSSatish Balay sizes[j] = 1; 84149b5e25fSSatish Balay i++; 84249b5e25fSSatish Balay } 84349b5e25fSSatish Balay } 84449b5e25fSSatish Balay } 84549b5e25fSSatish Balay *bs_max = j; 84649b5e25fSSatish Balay PetscFunctionReturn(0); 84749b5e25fSSatish Balay } 84849b5e25fSSatish Balay 8494a2ae208SSatish Balay #undef __FUNCT__ 8504a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 85187828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag) 85249b5e25fSSatish Balay { 85349b5e25fSSatish Balay PetscFunctionBegin; 854c0f24835SHong Zhang SETERRQ(PETSC_ERR_SUP,"No support for this function yet"); 85549b5e25fSSatish Balay } 85649b5e25fSSatish Balay 85749b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 85849b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 85949b5e25fSSatish Balay */ 86049b5e25fSSatish Balay 8614a2ae208SSatish Balay #undef __FUNCT__ 8624a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 86387828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 86449b5e25fSSatish Balay { 86549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 86649b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 86749b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 86849b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 86949b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 87049b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 87149b5e25fSSatish Balay 87249b5e25fSSatish Balay PetscFunctionBegin; 87349b5e25fSSatish Balay 87449b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 87549b5e25fSSatish Balay row = im[k]; /* row number */ 87649b5e25fSSatish Balay brow = row/bs; /* block row number */ 87749b5e25fSSatish Balay if (row < 0) continue; 87849b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 879c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 88049b5e25fSSatish Balay #endif 88149b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 88249b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 88349b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 88449b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 88549b5e25fSSatish Balay low = 0; 88649b5e25fSSatish Balay 88749b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 88849b5e25fSSatish Balay if (in[l] < 0) continue; 88949b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 890c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 89149b5e25fSSatish Balay #endif 89249b5e25fSSatish Balay col = in[l]; 89349b5e25fSSatish Balay bcol = col/bs; /* block col number */ 89449b5e25fSSatish Balay 89549b5e25fSSatish Balay if (brow <= bcol){ 89649b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8978549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 89849b5e25fSSatish Balay /* element value a(k,l) */ 89949b5e25fSSatish Balay if (roworiented) { 90049b5e25fSSatish Balay value = v[l + k*n]; 90149b5e25fSSatish Balay } else { 90249b5e25fSSatish Balay value = v[k + l*m]; 90349b5e25fSSatish Balay } 90449b5e25fSSatish Balay 90549b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 90649b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 90749b5e25fSSatish Balay while (high-low > 7) { 90849b5e25fSSatish Balay t = (low+high)/2; 90949b5e25fSSatish Balay if (rp[t] > bcol) high = t; 91049b5e25fSSatish Balay else low = t; 91149b5e25fSSatish Balay } 91249b5e25fSSatish Balay for (i=low; i<high; i++) { 91349b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 91449b5e25fSSatish Balay if (rp[i] > bcol) break; 91549b5e25fSSatish Balay if (rp[i] == bcol) { 91649b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 91749b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 91849b5e25fSSatish Balay else *bap = value; 9198549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 9208549e402SHong Zhang if (brow == bcol && ridx < cidx){ 9218549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 9228549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 9238549e402SHong Zhang else *bap = value; 9248549e402SHong Zhang } 92549b5e25fSSatish Balay goto noinsert1; 92649b5e25fSSatish Balay } 92749b5e25fSSatish Balay } 92849b5e25fSSatish Balay 92949b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 93029bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 93149b5e25fSSatish Balay if (nrow >= rmax) { 93249b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 93349b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 93449b5e25fSSatish Balay MatScalar *new_a; 93549b5e25fSSatish Balay 93629bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 93749b5e25fSSatish Balay 93849b5e25fSSatish Balay /* Malloc new storage space */ 93949b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 94082502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 94149b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 94249b5e25fSSatish Balay new_i = new_j + new_nz; 94349b5e25fSSatish Balay 94449b5e25fSSatish Balay /* copy over old data into new slots */ 94549b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 94649b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 94749b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 94849b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 94949b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 95049b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 95149b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 95249b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 95349b5e25fSSatish Balay /* free up old matrix storage */ 95449b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 95549b5e25fSSatish Balay if (!a->singlemalloc) { 95649b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 95749b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 95849b5e25fSSatish Balay } 95949b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 96049b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 96149b5e25fSSatish Balay 96249b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 96349b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 964b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 96549b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 96649b5e25fSSatish Balay a->reallocs++; 96749b5e25fSSatish Balay a->s_nz++; 96849b5e25fSSatish Balay } 96949b5e25fSSatish Balay 97049b5e25fSSatish Balay N = nrow++ - 1; 97149b5e25fSSatish Balay /* shift up all the later entries in this row */ 97249b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 97349b5e25fSSatish Balay rp[ii+1] = rp[ii]; 97449b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 97549b5e25fSSatish Balay } 97649b5e25fSSatish Balay if (N>=i) { 97749b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 97849b5e25fSSatish Balay } 97949b5e25fSSatish Balay rp[i] = bcol; 98049b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 98149b5e25fSSatish Balay noinsert1:; 98249b5e25fSSatish Balay low = i; 98349b5e25fSSatish Balay /* } */ 9848549e402SHong Zhang } 98549b5e25fSSatish Balay } /* end of if .. if.. */ 98649b5e25fSSatish Balay } /* end of loop over added columns */ 98749b5e25fSSatish Balay ailen[brow] = nrow; 98849b5e25fSSatish Balay } /* end of loop over added rows */ 98949b5e25fSSatish Balay 99049b5e25fSSatish Balay PetscFunctionReturn(0); 99149b5e25fSSatish Balay } 99249b5e25fSSatish Balay 993*15e8a5b3SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 994*15e8a5b3SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 99549b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 99649b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 99749b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 99849b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 99949b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 100087828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat); 100149b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 100249b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 100349b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 100449b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 100549b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 100649b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 100713a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 100849b5e25fSSatish Balay 100949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 101049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 101149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 101249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 101349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 101449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 101549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 101649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 101749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 101849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 101949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 102049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 102149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 102249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 102349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 102449b5e25fSSatish Balay 1025d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 1026d59c15a7SBarry Smith 102707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 102807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 102907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 103007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 103107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 103207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 103307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 103407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 103507f98182SSatish Balay 10365f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 10375f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 10385f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 10395f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 10405f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 10415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 10425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 104357d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 10443e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*); 104549b5e25fSSatish Balay 104649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 104749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 104849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 104949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 105049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 105149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 105249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 105349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 105449b5e25fSSatish Balay 105549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 105649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 105749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 105849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 105949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 106049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 106149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 106249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 106349b5e25fSSatish Balay 10644d101231SSatish Balay #ifdef HAVE_MatICCFactor 10651a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 10664a2ae208SSatish Balay #undef __FUNCT__ 10674d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 1068*15e8a5b3SHong Zhang int MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 106949b5e25fSSatish Balay { 10704ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 107149b5e25fSSatish Balay Mat outA; 107249b5e25fSSatish Balay int ierr; 107349b5e25fSSatish Balay PetscTruth row_identity,col_identity; 107449b5e25fSSatish Balay 107549b5e25fSSatish Balay PetscFunctionBegin; 10761a3463dfSHong Zhang /* 107729bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 107849b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 107949b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 108049b5e25fSSatish Balay if (!row_identity || !col_identity) { 108129bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 108249b5e25fSSatish Balay } 10831a3463dfSHong Zhang */ 108449b5e25fSSatish Balay 108549b5e25fSSatish Balay outA = inA; 10861260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 108749b5e25fSSatish Balay 108849b5e25fSSatish Balay if (!a->diag) { 10891a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 109049b5e25fSSatish Balay } 109149b5e25fSSatish Balay /* 109249b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 109349b5e25fSSatish Balay for ILU(0) factorization with natural ordering 109449b5e25fSSatish Balay */ 109549b5e25fSSatish Balay switch (a->bs) { 109649b5e25fSSatish Balay case 1: 10970fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 109807f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1099d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 11004d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 110149b5e25fSSatish Balay case 2: 11021a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 11031a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 11041a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 11054d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 110649b5e25fSSatish Balay break; 110749b5e25fSSatish Balay case 3: 11081a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 11091a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 11101a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 11114d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 111249b5e25fSSatish Balay break; 111349b5e25fSSatish Balay case 4: 11141a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 11151a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 11161a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 11174d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 111849b5e25fSSatish Balay break; 111949b5e25fSSatish Balay case 5: 11201a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 11211a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 11221a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 11234d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 112449b5e25fSSatish Balay break; 112549b5e25fSSatish Balay case 6: 11261a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 11271a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 11281a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 11294d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 113049b5e25fSSatish Balay break; 113149b5e25fSSatish Balay case 7: 11321a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 11331a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 11341a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 11354d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 113649b5e25fSSatish Balay break; 113749b5e25fSSatish Balay default: 113849b5e25fSSatish Balay a->row = row; 11391a3463dfSHong Zhang a->icol = col; 114049b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 114149b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 114249b5e25fSSatish Balay 114349b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 114449b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1145b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 114649b5e25fSSatish Balay 114749b5e25fSSatish Balay if (!a->solve_work) { 114887828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 114987828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 115049b5e25fSSatish Balay } 115149b5e25fSSatish Balay } 115249b5e25fSSatish Balay 11531a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 115449b5e25fSSatish Balay 115549b5e25fSSatish Balay PetscFunctionReturn(0); 115649b5e25fSSatish Balay } 1157950f1e5bSHong Zhang #endif 1158950f1e5bSHong Zhang 11594a2ae208SSatish Balay #undef __FUNCT__ 11604a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 116149b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 116249b5e25fSSatish Balay { 116349b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 116449b5e25fSSatish Balay MPI_Comm comm = A->comm; 116549b5e25fSSatish Balay int ierr; 116649b5e25fSSatish Balay 116749b5e25fSSatish Balay PetscFunctionBegin; 116849b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 116949b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 117049b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 117149b5e25fSSatish Balay PetscFunctionReturn(0); 117249b5e25fSSatish Balay } 117349b5e25fSSatish Balay 117449b5e25fSSatish Balay EXTERN_C_BEGIN 11754a2ae208SSatish Balay #undef __FUNCT__ 11764a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 117749b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 117849b5e25fSSatish Balay { 1179045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 118049b5e25fSSatish Balay int i,nz,n; 118149b5e25fSSatish Balay 118249b5e25fSSatish Balay PetscFunctionBegin; 1183045c9aa0SHong Zhang nz = baij->s_maxnz; 1184c464158bSHong Zhang n = mat->n; 118549b5e25fSSatish Balay for (i=0; i<nz; i++) { 118649b5e25fSSatish Balay baij->j[i] = indices[i]; 118749b5e25fSSatish Balay } 1188045c9aa0SHong Zhang baij->s_nz = nz; 118949b5e25fSSatish Balay for (i=0; i<n; i++) { 119049b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 119149b5e25fSSatish Balay } 119249b5e25fSSatish Balay 119349b5e25fSSatish Balay PetscFunctionReturn(0); 119449b5e25fSSatish Balay } 119549b5e25fSSatish Balay EXTERN_C_END 119649b5e25fSSatish Balay 11974a2ae208SSatish Balay #undef __FUNCT__ 11984a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 119949b5e25fSSatish Balay /*@ 120019585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 120149b5e25fSSatish Balay in the matrix. 120249b5e25fSSatish Balay 120349b5e25fSSatish Balay Input Parameters: 120419585528SSatish Balay + mat - the SeqSBAIJ matrix 120549b5e25fSSatish Balay - indices - the column indices 120649b5e25fSSatish Balay 120749b5e25fSSatish Balay Level: advanced 120849b5e25fSSatish Balay 120949b5e25fSSatish Balay Notes: 121049b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 121149b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 121249b5e25fSSatish Balay of the MatSetValues() operation. 121349b5e25fSSatish Balay 121449b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 121519585528SSatish Balay MatCreateSeqSBAIJ(). 121649b5e25fSSatish Balay 1217ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 121849b5e25fSSatish Balay 1219ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 122049b5e25fSSatish Balay @*/ 122149b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 122249b5e25fSSatish Balay { 122349b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 122449b5e25fSSatish Balay 122549b5e25fSSatish Balay PetscFunctionBegin; 122649b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1227c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 122849b5e25fSSatish Balay if (f) { 122949b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 123049b5e25fSSatish Balay } else { 123129bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 123249b5e25fSSatish Balay } 123349b5e25fSSatish Balay PetscFunctionReturn(0); 123449b5e25fSSatish Balay } 123549b5e25fSSatish Balay 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1238273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1239273d9f13SBarry Smith { 1240273d9f13SBarry Smith int ierr; 1241273d9f13SBarry Smith 1242273d9f13SBarry Smith PetscFunctionBegin; 1243273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1244273d9f13SBarry Smith PetscFunctionReturn(0); 1245273d9f13SBarry Smith } 1246273d9f13SBarry Smith 1247a6ece127SHong Zhang #undef __FUNCT__ 1248a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 1249a6ece127SHong Zhang int MatGetArray_SeqSBAIJ(Mat A,PetscScalar **array) 1250a6ece127SHong Zhang { 1251a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1252a6ece127SHong Zhang PetscFunctionBegin; 1253a6ece127SHong Zhang *array = a->a; 1254a6ece127SHong Zhang PetscFunctionReturn(0); 1255a6ece127SHong Zhang } 1256a6ece127SHong Zhang 1257a6ece127SHong Zhang #undef __FUNCT__ 1258a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 1259a6ece127SHong Zhang int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar **array) 1260a6ece127SHong Zhang { 1261a6ece127SHong Zhang PetscFunctionBegin; 1262a6ece127SHong Zhang PetscFunctionReturn(0); 1263a6ece127SHong Zhang } 1264a6ece127SHong Zhang 126542ee4b1aSHong Zhang #include "petscblaslapack.h" 126642ee4b1aSHong Zhang #undef __FUNCT__ 126742ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 126842ee4b1aSHong Zhang int MatAXPY_SeqSBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str) 126942ee4b1aSHong Zhang { 127042ee4b1aSHong Zhang int ierr,one=1; 127142ee4b1aSHong Zhang Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ *)X->data,*y = (Mat_SeqSBAIJ *)Y->data; 127242ee4b1aSHong Zhang 127342ee4b1aSHong Zhang PetscFunctionBegin; 127442ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 127542ee4b1aSHong Zhang BLaxpy_(&x->s_nz,a,x->a,&one,y->a,&one); 127642ee4b1aSHong Zhang } else { 127742ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 127842ee4b1aSHong Zhang } 127942ee4b1aSHong Zhang PetscFunctionReturn(0); 128042ee4b1aSHong Zhang } 128142ee4b1aSHong Zhang 128249b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 128349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 128449b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 128549b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 128649b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 128749b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 128849b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 128949b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 129049b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 129149b5e25fSSatish Balay 0, 129249b5e25fSSatish Balay 0, 129349b5e25fSSatish Balay 0, 129449b5e25fSSatish Balay 0, 12955f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1296d06b337dSHong Zhang MatRelax_SeqSBAIJ, 129749b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 129849b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 129949b5e25fSSatish Balay MatEqual_SeqSBAIJ, 130049b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 130149b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 130249b5e25fSSatish Balay MatNorm_SeqSBAIJ, 130349b5e25fSSatish Balay 0, 130449b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 130549b5e25fSSatish Balay 0, 130649b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 130749b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 130849b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 130949b5e25fSSatish Balay 0, 131049b5e25fSSatish Balay 0, 13115f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 13125f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1313273d9f13SBarry Smith MatSetUpPreallocation_SeqSBAIJ, 1314c464158bSHong Zhang 0, 13154d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1316a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1317a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 131849b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 131949b5e25fSSatish Balay 0, 132049b5e25fSSatish Balay 0, 132149b5e25fSSatish Balay 0, 1322950f1e5bSHong Zhang 0, 132342ee4b1aSHong Zhang MatAXPY_SeqSBAIJ, 132449b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 132549b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 132649b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 132749b5e25fSSatish Balay 0, 132849b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 132949b5e25fSSatish Balay MatScale_SeqSBAIJ, 133049b5e25fSSatish Balay 0, 133149b5e25fSSatish Balay 0, 133249b5e25fSSatish Balay 0, 133349b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 133449b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 133549b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 133649b5e25fSSatish Balay 0, 133749b5e25fSSatish Balay 0, 133849b5e25fSSatish Balay 0, 133949b5e25fSSatish Balay 0, 134049b5e25fSSatish Balay 0, 134149b5e25fSSatish Balay 0, 134249b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 134349b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 134449b5e25fSSatish Balay 0, 134549b5e25fSSatish Balay 0, 13468a124369SBarry Smith MatGetPetscMaps_Petsc, 1347d959ec07SHong Zhang 0, 1348d959ec07SHong Zhang 0, 1349d959ec07SHong Zhang 0, 1350d959ec07SHong Zhang 0, 1351d959ec07SHong Zhang 0, 1352d959ec07SHong Zhang 0, 13533e0d88b5SBarry Smith MatGetRowMax_SeqSBAIJ, 13543e0d88b5SBarry Smith 0, 13553e0d88b5SBarry Smith 0, 13563e0d88b5SBarry Smith 0, 13573e0d88b5SBarry Smith 0, 13583e0d88b5SBarry Smith 0, 13593e0d88b5SBarry Smith 0, 13603e0d88b5SBarry Smith 0, 13613e0d88b5SBarry Smith 0, 13623e0d88b5SBarry Smith 0, 13633e0d88b5SBarry Smith 0, 13643e0d88b5SBarry Smith 0, 13653e0d88b5SBarry Smith 0, 13663e0d88b5SBarry Smith 0, 13673e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 13683e0d88b5SBarry Smith MatGetInertia_SeqSBAIJ 13693e0d88b5SBarry Smith #else 13703e0d88b5SBarry Smith 0 13713e0d88b5SBarry Smith #endif 13723e0d88b5SBarry Smith }; 137349b5e25fSSatish Balay 137449b5e25fSSatish Balay EXTERN_C_BEGIN 13754a2ae208SSatish Balay #undef __FUNCT__ 13764a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 137749b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 137849b5e25fSSatish Balay { 13794afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1380c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 138149b5e25fSSatish Balay int ierr; 138249b5e25fSSatish Balay 138349b5e25fSSatish Balay PetscFunctionBegin; 138449b5e25fSSatish Balay if (aij->nonew != 1) { 138529bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 138649b5e25fSSatish Balay } 138749b5e25fSSatish Balay 138849b5e25fSSatish Balay /* allocate space for values if not already there */ 138949b5e25fSSatish Balay if (!aij->saved_values) { 139087828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 139149b5e25fSSatish Balay } 139249b5e25fSSatish Balay 139349b5e25fSSatish Balay /* copy values over */ 139487828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 139549b5e25fSSatish Balay PetscFunctionReturn(0); 139649b5e25fSSatish Balay } 139749b5e25fSSatish Balay EXTERN_C_END 139849b5e25fSSatish Balay 139949b5e25fSSatish Balay EXTERN_C_BEGIN 14004a2ae208SSatish Balay #undef __FUNCT__ 14014a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 140249b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 140349b5e25fSSatish Balay { 14044afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1405c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 140649b5e25fSSatish Balay 140749b5e25fSSatish Balay PetscFunctionBegin; 140849b5e25fSSatish Balay if (aij->nonew != 1) { 140929bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 141049b5e25fSSatish Balay } 141149b5e25fSSatish Balay if (!aij->saved_values) { 141229bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 141349b5e25fSSatish Balay } 141449b5e25fSSatish Balay 141549b5e25fSSatish Balay /* copy values over */ 141687828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 141749b5e25fSSatish Balay PetscFunctionReturn(0); 141849b5e25fSSatish Balay } 141949b5e25fSSatish Balay EXTERN_C_END 142049b5e25fSSatish Balay 14218549e402SHong Zhang EXTERN_C_BEGIN 14224a2ae208SSatish Balay #undef __FUNCT__ 14234a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ" 1424c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1425c464158bSHong Zhang { 1426c464158bSHong Zhang Mat_SeqSBAIJ *b; 14270c955e93SHong Zhang int ierr,size; 1428c464158bSHong Zhang 1429c464158bSHong Zhang PetscFunctionBegin; 1430c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1431c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1432273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1433273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1434c464158bSHong Zhang 1435b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1436b0a32e0cSBarry Smith B->data = (void*)b; 1437c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1438c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1439c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1440c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1441c464158bSHong Zhang B->factor = 0; 1442c464158bSHong Zhang B->lupivotthreshold = 1.0; 1443c464158bSHong Zhang B->mapping = 0; 1444c464158bSHong Zhang b->row = 0; 1445c464158bSHong Zhang b->icol = 0; 1446c464158bSHong Zhang b->reallocs = 0; 1447c464158bSHong Zhang b->saved_values = 0; 1448c464158bSHong Zhang 14498a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 14508a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1451c464158bSHong Zhang 1452c464158bSHong Zhang b->sorted = PETSC_FALSE; 1453c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1454c464158bSHong Zhang b->nonew = 0; 1455c464158bSHong Zhang b->diag = 0; 1456c464158bSHong Zhang b->solve_work = 0; 1457c464158bSHong Zhang b->mult_work = 0; 14582a1b7f2aSHong Zhang B->spptr = 0; 1459c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1460c464158bSHong Zhang 1461c464158bSHong Zhang b->inew = 0; 1462c464158bSHong Zhang b->jnew = 0; 1463c464158bSHong Zhang b->anew = 0; 1464c464158bSHong Zhang b->a2anew = 0; 1465c464158bSHong Zhang b->permute = PETSC_FALSE; 1466c464158bSHong Zhang 1467c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1468c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1469c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1470c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1471c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1472c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1473c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1474c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1475c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1476c464158bSHong Zhang PetscFunctionReturn(0); 1477c464158bSHong Zhang } 14788549e402SHong Zhang EXTERN_C_END 1479c464158bSHong Zhang 14804a2ae208SSatish Balay #undef __FUNCT__ 14814a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 148249b5e25fSSatish Balay /*@C 1483273d9f13SBarry Smith MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 148449b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 148549b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 148649b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 148749b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 148849b5e25fSSatish Balay 1489c464158bSHong Zhang Collective on Mat 149049b5e25fSSatish Balay 149149b5e25fSSatish Balay Input Parameters: 1492c464158bSHong Zhang + A - the symmetric matrix 149349b5e25fSSatish Balay . bs - size of block 149449b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 1495744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1496744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 149749b5e25fSSatish Balay 149849b5e25fSSatish Balay Options Database Keys: 149949b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 150049b5e25fSSatish Balay block calculations (much slower) 150149b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 150249b5e25fSSatish Balay 150349b5e25fSSatish Balay Level: intermediate 150449b5e25fSSatish Balay 150549b5e25fSSatish Balay Notes: 150649b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 150749b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 150849b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 150949b5e25fSSatish Balay matrices. 151049b5e25fSSatish Balay 1511a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 151249b5e25fSSatish Balay @*/ 1513273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 151449b5e25fSSatish Balay { 1515c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 15160c955e93SHong Zhang int i,len,ierr,mbs,bs2; 151749b5e25fSSatish Balay PetscTruth flg; 15184afc71dfSHong Zhang int s_nz; 151949b5e25fSSatish Balay 152049b5e25fSSatish Balay PetscFunctionBegin; 1521273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1522e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1523c464158bSHong Zhang mbs = B->m/bs; 152449b5e25fSSatish Balay bs2 = bs*bs; 152549b5e25fSSatish Balay 1526c464158bSHong Zhang if (mbs*bs != B->m) { 152729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 152849b5e25fSSatish Balay } 152949b5e25fSSatish Balay 1530435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1531435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 153249b5e25fSSatish Balay if (nnz) { 153349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 153429bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 153580fe4e49SBarry 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); 153649b5e25fSSatish Balay } 153749b5e25fSSatish Balay } 153849b5e25fSSatish Balay 1539e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 154049b5e25fSSatish Balay if (!flg) { 154149b5e25fSSatish Balay switch (bs) { 154249b5e25fSSatish Balay case 1: 15434ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 154449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1545d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 154649b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 154749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 154849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 154949b5e25fSSatish Balay break; 155049b5e25fSSatish Balay case 2: 15514ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 155249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 155349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 155449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 155549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 155649b5e25fSSatish Balay break; 155749b5e25fSSatish Balay case 3: 15585f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 155949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 156049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 156149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 156249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 156349b5e25fSSatish Balay break; 156449b5e25fSSatish Balay case 4: 15655f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 156649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 156749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 156849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 156949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 157049b5e25fSSatish Balay break; 157149b5e25fSSatish Balay case 5: 15725f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 157349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 157449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 157549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 157649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 157749b5e25fSSatish Balay break; 157849b5e25fSSatish Balay case 6: 15795f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 158049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 158149b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 158249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 158349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 158449b5e25fSSatish Balay break; 158549b5e25fSSatish Balay case 7: 1586de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 158749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 158849b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1589de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 159049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 159149b5e25fSSatish Balay break; 159249b5e25fSSatish Balay } 159349b5e25fSSatish Balay } 159449b5e25fSSatish Balay 159549b5e25fSSatish Balay b->mbs = mbs; 15964afc71dfSHong Zhang b->nbs = mbs; 1597b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 159849b5e25fSSatish Balay if (!nnz) { 1599435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 160049b5e25fSSatish Balay else if (nz <= 0) nz = 1; 160149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 16028cef66ccSHong Zhang b->imax[i] = nz; 160349b5e25fSSatish Balay } 1604153ea458SHong Zhang nz = nz*mbs; /* total nz */ 160549b5e25fSSatish Balay } else { 160649b5e25fSSatish Balay nz = 0; 16078cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 160849b5e25fSSatish Balay } 16098cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 16108cef66ccSHong Zhang s_nz = nz; 161149b5e25fSSatish Balay 161249b5e25fSSatish Balay /* allocate the matrix space */ 1613c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 161482502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 161549b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 161649b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 161749b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 161849b5e25fSSatish Balay b->i = b->j + s_nz; 161949b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 162049b5e25fSSatish Balay 162149b5e25fSSatish Balay /* pointer to beginning of each row */ 162249b5e25fSSatish Balay b->i[0] = 0; 162349b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 162449b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 162549b5e25fSSatish Balay } 162649b5e25fSSatish Balay 162749b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 16285033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1629b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 163049b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 163149b5e25fSSatish Balay 163249b5e25fSSatish Balay b->bs = bs; 163349b5e25fSSatish Balay b->bs2 = bs2; 163449b5e25fSSatish Balay b->s_nz = 0; 163549b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1636153ea458SHong Zhang 163716cdd363SHong Zhang b->inew = 0; 163816cdd363SHong Zhang b->jnew = 0; 163916cdd363SHong Zhang b->anew = 0; 164016cdd363SHong Zhang b->a2anew = 0; 16411a3463dfSHong Zhang b->permute = PETSC_FALSE; 1642c464158bSHong Zhang PetscFunctionReturn(0); 1643c464158bSHong Zhang } 1644153ea458SHong Zhang 164549b5e25fSSatish Balay 16464a2ae208SSatish Balay #undef __FUNCT__ 16474a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1648c464158bSHong Zhang /*@C 1649c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1650c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1651c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1652c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1653c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 165449b5e25fSSatish Balay 1655c464158bSHong Zhang Collective on MPI_Comm 1656c464158bSHong Zhang 1657c464158bSHong Zhang Input Parameters: 1658c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1659c464158bSHong Zhang . bs - size of block 1660c464158bSHong Zhang . m - number of rows, or number of columns 1661c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1662744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1663744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1664c464158bSHong Zhang 1665c464158bSHong Zhang Output Parameter: 1666c464158bSHong Zhang . A - the symmetric matrix 1667c464158bSHong Zhang 1668c464158bSHong Zhang Options Database Keys: 1669c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1670c464158bSHong Zhang block calculations (much slower) 1671c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1672c464158bSHong Zhang 1673c464158bSHong Zhang Level: intermediate 1674c464158bSHong Zhang 1675c464158bSHong Zhang Notes: 1676c464158bSHong Zhang 1677c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1678c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1679c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1680c464158bSHong Zhang matrices. 1681c464158bSHong Zhang 1682c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1683c464158bSHong Zhang @*/ 1684c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1685c464158bSHong Zhang { 1686c464158bSHong Zhang int ierr; 1687c464158bSHong Zhang 1688c464158bSHong Zhang PetscFunctionBegin; 1689c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1690c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1691273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 169249b5e25fSSatish Balay PetscFunctionReturn(0); 169349b5e25fSSatish Balay } 169449b5e25fSSatish Balay 16954a2ae208SSatish Balay #undef __FUNCT__ 16964a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 169749b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 169849b5e25fSSatish Balay { 169949b5e25fSSatish Balay Mat C; 170049b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 170149b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 170249b5e25fSSatish Balay 170349b5e25fSSatish Balay PetscFunctionBegin; 170429bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 170549b5e25fSSatish Balay 170649b5e25fSSatish Balay *B = 0; 1707692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1708692f9cbeSHong Zhang ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 1709692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1710692f9cbeSHong Zhang 171149b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1712273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 171349b5e25fSSatish Balay C->factor = A->factor; 171449b5e25fSSatish Balay c->row = 0; 171549b5e25fSSatish Balay c->icol = 0; 171649b5e25fSSatish Balay c->saved_values = 0; 171749b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 171849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 171949b5e25fSSatish Balay 172049b5e25fSSatish Balay c->bs = a->bs; 172149b5e25fSSatish Balay c->bs2 = a->bs2; 172249b5e25fSSatish Balay c->mbs = a->mbs; 172349b5e25fSSatish Balay c->nbs = a->nbs; 172449b5e25fSSatish Balay 1725b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1726b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 172749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 172849b5e25fSSatish Balay c->imax[i] = a->imax[i]; 172949b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 173049b5e25fSSatish Balay } 173149b5e25fSSatish Balay 173249b5e25fSSatish Balay /* allocate the matrix space */ 173349b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 173449b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 173582502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 173649b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 173749b5e25fSSatish Balay c->i = c->j + nz; 173849b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 173949b5e25fSSatish Balay if (mbs > 0) { 174049b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 174149b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 174249b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 174349b5e25fSSatish Balay } else { 174449b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 174549b5e25fSSatish Balay } 174649b5e25fSSatish Balay } 174749b5e25fSSatish Balay 1748b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 174949b5e25fSSatish Balay c->sorted = a->sorted; 175049b5e25fSSatish Balay c->roworiented = a->roworiented; 175149b5e25fSSatish Balay c->nonew = a->nonew; 175249b5e25fSSatish Balay 175349b5e25fSSatish Balay if (a->diag) { 1754b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1755b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 175649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 175749b5e25fSSatish Balay c->diag[i] = a->diag[i]; 175849b5e25fSSatish Balay } 175949b5e25fSSatish Balay } else c->diag = 0; 176049b5e25fSSatish Balay c->s_nz = a->s_nz; 176149b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 176249b5e25fSSatish Balay c->solve_work = 0; 17632a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 176449b5e25fSSatish Balay c->mult_work = 0; 176549b5e25fSSatish Balay *B = C; 1766b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 176749b5e25fSSatish Balay PetscFunctionReturn(0); 176849b5e25fSSatish Balay } 176949b5e25fSSatish Balay 17708549e402SHong Zhang EXTERN_C_BEGIN 17714a2ae208SSatish Balay #undef __FUNCT__ 17724a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1773b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A) 177449b5e25fSSatish Balay { 177549b5e25fSSatish Balay Mat_SeqSBAIJ *a; 177649b5e25fSSatish Balay Mat B; 177749b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 17783f7326a9SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 177949b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 178049b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 178187828ca2SBarry Smith PetscScalar *aa; 178249b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 178327fa2452SHong Zhang #if defined(PETSC_HAVE_SPOOLES) 1784ecc4ab6dSBarry Smith PetscTruth flag; 1785ecc4ab6dSBarry Smith #endif 178649b5e25fSSatish Balay 178749b5e25fSSatish Balay PetscFunctionBegin; 1788b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 178949b5e25fSSatish Balay bs2 = bs*bs; 179049b5e25fSSatish Balay 179149b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 179229bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1793b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 179449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1795552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 179649b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 179749b5e25fSSatish Balay 179849b5e25fSSatish Balay if (header[3] < 0) { 179929bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 180049b5e25fSSatish Balay } 180149b5e25fSSatish Balay 180229bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 180349b5e25fSSatish Balay 180449b5e25fSSatish Balay /* 180549b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 180649b5e25fSSatish Balay divisible by the blocksize 180749b5e25fSSatish Balay */ 180849b5e25fSSatish Balay mbs = M/bs; 180949b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 181049b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 181149b5e25fSSatish Balay else mbs++; 181249b5e25fSSatish Balay if (extra_rows) { 1813b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 181449b5e25fSSatish Balay } 181549b5e25fSSatish Balay 181649b5e25fSSatish Balay /* read in row lengths */ 1817b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 181849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 181949b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 182049b5e25fSSatish Balay 182149b5e25fSSatish Balay /* read in column indices */ 1822b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 182349b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 182449b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 182549b5e25fSSatish Balay 182649b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 182782502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 1828a91a614fSBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 182982502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 183049b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 183149b5e25fSSatish Balay masked = mask + mbs; 183249b5e25fSSatish Balay rowcount = 0; nzcount = 0; 183349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 183449b5e25fSSatish Balay nmask = 0; 183549b5e25fSSatish Balay for (j=0; j<bs; j++) { 183649b5e25fSSatish Balay kmax = rowlengths[rowcount]; 183749b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18382d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 183903630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 184049b5e25fSSatish Balay } 184149b5e25fSSatish Balay rowcount++; 184249b5e25fSSatish Balay } 1843574b2666SHong Zhang s_browlengths[i] += nmask; 1844574b2666SHong Zhang 184549b5e25fSSatish Balay /* zero out the mask elements we set */ 184649b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 184749b5e25fSSatish Balay } 184849b5e25fSSatish Balay 184949b5e25fSSatish Balay /* create our matrix */ 18507fe2be48SHong Zhang ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr); 185149b5e25fSSatish Balay B = *A; 185249b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 185349b5e25fSSatish Balay 185449b5e25fSSatish Balay /* set matrix "i" values */ 185549b5e25fSSatish Balay a->i[0] = 0; 185649b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1857574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1858574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 185949b5e25fSSatish Balay } 18607fe2be48SHong Zhang a->s_nz = a->i[mbs]; 186149b5e25fSSatish Balay 186249b5e25fSSatish Balay /* read in nonzero values */ 186387828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 186449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 186549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 186649b5e25fSSatish Balay 186749b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 186849b5e25fSSatish Balay nzcount = 0; jcount = 0; 186949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 187049b5e25fSSatish Balay nzcountb = nzcount; 187149b5e25fSSatish Balay nmask = 0; 187249b5e25fSSatish Balay for (j=0; j<bs; j++) { 187349b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 187449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18752d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 187603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 18772d703238SHong Zhang } 18782d703238SHong Zhang } 18792d703238SHong Zhang /* sort the masked values */ 18802d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 18812d703238SHong Zhang 18822d703238SHong Zhang /* set "j" values into matrix */ 18832d703238SHong Zhang maskcount = 1; 18842d703238SHong Zhang for (j=0; j<nmask; j++) { 188549b5e25fSSatish Balay a->j[jcount++] = masked[j]; 188649b5e25fSSatish Balay mask[masked[j]] = maskcount++; 188749b5e25fSSatish Balay } 1888574b2666SHong Zhang 188949b5e25fSSatish Balay /* set "a" values into matrix */ 189049b5e25fSSatish Balay ishift = bs2*a->i[i]; 189149b5e25fSSatish Balay for (j=0; j<bs; j++) { 189249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 189349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1894574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1895574b2666SHong Zhang if (tmp >= i){ 189649b5e25fSSatish Balay block = mask[tmp] - 1; 189749b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 189849b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1899574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1900574b2666SHong Zhang } 1901574b2666SHong Zhang nzcountb++; 190249b5e25fSSatish Balay } 190349b5e25fSSatish Balay } 190449b5e25fSSatish Balay /* zero out the mask elements we set */ 190549b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 190649b5e25fSSatish Balay } 190729bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 190849b5e25fSSatish Balay 190949b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1910574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 191149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 191249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 191349b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 191449b5e25fSSatish Balay 191549b5e25fSSatish Balay B->assembled = PETSC_TRUE; 191627fa2452SHong Zhang #if defined(PETSC_HAVE_SPOOLES) 191727fa2452SHong Zhang ierr = PetscOptionsHasName(B->prefix,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr); 191827fa2452SHong Zhang if (flag) { ierr = MatUseSpooles_SeqSBAIJ(B);CHKERRQ(ierr); } 1919ecc4ab6dSBarry Smith #endif 192049b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 192149b5e25fSSatish Balay PetscFunctionReturn(0); 192249b5e25fSSatish Balay } 19238549e402SHong Zhang EXTERN_C_END 1924574b2666SHong Zhang 1925d06b337dSHong Zhang #undef __FUNCT__ 1926d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1927c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx) 1928d06b337dSHong Zhang { 1929d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1930d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1931d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 1932d06b337dSHong Zhang int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr; 1933d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1934d06b337dSHong Zhang 1935d06b337dSHong Zhang PetscFunctionBegin; 1936b965ef7fSBarry Smith its = its*lits; 193791723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 1938d06b337dSHong Zhang 1939d06b337dSHong Zhang if (bs > 1) 1940d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1941d06b337dSHong Zhang 1942d06b337dSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1943d06b337dSHong Zhang if (xx != bb) { 1944d06b337dSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1945d06b337dSHong Zhang } else { 1946d06b337dSHong Zhang b = x; 1947d06b337dSHong Zhang } 1948d06b337dSHong Zhang 1949d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1950d06b337dSHong Zhang 1951d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19522798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1953d06b337dSHong Zhang for (i=0; i<m; i++) 1954d06b337dSHong Zhang t[i] = b[i]; 1955d06b337dSHong Zhang 1956d06b337dSHong Zhang for (i=0; i<m; i++){ 195744706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1958d06b337dSHong Zhang v = aa + ai[i] + 1; 1959d06b337dSHong Zhang vj = aj + ai[i] + 1; 1960d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1961d06b337dSHong Zhang x[i] = omega*t[i]/d; 1962d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 196344706875SHong Zhang PetscLogFlops(2*nz-1); 1964d06b337dSHong Zhang } 1965d06b337dSHong Zhang } 1966d06b337dSHong Zhang 19672798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1968d06b337dSHong Zhang for (i=0; i<m; i++) 1969d06b337dSHong Zhang t[i] = b[i]; 1970d06b337dSHong Zhang 1971d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1972d06b337dSHong Zhang v = aa + ai[i] + 1; 1973d06b337dSHong Zhang vj = aj + ai[i] + 1; 1974d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1975d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 197644706875SHong Zhang PetscLogFlops(2*nz-1); 1977d06b337dSHong Zhang } 1978d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1979d06b337dSHong Zhang d = *(aa + ai[i]); 1980d06b337dSHong Zhang v = aa + ai[i] + 1; 1981d06b337dSHong Zhang vj = aj + ai[i] + 1; 1982d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1983d06b337dSHong Zhang sum = t[i]; 1984d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 198544706875SHong Zhang PetscLogFlops(2*nz-1); 1986d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1987d06b337dSHong Zhang } 1988d06b337dSHong Zhang } 1989d06b337dSHong Zhang its--; 1990d06b337dSHong Zhang } 1991d06b337dSHong Zhang 1992d06b337dSHong Zhang while (its--) { 199344706875SHong Zhang /* 199444706875SHong Zhang forward sweep: 199544706875SHong Zhang for i=0,...,m-1: 199644706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 199744706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 199844706875SHong Zhang b = b - x[i]*U^T(i,:); 1999d06b337dSHong Zhang 200044706875SHong Zhang */ 20012798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2002d06b337dSHong Zhang for (i=0; i<m; i++) 2003d06b337dSHong Zhang t[i] = b[i]; 2004d06b337dSHong Zhang 2005d06b337dSHong Zhang for (i=0; i<m; i++){ 200644706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2007d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2008d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2009d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2010d06b337dSHong Zhang sum = t[i]; 2011d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2012d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2013d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 201444706875SHong Zhang PetscLogFlops(4*nz-2); 2015d06b337dSHong Zhang } 2016d06b337dSHong Zhang } 2017d06b337dSHong Zhang 20182798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 201944706875SHong Zhang /* 202044706875SHong Zhang backward sweep: 202144706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 202244706875SHong Zhang for i=m-1,...,0: 202344706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 202444706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 202544706875SHong Zhang */ 2026d06b337dSHong Zhang for (i=0; i<m; i++) 2027d06b337dSHong Zhang t[i] = b[i]; 2028d06b337dSHong Zhang 2029d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2030d06b337dSHong Zhang v = aa + ai[i] + 1; 2031d06b337dSHong Zhang vj = aj + ai[i] + 1; 2032d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2033d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 203444706875SHong Zhang PetscLogFlops(2*nz-1); 2035d06b337dSHong Zhang } 2036d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2037d06b337dSHong Zhang d = *(aa + ai[i]); 2038d06b337dSHong Zhang v = aa + ai[i] + 1; 2039d06b337dSHong Zhang vj = aj + ai[i] + 1; 2040d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2041d06b337dSHong Zhang sum = t[i]; 2042d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 204344706875SHong Zhang PetscLogFlops(2*nz-1); 2044d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2045d06b337dSHong Zhang } 2046d06b337dSHong Zhang } 2047d06b337dSHong Zhang } 2048d06b337dSHong Zhang 2049d06b337dSHong Zhang ierr = PetscFree(t); CHKERRQ(ierr); 2050d06b337dSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2051d06b337dSHong Zhang if (bb != xx) { 2052d06b337dSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2053d06b337dSHong Zhang } 20542798e883SHong Zhang 2055d06b337dSHong Zhang PetscFunctionReturn(0); 2056d06b337dSHong Zhang } 2057d06b337dSHong Zhang 2058d06b337dSHong Zhang 2059d06b337dSHong Zhang 2060d06b337dSHong Zhang 206149b5e25fSSatish Balay 206249b5e25fSSatish Balay 2063