173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 449b5e25fSSatish Balay Defines the basic matrix operations for the BAIJ (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" 11*cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 12*cf4441caSHong Zhang EXTERN int MatUseSpooles_SeqSBAIJ(Mat); 13*cf4441caSHong Zhang #endif 1449b5e25fSSatish Balay 1549b5e25fSSatish Balay #define CHUNKSIZE 10 1649b5e25fSSatish Balay 1749b5e25fSSatish Balay /* 1849b5e25fSSatish Balay Checks for missing diagonals 1949b5e25fSSatish Balay */ 204a2ae208SSatish Balay #undef __FUNCT__ 214a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 2249b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2349b5e25fSSatish Balay { 24045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2549b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2649b5e25fSSatish Balay 2749b5e25fSSatish Balay PetscFunctionBegin; 28045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2949b5e25fSSatish Balay diag = a->diag; 3049b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31e11b0e14SHong Zhang if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i); 3249b5e25fSSatish Balay } 3349b5e25fSSatish Balay PetscFunctionReturn(0); 3449b5e25fSSatish Balay } 3549b5e25fSSatish Balay 364a2ae208SSatish Balay #undef __FUNCT__ 374a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 3849b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3949b5e25fSSatish Balay { 40045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 41b424e231SHong Zhang int i,mbs = a->mbs,ierr; 4249b5e25fSSatish Balay 4349b5e25fSSatish Balay PetscFunctionBegin; 4449b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4549b5e25fSSatish Balay 46b424e231SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr); 47b424e231SHong Zhang PetscLogObjectMemory(A,(mbs+1)*sizeof(int)); 48b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4949b5e25fSSatish Balay PetscFunctionReturn(0); 5049b5e25fSSatish Balay } 5149b5e25fSSatish Balay 5249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 5349b5e25fSSatish Balay 544a2ae208SSatish Balay #undef __FUNCT__ 554a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 5649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 5749b5e25fSSatish Balay { 5849b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 5949b5e25fSSatish Balay 6049b5e25fSSatish Balay PetscFunctionBegin; 61045c9aa0SHong Zhang if (ia) { 6229bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 6349b5e25fSSatish Balay } 64045c9aa0SHong Zhang *nn = a->mbs; 6549b5e25fSSatish Balay PetscFunctionReturn(0); 6649b5e25fSSatish Balay } 6749b5e25fSSatish Balay 684a2ae208SSatish Balay #undef __FUNCT__ 694a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 70045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 7149b5e25fSSatish Balay { 7249b5e25fSSatish Balay PetscFunctionBegin; 7349b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 7429bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 7516cdd363SHong Zhang /* PetscFunctionReturn(0); */ 7649b5e25fSSatish Balay } 7749b5e25fSSatish Balay 784a2ae208SSatish Balay #undef __FUNCT__ 794a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 8049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 8149b5e25fSSatish Balay { 8249b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 8349b5e25fSSatish Balay 8449b5e25fSSatish Balay PetscFunctionBegin; 8549b5e25fSSatish Balay *bs = sbaij->bs; 8649b5e25fSSatish Balay PetscFunctionReturn(0); 8749b5e25fSSatish Balay } 8849b5e25fSSatish Balay 894a2ae208SSatish Balay #undef __FUNCT__ 904a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 9149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 9249b5e25fSSatish Balay { 9349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 9449b5e25fSSatish Balay int ierr; 9549b5e25fSSatish Balay 9649b5e25fSSatish Balay PetscFunctionBegin; 9749b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 98b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 9949b5e25fSSatish Balay #endif 10049b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 10149b5e25fSSatish Balay if (!a->singlemalloc) { 10249b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 10363c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 10449b5e25fSSatish Balay } 10549b5e25fSSatish Balay if (a->row) { 10649b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 10749b5e25fSSatish Balay } 10849b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 10949b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 11049b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 11149b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 11249b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 11349b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 11449b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1151a3463dfSHong Zhang 1161a3463dfSHong Zhang if (a->inew){ 1171a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1181a3463dfSHong Zhang a->inew = 0; 1191a3463dfSHong Zhang } 12049b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 12149b5e25fSSatish Balay PetscFunctionReturn(0); 12249b5e25fSSatish Balay } 12349b5e25fSSatish Balay 1244a2ae208SSatish Balay #undef __FUNCT__ 1254a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 12649b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 12749b5e25fSSatish Balay { 128045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 12949b5e25fSSatish Balay 13049b5e25fSSatish Balay PetscFunctionBegin; 1314d9d31abSKris Buschelman switch (op) { 1324d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1334d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1344d9d31abSKris Buschelman break; 1354d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1364d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1374d9d31abSKris Buschelman break; 1384d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1394d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1404d9d31abSKris Buschelman break; 1414d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1424d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1434d9d31abSKris Buschelman break; 1444d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1454d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1464d9d31abSKris Buschelman break; 1474d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1484d9d31abSKris Buschelman a->nonew = 1; 1494d9d31abSKris Buschelman break; 1504d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1514d9d31abSKris Buschelman a->nonew = -1; 1524d9d31abSKris Buschelman break; 1534d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1544d9d31abSKris Buschelman a->nonew = -2; 1554d9d31abSKris Buschelman break; 1564d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1574d9d31abSKris Buschelman a->nonew = 0; 1584d9d31abSKris Buschelman break; 1594d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1604d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1614d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1624d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1634d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 164d03495bdSKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 165b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1664d9d31abSKris Buschelman break; 1674d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 16829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1694d9d31abSKris Buschelman default: 17029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 17149b5e25fSSatish Balay } 17249b5e25fSSatish Balay PetscFunctionReturn(0); 17349b5e25fSSatish Balay } 17449b5e25fSSatish Balay 1754a2ae208SSatish Balay #undef __FUNCT__ 1764a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 17787828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v) 17849b5e25fSSatish Balay { 17949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 18082502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr; 18149b5e25fSSatish Balay MatScalar *aa,*aa_i; 18287828ca2SBarry Smith PetscScalar *v_i; 18349b5e25fSSatish Balay 18449b5e25fSSatish Balay PetscFunctionBegin; 18549b5e25fSSatish Balay bs = a->bs; 18649b5e25fSSatish Balay ai = a->i; 18749b5e25fSSatish Balay aj = a->j; 18849b5e25fSSatish Balay aa = a->a; 18949b5e25fSSatish Balay bs2 = a->bs2; 19049b5e25fSSatish Balay 191c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 19249b5e25fSSatish Balay 19349b5e25fSSatish Balay bn = row/bs; /* Block number */ 19449b5e25fSSatish Balay bp = row % bs; /* Block position */ 19549b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 19649b5e25fSSatish Balay *ncols = bs*M; 19749b5e25fSSatish Balay 19849b5e25fSSatish Balay if (v) { 19949b5e25fSSatish Balay *v = 0; 20049b5e25fSSatish Balay if (*ncols) { 20187828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 20249b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 20349b5e25fSSatish Balay v_i = *v + i*bs; 20449b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 20549b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 20649b5e25fSSatish Balay } 20749b5e25fSSatish Balay } 20849b5e25fSSatish Balay } 20949b5e25fSSatish Balay 21049b5e25fSSatish Balay if (cols) { 21149b5e25fSSatish Balay *cols = 0; 21249b5e25fSSatish Balay if (*ncols) { 21382502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 21449b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 21549b5e25fSSatish Balay cols_i = *cols + i*bs; 21649b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 21749b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 21849b5e25fSSatish Balay } 21949b5e25fSSatish Balay } 22049b5e25fSSatish Balay } 22149b5e25fSSatish Balay 22249b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2235ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2245ddb2528SHong Zhang #ifdef column_search 22549b5e25fSSatish Balay v_i = *v + M*bs; 22649b5e25fSSatish Balay cols_i = *cols + M*bs; 22749b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 22849b5e25fSSatish Balay M = ai[i+1] - ai[i]; 22949b5e25fSSatish Balay for (j=0; j<M; j++){ 23049b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 23149b5e25fSSatish Balay if (itmp == bn){ 23249b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 23349b5e25fSSatish Balay for (k=0; k<bs; k++) { 23449b5e25fSSatish Balay *cols_i++ = i*bs+k; 23549b5e25fSSatish Balay *v_i++ = aa_i[k]; 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay *ncols += bs; 23849b5e25fSSatish Balay break; 23949b5e25fSSatish Balay } 24049b5e25fSSatish Balay } 24149b5e25fSSatish Balay } 2425ddb2528SHong Zhang #endif 24349b5e25fSSatish Balay 24449b5e25fSSatish Balay PetscFunctionReturn(0); 24549b5e25fSSatish Balay } 24649b5e25fSSatish Balay 2474a2ae208SSatish Balay #undef __FUNCT__ 2484a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 24987828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 25049b5e25fSSatish Balay { 25149b5e25fSSatish Balay int ierr; 25249b5e25fSSatish Balay 25349b5e25fSSatish Balay PetscFunctionBegin; 25449b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 25549b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 25649b5e25fSSatish Balay PetscFunctionReturn(0); 25749b5e25fSSatish Balay } 25849b5e25fSSatish Balay 2594a2ae208SSatish Balay #undef __FUNCT__ 2604a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 26149b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 26249b5e25fSSatish Balay { 26349b5e25fSSatish Balay PetscFunctionBegin; 26429bbc08cSBarry Smith SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called"); 26516cdd363SHong Zhang /* PetscFunctionReturn(0); */ 26649b5e25fSSatish Balay } 26749b5e25fSSatish Balay 2684a2ae208SSatish Balay #undef __FUNCT__ 2694a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary" 270b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer) 27149b5e25fSSatish Balay { 27249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 27349b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 27487828ca2SBarry Smith PetscScalar *aa; 27549b5e25fSSatish Balay FILE *file; 27649b5e25fSSatish Balay 27749b5e25fSSatish Balay PetscFunctionBegin; 278b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 27982502324SSatish Balay ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 280552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 28149b5e25fSSatish Balay 282c464158bSHong Zhang col_lens[1] = A->m; 283c464158bSHong Zhang col_lens[2] = A->m; 28449b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 28549b5e25fSSatish Balay 28649b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 28749b5e25fSSatish Balay count = 0; 28849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 28949b5e25fSSatish Balay for (j=0; j<bs; j++) { 29049b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 29149b5e25fSSatish Balay } 29249b5e25fSSatish Balay } 293c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 29449b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 29549b5e25fSSatish Balay 29649b5e25fSSatish Balay /* store column indices (zero start index) */ 29782502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 29849b5e25fSSatish Balay count = 0; 29949b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30049b5e25fSSatish Balay for (j=0; j<bs; j++) { 30149b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 30249b5e25fSSatish Balay for (l=0; l<bs; l++) { 30349b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 30449b5e25fSSatish Balay } 30549b5e25fSSatish Balay } 30649b5e25fSSatish Balay } 30749b5e25fSSatish Balay } 30849b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 30949b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 31049b5e25fSSatish Balay 31149b5e25fSSatish Balay /* store nonzero values */ 31287828ca2SBarry Smith ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 31349b5e25fSSatish Balay count = 0; 31449b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31549b5e25fSSatish Balay for (j=0; j<bs; j++) { 31649b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31749b5e25fSSatish Balay for (l=0; l<bs; l++) { 31849b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 31949b5e25fSSatish Balay } 32049b5e25fSSatish Balay } 32149b5e25fSSatish Balay } 32249b5e25fSSatish Balay } 32349b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 32449b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 32549b5e25fSSatish Balay 326b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 32749b5e25fSSatish Balay if (file) { 32849b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 32949b5e25fSSatish Balay } 33049b5e25fSSatish Balay PetscFunctionReturn(0); 33149b5e25fSSatish Balay } 33249b5e25fSSatish Balay 3334a2ae208SSatish Balay #undef __FUNCT__ 3344a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 335b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 33649b5e25fSSatish Balay { 33749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 338fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 339fb9695e5SSatish Balay char *name; 340f3ef73ceSBarry Smith PetscViewerFormat format; 34149b5e25fSSatish Balay 34249b5e25fSSatish Balay PetscFunctionBegin; 34380fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 344b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 345fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) { 346b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 347fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 34829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 349fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 350b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 35149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 35249b5e25fSSatish Balay for (j=0; j<bs; j++) { 353b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 35449b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 35549b5e25fSSatish Balay for (l=0; l<bs; l++) { 35649b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 35749b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 35862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 35949b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36049b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 36162b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 36249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36349b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 36462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36549b5e25fSSatish Balay } 36649b5e25fSSatish Balay #else 36749b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 36862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 36949b5e25fSSatish Balay } 37049b5e25fSSatish Balay #endif 37149b5e25fSSatish Balay } 37249b5e25fSSatish Balay } 373b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 37449b5e25fSSatish Balay } 37549b5e25fSSatish Balay } 376b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 37749b5e25fSSatish Balay } else { 378b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 37949b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 38049b5e25fSSatish Balay for (j=0; j<bs; j++) { 381b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 38249b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 38349b5e25fSSatish Balay for (l=0; l<bs; l++) { 38449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 38549b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 38662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 38749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38849b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 38962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 39049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39149b5e25fSSatish Balay } else { 39262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39349b5e25fSSatish Balay } 39449b5e25fSSatish Balay #else 395b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 39649b5e25fSSatish Balay #endif 39749b5e25fSSatish Balay } 39849b5e25fSSatish Balay } 399b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 40049b5e25fSSatish Balay } 40149b5e25fSSatish Balay } 402b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 40349b5e25fSSatish Balay } 404b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 40549b5e25fSSatish Balay PetscFunctionReturn(0); 40649b5e25fSSatish Balay } 40749b5e25fSSatish Balay 4084a2ae208SSatish Balay #undef __FUNCT__ 4094a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 410b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 41149b5e25fSSatish Balay { 41249b5e25fSSatish Balay Mat A = (Mat) Aa; 41349b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 41449b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 41549b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 41649b5e25fSSatish Balay MatScalar *aa; 41749b5e25fSSatish Balay MPI_Comm comm; 418b0a32e0cSBarry Smith PetscViewer viewer; 41949b5e25fSSatish Balay 42049b5e25fSSatish Balay PetscFunctionBegin; 42149b5e25fSSatish Balay /* 42249b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 42349b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 42449b5e25fSSatish Balay rest should return immediately. 42549b5e25fSSatish Balay */ 42649b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 42749b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 42849b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 42949b5e25fSSatish Balay 43049b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 43149b5e25fSSatish Balay 432b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 433b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 43449b5e25fSSatish Balay 43549b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 436b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 43749b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 43849b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 439c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 44049b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 44149b5e25fSSatish Balay aa = a->a + j*bs2; 44249b5e25fSSatish Balay for (k=0; k<bs; k++) { 44349b5e25fSSatish Balay for (l=0; l<bs; l++) { 44449b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 445b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44649b5e25fSSatish Balay } 44749b5e25fSSatish Balay } 44849b5e25fSSatish Balay } 44949b5e25fSSatish Balay } 450b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 45149b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45249b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 453c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45449b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45549b5e25fSSatish Balay aa = a->a + j*bs2; 45649b5e25fSSatish Balay for (k=0; k<bs; k++) { 45749b5e25fSSatish Balay for (l=0; l<bs; l++) { 45849b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 459b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46049b5e25fSSatish Balay } 46149b5e25fSSatish Balay } 46249b5e25fSSatish Balay } 46349b5e25fSSatish Balay } 46449b5e25fSSatish Balay 465b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 46649b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46749b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 468c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 46949b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47049b5e25fSSatish Balay aa = a->a + j*bs2; 47149b5e25fSSatish Balay for (k=0; k<bs; k++) { 47249b5e25fSSatish Balay for (l=0; l<bs; l++) { 47349b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 474b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 47549b5e25fSSatish Balay } 47649b5e25fSSatish Balay } 47749b5e25fSSatish Balay } 47849b5e25fSSatish Balay } 47949b5e25fSSatish Balay PetscFunctionReturn(0); 48049b5e25fSSatish Balay } 48149b5e25fSSatish Balay 4824a2ae208SSatish Balay #undef __FUNCT__ 4834a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 484b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 48549b5e25fSSatish Balay { 48649b5e25fSSatish Balay int ierr; 48749b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 488b0a32e0cSBarry Smith PetscDraw draw; 48949b5e25fSSatish Balay PetscTruth isnull; 49049b5e25fSSatish Balay 49149b5e25fSSatish Balay PetscFunctionBegin; 49249b5e25fSSatish Balay 493b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 494b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 49549b5e25fSSatish Balay 49649b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 497c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 49849b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 499b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 500b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 50149b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 50249b5e25fSSatish Balay PetscFunctionReturn(0); 50349b5e25fSSatish Balay } 50449b5e25fSSatish Balay 5054a2ae208SSatish Balay #undef __FUNCT__ 5064a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 507b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 50849b5e25fSSatish Balay { 50949b5e25fSSatish Balay int ierr; 51049b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 51149b5e25fSSatish Balay 51249b5e25fSSatish Balay PetscFunctionBegin; 513b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 514b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 515fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 516fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 51749b5e25fSSatish Balay if (issocket) { 51829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 51949b5e25fSSatish Balay } else if (isascii){ 52049b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 52149b5e25fSSatish Balay } else if (isbinary) { 52249b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 52349b5e25fSSatish Balay } else if (isdraw) { 52449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 52549b5e25fSSatish Balay } else { 52629bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 52749b5e25fSSatish Balay } 52849b5e25fSSatish Balay PetscFunctionReturn(0); 52949b5e25fSSatish Balay } 53049b5e25fSSatish Balay 53149b5e25fSSatish Balay 5324a2ae208SSatish Balay #undef __FUNCT__ 5334a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 53487828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v) 53549b5e25fSSatish Balay { 536045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 53749b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 53849b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 53949b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 54049b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 54149b5e25fSSatish Balay 54249b5e25fSSatish Balay PetscFunctionBegin; 54349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 54449b5e25fSSatish Balay row = im[k]; brow = row/bs; 54529bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 546c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 54749b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 54849b5e25fSSatish Balay nrow = ailen[brow]; 54949b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 55029bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 551c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 55249b5e25fSSatish Balay col = in[l] ; 55349b5e25fSSatish Balay bcol = col/bs; 55449b5e25fSSatish Balay cidx = col%bs; 55549b5e25fSSatish Balay ridx = row%bs; 55649b5e25fSSatish Balay high = nrow; 55749b5e25fSSatish Balay low = 0; /* assume unsorted */ 55849b5e25fSSatish Balay while (high-low > 5) { 55949b5e25fSSatish Balay t = (low+high)/2; 56049b5e25fSSatish Balay if (rp[t] > bcol) high = t; 56149b5e25fSSatish Balay else low = t; 56249b5e25fSSatish Balay } 56349b5e25fSSatish Balay for (i=low; i<high; i++) { 56449b5e25fSSatish Balay if (rp[i] > bcol) break; 56549b5e25fSSatish Balay if (rp[i] == bcol) { 56649b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 56749b5e25fSSatish Balay goto finished; 56849b5e25fSSatish Balay } 56949b5e25fSSatish Balay } 57049b5e25fSSatish Balay *v++ = zero; 57149b5e25fSSatish Balay finished:; 57249b5e25fSSatish Balay } 57349b5e25fSSatish Balay } 57449b5e25fSSatish Balay PetscFunctionReturn(0); 57549b5e25fSSatish Balay } 57649b5e25fSSatish Balay 57749b5e25fSSatish Balay 5784a2ae208SSatish Balay #undef __FUNCT__ 5794a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 58087828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 58149b5e25fSSatish Balay { 58249b5e25fSSatish Balay PetscFunctionBegin; 58329bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format"); 58449b5e25fSSatish Balay } 58549b5e25fSSatish Balay 5864a2ae208SSatish Balay #undef __FUNCT__ 5874a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 58849b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 58949b5e25fSSatish Balay { 59049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 59149b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 592c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 59349b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 59449b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 595*cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 596*cf4441caSHong Zhang PetscTruth flag; 597*cf4441caSHong Zhang #endif 59849b5e25fSSatish Balay 59949b5e25fSSatish Balay PetscFunctionBegin; 60049b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 60149b5e25fSSatish Balay 60249b5e25fSSatish Balay if (m) rmax = ailen[0]; 60349b5e25fSSatish Balay for (i=1; i<mbs; i++) { 60449b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 60549b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 60649b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 60749b5e25fSSatish Balay if (fshift) { 60849b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 60949b5e25fSSatish Balay N = ailen[i]; 61049b5e25fSSatish Balay for (j=0; j<N; j++) { 61149b5e25fSSatish Balay ip[j-fshift] = ip[j]; 61249b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 61349b5e25fSSatish Balay } 61449b5e25fSSatish Balay } 61549b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 61649b5e25fSSatish Balay } 61749b5e25fSSatish Balay if (mbs) { 61849b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 61949b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 62049b5e25fSSatish Balay } 62149b5e25fSSatish Balay /* reset ilen and imax for each row */ 62249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 62349b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 62449b5e25fSSatish Balay } 62549b5e25fSSatish Balay a->s_nz = ai[mbs]; 62649b5e25fSSatish Balay 627b424e231SHong Zhang /* diagonals may have moved, reset it */ 628b424e231SHong Zhang if (a->diag) { 629b424e231SHong Zhang ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr); 63049b5e25fSSatish Balay } 631b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 632c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 633b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 63449b5e25fSSatish Balay a->reallocs); 635b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 63649b5e25fSSatish Balay a->reallocs = 0; 63749b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 63849b5e25fSSatish Balay 639*cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 640*cf4441caSHong Zhang ierr = PetscOptionsHasName(PETSC_NULL,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr); 641*cf4441caSHong Zhang if (flag) { ierr = MatUseSpooles_SeqSBAIJ(A);CHKERRQ(ierr); } 642*cf4441caSHong Zhang #endif 643*cf4441caSHong Zhang 64449b5e25fSSatish Balay PetscFunctionReturn(0); 64549b5e25fSSatish Balay } 64649b5e25fSSatish Balay 64749b5e25fSSatish Balay /* 64849b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 64949b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 65049b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 65149b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 65249b5e25fSSatish Balay */ 6534a2ae208SSatish Balay #undef __FUNCT__ 6544a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 65549b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 65649b5e25fSSatish Balay { 65749b5e25fSSatish Balay int i,j,k,row; 65849b5e25fSSatish Balay PetscTruth flg; 65949b5e25fSSatish Balay 66049b5e25fSSatish Balay PetscFunctionBegin; 66149b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 66249b5e25fSSatish Balay row = idx[i]; 66349b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 66449b5e25fSSatish Balay sizes[j] = 1; 66549b5e25fSSatish Balay i++; 66649b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 66749b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 66849b5e25fSSatish Balay i++; 66949b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 67049b5e25fSSatish Balay flg = PETSC_TRUE; 67149b5e25fSSatish Balay for (k=1; k<bs; k++) { 67249b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 67349b5e25fSSatish Balay flg = PETSC_FALSE; 67449b5e25fSSatish Balay break; 67549b5e25fSSatish Balay } 67649b5e25fSSatish Balay } 67749b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 67849b5e25fSSatish Balay sizes[j] = bs; 67949b5e25fSSatish Balay i+= bs; 68049b5e25fSSatish Balay } else { 68149b5e25fSSatish Balay sizes[j] = 1; 68249b5e25fSSatish Balay i++; 68349b5e25fSSatish Balay } 68449b5e25fSSatish Balay } 68549b5e25fSSatish Balay } 68649b5e25fSSatish Balay *bs_max = j; 68749b5e25fSSatish Balay PetscFunctionReturn(0); 68849b5e25fSSatish Balay } 68949b5e25fSSatish Balay 6904a2ae208SSatish Balay #undef __FUNCT__ 6914a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 69287828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag) 69349b5e25fSSatish Balay { 69449b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data; 69549b5e25fSSatish Balay int ierr,i,j,k,count,is_n,*is_idx,*rows; 69649b5e25fSSatish Balay int bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max; 69787828ca2SBarry Smith PetscScalar zero = 0.0; 69849b5e25fSSatish Balay MatScalar *aa; 69949b5e25fSSatish Balay 70049b5e25fSSatish Balay PetscFunctionBegin; 70149b5e25fSSatish Balay /* Make a copy of the IS and sort it */ 70249b5e25fSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 70349b5e25fSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 70449b5e25fSSatish Balay 70549b5e25fSSatish Balay /* allocate memory for rows,sizes */ 706b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 70749b5e25fSSatish Balay sizes = rows + is_n; 70849b5e25fSSatish Balay 70949b5e25fSSatish Balay /* initialize copy IS values to rows, and sort them */ 71049b5e25fSSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 71149b5e25fSSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 71249b5e25fSSatish Balay if (sbaij->keepzeroedrows) { /* do not change nonzero structure */ 71349b5e25fSSatish Balay for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */ 71449b5e25fSSatish Balay bs_max = is_n; /* bs_max: num. of contiguous block row in the row */ 71549b5e25fSSatish Balay } else { 71649b5e25fSSatish Balay ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 71749b5e25fSSatish Balay } 71849b5e25fSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 71949b5e25fSSatish Balay 72049b5e25fSSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 72149b5e25fSSatish Balay row = rows[j]; /* row to be zeroed */ 722c464158bSHong Zhang if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 72349b5e25fSSatish Balay count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */ 72449b5e25fSSatish Balay aa = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs); 72549b5e25fSSatish Balay if (sizes[i] == bs && !sbaij->keepzeroedrows) { 72649b5e25fSSatish Balay if (diag) { 72749b5e25fSSatish Balay if (sbaij->ilen[row/bs] > 0) { 72849b5e25fSSatish Balay sbaij->ilen[row/bs] = 1; 72949b5e25fSSatish Balay sbaij->j[sbaij->i[row/bs]] = row/bs; 73049b5e25fSSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 73149b5e25fSSatish Balay } 73249b5e25fSSatish Balay /* Now insert all the diagoanl values for this bs */ 73349b5e25fSSatish Balay for (k=0; k<bs; k++) { 73449b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 73549b5e25fSSatish Balay } 73649b5e25fSSatish Balay } else { /* (!diag) */ 73749b5e25fSSatish Balay sbaij->ilen[row/bs] = 0; 73849b5e25fSSatish Balay } /* end (!diag) */ 73949b5e25fSSatish Balay } else { /* (sizes[i] != bs), broken block */ 74049b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG) 74129bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 74249b5e25fSSatish Balay #endif 74349b5e25fSSatish Balay for (k=0; k<count; k++) { 74449b5e25fSSatish Balay aa[0] = zero; 74549b5e25fSSatish Balay aa+=bs; 74649b5e25fSSatish Balay } 74749b5e25fSSatish Balay if (diag) { 74849b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 74949b5e25fSSatish Balay } 75049b5e25fSSatish Balay } 75149b5e25fSSatish Balay } 75249b5e25fSSatish Balay 75349b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 75449b5e25fSSatish Balay ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75549b5e25fSSatish Balay ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75649b5e25fSSatish Balay PetscFunctionReturn(0); 75749b5e25fSSatish Balay } 75849b5e25fSSatish Balay 75949b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 76049b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 76149b5e25fSSatish Balay */ 76249b5e25fSSatish Balay 7634a2ae208SSatish Balay #undef __FUNCT__ 7644a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 76587828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 76649b5e25fSSatish Balay { 76749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 76849b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 76949b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 77049b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 77149b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 77249b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 77349b5e25fSSatish Balay 77449b5e25fSSatish Balay PetscFunctionBegin; 77549b5e25fSSatish Balay 77649b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 77749b5e25fSSatish Balay row = im[k]; /* row number */ 77849b5e25fSSatish Balay brow = row/bs; /* block row number */ 77949b5e25fSSatish Balay if (row < 0) continue; 78049b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 781c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 78249b5e25fSSatish Balay #endif 78349b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 78449b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 78549b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 78649b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 78749b5e25fSSatish Balay low = 0; 78849b5e25fSSatish Balay 78949b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 79049b5e25fSSatish Balay if (in[l] < 0) continue; 79149b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 792c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 79349b5e25fSSatish Balay #endif 79449b5e25fSSatish Balay col = in[l]; 79549b5e25fSSatish Balay bcol = col/bs; /* block col number */ 79649b5e25fSSatish Balay 79749b5e25fSSatish Balay if (brow <= bcol){ 79849b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 7998549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 80049b5e25fSSatish Balay /* element value a(k,l) */ 80149b5e25fSSatish Balay if (roworiented) { 80249b5e25fSSatish Balay value = v[l + k*n]; 80349b5e25fSSatish Balay } else { 80449b5e25fSSatish Balay value = v[k + l*m]; 80549b5e25fSSatish Balay } 80649b5e25fSSatish Balay 80749b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 80849b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 80949b5e25fSSatish Balay while (high-low > 7) { 81049b5e25fSSatish Balay t = (low+high)/2; 81149b5e25fSSatish Balay if (rp[t] > bcol) high = t; 81249b5e25fSSatish Balay else low = t; 81349b5e25fSSatish Balay } 81449b5e25fSSatish Balay for (i=low; i<high; i++) { 81549b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 81649b5e25fSSatish Balay if (rp[i] > bcol) break; 81749b5e25fSSatish Balay if (rp[i] == bcol) { 81849b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 81949b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 82049b5e25fSSatish Balay else *bap = value; 8218549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8228549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8238549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8248549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8258549e402SHong Zhang else *bap = value; 8268549e402SHong Zhang } 82749b5e25fSSatish Balay goto noinsert1; 82849b5e25fSSatish Balay } 82949b5e25fSSatish Balay } 83049b5e25fSSatish Balay 83149b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 83229bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 83349b5e25fSSatish Balay if (nrow >= rmax) { 83449b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 83549b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 83649b5e25fSSatish Balay MatScalar *new_a; 83749b5e25fSSatish Balay 83829bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 83949b5e25fSSatish Balay 84049b5e25fSSatish Balay /* Malloc new storage space */ 84149b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 84282502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 84349b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 84449b5e25fSSatish Balay new_i = new_j + new_nz; 84549b5e25fSSatish Balay 84649b5e25fSSatish Balay /* copy over old data into new slots */ 84749b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 84849b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 84949b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 85049b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 85149b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 85249b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 85349b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 85449b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 85549b5e25fSSatish Balay /* free up old matrix storage */ 85649b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 85749b5e25fSSatish Balay if (!a->singlemalloc) { 85849b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 85949b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 86049b5e25fSSatish Balay } 86149b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 86249b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 86349b5e25fSSatish Balay 86449b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 86549b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 866b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 86749b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 86849b5e25fSSatish Balay a->reallocs++; 86949b5e25fSSatish Balay a->s_nz++; 87049b5e25fSSatish Balay } 87149b5e25fSSatish Balay 87249b5e25fSSatish Balay N = nrow++ - 1; 87349b5e25fSSatish Balay /* shift up all the later entries in this row */ 87449b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 87549b5e25fSSatish Balay rp[ii+1] = rp[ii]; 87649b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 87749b5e25fSSatish Balay } 87849b5e25fSSatish Balay if (N>=i) { 87949b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 88049b5e25fSSatish Balay } 88149b5e25fSSatish Balay rp[i] = bcol; 88249b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 88349b5e25fSSatish Balay noinsert1:; 88449b5e25fSSatish Balay low = i; 88549b5e25fSSatish Balay /* } */ 8868549e402SHong Zhang } 88749b5e25fSSatish Balay } /* end of if .. if.. */ 88849b5e25fSSatish Balay } /* end of loop over added columns */ 88949b5e25fSSatish Balay ailen[brow] = nrow; 89049b5e25fSSatish Balay } /* end of loop over added rows */ 89149b5e25fSSatish Balay 89249b5e25fSSatish Balay PetscFunctionReturn(0); 89349b5e25fSSatish Balay } 89449b5e25fSSatish Balay 895915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*); 896915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal); 89749b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 89849b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 89949b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 90049b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 90149b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 90287828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat); 90349b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 90449b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 90549b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 90649b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 90749b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 90849b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 90913a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 91049b5e25fSSatish Balay 91149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 91249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 91349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 91449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 91549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 91649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 91749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 91849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 91949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 92049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 92149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 92249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 92349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 92449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 92549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 92649b5e25fSSatish Balay 92707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 92807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 92907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 93007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 93107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 93207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 93307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 93407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 93507f98182SSatish Balay 9365f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9375f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9385f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9395f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9405f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 94357d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 94449b5e25fSSatish Balay 94549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 94649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 94749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 94849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 94949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 95049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 95149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 95249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 95349b5e25fSSatish Balay 95449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 95549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 95649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 95749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 95849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 95949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 96049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 96149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 96249b5e25fSSatish Balay 9634d101231SSatish Balay #ifdef HAVE_MatICCFactor 9641a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 9654a2ae208SSatish Balay #undef __FUNCT__ 9664d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 9674d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level) 96849b5e25fSSatish Balay { 9694ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 97049b5e25fSSatish Balay Mat outA; 97149b5e25fSSatish Balay int ierr; 97249b5e25fSSatish Balay PetscTruth row_identity,col_identity; 97349b5e25fSSatish Balay 97449b5e25fSSatish Balay PetscFunctionBegin; 9751a3463dfSHong Zhang /* 97629bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 97749b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 97849b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 97949b5e25fSSatish Balay if (!row_identity || !col_identity) { 98029bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 98149b5e25fSSatish Balay } 9821a3463dfSHong Zhang */ 98349b5e25fSSatish Balay 98449b5e25fSSatish Balay outA = inA; 9851260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 98649b5e25fSSatish Balay 98749b5e25fSSatish Balay if (!a->diag) { 9881a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 98949b5e25fSSatish Balay } 99049b5e25fSSatish Balay /* 99149b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 99249b5e25fSSatish Balay for ILU(0) factorization with natural ordering 99349b5e25fSSatish Balay */ 99449b5e25fSSatish Balay switch (a->bs) { 99549b5e25fSSatish Balay case 1: 9960fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 99707f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 9984d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 99949b5e25fSSatish Balay case 2: 10001a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10011a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10021a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10034d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 100449b5e25fSSatish Balay break; 100549b5e25fSSatish Balay case 3: 10061a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10071a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10081a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10094d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 101049b5e25fSSatish Balay break; 101149b5e25fSSatish Balay case 4: 10121a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10131a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10141a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10154d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 101649b5e25fSSatish Balay break; 101749b5e25fSSatish Balay case 5: 10181a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10191a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10201a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10214d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 102249b5e25fSSatish Balay break; 102349b5e25fSSatish Balay case 6: 10241a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10251a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10261a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10274d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 102849b5e25fSSatish Balay break; 102949b5e25fSSatish Balay case 7: 10301a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10311a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10321a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10334d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 103449b5e25fSSatish Balay break; 103549b5e25fSSatish Balay default: 103649b5e25fSSatish Balay a->row = row; 10371a3463dfSHong Zhang a->icol = col; 103849b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 103949b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 104049b5e25fSSatish Balay 104149b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 104249b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1043b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 104449b5e25fSSatish Balay 104549b5e25fSSatish Balay if (!a->solve_work) { 104687828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 104787828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 104849b5e25fSSatish Balay } 104949b5e25fSSatish Balay } 105049b5e25fSSatish Balay 10511a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 105249b5e25fSSatish Balay 105349b5e25fSSatish Balay PetscFunctionReturn(0); 105449b5e25fSSatish Balay } 1055950f1e5bSHong Zhang #endif 1056950f1e5bSHong Zhang 10574a2ae208SSatish Balay #undef __FUNCT__ 10584a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 105949b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 106049b5e25fSSatish Balay { 106149b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 106249b5e25fSSatish Balay MPI_Comm comm = A->comm; 106349b5e25fSSatish Balay int ierr; 106449b5e25fSSatish Balay 106549b5e25fSSatish Balay PetscFunctionBegin; 106649b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 106749b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 106849b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 106949b5e25fSSatish Balay PetscFunctionReturn(0); 107049b5e25fSSatish Balay } 107149b5e25fSSatish Balay 107249b5e25fSSatish Balay EXTERN_C_BEGIN 10734a2ae208SSatish Balay #undef __FUNCT__ 10744a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 107549b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 107649b5e25fSSatish Balay { 1077045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 107849b5e25fSSatish Balay int i,nz,n; 107949b5e25fSSatish Balay 108049b5e25fSSatish Balay PetscFunctionBegin; 1081045c9aa0SHong Zhang nz = baij->s_maxnz; 1082c464158bSHong Zhang n = mat->n; 108349b5e25fSSatish Balay for (i=0; i<nz; i++) { 108449b5e25fSSatish Balay baij->j[i] = indices[i]; 108549b5e25fSSatish Balay } 1086045c9aa0SHong Zhang baij->s_nz = nz; 108749b5e25fSSatish Balay for (i=0; i<n; i++) { 108849b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 108949b5e25fSSatish Balay } 109049b5e25fSSatish Balay 109149b5e25fSSatish Balay PetscFunctionReturn(0); 109249b5e25fSSatish Balay } 109349b5e25fSSatish Balay EXTERN_C_END 109449b5e25fSSatish Balay 10954a2ae208SSatish Balay #undef __FUNCT__ 10964a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 109749b5e25fSSatish Balay /*@ 109819585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 109949b5e25fSSatish Balay in the matrix. 110049b5e25fSSatish Balay 110149b5e25fSSatish Balay Input Parameters: 110219585528SSatish Balay + mat - the SeqSBAIJ matrix 110349b5e25fSSatish Balay - indices - the column indices 110449b5e25fSSatish Balay 110549b5e25fSSatish Balay Level: advanced 110649b5e25fSSatish Balay 110749b5e25fSSatish Balay Notes: 110849b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 110949b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 111049b5e25fSSatish Balay of the MatSetValues() operation. 111149b5e25fSSatish Balay 111249b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 111319585528SSatish Balay MatCreateSeqSBAIJ(). 111449b5e25fSSatish Balay 1115ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 111649b5e25fSSatish Balay 1117ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 111849b5e25fSSatish Balay @*/ 111949b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 112049b5e25fSSatish Balay { 112149b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 112249b5e25fSSatish Balay 112349b5e25fSSatish Balay PetscFunctionBegin; 112449b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1125c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 112649b5e25fSSatish Balay if (f) { 112749b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 112849b5e25fSSatish Balay } else { 112929bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 113049b5e25fSSatish Balay } 113149b5e25fSSatish Balay PetscFunctionReturn(0); 113249b5e25fSSatish Balay } 113349b5e25fSSatish Balay 11344a2ae208SSatish Balay #undef __FUNCT__ 11354a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1136273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1137273d9f13SBarry Smith { 1138273d9f13SBarry Smith int ierr; 1139273d9f13SBarry Smith 1140273d9f13SBarry Smith PetscFunctionBegin; 1141273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1142273d9f13SBarry Smith PetscFunctionReturn(0); 1143273d9f13SBarry Smith } 1144273d9f13SBarry Smith 114549b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 114649b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 114749b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 114849b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 114949b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 115049b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 115149b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 115249b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 115349b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 115449b5e25fSSatish Balay 0, 115549b5e25fSSatish Balay 0, 115649b5e25fSSatish Balay 0, 115749b5e25fSSatish Balay 0, 11585f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1159d06b337dSHong Zhang MatRelax_SeqSBAIJ, 116049b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 116149b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 116249b5e25fSSatish Balay MatEqual_SeqSBAIJ, 116349b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 116449b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 116549b5e25fSSatish Balay MatNorm_SeqSBAIJ, 116649b5e25fSSatish Balay 0, 116749b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 116849b5e25fSSatish Balay 0, 116949b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 117049b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 117149b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 117249b5e25fSSatish Balay 0, 117349b5e25fSSatish Balay 0, 11745f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 11755f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1176273d9f13SBarry Smith MatSetUpPreallocation_SeqSBAIJ, 1177c464158bSHong Zhang 0, 11784d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 117949b5e25fSSatish Balay 0, 118049b5e25fSSatish Balay 0, 118149b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 118249b5e25fSSatish Balay 0, 118349b5e25fSSatish Balay 0, 118449b5e25fSSatish Balay 0, 1185950f1e5bSHong Zhang 0, 118649b5e25fSSatish Balay 0, 118749b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 118849b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 118949b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 119049b5e25fSSatish Balay 0, 119149b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 119249b5e25fSSatish Balay MatScale_SeqSBAIJ, 119349b5e25fSSatish Balay 0, 119449b5e25fSSatish Balay 0, 119549b5e25fSSatish Balay 0, 119649b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 119749b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 119849b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 119949b5e25fSSatish Balay 0, 120049b5e25fSSatish Balay 0, 120149b5e25fSSatish Balay 0, 120249b5e25fSSatish Balay 0, 120349b5e25fSSatish Balay 0, 120449b5e25fSSatish Balay 0, 120549b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 120649b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 120749b5e25fSSatish Balay 0, 120849b5e25fSSatish Balay 0, 12098a124369SBarry Smith MatGetPetscMaps_Petsc, 1210d959ec07SHong Zhang 0, 1211d959ec07SHong Zhang 0, 1212d959ec07SHong Zhang 0, 1213d959ec07SHong Zhang 0, 1214d959ec07SHong Zhang 0, 1215d959ec07SHong Zhang 0, 121624d5174aSHong Zhang MatGetRowMax_SeqSBAIJ}; 121749b5e25fSSatish Balay 121849b5e25fSSatish Balay EXTERN_C_BEGIN 12194a2ae208SSatish Balay #undef __FUNCT__ 12204a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 122149b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 122249b5e25fSSatish Balay { 12234afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1224c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 122549b5e25fSSatish Balay int ierr; 122649b5e25fSSatish Balay 122749b5e25fSSatish Balay PetscFunctionBegin; 122849b5e25fSSatish Balay if (aij->nonew != 1) { 122929bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 123049b5e25fSSatish Balay } 123149b5e25fSSatish Balay 123249b5e25fSSatish Balay /* allocate space for values if not already there */ 123349b5e25fSSatish Balay if (!aij->saved_values) { 123487828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 123549b5e25fSSatish Balay } 123649b5e25fSSatish Balay 123749b5e25fSSatish Balay /* copy values over */ 123887828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 123949b5e25fSSatish Balay PetscFunctionReturn(0); 124049b5e25fSSatish Balay } 124149b5e25fSSatish Balay EXTERN_C_END 124249b5e25fSSatish Balay 124349b5e25fSSatish Balay EXTERN_C_BEGIN 12444a2ae208SSatish Balay #undef __FUNCT__ 12454a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 124649b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 124749b5e25fSSatish Balay { 12484afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1249c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 125049b5e25fSSatish Balay 125149b5e25fSSatish Balay PetscFunctionBegin; 125249b5e25fSSatish Balay if (aij->nonew != 1) { 125329bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 125449b5e25fSSatish Balay } 125549b5e25fSSatish Balay if (!aij->saved_values) { 125629bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 125749b5e25fSSatish Balay } 125849b5e25fSSatish Balay 125949b5e25fSSatish Balay /* copy values over */ 126087828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 126149b5e25fSSatish Balay PetscFunctionReturn(0); 126249b5e25fSSatish Balay } 126349b5e25fSSatish Balay EXTERN_C_END 126449b5e25fSSatish Balay 12658549e402SHong Zhang EXTERN_C_BEGIN 12664a2ae208SSatish Balay #undef __FUNCT__ 12674a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ" 1268c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1269c464158bSHong Zhang { 1270c464158bSHong Zhang Mat_SeqSBAIJ *b; 12710c955e93SHong Zhang int ierr,size; 1272c464158bSHong Zhang 1273c464158bSHong Zhang PetscFunctionBegin; 1274c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1275c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1276273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1277273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1278c464158bSHong Zhang 1279b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1280b0a32e0cSBarry Smith B->data = (void*)b; 1281c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1282c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1283c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1284c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1285c464158bSHong Zhang B->factor = 0; 1286c464158bSHong Zhang B->lupivotthreshold = 1.0; 1287c464158bSHong Zhang B->mapping = 0; 1288c464158bSHong Zhang b->row = 0; 1289c464158bSHong Zhang b->icol = 0; 1290c464158bSHong Zhang b->reallocs = 0; 1291c464158bSHong Zhang b->saved_values = 0; 1292c464158bSHong Zhang 12938a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 12948a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1295c464158bSHong Zhang 1296c464158bSHong Zhang b->sorted = PETSC_FALSE; 1297c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1298c464158bSHong Zhang b->nonew = 0; 1299c464158bSHong Zhang b->diag = 0; 1300c464158bSHong Zhang b->solve_work = 0; 1301c464158bSHong Zhang b->mult_work = 0; 13022a1b7f2aSHong Zhang B->spptr = 0; 1303c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1304c464158bSHong Zhang 1305c464158bSHong Zhang b->inew = 0; 1306c464158bSHong Zhang b->jnew = 0; 1307c464158bSHong Zhang b->anew = 0; 1308c464158bSHong Zhang b->a2anew = 0; 1309c464158bSHong Zhang b->permute = PETSC_FALSE; 1310c464158bSHong Zhang 1311c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1312c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1313c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1314c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1315c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1316c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1317c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1318c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1319c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1320c464158bSHong Zhang PetscFunctionReturn(0); 1321c464158bSHong Zhang } 13228549e402SHong Zhang EXTERN_C_END 1323c464158bSHong Zhang 13244a2ae208SSatish Balay #undef __FUNCT__ 13254a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 132649b5e25fSSatish Balay /*@C 1327273d9f13SBarry Smith MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 132849b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 132949b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 133049b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 133149b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 133249b5e25fSSatish Balay 1333c464158bSHong Zhang Collective on Mat 133449b5e25fSSatish Balay 133549b5e25fSSatish Balay Input Parameters: 1336c464158bSHong Zhang + A - the symmetric matrix 133749b5e25fSSatish Balay . bs - size of block 133849b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 133949b5e25fSSatish Balay - nnz - array containing the number of block nonzeros in the various block rows 134049b5e25fSSatish Balay (possibly different for each block row) or PETSC_NULL 134149b5e25fSSatish Balay 134249b5e25fSSatish Balay Options Database Keys: 134349b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 134449b5e25fSSatish Balay block calculations (much slower) 134549b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 134649b5e25fSSatish Balay 134749b5e25fSSatish Balay Level: intermediate 134849b5e25fSSatish Balay 134949b5e25fSSatish Balay Notes: 1350c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 135149b5e25fSSatish Balay storage. That is, the stored row and column indices can begin at 135249b5e25fSSatish Balay either one (as in Fortran) or zero. See the users' manual for details. 135349b5e25fSSatish Balay 135449b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 135549b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 135649b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 135749b5e25fSSatish Balay matrices. 135849b5e25fSSatish Balay 1359a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 136049b5e25fSSatish Balay @*/ 1361273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 136249b5e25fSSatish Balay { 1363c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13640c955e93SHong Zhang int i,len,ierr,mbs,bs2; 136549b5e25fSSatish Balay PetscTruth flg; 13664afc71dfSHong Zhang int s_nz; 136749b5e25fSSatish Balay 136849b5e25fSSatish Balay PetscFunctionBegin; 1369273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1370b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1371c464158bSHong Zhang mbs = B->m/bs; 137249b5e25fSSatish Balay bs2 = bs*bs; 137349b5e25fSSatish Balay 1374c464158bSHong Zhang if (mbs*bs != B->m) { 137529bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 137649b5e25fSSatish Balay } 137749b5e25fSSatish Balay 1378435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1379435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 138049b5e25fSSatish Balay if (nnz) { 138149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 138229bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 138380fe4e49SBarry 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); 138449b5e25fSSatish Balay } 138549b5e25fSSatish Balay } 138649b5e25fSSatish Balay 1387b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 138849b5e25fSSatish Balay if (!flg) { 138949b5e25fSSatish Balay switch (bs) { 139049b5e25fSSatish Balay case 1: 13914ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 139249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 139349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 139449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 139549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 139649b5e25fSSatish Balay break; 139749b5e25fSSatish Balay case 2: 13984ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 139949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 140049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 140149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 140249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 140349b5e25fSSatish Balay break; 140449b5e25fSSatish Balay case 3: 14055f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 140649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 140749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 140849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 140949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 141049b5e25fSSatish Balay break; 141149b5e25fSSatish Balay case 4: 14125f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 141349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 141449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 141549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 141649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 141749b5e25fSSatish Balay break; 141849b5e25fSSatish Balay case 5: 14195f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 142049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 142149b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 142249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 142349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 142449b5e25fSSatish Balay break; 142549b5e25fSSatish Balay case 6: 14265f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 142749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 142849b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 142949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 143049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 143149b5e25fSSatish Balay break; 143249b5e25fSSatish Balay case 7: 1433de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 143449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 143549b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1436de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 143749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 143849b5e25fSSatish Balay break; 143949b5e25fSSatish Balay } 144049b5e25fSSatish Balay } 144149b5e25fSSatish Balay 144249b5e25fSSatish Balay b->mbs = mbs; 14434afc71dfSHong Zhang b->nbs = mbs; 1444b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 144549b5e25fSSatish Balay if (!nnz) { 1446435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 144749b5e25fSSatish Balay else if (nz <= 0) nz = 1; 144849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14498cef66ccSHong Zhang b->imax[i] = nz; 145049b5e25fSSatish Balay } 1451153ea458SHong Zhang nz = nz*mbs; /* total nz */ 145249b5e25fSSatish Balay } else { 145349b5e25fSSatish Balay nz = 0; 14548cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 145549b5e25fSSatish Balay } 14568cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 14578cef66ccSHong Zhang s_nz = nz; 145849b5e25fSSatish Balay 145949b5e25fSSatish Balay /* allocate the matrix space */ 1460c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 146182502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 146249b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 146349b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 146449b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 146549b5e25fSSatish Balay b->i = b->j + s_nz; 146649b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 146749b5e25fSSatish Balay 146849b5e25fSSatish Balay /* pointer to beginning of each row */ 146949b5e25fSSatish Balay b->i[0] = 0; 147049b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 147149b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 147249b5e25fSSatish Balay } 147349b5e25fSSatish Balay 147449b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 14755033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1476b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 147749b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 147849b5e25fSSatish Balay 147949b5e25fSSatish Balay b->bs = bs; 148049b5e25fSSatish Balay b->bs2 = bs2; 148149b5e25fSSatish Balay b->s_nz = 0; 148249b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1483153ea458SHong Zhang 148416cdd363SHong Zhang b->inew = 0; 148516cdd363SHong Zhang b->jnew = 0; 148616cdd363SHong Zhang b->anew = 0; 148716cdd363SHong Zhang b->a2anew = 0; 14881a3463dfSHong Zhang b->permute = PETSC_FALSE; 1489c464158bSHong Zhang PetscFunctionReturn(0); 1490c464158bSHong Zhang } 1491153ea458SHong Zhang 149249b5e25fSSatish Balay 14934a2ae208SSatish Balay #undef __FUNCT__ 14944a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1495c464158bSHong Zhang /*@C 1496c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1497c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1498c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1499c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1500c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 150149b5e25fSSatish Balay 1502c464158bSHong Zhang Collective on MPI_Comm 1503c464158bSHong Zhang 1504c464158bSHong Zhang Input Parameters: 1505c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1506c464158bSHong Zhang . bs - size of block 1507c464158bSHong Zhang . m - number of rows, or number of columns 1508c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1509c464158bSHong Zhang - nnz - array containing the number of block nonzeros in the various block rows 1510c464158bSHong Zhang (possibly different for each block row) or PETSC_NULL 1511c464158bSHong Zhang 1512c464158bSHong Zhang Output Parameter: 1513c464158bSHong Zhang . A - the symmetric matrix 1514c464158bSHong Zhang 1515c464158bSHong Zhang Options Database Keys: 1516c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1517c464158bSHong Zhang block calculations (much slower) 1518c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1519c464158bSHong Zhang 1520c464158bSHong Zhang Level: intermediate 1521c464158bSHong Zhang 1522c464158bSHong Zhang Notes: 1523c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 1524c464158bSHong Zhang storage. That is, the stored row and column indices can begin at 1525c464158bSHong Zhang either one (as in Fortran) or zero. See the users' manual for details. 1526c464158bSHong Zhang 1527c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1528c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1529c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1530c464158bSHong Zhang matrices. 1531c464158bSHong Zhang 1532c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1533c464158bSHong Zhang @*/ 1534c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1535c464158bSHong Zhang { 1536c464158bSHong Zhang int ierr; 1537c464158bSHong Zhang 1538c464158bSHong Zhang PetscFunctionBegin; 1539c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1540c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1541273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 154249b5e25fSSatish Balay PetscFunctionReturn(0); 154349b5e25fSSatish Balay } 154449b5e25fSSatish Balay 15454a2ae208SSatish Balay #undef __FUNCT__ 15464a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 154749b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 154849b5e25fSSatish Balay { 154949b5e25fSSatish Balay Mat C; 155049b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 155149b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 155249b5e25fSSatish Balay 155349b5e25fSSatish Balay PetscFunctionBegin; 155429bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 155549b5e25fSSatish Balay 155649b5e25fSSatish Balay *B = 0; 1557692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1558692f9cbeSHong Zhang ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 1559692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1560692f9cbeSHong Zhang 156149b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1562273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 156349b5e25fSSatish Balay C->factor = A->factor; 156449b5e25fSSatish Balay c->row = 0; 156549b5e25fSSatish Balay c->icol = 0; 156649b5e25fSSatish Balay c->saved_values = 0; 156749b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 156849b5e25fSSatish Balay C->assembled = PETSC_TRUE; 156949b5e25fSSatish Balay 157049b5e25fSSatish Balay c->bs = a->bs; 157149b5e25fSSatish Balay c->bs2 = a->bs2; 157249b5e25fSSatish Balay c->mbs = a->mbs; 157349b5e25fSSatish Balay c->nbs = a->nbs; 157449b5e25fSSatish Balay 1575b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1576b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 157749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 157849b5e25fSSatish Balay c->imax[i] = a->imax[i]; 157949b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 158049b5e25fSSatish Balay } 158149b5e25fSSatish Balay 158249b5e25fSSatish Balay /* allocate the matrix space */ 158349b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 158449b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 158582502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 158649b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 158749b5e25fSSatish Balay c->i = c->j + nz; 158849b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 158949b5e25fSSatish Balay if (mbs > 0) { 159049b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 159149b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 159249b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 159349b5e25fSSatish Balay } else { 159449b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 159549b5e25fSSatish Balay } 159649b5e25fSSatish Balay } 159749b5e25fSSatish Balay 1598b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 159949b5e25fSSatish Balay c->sorted = a->sorted; 160049b5e25fSSatish Balay c->roworiented = a->roworiented; 160149b5e25fSSatish Balay c->nonew = a->nonew; 160249b5e25fSSatish Balay 160349b5e25fSSatish Balay if (a->diag) { 1604b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1605b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 160649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 160749b5e25fSSatish Balay c->diag[i] = a->diag[i]; 160849b5e25fSSatish Balay } 160949b5e25fSSatish Balay } else c->diag = 0; 161049b5e25fSSatish Balay c->s_nz = a->s_nz; 161149b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 161249b5e25fSSatish Balay c->solve_work = 0; 16132a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 161449b5e25fSSatish Balay c->mult_work = 0; 161549b5e25fSSatish Balay *B = C; 1616b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 161749b5e25fSSatish Balay PetscFunctionReturn(0); 161849b5e25fSSatish Balay } 161949b5e25fSSatish Balay 16208549e402SHong Zhang EXTERN_C_BEGIN 16214a2ae208SSatish Balay #undef __FUNCT__ 16224a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1623b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A) 162449b5e25fSSatish Balay { 162549b5e25fSSatish Balay Mat_SeqSBAIJ *a; 162649b5e25fSSatish Balay Mat B; 162749b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 16283f7326a9SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 162949b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 163049b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 163187828ca2SBarry Smith PetscScalar *aa; 163249b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 163349b5e25fSSatish Balay 163449b5e25fSSatish Balay PetscFunctionBegin; 1635b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 163649b5e25fSSatish Balay bs2 = bs*bs; 163749b5e25fSSatish Balay 163849b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 163929bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1640b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 164149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1642552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 164349b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 164449b5e25fSSatish Balay 164549b5e25fSSatish Balay if (header[3] < 0) { 164629bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 164749b5e25fSSatish Balay } 164849b5e25fSSatish Balay 164929bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 165049b5e25fSSatish Balay 165149b5e25fSSatish Balay /* 165249b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 165349b5e25fSSatish Balay divisible by the blocksize 165449b5e25fSSatish Balay */ 165549b5e25fSSatish Balay mbs = M/bs; 165649b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 165749b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 165849b5e25fSSatish Balay else mbs++; 165949b5e25fSSatish Balay if (extra_rows) { 1660b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 166149b5e25fSSatish Balay } 166249b5e25fSSatish Balay 166349b5e25fSSatish Balay /* read in row lengths */ 1664b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 166549b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 166649b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 166749b5e25fSSatish Balay 166849b5e25fSSatish Balay /* read in column indices */ 1669b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 167049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 167149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 167249b5e25fSSatish Balay 167349b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 167482502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 167582502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 167649b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 167749b5e25fSSatish Balay masked = mask + mbs; 167849b5e25fSSatish Balay rowcount = 0; nzcount = 0; 167949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 168049b5e25fSSatish Balay nmask = 0; 168149b5e25fSSatish Balay for (j=0; j<bs; j++) { 168249b5e25fSSatish Balay kmax = rowlengths[rowcount]; 168349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 16842d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 168503630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 168649b5e25fSSatish Balay } 168749b5e25fSSatish Balay rowcount++; 168849b5e25fSSatish Balay } 1689574b2666SHong Zhang s_browlengths[i] += nmask; 1690574b2666SHong Zhang 169149b5e25fSSatish Balay /* zero out the mask elements we set */ 169249b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 169349b5e25fSSatish Balay } 169449b5e25fSSatish Balay 169549b5e25fSSatish Balay /* create our matrix */ 16967fe2be48SHong Zhang ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr); 169749b5e25fSSatish Balay B = *A; 169849b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 169949b5e25fSSatish Balay 170049b5e25fSSatish Balay /* set matrix "i" values */ 170149b5e25fSSatish Balay a->i[0] = 0; 170249b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1703574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1704574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 170549b5e25fSSatish Balay } 17067fe2be48SHong Zhang a->s_nz = a->i[mbs]; 170749b5e25fSSatish Balay 170849b5e25fSSatish Balay /* read in nonzero values */ 170987828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 171049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 171149b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 171249b5e25fSSatish Balay 171349b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 171449b5e25fSSatish Balay nzcount = 0; jcount = 0; 171549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 171649b5e25fSSatish Balay nzcountb = nzcount; 171749b5e25fSSatish Balay nmask = 0; 171849b5e25fSSatish Balay for (j=0; j<bs; j++) { 171949b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 172049b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17212d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 172203630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 17232d703238SHong Zhang } 17242d703238SHong Zhang } 17252d703238SHong Zhang /* sort the masked values */ 17262d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 17272d703238SHong Zhang 17282d703238SHong Zhang /* set "j" values into matrix */ 17292d703238SHong Zhang maskcount = 1; 17302d703238SHong Zhang for (j=0; j<nmask; j++) { 173149b5e25fSSatish Balay a->j[jcount++] = masked[j]; 173249b5e25fSSatish Balay mask[masked[j]] = maskcount++; 173349b5e25fSSatish Balay } 1734574b2666SHong Zhang 173549b5e25fSSatish Balay /* set "a" values into matrix */ 173649b5e25fSSatish Balay ishift = bs2*a->i[i]; 173749b5e25fSSatish Balay for (j=0; j<bs; j++) { 173849b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 173949b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1740574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1741574b2666SHong Zhang if (tmp >= i){ 174249b5e25fSSatish Balay block = mask[tmp] - 1; 174349b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 174449b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1745574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1746574b2666SHong Zhang } 1747574b2666SHong Zhang nzcountb++; 174849b5e25fSSatish Balay } 174949b5e25fSSatish Balay } 175049b5e25fSSatish Balay /* zero out the mask elements we set */ 175149b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 175249b5e25fSSatish Balay } 175329bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 175449b5e25fSSatish Balay 175549b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1756574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 175749b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 175849b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 175949b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 176049b5e25fSSatish Balay 176149b5e25fSSatish Balay B->assembled = PETSC_TRUE; 176249b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 176349b5e25fSSatish Balay PetscFunctionReturn(0); 176449b5e25fSSatish Balay } 17658549e402SHong Zhang EXTERN_C_END 1766574b2666SHong Zhang 1767d06b337dSHong Zhang #undef __FUNCT__ 1768d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1769c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx) 1770d06b337dSHong Zhang { 1771d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1772d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1773d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 1774d06b337dSHong Zhang int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr; 1775d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1776d06b337dSHong Zhang 1777d06b337dSHong Zhang PetscFunctionBegin; 1778b965ef7fSBarry Smith its = its*lits; 177991723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 1780d06b337dSHong Zhang 1781d06b337dSHong Zhang if (bs > 1) 1782d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1783d06b337dSHong Zhang 1784d06b337dSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1785d06b337dSHong Zhang if (xx != bb) { 1786d06b337dSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1787d06b337dSHong Zhang } else { 1788d06b337dSHong Zhang b = x; 1789d06b337dSHong Zhang } 1790d06b337dSHong Zhang 1791d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1792d06b337dSHong Zhang 1793d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 17942798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1795d06b337dSHong Zhang for (i=0; i<m; i++) 1796d06b337dSHong Zhang t[i] = b[i]; 1797d06b337dSHong Zhang 1798d06b337dSHong Zhang for (i=0; i<m; i++){ 179944706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1800d06b337dSHong Zhang v = aa + ai[i] + 1; 1801d06b337dSHong Zhang vj = aj + ai[i] + 1; 1802d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1803d06b337dSHong Zhang x[i] = omega*t[i]/d; 1804d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 180544706875SHong Zhang PetscLogFlops(2*nz-1); 1806d06b337dSHong Zhang } 1807d06b337dSHong Zhang } 1808d06b337dSHong Zhang 18092798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1810d06b337dSHong Zhang for (i=0; i<m; i++) 1811d06b337dSHong Zhang t[i] = b[i]; 1812d06b337dSHong Zhang 1813d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1814d06b337dSHong Zhang v = aa + ai[i] + 1; 1815d06b337dSHong Zhang vj = aj + ai[i] + 1; 1816d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1817d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 181844706875SHong Zhang PetscLogFlops(2*nz-1); 1819d06b337dSHong Zhang } 1820d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1821d06b337dSHong Zhang d = *(aa + ai[i]); 1822d06b337dSHong Zhang v = aa + ai[i] + 1; 1823d06b337dSHong Zhang vj = aj + ai[i] + 1; 1824d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1825d06b337dSHong Zhang sum = t[i]; 1826d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 182744706875SHong Zhang PetscLogFlops(2*nz-1); 1828d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1829d06b337dSHong Zhang } 1830d06b337dSHong Zhang } 1831d06b337dSHong Zhang its--; 1832d06b337dSHong Zhang } 1833d06b337dSHong Zhang 1834d06b337dSHong Zhang while (its--) { 183544706875SHong Zhang /* 183644706875SHong Zhang forward sweep: 183744706875SHong Zhang for i=0,...,m-1: 183844706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 183944706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 184044706875SHong Zhang b = b - x[i]*U^T(i,:); 1841d06b337dSHong Zhang 184244706875SHong Zhang */ 18432798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1844d06b337dSHong Zhang for (i=0; i<m; i++) 1845d06b337dSHong Zhang t[i] = b[i]; 1846d06b337dSHong Zhang 1847d06b337dSHong Zhang for (i=0; i<m; i++){ 184844706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1849d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 1850d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 1851d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 1852d06b337dSHong Zhang sum = t[i]; 1853d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 1854d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1855d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 185644706875SHong Zhang PetscLogFlops(4*nz-2); 1857d06b337dSHong Zhang } 1858d06b337dSHong Zhang } 1859d06b337dSHong Zhang 18602798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 186144706875SHong Zhang /* 186244706875SHong Zhang backward sweep: 186344706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 186444706875SHong Zhang for i=m-1,...,0: 186544706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 186644706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 186744706875SHong Zhang */ 1868d06b337dSHong Zhang for (i=0; i<m; i++) 1869d06b337dSHong Zhang t[i] = b[i]; 1870d06b337dSHong Zhang 1871d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1872d06b337dSHong Zhang v = aa + ai[i] + 1; 1873d06b337dSHong Zhang vj = aj + ai[i] + 1; 1874d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1875d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 187644706875SHong Zhang PetscLogFlops(2*nz-1); 1877d06b337dSHong Zhang } 1878d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1879d06b337dSHong Zhang d = *(aa + ai[i]); 1880d06b337dSHong Zhang v = aa + ai[i] + 1; 1881d06b337dSHong Zhang vj = aj + ai[i] + 1; 1882d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1883d06b337dSHong Zhang sum = t[i]; 1884d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 188544706875SHong Zhang PetscLogFlops(2*nz-1); 1886d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1887d06b337dSHong Zhang } 1888d06b337dSHong Zhang } 1889d06b337dSHong Zhang } 1890d06b337dSHong Zhang 1891d06b337dSHong Zhang ierr = PetscFree(t); CHKERRQ(ierr); 1892d06b337dSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1893d06b337dSHong Zhang if (bb != xx) { 1894d06b337dSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1895d06b337dSHong Zhang } 18962798e883SHong Zhang 1897d06b337dSHong Zhang PetscFunctionReturn(0); 1898d06b337dSHong Zhang } 1899d06b337dSHong Zhang 1900d06b337dSHong Zhang 1901d06b337dSHong Zhang 1902d06b337dSHong Zhang 190349b5e25fSSatish Balay 190449b5e25fSSatish Balay 1905