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 */ 7ab9f2c04SSatish Balay #include "petscsys.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h" 949b5e25fSSatish Balay #include "src/vec/vecimpl.h" 1049b5e25fSSatish Balay #include "src/inline/spops.h" 11aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay #define CHUNKSIZE 10 1449b5e25fSSatish Balay 1549b5e25fSSatish Balay /* 1649b5e25fSSatish Balay Checks for missing diagonals 1749b5e25fSSatish Balay */ 184a2ae208SSatish Balay #undef __FUNCT__ 194a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 2049b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2149b5e25fSSatish Balay { 22045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2349b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2449b5e25fSSatish Balay 2549b5e25fSSatish Balay PetscFunctionBegin; 26045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2749b5e25fSSatish Balay diag = a->diag; 2849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2949b5e25fSSatish Balay if (jj[diag[i]] != i) { 3029bbc08cSBarry Smith SETERRQ1(1,"Matrix is missing diagonal number %d",i); 3149b5e25fSSatish Balay } 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; 4182502324SSatish Balay int i,j,*diag,m = a->mbs,ierr; 4249b5e25fSSatish Balay 4349b5e25fSSatish Balay PetscFunctionBegin; 4449b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4549b5e25fSSatish Balay 46b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 47b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 4849b5e25fSSatish Balay for (i=0; i<m; i++) { 4949b5e25fSSatish Balay diag[i] = a->i[i+1]; 5049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 5149b5e25fSSatish Balay if (a->j[j] == i) { 5249b5e25fSSatish Balay diag[i] = j; 5349b5e25fSSatish Balay break; 5449b5e25fSSatish Balay } 5549b5e25fSSatish Balay } 5649b5e25fSSatish Balay } 5749b5e25fSSatish Balay a->diag = diag; 5849b5e25fSSatish Balay PetscFunctionReturn(0); 5949b5e25fSSatish Balay } 6049b5e25fSSatish Balay 6149b5e25fSSatish Balay 6249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 6349b5e25fSSatish Balay 644a2ae208SSatish Balay #undef __FUNCT__ 654a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 6649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 6749b5e25fSSatish Balay { 6849b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6949b5e25fSSatish Balay 7049b5e25fSSatish Balay PetscFunctionBegin; 71045c9aa0SHong Zhang if (ia) { 7229bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 7349b5e25fSSatish Balay } 74045c9aa0SHong Zhang *nn = a->mbs; 7549b5e25fSSatish Balay PetscFunctionReturn(0); 7649b5e25fSSatish Balay } 7749b5e25fSSatish Balay 784a2ae208SSatish Balay #undef __FUNCT__ 794a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 80045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 8149b5e25fSSatish Balay { 8249b5e25fSSatish Balay PetscFunctionBegin; 8349b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 8429bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 8516cdd363SHong Zhang /* PetscFunctionReturn(0); */ 8649b5e25fSSatish Balay } 8749b5e25fSSatish Balay 884a2ae208SSatish Balay #undef __FUNCT__ 894a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 9049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9149b5e25fSSatish Balay { 9249b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 9349b5e25fSSatish Balay 9449b5e25fSSatish Balay PetscFunctionBegin; 9549b5e25fSSatish Balay *bs = sbaij->bs; 9649b5e25fSSatish Balay PetscFunctionReturn(0); 9749b5e25fSSatish Balay } 9849b5e25fSSatish Balay 994a2ae208SSatish Balay #undef __FUNCT__ 1004a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 10149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 10249b5e25fSSatish Balay { 10349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10449b5e25fSSatish Balay int ierr; 10549b5e25fSSatish Balay 10649b5e25fSSatish Balay PetscFunctionBegin; 10749b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 108b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 10949b5e25fSSatish Balay #endif 11049b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 11149b5e25fSSatish Balay if (!a->singlemalloc) { 11249b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 11363c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 11449b5e25fSSatish Balay } 11549b5e25fSSatish Balay if (a->row) { 11649b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 11749b5e25fSSatish Balay } 11849b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 11949b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 12049b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 12149b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 12249b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 12349b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 12449b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1251a3463dfSHong Zhang 1261a3463dfSHong Zhang if (a->inew){ 1271a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1281a3463dfSHong Zhang a->inew = 0; 1291a3463dfSHong Zhang } 13049b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 13149b5e25fSSatish Balay PetscFunctionReturn(0); 13249b5e25fSSatish Balay } 13349b5e25fSSatish Balay 1344a2ae208SSatish Balay #undef __FUNCT__ 1354a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 13649b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 13749b5e25fSSatish Balay { 138045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 13949b5e25fSSatish Balay 14049b5e25fSSatish Balay PetscFunctionBegin; 1414d9d31abSKris Buschelman switch (op) { 1424d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1434d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1444d9d31abSKris Buschelman break; 1454d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1464d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1474d9d31abSKris Buschelman break; 1484d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1494d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1504d9d31abSKris Buschelman break; 1514d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1524d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1554d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1584d9d31abSKris Buschelman a->nonew = 1; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1614d9d31abSKris Buschelman a->nonew = -1; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1644d9d31abSKris Buschelman a->nonew = -2; 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1674d9d31abSKris Buschelman a->nonew = 0; 1684d9d31abSKris Buschelman break; 1694d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1704d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1714d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1724d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1734d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 174d03495bdSKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 175b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1764d9d31abSKris Buschelman break; 1774d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 17829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1794d9d31abSKris Buschelman break; 1804d9d31abSKris Buschelman default: 18129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 1824d9d31abSKris Buschelman break; 18349b5e25fSSatish Balay } 18449b5e25fSSatish Balay PetscFunctionReturn(0); 18549b5e25fSSatish Balay } 18649b5e25fSSatish Balay 1874a2ae208SSatish Balay #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 18987828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v) 19049b5e25fSSatish Balay { 19149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 19282502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr; 19349b5e25fSSatish Balay MatScalar *aa,*aa_i; 19487828ca2SBarry Smith PetscScalar *v_i; 19549b5e25fSSatish Balay 19649b5e25fSSatish Balay PetscFunctionBegin; 19749b5e25fSSatish Balay bs = a->bs; 19849b5e25fSSatish Balay ai = a->i; 19949b5e25fSSatish Balay aj = a->j; 20049b5e25fSSatish Balay aa = a->a; 20149b5e25fSSatish Balay bs2 = a->bs2; 20249b5e25fSSatish Balay 203c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 20449b5e25fSSatish Balay 20549b5e25fSSatish Balay bn = row/bs; /* Block number */ 20649b5e25fSSatish Balay bp = row % bs; /* Block position */ 20749b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 20849b5e25fSSatish Balay *ncols = bs*M; 20949b5e25fSSatish Balay 21049b5e25fSSatish Balay if (v) { 21149b5e25fSSatish Balay *v = 0; 21249b5e25fSSatish Balay if (*ncols) { 21387828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 21449b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 21549b5e25fSSatish Balay v_i = *v + i*bs; 21649b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 21749b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 21849b5e25fSSatish Balay } 21949b5e25fSSatish Balay } 22049b5e25fSSatish Balay } 22149b5e25fSSatish Balay 22249b5e25fSSatish Balay if (cols) { 22349b5e25fSSatish Balay *cols = 0; 22449b5e25fSSatish Balay if (*ncols) { 22582502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 22649b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22749b5e25fSSatish Balay cols_i = *cols + i*bs; 22849b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 22949b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 23049b5e25fSSatish Balay } 23149b5e25fSSatish Balay } 23249b5e25fSSatish Balay } 23349b5e25fSSatish Balay 23449b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2355ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2365ddb2528SHong Zhang #ifdef column_search 23749b5e25fSSatish Balay v_i = *v + M*bs; 23849b5e25fSSatish Balay cols_i = *cols + M*bs; 23949b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 24049b5e25fSSatish Balay M = ai[i+1] - ai[i]; 24149b5e25fSSatish Balay for (j=0; j<M; j++){ 24249b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 24349b5e25fSSatish Balay if (itmp == bn){ 24449b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 24549b5e25fSSatish Balay for (k=0; k<bs; k++) { 24649b5e25fSSatish Balay *cols_i++ = i*bs+k; 24749b5e25fSSatish Balay *v_i++ = aa_i[k]; 24849b5e25fSSatish Balay } 24949b5e25fSSatish Balay *ncols += bs; 25049b5e25fSSatish Balay break; 25149b5e25fSSatish Balay } 25249b5e25fSSatish Balay } 25349b5e25fSSatish Balay } 2545ddb2528SHong Zhang #endif 25549b5e25fSSatish Balay 25649b5e25fSSatish Balay PetscFunctionReturn(0); 25749b5e25fSSatish Balay } 25849b5e25fSSatish Balay 2594a2ae208SSatish Balay #undef __FUNCT__ 2604a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 26187828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 26249b5e25fSSatish Balay { 26349b5e25fSSatish Balay int ierr; 26449b5e25fSSatish Balay 26549b5e25fSSatish Balay PetscFunctionBegin; 26649b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 26749b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 26849b5e25fSSatish Balay PetscFunctionReturn(0); 26949b5e25fSSatish Balay } 27049b5e25fSSatish Balay 2714a2ae208SSatish Balay #undef __FUNCT__ 2724a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 27349b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 27449b5e25fSSatish Balay { 27549b5e25fSSatish Balay PetscFunctionBegin; 27629bbc08cSBarry Smith SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called"); 27716cdd363SHong Zhang /* PetscFunctionReturn(0); */ 27849b5e25fSSatish Balay } 27949b5e25fSSatish Balay 2804a2ae208SSatish Balay #undef __FUNCT__ 2814a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary" 282b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer) 28349b5e25fSSatish Balay { 28449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 28549b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 28687828ca2SBarry Smith PetscScalar *aa; 28749b5e25fSSatish Balay FILE *file; 28849b5e25fSSatish Balay 28949b5e25fSSatish Balay PetscFunctionBegin; 290b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 29182502324SSatish Balay ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 29249b5e25fSSatish Balay col_lens[0] = MAT_COOKIE; 29349b5e25fSSatish Balay 294c464158bSHong Zhang col_lens[1] = A->m; 295c464158bSHong Zhang col_lens[2] = A->m; 29649b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 29749b5e25fSSatish Balay 29849b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 29949b5e25fSSatish Balay count = 0; 30049b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30149b5e25fSSatish Balay for (j=0; j<bs; j++) { 30249b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 30349b5e25fSSatish Balay } 30449b5e25fSSatish Balay } 305c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 30649b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 30749b5e25fSSatish Balay 30849b5e25fSSatish Balay /* store column indices (zero start index) */ 30982502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 31049b5e25fSSatish Balay count = 0; 31149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31249b5e25fSSatish Balay for (j=0; j<bs; j++) { 31349b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31449b5e25fSSatish Balay for (l=0; l<bs; l++) { 31549b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 31649b5e25fSSatish Balay } 31749b5e25fSSatish Balay } 31849b5e25fSSatish Balay } 31949b5e25fSSatish Balay } 32049b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 32149b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 32249b5e25fSSatish Balay 32349b5e25fSSatish Balay /* store nonzero values */ 32487828ca2SBarry Smith ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 32549b5e25fSSatish Balay count = 0; 32649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 32749b5e25fSSatish Balay for (j=0; j<bs; j++) { 32849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 32949b5e25fSSatish Balay for (l=0; l<bs; l++) { 33049b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 33149b5e25fSSatish Balay } 33249b5e25fSSatish Balay } 33349b5e25fSSatish Balay } 33449b5e25fSSatish Balay } 33549b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 33649b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 33749b5e25fSSatish Balay 338b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 33949b5e25fSSatish Balay if (file) { 34049b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 34149b5e25fSSatish Balay } 34249b5e25fSSatish Balay PetscFunctionReturn(0); 34349b5e25fSSatish Balay } 34449b5e25fSSatish Balay 3454a2ae208SSatish Balay #undef __FUNCT__ 3464a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 347b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 34849b5e25fSSatish Balay { 34949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 350fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 351fb9695e5SSatish Balay char *name; 352f3ef73ceSBarry Smith PetscViewerFormat format; 35349b5e25fSSatish Balay 35449b5e25fSSatish Balay PetscFunctionBegin; 35580fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 356b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 357fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) { 358b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 359fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 36029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 361fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 362b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36349b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 36449b5e25fSSatish Balay for (j=0; j<bs; j++) { 365b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 36649b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 36749b5e25fSSatish Balay for (l=0; l<bs; l++) { 36849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 36949b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 370b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 37149b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37249b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 373b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 37449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37549b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 376b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37749b5e25fSSatish Balay } 37849b5e25fSSatish Balay #else 37949b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 380b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38149b5e25fSSatish Balay } 38249b5e25fSSatish Balay #endif 38349b5e25fSSatish Balay } 38449b5e25fSSatish Balay } 385b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 38649b5e25fSSatish Balay } 38749b5e25fSSatish Balay } 388b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 38949b5e25fSSatish Balay } else { 390b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 39149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 39249b5e25fSSatish Balay for (j=0; j<bs; j++) { 393b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 39449b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 39549b5e25fSSatish Balay for (l=0; l<bs; l++) { 39649b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 39749b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 398b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 39949b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40049b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 401b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 40249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40349b5e25fSSatish Balay } else { 404b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40549b5e25fSSatish Balay } 40649b5e25fSSatish Balay #else 407b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 40849b5e25fSSatish Balay #endif 40949b5e25fSSatish Balay } 41049b5e25fSSatish Balay } 411b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 41249b5e25fSSatish Balay } 41349b5e25fSSatish Balay } 414b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 41549b5e25fSSatish Balay } 416b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 41749b5e25fSSatish Balay PetscFunctionReturn(0); 41849b5e25fSSatish Balay } 41949b5e25fSSatish Balay 4204a2ae208SSatish Balay #undef __FUNCT__ 4214a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 422b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 42349b5e25fSSatish Balay { 42449b5e25fSSatish Balay Mat A = (Mat) Aa; 42549b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 42649b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 42749b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 42849b5e25fSSatish Balay MatScalar *aa; 42949b5e25fSSatish Balay MPI_Comm comm; 430b0a32e0cSBarry Smith PetscViewer viewer; 43149b5e25fSSatish Balay 43249b5e25fSSatish Balay PetscFunctionBegin; 43349b5e25fSSatish Balay /* 43449b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 43549b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 43649b5e25fSSatish Balay rest should return immediately. 43749b5e25fSSatish Balay */ 43849b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 43949b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 44049b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 44149b5e25fSSatish Balay 44249b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 44349b5e25fSSatish Balay 444b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 445b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 44649b5e25fSSatish Balay 44749b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 448b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 44949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 451c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45349b5e25fSSatish Balay aa = a->a + j*bs2; 45449b5e25fSSatish Balay for (k=0; k<bs; k++) { 45549b5e25fSSatish Balay for (l=0; l<bs; l++) { 45649b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 457b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 45849b5e25fSSatish Balay } 45949b5e25fSSatish Balay } 46049b5e25fSSatish Balay } 46149b5e25fSSatish Balay } 462b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 46349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 465c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 46649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46749b5e25fSSatish Balay aa = a->a + j*bs2; 46849b5e25fSSatish Balay for (k=0; k<bs; k++) { 46949b5e25fSSatish Balay for (l=0; l<bs; l++) { 47049b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 471b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 47249b5e25fSSatish Balay } 47349b5e25fSSatish Balay } 47449b5e25fSSatish Balay } 47549b5e25fSSatish Balay } 47649b5e25fSSatish Balay 477b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 47849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 47949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 480c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 48149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 48249b5e25fSSatish Balay aa = a->a + j*bs2; 48349b5e25fSSatish Balay for (k=0; k<bs; k++) { 48449b5e25fSSatish Balay for (l=0; l<bs; l++) { 48549b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 486b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 48749b5e25fSSatish Balay } 48849b5e25fSSatish Balay } 48949b5e25fSSatish Balay } 49049b5e25fSSatish Balay } 49149b5e25fSSatish Balay PetscFunctionReturn(0); 49249b5e25fSSatish Balay } 49349b5e25fSSatish Balay 4944a2ae208SSatish Balay #undef __FUNCT__ 4954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 496b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 49749b5e25fSSatish Balay { 49849b5e25fSSatish Balay int ierr; 49949b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 500b0a32e0cSBarry Smith PetscDraw draw; 50149b5e25fSSatish Balay PetscTruth isnull; 50249b5e25fSSatish Balay 50349b5e25fSSatish Balay PetscFunctionBegin; 50449b5e25fSSatish Balay 505b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 506b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 50749b5e25fSSatish Balay 50849b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 509c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 51049b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 511b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 512b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 51349b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 51449b5e25fSSatish Balay PetscFunctionReturn(0); 51549b5e25fSSatish Balay } 51649b5e25fSSatish Balay 5174a2ae208SSatish Balay #undef __FUNCT__ 5184a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 519b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 52049b5e25fSSatish Balay { 52149b5e25fSSatish Balay int ierr; 52249b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 52349b5e25fSSatish Balay 52449b5e25fSSatish Balay PetscFunctionBegin; 525b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 526b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 527fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 528fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 52949b5e25fSSatish Balay if (issocket) { 53029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 53149b5e25fSSatish Balay } else if (isascii){ 53249b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 53349b5e25fSSatish Balay } else if (isbinary) { 53449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 53549b5e25fSSatish Balay } else if (isdraw) { 53649b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 53749b5e25fSSatish Balay } else { 53829bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 53949b5e25fSSatish Balay } 54049b5e25fSSatish Balay PetscFunctionReturn(0); 54149b5e25fSSatish Balay } 54249b5e25fSSatish Balay 54349b5e25fSSatish Balay 5444a2ae208SSatish Balay #undef __FUNCT__ 5454a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 54687828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v) 54749b5e25fSSatish Balay { 548045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 54949b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 55049b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 55149b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 55249b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 55349b5e25fSSatish Balay 55449b5e25fSSatish Balay PetscFunctionBegin; 55549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 55649b5e25fSSatish Balay row = im[k]; brow = row/bs; 55729bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 558c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 55949b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 56049b5e25fSSatish Balay nrow = ailen[brow]; 56149b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 56229bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 563c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 56449b5e25fSSatish Balay col = in[l] ; 56549b5e25fSSatish Balay bcol = col/bs; 56649b5e25fSSatish Balay cidx = col%bs; 56749b5e25fSSatish Balay ridx = row%bs; 56849b5e25fSSatish Balay high = nrow; 56949b5e25fSSatish Balay low = 0; /* assume unsorted */ 57049b5e25fSSatish Balay while (high-low > 5) { 57149b5e25fSSatish Balay t = (low+high)/2; 57249b5e25fSSatish Balay if (rp[t] > bcol) high = t; 57349b5e25fSSatish Balay else low = t; 57449b5e25fSSatish Balay } 57549b5e25fSSatish Balay for (i=low; i<high; i++) { 57649b5e25fSSatish Balay if (rp[i] > bcol) break; 57749b5e25fSSatish Balay if (rp[i] == bcol) { 57849b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 57949b5e25fSSatish Balay goto finished; 58049b5e25fSSatish Balay } 58149b5e25fSSatish Balay } 58249b5e25fSSatish Balay *v++ = zero; 58349b5e25fSSatish Balay finished:; 58449b5e25fSSatish Balay } 58549b5e25fSSatish Balay } 58649b5e25fSSatish Balay PetscFunctionReturn(0); 58749b5e25fSSatish Balay } 58849b5e25fSSatish Balay 58949b5e25fSSatish Balay 5904a2ae208SSatish Balay #undef __FUNCT__ 5914a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 59287828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 59349b5e25fSSatish Balay { 59449b5e25fSSatish Balay PetscFunctionBegin; 59529bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format"); 59616cdd363SHong Zhang /* PetscFunctionReturn(0); */ 59749b5e25fSSatish Balay } 59849b5e25fSSatish Balay 5994a2ae208SSatish Balay #undef __FUNCT__ 6004a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 60149b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 60249b5e25fSSatish Balay { 60349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 60449b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 605c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 60649b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 60749b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 60849b5e25fSSatish Balay 60949b5e25fSSatish Balay PetscFunctionBegin; 61049b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 61149b5e25fSSatish Balay 61249b5e25fSSatish Balay if (m) rmax = ailen[0]; 61349b5e25fSSatish Balay for (i=1; i<mbs; i++) { 61449b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 61549b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 61649b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 61749b5e25fSSatish Balay if (fshift) { 61849b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 61949b5e25fSSatish Balay N = ailen[i]; 62049b5e25fSSatish Balay for (j=0; j<N; j++) { 62149b5e25fSSatish Balay ip[j-fshift] = ip[j]; 62249b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 62349b5e25fSSatish Balay } 62449b5e25fSSatish Balay } 62549b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 62649b5e25fSSatish Balay } 62749b5e25fSSatish Balay if (mbs) { 62849b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 62949b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 63049b5e25fSSatish Balay } 63149b5e25fSSatish Balay /* reset ilen and imax for each row */ 63249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 63349b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 63449b5e25fSSatish Balay } 63549b5e25fSSatish Balay a->s_nz = ai[mbs]; 63649b5e25fSSatish Balay 63749b5e25fSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 63849b5e25fSSatish Balay if (fshift && a->diag) { 63949b5e25fSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 640b0a32e0cSBarry Smith PetscLogObjectMemory(A,-(m+1)*sizeof(int)); 64149b5e25fSSatish Balay a->diag = 0; 64249b5e25fSSatish Balay } 643b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 644c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 645b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 64649b5e25fSSatish Balay a->reallocs); 647b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 64849b5e25fSSatish Balay a->reallocs = 0; 64949b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 65049b5e25fSSatish Balay 65149b5e25fSSatish Balay PetscFunctionReturn(0); 65249b5e25fSSatish Balay } 65349b5e25fSSatish Balay 65449b5e25fSSatish Balay /* 65549b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 65649b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 65749b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 65849b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 65949b5e25fSSatish Balay */ 6604a2ae208SSatish Balay #undef __FUNCT__ 6614a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 66249b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 66349b5e25fSSatish Balay { 66449b5e25fSSatish Balay int i,j,k,row; 66549b5e25fSSatish Balay PetscTruth flg; 66649b5e25fSSatish Balay 66749b5e25fSSatish Balay PetscFunctionBegin; 66849b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 66949b5e25fSSatish Balay row = idx[i]; 67049b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 67149b5e25fSSatish Balay sizes[j] = 1; 67249b5e25fSSatish Balay i++; 67349b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 67449b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 67549b5e25fSSatish Balay i++; 67649b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 67749b5e25fSSatish Balay flg = PETSC_TRUE; 67849b5e25fSSatish Balay for (k=1; k<bs; k++) { 67949b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 68049b5e25fSSatish Balay flg = PETSC_FALSE; 68149b5e25fSSatish Balay break; 68249b5e25fSSatish Balay } 68349b5e25fSSatish Balay } 68449b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 68549b5e25fSSatish Balay sizes[j] = bs; 68649b5e25fSSatish Balay i+= bs; 68749b5e25fSSatish Balay } else { 68849b5e25fSSatish Balay sizes[j] = 1; 68949b5e25fSSatish Balay i++; 69049b5e25fSSatish Balay } 69149b5e25fSSatish Balay } 69249b5e25fSSatish Balay } 69349b5e25fSSatish Balay *bs_max = j; 69449b5e25fSSatish Balay PetscFunctionReturn(0); 69549b5e25fSSatish Balay } 69649b5e25fSSatish Balay 6974a2ae208SSatish Balay #undef __FUNCT__ 6984a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 69987828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag) 70049b5e25fSSatish Balay { 70149b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data; 70249b5e25fSSatish Balay int ierr,i,j,k,count,is_n,*is_idx,*rows; 70349b5e25fSSatish Balay int bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max; 70487828ca2SBarry Smith PetscScalar zero = 0.0; 70549b5e25fSSatish Balay MatScalar *aa; 70649b5e25fSSatish Balay 70749b5e25fSSatish Balay PetscFunctionBegin; 70849b5e25fSSatish Balay /* Make a copy of the IS and sort it */ 70949b5e25fSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 71049b5e25fSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 71149b5e25fSSatish Balay 71249b5e25fSSatish Balay /* allocate memory for rows,sizes */ 713b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 71449b5e25fSSatish Balay sizes = rows + is_n; 71549b5e25fSSatish Balay 71649b5e25fSSatish Balay /* initialize copy IS values to rows, and sort them */ 71749b5e25fSSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 71849b5e25fSSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 71949b5e25fSSatish Balay if (sbaij->keepzeroedrows) { /* do not change nonzero structure */ 72049b5e25fSSatish Balay for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */ 72149b5e25fSSatish Balay bs_max = is_n; /* bs_max: num. of contiguous block row in the row */ 72249b5e25fSSatish Balay } else { 72349b5e25fSSatish Balay ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 72449b5e25fSSatish Balay } 72549b5e25fSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 72649b5e25fSSatish Balay 72749b5e25fSSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 72849b5e25fSSatish Balay row = rows[j]; /* row to be zeroed */ 729c464158bSHong Zhang if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 73049b5e25fSSatish Balay count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */ 73149b5e25fSSatish Balay aa = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs); 73249b5e25fSSatish Balay if (sizes[i] == bs && !sbaij->keepzeroedrows) { 73349b5e25fSSatish Balay if (diag) { 73449b5e25fSSatish Balay if (sbaij->ilen[row/bs] > 0) { 73549b5e25fSSatish Balay sbaij->ilen[row/bs] = 1; 73649b5e25fSSatish Balay sbaij->j[sbaij->i[row/bs]] = row/bs; 73749b5e25fSSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 73849b5e25fSSatish Balay } 73949b5e25fSSatish Balay /* Now insert all the diagoanl values for this bs */ 74049b5e25fSSatish Balay for (k=0; k<bs; k++) { 74149b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 74249b5e25fSSatish Balay } 74349b5e25fSSatish Balay } else { /* (!diag) */ 74449b5e25fSSatish Balay sbaij->ilen[row/bs] = 0; 74549b5e25fSSatish Balay } /* end (!diag) */ 74649b5e25fSSatish Balay } else { /* (sizes[i] != bs), broken block */ 74749b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG) 74829bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 74949b5e25fSSatish Balay #endif 75049b5e25fSSatish Balay for (k=0; k<count; k++) { 75149b5e25fSSatish Balay aa[0] = zero; 75249b5e25fSSatish Balay aa+=bs; 75349b5e25fSSatish Balay } 75449b5e25fSSatish Balay if (diag) { 75549b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 75649b5e25fSSatish Balay } 75749b5e25fSSatish Balay } 75849b5e25fSSatish Balay } 75949b5e25fSSatish Balay 76049b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 76149b5e25fSSatish Balay ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76249b5e25fSSatish Balay ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76349b5e25fSSatish Balay PetscFunctionReturn(0); 76449b5e25fSSatish Balay } 76549b5e25fSSatish Balay 76649b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 76749b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 76849b5e25fSSatish Balay */ 76949b5e25fSSatish Balay 7704a2ae208SSatish Balay #undef __FUNCT__ 7714a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 77287828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 77349b5e25fSSatish Balay { 77449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 77549b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 77649b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 77749b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 77849b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 77949b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 78049b5e25fSSatish Balay 78149b5e25fSSatish Balay PetscFunctionBegin; 78249b5e25fSSatish Balay 78349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 78449b5e25fSSatish Balay row = im[k]; /* row number */ 78549b5e25fSSatish Balay brow = row/bs; /* block row number */ 78649b5e25fSSatish Balay if (row < 0) continue; 78749b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 788c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 78949b5e25fSSatish Balay #endif 79049b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 79149b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 79249b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 79349b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 79449b5e25fSSatish Balay low = 0; 79549b5e25fSSatish Balay 79649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 79749b5e25fSSatish Balay if (in[l] < 0) continue; 79849b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 799c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 80049b5e25fSSatish Balay #endif 80149b5e25fSSatish Balay col = in[l]; 80249b5e25fSSatish Balay bcol = col/bs; /* block col number */ 80349b5e25fSSatish Balay 80449b5e25fSSatish Balay if (brow <= bcol){ 80549b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8068549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 80749b5e25fSSatish Balay /* element value a(k,l) */ 80849b5e25fSSatish Balay if (roworiented) { 80949b5e25fSSatish Balay value = v[l + k*n]; 81049b5e25fSSatish Balay } else { 81149b5e25fSSatish Balay value = v[k + l*m]; 81249b5e25fSSatish Balay } 81349b5e25fSSatish Balay 81449b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 81549b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 81649b5e25fSSatish Balay while (high-low > 7) { 81749b5e25fSSatish Balay t = (low+high)/2; 81849b5e25fSSatish Balay if (rp[t] > bcol) high = t; 81949b5e25fSSatish Balay else low = t; 82049b5e25fSSatish Balay } 82149b5e25fSSatish Balay for (i=low; i<high; i++) { 82249b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 82349b5e25fSSatish Balay if (rp[i] > bcol) break; 82449b5e25fSSatish Balay if (rp[i] == bcol) { 82549b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 82649b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 82749b5e25fSSatish Balay else *bap = value; 8288549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8298549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8308549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8318549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8328549e402SHong Zhang else *bap = value; 8338549e402SHong Zhang } 83449b5e25fSSatish Balay goto noinsert1; 83549b5e25fSSatish Balay } 83649b5e25fSSatish Balay } 83749b5e25fSSatish Balay 83849b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 83929bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 84049b5e25fSSatish Balay if (nrow >= rmax) { 84149b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 84249b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 84349b5e25fSSatish Balay MatScalar *new_a; 84449b5e25fSSatish Balay 84529bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 84649b5e25fSSatish Balay 84749b5e25fSSatish Balay /* Malloc new storage space */ 84849b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 84982502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 85049b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 85149b5e25fSSatish Balay new_i = new_j + new_nz; 85249b5e25fSSatish Balay 85349b5e25fSSatish Balay /* copy over old data into new slots */ 85449b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 85549b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 85649b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 85749b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 85849b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 85949b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 86049b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 86149b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 86249b5e25fSSatish Balay /* free up old matrix storage */ 86349b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 86449b5e25fSSatish Balay if (!a->singlemalloc) { 86549b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 86649b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 86749b5e25fSSatish Balay } 86849b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 86949b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 87049b5e25fSSatish Balay 87149b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 87249b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 873b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 87449b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 87549b5e25fSSatish Balay a->reallocs++; 87649b5e25fSSatish Balay a->s_nz++; 87749b5e25fSSatish Balay } 87849b5e25fSSatish Balay 87949b5e25fSSatish Balay N = nrow++ - 1; 88049b5e25fSSatish Balay /* shift up all the later entries in this row */ 88149b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 88249b5e25fSSatish Balay rp[ii+1] = rp[ii]; 88349b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 88449b5e25fSSatish Balay } 88549b5e25fSSatish Balay if (N>=i) { 88649b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 88749b5e25fSSatish Balay } 88849b5e25fSSatish Balay rp[i] = bcol; 88949b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 89049b5e25fSSatish Balay noinsert1:; 89149b5e25fSSatish Balay low = i; 89249b5e25fSSatish Balay /* } */ 8938549e402SHong Zhang } 89449b5e25fSSatish Balay } /* end of if .. if.. */ 89549b5e25fSSatish Balay } /* end of loop over added columns */ 89649b5e25fSSatish Balay ailen[brow] = nrow; 89749b5e25fSSatish Balay } /* end of loop over added rows */ 89849b5e25fSSatish Balay 89949b5e25fSSatish Balay PetscFunctionReturn(0); 90049b5e25fSSatish Balay } 90149b5e25fSSatish Balay 902915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*); 903915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal); 90449b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 90549b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 90649b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 90749b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 90849b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 90987828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat); 91049b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 91149b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 91249b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 91349b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 91449b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 91549b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 91613a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 91749b5e25fSSatish Balay 91849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 91949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 92049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 92149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 92249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 92349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 92449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 92549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 92649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 92749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 92849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 92949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 93049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 93149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 93249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 93349b5e25fSSatish Balay 93407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 93507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 93607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 93707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 93807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 93907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 94007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 94107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 94207f98182SSatish Balay 9435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9445f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9455f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9465f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9475f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9485f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9495f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 95057d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 95149b5e25fSSatish Balay 95249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 95349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 95449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 95549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 95649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 95749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 95849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 95949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 96049b5e25fSSatish Balay 96149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 96249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 96349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 96449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 96549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 96649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 96749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 96849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 96949b5e25fSSatish Balay 9704d101231SSatish Balay #ifdef HAVE_MatICCFactor 9711a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 9724a2ae208SSatish Balay #undef __FUNCT__ 9734d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 9744d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level) 97549b5e25fSSatish Balay { 9764ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 97749b5e25fSSatish Balay Mat outA; 97849b5e25fSSatish Balay int ierr; 97949b5e25fSSatish Balay PetscTruth row_identity,col_identity; 98049b5e25fSSatish Balay 98149b5e25fSSatish Balay PetscFunctionBegin; 9821a3463dfSHong Zhang /* 98329bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 98449b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 98549b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 98649b5e25fSSatish Balay if (!row_identity || !col_identity) { 98729bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 98849b5e25fSSatish Balay } 9891a3463dfSHong Zhang */ 99049b5e25fSSatish Balay 99149b5e25fSSatish Balay outA = inA; 9921260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 99349b5e25fSSatish Balay 99449b5e25fSSatish Balay if (!a->diag) { 9951a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 99649b5e25fSSatish Balay } 99749b5e25fSSatish Balay /* 99849b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 99949b5e25fSSatish Balay for ILU(0) factorization with natural ordering 100049b5e25fSSatish Balay */ 100149b5e25fSSatish Balay switch (a->bs) { 100249b5e25fSSatish Balay case 1: 10030fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 100407f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 10054d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 100649b5e25fSSatish Balay case 2: 10071a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10081a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10091a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10104d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 101149b5e25fSSatish Balay break; 101249b5e25fSSatish Balay case 3: 10131a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10141a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10151a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10164d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 101749b5e25fSSatish Balay break; 101849b5e25fSSatish Balay case 4: 10191a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10201a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10211a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10224d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 102349b5e25fSSatish Balay break; 102449b5e25fSSatish Balay case 5: 10251a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10261a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10271a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10284d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 102949b5e25fSSatish Balay break; 103049b5e25fSSatish Balay case 6: 10311a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10321a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10331a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10344d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 103549b5e25fSSatish Balay break; 103649b5e25fSSatish Balay case 7: 10371a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10381a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10391a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10404d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 104149b5e25fSSatish Balay break; 104249b5e25fSSatish Balay default: 104349b5e25fSSatish Balay a->row = row; 10441a3463dfSHong Zhang a->icol = col; 104549b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 104649b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 104749b5e25fSSatish Balay 104849b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 104949b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1050b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 105149b5e25fSSatish Balay 105249b5e25fSSatish Balay if (!a->solve_work) { 105387828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 105487828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 105549b5e25fSSatish Balay } 105649b5e25fSSatish Balay } 105749b5e25fSSatish Balay 10581a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 105949b5e25fSSatish Balay 106049b5e25fSSatish Balay PetscFunctionReturn(0); 106149b5e25fSSatish Balay } 1062950f1e5bSHong Zhang #endif 1063950f1e5bSHong Zhang 10644a2ae208SSatish Balay #undef __FUNCT__ 10654a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 106649b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 106749b5e25fSSatish Balay { 106849b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 106949b5e25fSSatish Balay MPI_Comm comm = A->comm; 107049b5e25fSSatish Balay int ierr; 107149b5e25fSSatish Balay 107249b5e25fSSatish Balay PetscFunctionBegin; 107349b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 107449b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 107549b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 107649b5e25fSSatish Balay PetscFunctionReturn(0); 107749b5e25fSSatish Balay } 107849b5e25fSSatish Balay 107949b5e25fSSatish Balay EXTERN_C_BEGIN 10804a2ae208SSatish Balay #undef __FUNCT__ 10814a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 108249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 108349b5e25fSSatish Balay { 1084045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 108549b5e25fSSatish Balay int i,nz,n; 108649b5e25fSSatish Balay 108749b5e25fSSatish Balay PetscFunctionBegin; 1088045c9aa0SHong Zhang nz = baij->s_maxnz; 1089c464158bSHong Zhang n = mat->n; 109049b5e25fSSatish Balay for (i=0; i<nz; i++) { 109149b5e25fSSatish Balay baij->j[i] = indices[i]; 109249b5e25fSSatish Balay } 1093045c9aa0SHong Zhang baij->s_nz = nz; 109449b5e25fSSatish Balay for (i=0; i<n; i++) { 109549b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 109649b5e25fSSatish Balay } 109749b5e25fSSatish Balay 109849b5e25fSSatish Balay PetscFunctionReturn(0); 109949b5e25fSSatish Balay } 110049b5e25fSSatish Balay EXTERN_C_END 110149b5e25fSSatish Balay 11024a2ae208SSatish Balay #undef __FUNCT__ 11034a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 110449b5e25fSSatish Balay /*@ 110519585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 110649b5e25fSSatish Balay in the matrix. 110749b5e25fSSatish Balay 110849b5e25fSSatish Balay Input Parameters: 110919585528SSatish Balay + mat - the SeqSBAIJ matrix 111049b5e25fSSatish Balay - indices - the column indices 111149b5e25fSSatish Balay 111249b5e25fSSatish Balay Level: advanced 111349b5e25fSSatish Balay 111449b5e25fSSatish Balay Notes: 111549b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 111649b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 111749b5e25fSSatish Balay of the MatSetValues() operation. 111849b5e25fSSatish Balay 111949b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 112019585528SSatish Balay MatCreateSeqSBAIJ(). 112149b5e25fSSatish Balay 1122ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 112349b5e25fSSatish Balay 1124ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 112549b5e25fSSatish Balay @*/ 112649b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 112749b5e25fSSatish Balay { 112849b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 112949b5e25fSSatish Balay 113049b5e25fSSatish Balay PetscFunctionBegin; 113149b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1132c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 113349b5e25fSSatish Balay if (f) { 113449b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 113549b5e25fSSatish Balay } else { 113629bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 113749b5e25fSSatish Balay } 113849b5e25fSSatish Balay PetscFunctionReturn(0); 113949b5e25fSSatish Balay } 114049b5e25fSSatish Balay 11414a2ae208SSatish Balay #undef __FUNCT__ 11424a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1143273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1144273d9f13SBarry Smith { 1145273d9f13SBarry Smith int ierr; 1146273d9f13SBarry Smith 1147273d9f13SBarry Smith PetscFunctionBegin; 1148273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1149273d9f13SBarry Smith PetscFunctionReturn(0); 1150273d9f13SBarry Smith } 1151273d9f13SBarry Smith 115249b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 115349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 115449b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 115549b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 115649b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 115749b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 115849b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 115949b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 116049b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 116149b5e25fSSatish Balay 0, 116249b5e25fSSatish Balay 0, 116349b5e25fSSatish Balay 0, 116449b5e25fSSatish Balay 0, 11655f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1166*d06b337dSHong Zhang MatRelax_SeqSBAIJ, 116749b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 116849b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 116949b5e25fSSatish Balay MatEqual_SeqSBAIJ, 117049b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 117149b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 117249b5e25fSSatish Balay MatNorm_SeqSBAIJ, 117349b5e25fSSatish Balay 0, 117449b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 117549b5e25fSSatish Balay 0, 117649b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 117749b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 117849b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 117949b5e25fSSatish Balay 0, 118049b5e25fSSatish Balay 0, 11815f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 11825f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1183273d9f13SBarry Smith MatSetUpPreallocation_SeqSBAIJ, 1184c464158bSHong Zhang 0, 11854d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 118649b5e25fSSatish Balay 0, 118749b5e25fSSatish Balay 0, 118849b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 118949b5e25fSSatish Balay 0, 119049b5e25fSSatish Balay 0, 119149b5e25fSSatish Balay 0, 1192950f1e5bSHong Zhang 0, 119349b5e25fSSatish Balay 0, 119449b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 119549b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 119649b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 119749b5e25fSSatish Balay 0, 119849b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 119949b5e25fSSatish Balay MatScale_SeqSBAIJ, 120049b5e25fSSatish Balay 0, 120149b5e25fSSatish Balay 0, 120249b5e25fSSatish Balay 0, 120349b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 120449b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 120549b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 120649b5e25fSSatish Balay 0, 120749b5e25fSSatish Balay 0, 120849b5e25fSSatish Balay 0, 120949b5e25fSSatish Balay 0, 121049b5e25fSSatish Balay 0, 121149b5e25fSSatish Balay 0, 121249b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 121349b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 121449b5e25fSSatish Balay 0, 121549b5e25fSSatish Balay 0, 12168a124369SBarry Smith MatGetPetscMaps_Petsc, 1217d959ec07SHong Zhang 0, 1218d959ec07SHong Zhang 0, 1219d959ec07SHong Zhang 0, 1220d959ec07SHong Zhang 0, 1221d959ec07SHong Zhang 0, 1222d959ec07SHong Zhang 0, 122324d5174aSHong Zhang MatGetRowMax_SeqSBAIJ}; 122449b5e25fSSatish Balay 122549b5e25fSSatish Balay EXTERN_C_BEGIN 12264a2ae208SSatish Balay #undef __FUNCT__ 12274a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 122849b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 122949b5e25fSSatish Balay { 12304afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1231c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 123249b5e25fSSatish Balay int ierr; 123349b5e25fSSatish Balay 123449b5e25fSSatish Balay PetscFunctionBegin; 123549b5e25fSSatish Balay if (aij->nonew != 1) { 123629bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 123749b5e25fSSatish Balay } 123849b5e25fSSatish Balay 123949b5e25fSSatish Balay /* allocate space for values if not already there */ 124049b5e25fSSatish Balay if (!aij->saved_values) { 124187828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 124249b5e25fSSatish Balay } 124349b5e25fSSatish Balay 124449b5e25fSSatish Balay /* copy values over */ 124587828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 124649b5e25fSSatish Balay PetscFunctionReturn(0); 124749b5e25fSSatish Balay } 124849b5e25fSSatish Balay EXTERN_C_END 124949b5e25fSSatish Balay 125049b5e25fSSatish Balay EXTERN_C_BEGIN 12514a2ae208SSatish Balay #undef __FUNCT__ 12524a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 125349b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 125449b5e25fSSatish Balay { 12554afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1256c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 125749b5e25fSSatish Balay 125849b5e25fSSatish Balay PetscFunctionBegin; 125949b5e25fSSatish Balay if (aij->nonew != 1) { 126029bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 126149b5e25fSSatish Balay } 126249b5e25fSSatish Balay if (!aij->saved_values) { 126329bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 126449b5e25fSSatish Balay } 126549b5e25fSSatish Balay 126649b5e25fSSatish Balay /* copy values over */ 126787828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 126849b5e25fSSatish Balay PetscFunctionReturn(0); 126949b5e25fSSatish Balay } 127049b5e25fSSatish Balay EXTERN_C_END 127149b5e25fSSatish Balay 12728549e402SHong Zhang EXTERN_C_BEGIN 12734a2ae208SSatish Balay #undef __FUNCT__ 12744a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ" 1275c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1276c464158bSHong Zhang { 1277c464158bSHong Zhang Mat_SeqSBAIJ *b; 12780c955e93SHong Zhang int ierr,size; 1279c464158bSHong Zhang 1280c464158bSHong Zhang PetscFunctionBegin; 1281c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1282c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1283273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1284273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1285c464158bSHong Zhang 1286b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1287b0a32e0cSBarry Smith B->data = (void*)b; 1288c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1289c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1290c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1291c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1292c464158bSHong Zhang B->factor = 0; 1293c464158bSHong Zhang B->lupivotthreshold = 1.0; 1294c464158bSHong Zhang B->mapping = 0; 1295c464158bSHong Zhang b->row = 0; 1296c464158bSHong Zhang b->icol = 0; 1297c464158bSHong Zhang b->reallocs = 0; 1298c464158bSHong Zhang b->saved_values = 0; 1299c464158bSHong Zhang 13008a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 13018a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1302c464158bSHong Zhang 1303c464158bSHong Zhang b->sorted = PETSC_FALSE; 1304c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1305c464158bSHong Zhang b->nonew = 0; 1306c464158bSHong Zhang b->diag = 0; 1307c464158bSHong Zhang b->solve_work = 0; 1308c464158bSHong Zhang b->mult_work = 0; 1309c464158bSHong Zhang b->spptr = 0; 1310c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1311c464158bSHong Zhang 1312c464158bSHong Zhang b->inew = 0; 1313c464158bSHong Zhang b->jnew = 0; 1314c464158bSHong Zhang b->anew = 0; 1315c464158bSHong Zhang b->a2anew = 0; 1316c464158bSHong Zhang b->permute = PETSC_FALSE; 1317c464158bSHong Zhang 1318c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1319c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1320c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1321c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1322c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1323c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1324c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1325c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1326c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1327c464158bSHong Zhang PetscFunctionReturn(0); 1328c464158bSHong Zhang } 13298549e402SHong Zhang EXTERN_C_END 1330c464158bSHong Zhang 13314a2ae208SSatish Balay #undef __FUNCT__ 13324a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 133349b5e25fSSatish Balay /*@C 1334273d9f13SBarry Smith MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 133549b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 133649b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 133749b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 133849b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 133949b5e25fSSatish Balay 1340c464158bSHong Zhang Collective on Mat 134149b5e25fSSatish Balay 134249b5e25fSSatish Balay Input Parameters: 1343c464158bSHong Zhang + A - the symmetric matrix 134449b5e25fSSatish Balay . bs - size of block 134549b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 134649b5e25fSSatish Balay - nnz - array containing the number of block nonzeros in the various block rows 134749b5e25fSSatish Balay (possibly different for each block row) or PETSC_NULL 134849b5e25fSSatish Balay 134949b5e25fSSatish Balay Options Database Keys: 135049b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 135149b5e25fSSatish Balay block calculations (much slower) 135249b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 135349b5e25fSSatish Balay 135449b5e25fSSatish Balay Level: intermediate 135549b5e25fSSatish Balay 135649b5e25fSSatish Balay Notes: 1357c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 135849b5e25fSSatish Balay storage. That is, the stored row and column indices can begin at 135949b5e25fSSatish Balay either one (as in Fortran) or zero. See the users' manual for details. 136049b5e25fSSatish Balay 136149b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 136249b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 136349b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 136449b5e25fSSatish Balay matrices. 136549b5e25fSSatish Balay 1366a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 136749b5e25fSSatish Balay @*/ 1368273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 136949b5e25fSSatish Balay { 1370c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13710c955e93SHong Zhang int i,len,ierr,mbs,bs2; 137249b5e25fSSatish Balay PetscTruth flg; 13734afc71dfSHong Zhang int s_nz; 137449b5e25fSSatish Balay 137549b5e25fSSatish Balay PetscFunctionBegin; 1376273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1377b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1378c464158bSHong Zhang mbs = B->m/bs; 137949b5e25fSSatish Balay bs2 = bs*bs; 138049b5e25fSSatish Balay 1381c464158bSHong Zhang if (mbs*bs != B->m) { 138229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 138349b5e25fSSatish Balay } 138449b5e25fSSatish Balay 1385435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1386435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 138749b5e25fSSatish Balay if (nnz) { 138849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 138929bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 139080fe4e49SBarry 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); 139149b5e25fSSatish Balay } 139249b5e25fSSatish Balay } 139349b5e25fSSatish Balay 1394b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 139549b5e25fSSatish Balay if (!flg) { 139649b5e25fSSatish Balay switch (bs) { 139749b5e25fSSatish Balay case 1: 13984ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 139949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 140049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 140149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 140249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 140349b5e25fSSatish Balay break; 140449b5e25fSSatish Balay case 2: 14054ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 140649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 140749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 140849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 140949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 141049b5e25fSSatish Balay break; 141149b5e25fSSatish Balay case 3: 14125f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 141349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 141449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 141549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 141649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 141749b5e25fSSatish Balay break; 141849b5e25fSSatish Balay case 4: 14195f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 142049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 142149b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 142249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 142349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 142449b5e25fSSatish Balay break; 142549b5e25fSSatish Balay case 5: 14265f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 142749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 142849b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 142949b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 143049b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 143149b5e25fSSatish Balay break; 143249b5e25fSSatish Balay case 6: 14335f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 143449b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 143549b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 143649b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 143749b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 143849b5e25fSSatish Balay break; 143949b5e25fSSatish Balay case 7: 1440de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 144149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 144249b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1443de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 144449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 144549b5e25fSSatish Balay break; 144649b5e25fSSatish Balay } 144749b5e25fSSatish Balay } 144849b5e25fSSatish Balay 144949b5e25fSSatish Balay b->mbs = mbs; 14504afc71dfSHong Zhang b->nbs = mbs; 1451b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 145249b5e25fSSatish Balay if (!nnz) { 1453435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 145449b5e25fSSatish Balay else if (nz <= 0) nz = 1; 145549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14568cef66ccSHong Zhang b->imax[i] = nz; 145749b5e25fSSatish Balay } 1458153ea458SHong Zhang nz = nz*mbs; /* total nz */ 145949b5e25fSSatish Balay } else { 146049b5e25fSSatish Balay nz = 0; 14618cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 146249b5e25fSSatish Balay } 14638cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 14648cef66ccSHong Zhang s_nz = nz; 146549b5e25fSSatish Balay 146649b5e25fSSatish Balay /* allocate the matrix space */ 1467c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 146882502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 146949b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 147049b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 147149b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 147249b5e25fSSatish Balay b->i = b->j + s_nz; 147349b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 147449b5e25fSSatish Balay 147549b5e25fSSatish Balay /* pointer to beginning of each row */ 147649b5e25fSSatish Balay b->i[0] = 0; 147749b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 147849b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 147949b5e25fSSatish Balay } 148049b5e25fSSatish Balay 148149b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 14825033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1483b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 148449b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 148549b5e25fSSatish Balay 148649b5e25fSSatish Balay b->bs = bs; 148749b5e25fSSatish Balay b->bs2 = bs2; 148849b5e25fSSatish Balay b->s_nz = 0; 148949b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1490153ea458SHong Zhang 149116cdd363SHong Zhang b->inew = 0; 149216cdd363SHong Zhang b->jnew = 0; 149316cdd363SHong Zhang b->anew = 0; 149416cdd363SHong Zhang b->a2anew = 0; 14951a3463dfSHong Zhang b->permute = PETSC_FALSE; 1496c464158bSHong Zhang PetscFunctionReturn(0); 1497c464158bSHong Zhang } 1498153ea458SHong Zhang 149949b5e25fSSatish Balay 15004a2ae208SSatish Balay #undef __FUNCT__ 15014a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1502c464158bSHong Zhang /*@C 1503c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1504c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1505c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1506c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1507c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 150849b5e25fSSatish Balay 1509c464158bSHong Zhang Collective on MPI_Comm 1510c464158bSHong Zhang 1511c464158bSHong Zhang Input Parameters: 1512c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1513c464158bSHong Zhang . bs - size of block 1514c464158bSHong Zhang . m - number of rows, or number of columns 1515c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1516c464158bSHong Zhang - nnz - array containing the number of block nonzeros in the various block rows 1517c464158bSHong Zhang (possibly different for each block row) or PETSC_NULL 1518c464158bSHong Zhang 1519c464158bSHong Zhang Output Parameter: 1520c464158bSHong Zhang . A - the symmetric matrix 1521c464158bSHong Zhang 1522c464158bSHong Zhang Options Database Keys: 1523c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1524c464158bSHong Zhang block calculations (much slower) 1525c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1526c464158bSHong Zhang 1527c464158bSHong Zhang Level: intermediate 1528c464158bSHong Zhang 1529c464158bSHong Zhang Notes: 1530c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 1531c464158bSHong Zhang storage. That is, the stored row and column indices can begin at 1532c464158bSHong Zhang either one (as in Fortran) or zero. See the users' manual for details. 1533c464158bSHong Zhang 1534c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1535c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1536c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1537c464158bSHong Zhang matrices. 1538c464158bSHong Zhang 1539c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1540c464158bSHong Zhang @*/ 1541c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1542c464158bSHong Zhang { 1543c464158bSHong Zhang int ierr; 1544c464158bSHong Zhang 1545c464158bSHong Zhang PetscFunctionBegin; 1546c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1547c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1548273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 154949b5e25fSSatish Balay PetscFunctionReturn(0); 155049b5e25fSSatish Balay } 155149b5e25fSSatish Balay 15524a2ae208SSatish Balay #undef __FUNCT__ 15534a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 155449b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 155549b5e25fSSatish Balay { 155649b5e25fSSatish Balay Mat C; 155749b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 155849b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 155949b5e25fSSatish Balay 156049b5e25fSSatish Balay PetscFunctionBegin; 156129bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 156249b5e25fSSatish Balay 156349b5e25fSSatish Balay *B = 0; 1564692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1565692f9cbeSHong Zhang ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 1566692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1567692f9cbeSHong Zhang 156849b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1569273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 157049b5e25fSSatish Balay C->factor = A->factor; 157149b5e25fSSatish Balay c->row = 0; 157249b5e25fSSatish Balay c->icol = 0; 157349b5e25fSSatish Balay c->saved_values = 0; 157449b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 157549b5e25fSSatish Balay C->assembled = PETSC_TRUE; 157649b5e25fSSatish Balay 157749b5e25fSSatish Balay c->bs = a->bs; 157849b5e25fSSatish Balay c->bs2 = a->bs2; 157949b5e25fSSatish Balay c->mbs = a->mbs; 158049b5e25fSSatish Balay c->nbs = a->nbs; 158149b5e25fSSatish Balay 1582b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1583b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 158449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 158549b5e25fSSatish Balay c->imax[i] = a->imax[i]; 158649b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 158749b5e25fSSatish Balay } 158849b5e25fSSatish Balay 158949b5e25fSSatish Balay /* allocate the matrix space */ 159049b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 159149b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 159282502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 159349b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 159449b5e25fSSatish Balay c->i = c->j + nz; 159549b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 159649b5e25fSSatish Balay if (mbs > 0) { 159749b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 159849b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 159949b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 160049b5e25fSSatish Balay } else { 160149b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 160249b5e25fSSatish Balay } 160349b5e25fSSatish Balay } 160449b5e25fSSatish Balay 1605b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 160649b5e25fSSatish Balay c->sorted = a->sorted; 160749b5e25fSSatish Balay c->roworiented = a->roworiented; 160849b5e25fSSatish Balay c->nonew = a->nonew; 160949b5e25fSSatish Balay 161049b5e25fSSatish Balay if (a->diag) { 1611b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1612b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 161349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 161449b5e25fSSatish Balay c->diag[i] = a->diag[i]; 161549b5e25fSSatish Balay } 161649b5e25fSSatish Balay } else c->diag = 0; 161749b5e25fSSatish Balay c->s_nz = a->s_nz; 161849b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 161949b5e25fSSatish Balay c->solve_work = 0; 162049b5e25fSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 162149b5e25fSSatish Balay c->mult_work = 0; 162249b5e25fSSatish Balay *B = C; 1623b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 162449b5e25fSSatish Balay PetscFunctionReturn(0); 162549b5e25fSSatish Balay } 162649b5e25fSSatish Balay 16278549e402SHong Zhang EXTERN_C_BEGIN 16284a2ae208SSatish Balay #undef __FUNCT__ 16294a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1630b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A) 163149b5e25fSSatish Balay { 163249b5e25fSSatish Balay Mat_SeqSBAIJ *a; 163349b5e25fSSatish Balay Mat B; 163449b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1635574b2666SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount; 163649b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 163749b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 163887828ca2SBarry Smith PetscScalar *aa; 163949b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 164049b5e25fSSatish Balay 164149b5e25fSSatish Balay PetscFunctionBegin; 1642b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 164349b5e25fSSatish Balay bs2 = bs*bs; 164449b5e25fSSatish Balay 164549b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 164629bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1647b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 164849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 164929bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 165049b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 165149b5e25fSSatish Balay 165249b5e25fSSatish Balay if (header[3] < 0) { 165329bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 165449b5e25fSSatish Balay } 165549b5e25fSSatish Balay 165629bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 165749b5e25fSSatish Balay 165849b5e25fSSatish Balay /* 165949b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 166049b5e25fSSatish Balay divisible by the blocksize 166149b5e25fSSatish Balay */ 166249b5e25fSSatish Balay mbs = M/bs; 166349b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 166449b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 166549b5e25fSSatish Balay else mbs++; 166649b5e25fSSatish Balay if (extra_rows) { 1667b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 166849b5e25fSSatish Balay } 166949b5e25fSSatish Balay 167049b5e25fSSatish Balay /* read in row lengths */ 1671b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 167249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 167349b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 167449b5e25fSSatish Balay 167549b5e25fSSatish Balay /* read in column indices */ 1676b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 167749b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 167849b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 167949b5e25fSSatish Balay 168049b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 168182502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 168282502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 1683574b2666SHong Zhang ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 168482502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 168549b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 168649b5e25fSSatish Balay masked = mask + mbs; 168749b5e25fSSatish Balay rowcount = 0; nzcount = 0; 168849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 168949b5e25fSSatish Balay nmask = 0; 169049b5e25fSSatish Balay for (j=0; j<bs; j++) { 169149b5e25fSSatish Balay kmax = rowlengths[rowcount]; 169249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 16932d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 169403630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 169549b5e25fSSatish Balay } 169649b5e25fSSatish Balay rowcount++; 169749b5e25fSSatish Balay } 1698574b2666SHong Zhang s_browlengths[i] += nmask; 1699574b2666SHong Zhang browlengths[i] = 2*s_browlengths[i]; 1700574b2666SHong Zhang 170149b5e25fSSatish Balay /* zero out the mask elements we set */ 170249b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 170349b5e25fSSatish Balay } 170449b5e25fSSatish Balay 170549b5e25fSSatish Balay /* create our matrix */ 170649b5e25fSSatish Balay ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 170749b5e25fSSatish Balay B = *A; 170849b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 170949b5e25fSSatish Balay 171049b5e25fSSatish Balay /* set matrix "i" values */ 171149b5e25fSSatish Balay a->i[0] = 0; 171249b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1713574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1714574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 171549b5e25fSSatish Balay } 171649b5e25fSSatish Balay a->s_nz = 0; 1717574b2666SHong Zhang for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i]; 171849b5e25fSSatish Balay 171949b5e25fSSatish Balay /* read in nonzero values */ 172087828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 172149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 172249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 172349b5e25fSSatish Balay 172449b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 172549b5e25fSSatish Balay nzcount = 0; jcount = 0; 172649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 172749b5e25fSSatish Balay nzcountb = nzcount; 172849b5e25fSSatish Balay nmask = 0; 172949b5e25fSSatish Balay for (j=0; j<bs; j++) { 173049b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 173149b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17322d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 173303630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 17342d703238SHong Zhang } 17352d703238SHong Zhang } 17362d703238SHong Zhang /* sort the masked values */ 17372d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 17382d703238SHong Zhang 17392d703238SHong Zhang /* set "j" values into matrix */ 17402d703238SHong Zhang maskcount = 1; 17412d703238SHong Zhang for (j=0; j<nmask; j++) { 174249b5e25fSSatish Balay a->j[jcount++] = masked[j]; 174349b5e25fSSatish Balay mask[masked[j]] = maskcount++; 174449b5e25fSSatish Balay } 1745574b2666SHong Zhang 174649b5e25fSSatish Balay /* set "a" values into matrix */ 174749b5e25fSSatish Balay ishift = bs2*a->i[i]; 174849b5e25fSSatish Balay for (j=0; j<bs; j++) { 174949b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 175049b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1751574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1752574b2666SHong Zhang if (tmp >= i){ 175349b5e25fSSatish Balay block = mask[tmp] - 1; 175449b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 175549b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1756574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1757574b2666SHong Zhang } 1758574b2666SHong Zhang nzcountb++; 175949b5e25fSSatish Balay } 176049b5e25fSSatish Balay } 176149b5e25fSSatish Balay /* zero out the mask elements we set */ 176249b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 176349b5e25fSSatish Balay } 176429bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 176549b5e25fSSatish Balay 176649b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 176749b5e25fSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1768574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 176949b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 177049b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 177149b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 177249b5e25fSSatish Balay 177349b5e25fSSatish Balay B->assembled = PETSC_TRUE; 177449b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 177549b5e25fSSatish Balay PetscFunctionReturn(0); 177649b5e25fSSatish Balay } 17778549e402SHong Zhang EXTERN_C_END 1778574b2666SHong Zhang 1779*d06b337dSHong Zhang #undef __FUNCT__ 1780*d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1781*d06b337dSHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,Vec xx) 1782*d06b337dSHong Zhang { 1783*d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1784*d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1785*d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 1786*d06b337dSHong Zhang int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr; 1787*d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1788*d06b337dSHong Zhang 1789*d06b337dSHong Zhang PetscFunctionBegin; 1790*d06b337dSHong Zhang 1791*d06b337dSHong Zhang if (bs > 1) 1792*d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1793*d06b337dSHong Zhang 1794*d06b337dSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1795*d06b337dSHong Zhang if (xx != bb) { 1796*d06b337dSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1797*d06b337dSHong Zhang } else { 1798*d06b337dSHong Zhang b = x; 1799*d06b337dSHong Zhang } 1800*d06b337dSHong Zhang 1801*d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1802*d06b337dSHong Zhang 1803*d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 1804*d06b337dSHong Zhang if (flag & SOR_FORWARD_SWEEP){ 1805*d06b337dSHong Zhang for (i=0; i<m; i++) 1806*d06b337dSHong Zhang t[i] = b[i]; 1807*d06b337dSHong Zhang 1808*d06b337dSHong Zhang for (i=0; i<m; i++){ 1809*d06b337dSHong Zhang d = *(aa + ai[i]); /* diag */ 1810*d06b337dSHong Zhang v = aa + ai[i] + 1; 1811*d06b337dSHong Zhang vj = aj + ai[i] + 1; 1812*d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1813*d06b337dSHong Zhang x[i] = omega*t[i]/d; 1814*d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 1815*d06b337dSHong Zhang } 1816*d06b337dSHong Zhang } 1817*d06b337dSHong Zhang 1818*d06b337dSHong Zhang if (flag & SOR_BACKWARD_SWEEP){ 1819*d06b337dSHong Zhang for (i=0; i<m; i++) 1820*d06b337dSHong Zhang t[i] = b[i]; 1821*d06b337dSHong Zhang 1822*d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1823*d06b337dSHong Zhang v = aa + ai[i] + 1; 1824*d06b337dSHong Zhang vj = aj + ai[i] + 1; 1825*d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1826*d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 1827*d06b337dSHong Zhang } 1828*d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1829*d06b337dSHong Zhang d = *(aa + ai[i]); 1830*d06b337dSHong Zhang v = aa + ai[i] + 1; 1831*d06b337dSHong Zhang vj = aj + ai[i] + 1; 1832*d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1833*d06b337dSHong Zhang sum = t[i]; 1834*d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 1835*d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1836*d06b337dSHong Zhang } 1837*d06b337dSHong Zhang } 1838*d06b337dSHong Zhang its--; 1839*d06b337dSHong Zhang } 1840*d06b337dSHong Zhang 1841*d06b337dSHong Zhang while (its--) { 1842*d06b337dSHong Zhang 1843*d06b337dSHong Zhang if (flag & SOR_FORWARD_SWEEP ){ 1844*d06b337dSHong Zhang for (i=0; i<m; i++) 1845*d06b337dSHong Zhang t[i] = b[i]; 1846*d06b337dSHong Zhang 1847*d06b337dSHong Zhang for (i=0; i<m; i++){ 1848*d06b337dSHong Zhang d = *(aa + ai[i]); /* diag */ 1849*d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 1850*d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 1851*d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 1852*d06b337dSHong Zhang sum = t[i]; 1853*d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 1854*d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1855*d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 1856*d06b337dSHong Zhang } 1857*d06b337dSHong Zhang } 1858*d06b337dSHong Zhang 1859*d06b337dSHong Zhang if (flag & SOR_BACKWARD_SWEEP){ 1860*d06b337dSHong Zhang for (i=0; i<m; i++) 1861*d06b337dSHong Zhang t[i] = b[i]; 1862*d06b337dSHong Zhang 1863*d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1864*d06b337dSHong Zhang v = aa + ai[i] + 1; 1865*d06b337dSHong Zhang vj = aj + ai[i] + 1; 1866*d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1867*d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 1868*d06b337dSHong Zhang } 1869*d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1870*d06b337dSHong Zhang d = *(aa + ai[i]); 1871*d06b337dSHong Zhang v = aa + ai[i] + 1; 1872*d06b337dSHong Zhang vj = aj + ai[i] + 1; 1873*d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1874*d06b337dSHong Zhang sum = t[i]; 1875*d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 1876*d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1877*d06b337dSHong Zhang } 1878*d06b337dSHong Zhang } 1879*d06b337dSHong Zhang } 1880*d06b337dSHong Zhang 1881*d06b337dSHong Zhang /* PetscLogFlops(m); -- to be put later */ 1882*d06b337dSHong Zhang ierr = PetscFree(t); CHKERRQ(ierr); 1883*d06b337dSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1884*d06b337dSHong Zhang if (bb != xx) { 1885*d06b337dSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1886*d06b337dSHong Zhang } 1887*d06b337dSHong Zhang PetscFunctionReturn(0); 1888*d06b337dSHong Zhang } 1889*d06b337dSHong Zhang 1890*d06b337dSHong Zhang 1891*d06b337dSHong Zhang 1892*d06b337dSHong Zhang 189349b5e25fSSatish Balay 189449b5e25fSSatish Balay 1895