173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 449b5e25fSSatish Balay Defines the basic matrix operations for the BAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/vec/vecimpl.h" 949b5e25fSSatish Balay #include "src/inline/spops.h" 10aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1149b5e25fSSatish Balay 1249b5e25fSSatish Balay #define CHUNKSIZE 10 1349b5e25fSSatish Balay 1449b5e25fSSatish Balay /* 1549b5e25fSSatish Balay Checks for missing diagonals 1649b5e25fSSatish Balay */ 174a2ae208SSatish Balay #undef __FUNCT__ 184a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 1949b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2049b5e25fSSatish Balay { 21045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2249b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2349b5e25fSSatish Balay 2449b5e25fSSatish Balay PetscFunctionBegin; 25045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2649b5e25fSSatish Balay diag = a->diag; 2749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 28e11b0e14SHong Zhang if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i); 2949b5e25fSSatish Balay } 3049b5e25fSSatish Balay PetscFunctionReturn(0); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay 334a2ae208SSatish Balay #undef __FUNCT__ 344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 3549b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3649b5e25fSSatish Balay { 37045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 38b424e231SHong Zhang int i,mbs = a->mbs,ierr; 3949b5e25fSSatish Balay 4049b5e25fSSatish Balay PetscFunctionBegin; 4149b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4249b5e25fSSatish Balay 43b424e231SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr); 44b424e231SHong Zhang PetscLogObjectMemory(A,(mbs+1)*sizeof(int)); 45b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4649b5e25fSSatish Balay PetscFunctionReturn(0); 4749b5e25fSSatish Balay } 4849b5e25fSSatish Balay 4949b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 5049b5e25fSSatish Balay 514a2ae208SSatish Balay #undef __FUNCT__ 524a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 5349b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 5449b5e25fSSatish Balay { 5549b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 5649b5e25fSSatish Balay 5749b5e25fSSatish Balay PetscFunctionBegin; 58045c9aa0SHong Zhang if (ia) { 5929bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 6049b5e25fSSatish Balay } 61045c9aa0SHong Zhang *nn = a->mbs; 6249b5e25fSSatish Balay PetscFunctionReturn(0); 6349b5e25fSSatish Balay } 6449b5e25fSSatish Balay 654a2ae208SSatish Balay #undef __FUNCT__ 664a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 67045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 6849b5e25fSSatish Balay { 6949b5e25fSSatish Balay PetscFunctionBegin; 7049b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 7129bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 7216cdd363SHong Zhang /* PetscFunctionReturn(0); */ 7349b5e25fSSatish Balay } 7449b5e25fSSatish Balay 754a2ae208SSatish Balay #undef __FUNCT__ 764a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 7749b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 7849b5e25fSSatish Balay { 7949b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 8049b5e25fSSatish Balay 8149b5e25fSSatish Balay PetscFunctionBegin; 8249b5e25fSSatish Balay *bs = sbaij->bs; 8349b5e25fSSatish Balay PetscFunctionReturn(0); 8449b5e25fSSatish Balay } 8549b5e25fSSatish Balay 864a2ae208SSatish Balay #undef __FUNCT__ 874a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 8849b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 8949b5e25fSSatish Balay { 9049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 9149b5e25fSSatish Balay int ierr; 9249b5e25fSSatish Balay 9349b5e25fSSatish Balay PetscFunctionBegin; 9449b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 95b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 9649b5e25fSSatish Balay #endif 9749b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 9849b5e25fSSatish Balay if (!a->singlemalloc) { 9949b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 10063c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 10149b5e25fSSatish Balay } 10249b5e25fSSatish Balay if (a->row) { 10349b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 10449b5e25fSSatish Balay } 10549b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 10649b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 10749b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 10849b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 10949b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 11049b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 11149b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1121a3463dfSHong Zhang 1131a3463dfSHong Zhang if (a->inew){ 1141a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1151a3463dfSHong Zhang a->inew = 0; 1161a3463dfSHong Zhang } 11749b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 11849b5e25fSSatish Balay PetscFunctionReturn(0); 11949b5e25fSSatish Balay } 12049b5e25fSSatish Balay 1214a2ae208SSatish Balay #undef __FUNCT__ 1224a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 12349b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 12449b5e25fSSatish Balay { 125045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 12649b5e25fSSatish Balay 12749b5e25fSSatish Balay PetscFunctionBegin; 1284d9d31abSKris Buschelman switch (op) { 1294d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1304d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1314d9d31abSKris Buschelman break; 1324d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1334d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1344d9d31abSKris Buschelman break; 1354d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1364d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1374d9d31abSKris Buschelman break; 1384d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1394d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1404d9d31abSKris Buschelman break; 1414d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1424d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1434d9d31abSKris Buschelman break; 1444d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1454d9d31abSKris Buschelman a->nonew = 1; 1464d9d31abSKris Buschelman break; 1474d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1484d9d31abSKris Buschelman a->nonew = -1; 1494d9d31abSKris Buschelman break; 1504d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1514d9d31abSKris Buschelman a->nonew = -2; 1524d9d31abSKris Buschelman break; 1534d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1544d9d31abSKris Buschelman a->nonew = 0; 1554d9d31abSKris Buschelman break; 1564d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1574d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1584d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1594d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1604d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 161d03495bdSKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 162b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1634d9d31abSKris Buschelman break; 1644d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 16529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1664d9d31abSKris Buschelman default: 16729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 16849b5e25fSSatish Balay } 16949b5e25fSSatish Balay PetscFunctionReturn(0); 17049b5e25fSSatish Balay } 17149b5e25fSSatish Balay 1724a2ae208SSatish Balay #undef __FUNCT__ 1734a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 17487828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v) 17549b5e25fSSatish Balay { 17649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 17782502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr; 17849b5e25fSSatish Balay MatScalar *aa,*aa_i; 17987828ca2SBarry Smith PetscScalar *v_i; 18049b5e25fSSatish Balay 18149b5e25fSSatish Balay PetscFunctionBegin; 18249b5e25fSSatish Balay bs = a->bs; 18349b5e25fSSatish Balay ai = a->i; 18449b5e25fSSatish Balay aj = a->j; 18549b5e25fSSatish Balay aa = a->a; 18649b5e25fSSatish Balay bs2 = a->bs2; 18749b5e25fSSatish Balay 188c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 18949b5e25fSSatish Balay 19049b5e25fSSatish Balay bn = row/bs; /* Block number */ 19149b5e25fSSatish Balay bp = row % bs; /* Block position */ 19249b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 19349b5e25fSSatish Balay *ncols = bs*M; 19449b5e25fSSatish Balay 19549b5e25fSSatish Balay if (v) { 19649b5e25fSSatish Balay *v = 0; 19749b5e25fSSatish Balay if (*ncols) { 19887828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 19949b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 20049b5e25fSSatish Balay v_i = *v + i*bs; 20149b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 20249b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 20349b5e25fSSatish Balay } 20449b5e25fSSatish Balay } 20549b5e25fSSatish Balay } 20649b5e25fSSatish Balay 20749b5e25fSSatish Balay if (cols) { 20849b5e25fSSatish Balay *cols = 0; 20949b5e25fSSatish Balay if (*ncols) { 21082502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 21149b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 21249b5e25fSSatish Balay cols_i = *cols + i*bs; 21349b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 21449b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 21549b5e25fSSatish Balay } 21649b5e25fSSatish Balay } 21749b5e25fSSatish Balay } 21849b5e25fSSatish Balay 21949b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2205ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2215ddb2528SHong Zhang #ifdef column_search 22249b5e25fSSatish Balay v_i = *v + M*bs; 22349b5e25fSSatish Balay cols_i = *cols + M*bs; 22449b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 22549b5e25fSSatish Balay M = ai[i+1] - ai[i]; 22649b5e25fSSatish Balay for (j=0; j<M; j++){ 22749b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 22849b5e25fSSatish Balay if (itmp == bn){ 22949b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 23049b5e25fSSatish Balay for (k=0; k<bs; k++) { 23149b5e25fSSatish Balay *cols_i++ = i*bs+k; 23249b5e25fSSatish Balay *v_i++ = aa_i[k]; 23349b5e25fSSatish Balay } 23449b5e25fSSatish Balay *ncols += bs; 23549b5e25fSSatish Balay break; 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay } 23849b5e25fSSatish Balay } 2395ddb2528SHong Zhang #endif 24049b5e25fSSatish Balay 24149b5e25fSSatish Balay PetscFunctionReturn(0); 24249b5e25fSSatish Balay } 24349b5e25fSSatish Balay 2444a2ae208SSatish Balay #undef __FUNCT__ 2454a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 24687828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 24749b5e25fSSatish Balay { 24849b5e25fSSatish Balay int ierr; 24949b5e25fSSatish Balay 25049b5e25fSSatish Balay PetscFunctionBegin; 25149b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 25249b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 25349b5e25fSSatish Balay PetscFunctionReturn(0); 25449b5e25fSSatish Balay } 25549b5e25fSSatish Balay 2564a2ae208SSatish Balay #undef __FUNCT__ 2574a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 25849b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 25949b5e25fSSatish Balay { 26049b5e25fSSatish Balay PetscFunctionBegin; 26129bbc08cSBarry Smith SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called"); 26216cdd363SHong Zhang /* PetscFunctionReturn(0); */ 26349b5e25fSSatish Balay } 26449b5e25fSSatish Balay 2654a2ae208SSatish Balay #undef __FUNCT__ 2664a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary" 267b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer) 26849b5e25fSSatish Balay { 26949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 27049b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 27187828ca2SBarry Smith PetscScalar *aa; 27249b5e25fSSatish Balay FILE *file; 27349b5e25fSSatish Balay 27449b5e25fSSatish Balay PetscFunctionBegin; 275b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 27682502324SSatish Balay ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 277552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 27849b5e25fSSatish Balay 279c464158bSHong Zhang col_lens[1] = A->m; 280c464158bSHong Zhang col_lens[2] = A->m; 28149b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 28249b5e25fSSatish Balay 28349b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 28449b5e25fSSatish Balay count = 0; 28549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 28649b5e25fSSatish Balay for (j=0; j<bs; j++) { 28749b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 28849b5e25fSSatish Balay } 28949b5e25fSSatish Balay } 290c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 29149b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 29249b5e25fSSatish Balay 29349b5e25fSSatish Balay /* store column indices (zero start index) */ 29482502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 29549b5e25fSSatish Balay count = 0; 29649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 29749b5e25fSSatish Balay for (j=0; j<bs; j++) { 29849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 29949b5e25fSSatish Balay for (l=0; l<bs; l++) { 30049b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 30149b5e25fSSatish Balay } 30249b5e25fSSatish Balay } 30349b5e25fSSatish Balay } 30449b5e25fSSatish Balay } 30549b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 30649b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 30749b5e25fSSatish Balay 30849b5e25fSSatish Balay /* store nonzero values */ 30987828ca2SBarry Smith ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);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 aa[count++] = a->a[bs2*k + l*bs + j]; 31649b5e25fSSatish Balay } 31749b5e25fSSatish Balay } 31849b5e25fSSatish Balay } 31949b5e25fSSatish Balay } 32049b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 32149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 32249b5e25fSSatish Balay 323b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 32449b5e25fSSatish Balay if (file) { 32549b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 32649b5e25fSSatish Balay } 32749b5e25fSSatish Balay PetscFunctionReturn(0); 32849b5e25fSSatish Balay } 32949b5e25fSSatish Balay 3304a2ae208SSatish Balay #undef __FUNCT__ 3314a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 332b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 33349b5e25fSSatish Balay { 33449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 335fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 336fb9695e5SSatish Balay char *name; 337f3ef73ceSBarry Smith PetscViewerFormat format; 33849b5e25fSSatish Balay 33949b5e25fSSatish Balay PetscFunctionBegin; 34080fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 341b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 342fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) { 343b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 344fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 34529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 346fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 347b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 34849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 34949b5e25fSSatish Balay for (j=0; j<bs; j++) { 350b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 35149b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 35249b5e25fSSatish Balay for (l=0; l<bs; l++) { 35349b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 35449b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 35562b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 35649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35749b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 35862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 35949b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36049b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 36162b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 36249b5e25fSSatish Balay } 36349b5e25fSSatish Balay #else 36449b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 36562b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 36649b5e25fSSatish Balay } 36749b5e25fSSatish Balay #endif 36849b5e25fSSatish Balay } 36949b5e25fSSatish Balay } 370b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 37149b5e25fSSatish Balay } 37249b5e25fSSatish Balay } 373b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 37449b5e25fSSatish Balay } else { 375b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 37649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 37749b5e25fSSatish Balay for (j=0; j<bs; j++) { 378b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 37949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 38049b5e25fSSatish Balay for (l=0; l<bs; l++) { 38149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 38249b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 38362b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 38449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38549b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 38662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 38749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38849b5e25fSSatish Balay } else { 38962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 39049b5e25fSSatish Balay } 39149b5e25fSSatish Balay #else 392b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 39349b5e25fSSatish Balay #endif 39449b5e25fSSatish Balay } 39549b5e25fSSatish Balay } 396b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 39749b5e25fSSatish Balay } 39849b5e25fSSatish Balay } 399b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 40049b5e25fSSatish Balay } 401b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 40249b5e25fSSatish Balay PetscFunctionReturn(0); 40349b5e25fSSatish Balay } 40449b5e25fSSatish Balay 4054a2ae208SSatish Balay #undef __FUNCT__ 4064a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 407b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 40849b5e25fSSatish Balay { 40949b5e25fSSatish Balay Mat A = (Mat) Aa; 41049b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 41149b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 41249b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 41349b5e25fSSatish Balay MatScalar *aa; 41449b5e25fSSatish Balay MPI_Comm comm; 415b0a32e0cSBarry Smith PetscViewer viewer; 41649b5e25fSSatish Balay 41749b5e25fSSatish Balay PetscFunctionBegin; 41849b5e25fSSatish Balay /* 41949b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 42049b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 42149b5e25fSSatish Balay rest should return immediately. 42249b5e25fSSatish Balay */ 42349b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 42449b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 42549b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 42649b5e25fSSatish Balay 42749b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 42849b5e25fSSatish Balay 429b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 430b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 43149b5e25fSSatish Balay 43249b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 433b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 43449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 43549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 436c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 43749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 43849b5e25fSSatish Balay aa = a->a + j*bs2; 43949b5e25fSSatish Balay for (k=0; k<bs; k++) { 44049b5e25fSSatish Balay for (l=0; l<bs; l++) { 44149b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 442b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 44349b5e25fSSatish Balay } 44449b5e25fSSatish Balay } 44549b5e25fSSatish Balay } 44649b5e25fSSatish Balay } 447b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 44849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 44949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 450c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45249b5e25fSSatish Balay aa = a->a + j*bs2; 45349b5e25fSSatish Balay for (k=0; k<bs; k++) { 45449b5e25fSSatish Balay for (l=0; l<bs; l++) { 45549b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 456b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 45749b5e25fSSatish Balay } 45849b5e25fSSatish Balay } 45949b5e25fSSatish Balay } 46049b5e25fSSatish Balay } 46149b5e25fSSatish Balay 462b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 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 PetscFunctionReturn(0); 47749b5e25fSSatish Balay } 47849b5e25fSSatish Balay 4794a2ae208SSatish Balay #undef __FUNCT__ 4804a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 481b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 48249b5e25fSSatish Balay { 48349b5e25fSSatish Balay int ierr; 48449b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 485b0a32e0cSBarry Smith PetscDraw draw; 48649b5e25fSSatish Balay PetscTruth isnull; 48749b5e25fSSatish Balay 48849b5e25fSSatish Balay PetscFunctionBegin; 48949b5e25fSSatish Balay 490b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 491b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 49249b5e25fSSatish Balay 49349b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 494c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 49549b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 496b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 497b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 49849b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 49949b5e25fSSatish Balay PetscFunctionReturn(0); 50049b5e25fSSatish Balay } 50149b5e25fSSatish Balay 5024a2ae208SSatish Balay #undef __FUNCT__ 5034a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 504b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 50549b5e25fSSatish Balay { 50649b5e25fSSatish Balay int ierr; 50749b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 50849b5e25fSSatish Balay 50949b5e25fSSatish Balay PetscFunctionBegin; 510b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 511b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 512fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 513fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 51449b5e25fSSatish Balay if (issocket) { 51529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 51649b5e25fSSatish Balay } else if (isascii){ 51749b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 51849b5e25fSSatish Balay } else if (isbinary) { 51949b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 52049b5e25fSSatish Balay } else if (isdraw) { 52149b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 52249b5e25fSSatish Balay } else { 52329bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 52449b5e25fSSatish Balay } 52549b5e25fSSatish Balay PetscFunctionReturn(0); 52649b5e25fSSatish Balay } 52749b5e25fSSatish Balay 52849b5e25fSSatish Balay 5294a2ae208SSatish Balay #undef __FUNCT__ 5304a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 53187828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v) 53249b5e25fSSatish Balay { 533045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 53449b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 53549b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 53649b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 53749b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 53849b5e25fSSatish Balay 53949b5e25fSSatish Balay PetscFunctionBegin; 54049b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 54149b5e25fSSatish Balay row = im[k]; brow = row/bs; 54229bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 543c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 54449b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 54549b5e25fSSatish Balay nrow = ailen[brow]; 54649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 54729bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 548c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 54949b5e25fSSatish Balay col = in[l] ; 55049b5e25fSSatish Balay bcol = col/bs; 55149b5e25fSSatish Balay cidx = col%bs; 55249b5e25fSSatish Balay ridx = row%bs; 55349b5e25fSSatish Balay high = nrow; 55449b5e25fSSatish Balay low = 0; /* assume unsorted */ 55549b5e25fSSatish Balay while (high-low > 5) { 55649b5e25fSSatish Balay t = (low+high)/2; 55749b5e25fSSatish Balay if (rp[t] > bcol) high = t; 55849b5e25fSSatish Balay else low = t; 55949b5e25fSSatish Balay } 56049b5e25fSSatish Balay for (i=low; i<high; i++) { 56149b5e25fSSatish Balay if (rp[i] > bcol) break; 56249b5e25fSSatish Balay if (rp[i] == bcol) { 56349b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 56449b5e25fSSatish Balay goto finished; 56549b5e25fSSatish Balay } 56649b5e25fSSatish Balay } 56749b5e25fSSatish Balay *v++ = zero; 56849b5e25fSSatish Balay finished:; 56949b5e25fSSatish Balay } 57049b5e25fSSatish Balay } 57149b5e25fSSatish Balay PetscFunctionReturn(0); 57249b5e25fSSatish Balay } 57349b5e25fSSatish Balay 57449b5e25fSSatish Balay 5754a2ae208SSatish Balay #undef __FUNCT__ 5764a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 57787828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 57849b5e25fSSatish Balay { 57949b5e25fSSatish Balay PetscFunctionBegin; 58029bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format"); 58149b5e25fSSatish Balay } 58249b5e25fSSatish Balay 5834a2ae208SSatish Balay #undef __FUNCT__ 5844a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 58549b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 58649b5e25fSSatish Balay { 58749b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 58849b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 589c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 59049b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 59149b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 59249b5e25fSSatish Balay 59349b5e25fSSatish Balay PetscFunctionBegin; 59449b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 59549b5e25fSSatish Balay 59649b5e25fSSatish Balay if (m) rmax = ailen[0]; 59749b5e25fSSatish Balay for (i=1; i<mbs; i++) { 59849b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 59949b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 60049b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 60149b5e25fSSatish Balay if (fshift) { 60249b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 60349b5e25fSSatish Balay N = ailen[i]; 60449b5e25fSSatish Balay for (j=0; j<N; j++) { 60549b5e25fSSatish Balay ip[j-fshift] = ip[j]; 60649b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 60749b5e25fSSatish Balay } 60849b5e25fSSatish Balay } 60949b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 61049b5e25fSSatish Balay } 61149b5e25fSSatish Balay if (mbs) { 61249b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 61349b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 61449b5e25fSSatish Balay } 61549b5e25fSSatish Balay /* reset ilen and imax for each row */ 61649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 61749b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 61849b5e25fSSatish Balay } 61949b5e25fSSatish Balay a->s_nz = ai[mbs]; 62049b5e25fSSatish Balay 621b424e231SHong Zhang /* diagonals may have moved, reset it */ 622b424e231SHong Zhang if (a->diag) { 623b424e231SHong Zhang ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr); 62449b5e25fSSatish Balay } 625b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 626c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 627b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 62849b5e25fSSatish Balay a->reallocs); 629b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 63049b5e25fSSatish Balay a->reallocs = 0; 63149b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 63249b5e25fSSatish Balay 63349b5e25fSSatish Balay PetscFunctionReturn(0); 63449b5e25fSSatish Balay } 63549b5e25fSSatish Balay 63649b5e25fSSatish Balay /* 63749b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 63849b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 63949b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 64049b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 64149b5e25fSSatish Balay */ 6424a2ae208SSatish Balay #undef __FUNCT__ 6434a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 64449b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 64549b5e25fSSatish Balay { 64649b5e25fSSatish Balay int i,j,k,row; 64749b5e25fSSatish Balay PetscTruth flg; 64849b5e25fSSatish Balay 64949b5e25fSSatish Balay PetscFunctionBegin; 65049b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 65149b5e25fSSatish Balay row = idx[i]; 65249b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 65349b5e25fSSatish Balay sizes[j] = 1; 65449b5e25fSSatish Balay i++; 65549b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 65649b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 65749b5e25fSSatish Balay i++; 65849b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 65949b5e25fSSatish Balay flg = PETSC_TRUE; 66049b5e25fSSatish Balay for (k=1; k<bs; k++) { 66149b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 66249b5e25fSSatish Balay flg = PETSC_FALSE; 66349b5e25fSSatish Balay break; 66449b5e25fSSatish Balay } 66549b5e25fSSatish Balay } 66649b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 66749b5e25fSSatish Balay sizes[j] = bs; 66849b5e25fSSatish Balay i+= bs; 66949b5e25fSSatish Balay } else { 67049b5e25fSSatish Balay sizes[j] = 1; 67149b5e25fSSatish Balay i++; 67249b5e25fSSatish Balay } 67349b5e25fSSatish Balay } 67449b5e25fSSatish Balay } 67549b5e25fSSatish Balay *bs_max = j; 67649b5e25fSSatish Balay PetscFunctionReturn(0); 67749b5e25fSSatish Balay } 67849b5e25fSSatish Balay 6794a2ae208SSatish Balay #undef __FUNCT__ 6804a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 68187828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag) 68249b5e25fSSatish Balay { 68349b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data; 68449b5e25fSSatish Balay int ierr,i,j,k,count,is_n,*is_idx,*rows; 68549b5e25fSSatish Balay int bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max; 68687828ca2SBarry Smith PetscScalar zero = 0.0; 68749b5e25fSSatish Balay MatScalar *aa; 68849b5e25fSSatish Balay 68949b5e25fSSatish Balay PetscFunctionBegin; 69049b5e25fSSatish Balay /* Make a copy of the IS and sort it */ 69149b5e25fSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 69249b5e25fSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 69349b5e25fSSatish Balay 69449b5e25fSSatish Balay /* allocate memory for rows,sizes */ 695b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 69649b5e25fSSatish Balay sizes = rows + is_n; 69749b5e25fSSatish Balay 69849b5e25fSSatish Balay /* initialize copy IS values to rows, and sort them */ 69949b5e25fSSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 70049b5e25fSSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 70149b5e25fSSatish Balay if (sbaij->keepzeroedrows) { /* do not change nonzero structure */ 70249b5e25fSSatish Balay for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */ 70349b5e25fSSatish Balay bs_max = is_n; /* bs_max: num. of contiguous block row in the row */ 70449b5e25fSSatish Balay } else { 70549b5e25fSSatish Balay ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 70649b5e25fSSatish Balay } 70749b5e25fSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 70849b5e25fSSatish Balay 70949b5e25fSSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 71049b5e25fSSatish Balay row = rows[j]; /* row to be zeroed */ 711c464158bSHong Zhang if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 71249b5e25fSSatish Balay count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */ 71349b5e25fSSatish Balay aa = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs); 71449b5e25fSSatish Balay if (sizes[i] == bs && !sbaij->keepzeroedrows) { 71549b5e25fSSatish Balay if (diag) { 71649b5e25fSSatish Balay if (sbaij->ilen[row/bs] > 0) { 71749b5e25fSSatish Balay sbaij->ilen[row/bs] = 1; 71849b5e25fSSatish Balay sbaij->j[sbaij->i[row/bs]] = row/bs; 71949b5e25fSSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 72049b5e25fSSatish Balay } 72149b5e25fSSatish Balay /* Now insert all the diagoanl values for this bs */ 72249b5e25fSSatish Balay for (k=0; k<bs; k++) { 72349b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 72449b5e25fSSatish Balay } 72549b5e25fSSatish Balay } else { /* (!diag) */ 72649b5e25fSSatish Balay sbaij->ilen[row/bs] = 0; 72749b5e25fSSatish Balay } /* end (!diag) */ 72849b5e25fSSatish Balay } else { /* (sizes[i] != bs), broken block */ 72949b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG) 73029bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 73149b5e25fSSatish Balay #endif 73249b5e25fSSatish Balay for (k=0; k<count; k++) { 73349b5e25fSSatish Balay aa[0] = zero; 73449b5e25fSSatish Balay aa+=bs; 73549b5e25fSSatish Balay } 73649b5e25fSSatish Balay if (diag) { 73749b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 73849b5e25fSSatish Balay } 73949b5e25fSSatish Balay } 74049b5e25fSSatish Balay } 74149b5e25fSSatish Balay 74249b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 74349b5e25fSSatish Balay ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74449b5e25fSSatish Balay ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74549b5e25fSSatish Balay PetscFunctionReturn(0); 74649b5e25fSSatish Balay } 74749b5e25fSSatish Balay 74849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 74949b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 75049b5e25fSSatish Balay */ 75149b5e25fSSatish Balay 7524a2ae208SSatish Balay #undef __FUNCT__ 7534a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 75487828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 75549b5e25fSSatish Balay { 75649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 75749b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 75849b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 75949b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 76049b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 76149b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 76249b5e25fSSatish Balay 76349b5e25fSSatish Balay PetscFunctionBegin; 76449b5e25fSSatish Balay 76549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 76649b5e25fSSatish Balay row = im[k]; /* row number */ 76749b5e25fSSatish Balay brow = row/bs; /* block row number */ 76849b5e25fSSatish Balay if (row < 0) continue; 76949b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 770c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 77149b5e25fSSatish Balay #endif 77249b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 77349b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 77449b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 77549b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 77649b5e25fSSatish Balay low = 0; 77749b5e25fSSatish Balay 77849b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 77949b5e25fSSatish Balay if (in[l] < 0) continue; 78049b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 781c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 78249b5e25fSSatish Balay #endif 78349b5e25fSSatish Balay col = in[l]; 78449b5e25fSSatish Balay bcol = col/bs; /* block col number */ 78549b5e25fSSatish Balay 78649b5e25fSSatish Balay if (brow <= bcol){ 78749b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 7888549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 78949b5e25fSSatish Balay /* element value a(k,l) */ 79049b5e25fSSatish Balay if (roworiented) { 79149b5e25fSSatish Balay value = v[l + k*n]; 79249b5e25fSSatish Balay } else { 79349b5e25fSSatish Balay value = v[k + l*m]; 79449b5e25fSSatish Balay } 79549b5e25fSSatish Balay 79649b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 79749b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 79849b5e25fSSatish Balay while (high-low > 7) { 79949b5e25fSSatish Balay t = (low+high)/2; 80049b5e25fSSatish Balay if (rp[t] > bcol) high = t; 80149b5e25fSSatish Balay else low = t; 80249b5e25fSSatish Balay } 80349b5e25fSSatish Balay for (i=low; i<high; i++) { 80449b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 80549b5e25fSSatish Balay if (rp[i] > bcol) break; 80649b5e25fSSatish Balay if (rp[i] == bcol) { 80749b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 80849b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 80949b5e25fSSatish Balay else *bap = value; 8108549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8118549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8128549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8138549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8148549e402SHong Zhang else *bap = value; 8158549e402SHong Zhang } 81649b5e25fSSatish Balay goto noinsert1; 81749b5e25fSSatish Balay } 81849b5e25fSSatish Balay } 81949b5e25fSSatish Balay 82049b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 82129bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 82249b5e25fSSatish Balay if (nrow >= rmax) { 82349b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 82449b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 82549b5e25fSSatish Balay MatScalar *new_a; 82649b5e25fSSatish Balay 82729bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 82849b5e25fSSatish Balay 82949b5e25fSSatish Balay /* Malloc new storage space */ 83049b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 83182502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 83249b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 83349b5e25fSSatish Balay new_i = new_j + new_nz; 83449b5e25fSSatish Balay 83549b5e25fSSatish Balay /* copy over old data into new slots */ 83649b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 83749b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 83849b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 83949b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 84049b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 84149b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 84249b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 84349b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 84449b5e25fSSatish Balay /* free up old matrix storage */ 84549b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 84649b5e25fSSatish Balay if (!a->singlemalloc) { 84749b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 84849b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 84949b5e25fSSatish Balay } 85049b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 85149b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 85249b5e25fSSatish Balay 85349b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 85449b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 855b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 85649b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 85749b5e25fSSatish Balay a->reallocs++; 85849b5e25fSSatish Balay a->s_nz++; 85949b5e25fSSatish Balay } 86049b5e25fSSatish Balay 86149b5e25fSSatish Balay N = nrow++ - 1; 86249b5e25fSSatish Balay /* shift up all the later entries in this row */ 86349b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 86449b5e25fSSatish Balay rp[ii+1] = rp[ii]; 86549b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 86649b5e25fSSatish Balay } 86749b5e25fSSatish Balay if (N>=i) { 86849b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 86949b5e25fSSatish Balay } 87049b5e25fSSatish Balay rp[i] = bcol; 87149b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 87249b5e25fSSatish Balay noinsert1:; 87349b5e25fSSatish Balay low = i; 87449b5e25fSSatish Balay /* } */ 8758549e402SHong Zhang } 87649b5e25fSSatish Balay } /* end of if .. if.. */ 87749b5e25fSSatish Balay } /* end of loop over added columns */ 87849b5e25fSSatish Balay ailen[brow] = nrow; 87949b5e25fSSatish Balay } /* end of loop over added rows */ 88049b5e25fSSatish Balay 88149b5e25fSSatish Balay PetscFunctionReturn(0); 88249b5e25fSSatish Balay } 88349b5e25fSSatish Balay 884915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*); 885915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal); 88649b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 88749b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 88849b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 88949b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 89049b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 89187828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat); 89249b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 89349b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 89449b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 89549b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 89649b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 89749b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 89813a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 89949b5e25fSSatish Balay 90049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 90149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 90249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 90349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 90449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 90549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 90649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 90749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 90849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 90949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 91049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 91149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 91249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 91349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 91449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 91549b5e25fSSatish Balay 91607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 91707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 91807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 91907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 92007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 92107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 92207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 92307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 92407f98182SSatish Balay 9255f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9265f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9275f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9285f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9295f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9305f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9315f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 93257d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 93349b5e25fSSatish Balay 93449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 93549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 93649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 93749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 93849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 93949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 94049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 94149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 94249b5e25fSSatish Balay 94349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 94449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 94549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 94649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 94749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 94849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 94949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 95049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 95149b5e25fSSatish Balay 9524d101231SSatish Balay #ifdef HAVE_MatICCFactor 9531a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 9544a2ae208SSatish Balay #undef __FUNCT__ 9554d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 9564d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level) 95749b5e25fSSatish Balay { 9584ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 95949b5e25fSSatish Balay Mat outA; 96049b5e25fSSatish Balay int ierr; 96149b5e25fSSatish Balay PetscTruth row_identity,col_identity; 96249b5e25fSSatish Balay 96349b5e25fSSatish Balay PetscFunctionBegin; 9641a3463dfSHong Zhang /* 96529bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 96649b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 96749b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 96849b5e25fSSatish Balay if (!row_identity || !col_identity) { 96929bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 97049b5e25fSSatish Balay } 9711a3463dfSHong Zhang */ 97249b5e25fSSatish Balay 97349b5e25fSSatish Balay outA = inA; 9741260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 97549b5e25fSSatish Balay 97649b5e25fSSatish Balay if (!a->diag) { 9771a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 97849b5e25fSSatish Balay } 97949b5e25fSSatish Balay /* 98049b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 98149b5e25fSSatish Balay for ILU(0) factorization with natural ordering 98249b5e25fSSatish Balay */ 98349b5e25fSSatish Balay switch (a->bs) { 98449b5e25fSSatish Balay case 1: 9850fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 98607f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 9874d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 98849b5e25fSSatish Balay case 2: 9891a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9901a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9911a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9924d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 99349b5e25fSSatish Balay break; 99449b5e25fSSatish Balay case 3: 9951a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 9961a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9971a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 9984d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 99949b5e25fSSatish Balay break; 100049b5e25fSSatish Balay case 4: 10011a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10021a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10031a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10044d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 100549b5e25fSSatish Balay break; 100649b5e25fSSatish Balay case 5: 10071a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10081a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10091a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10104d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 101149b5e25fSSatish Balay break; 101249b5e25fSSatish Balay case 6: 10131a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10141a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10151a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10164d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 101749b5e25fSSatish Balay break; 101849b5e25fSSatish Balay case 7: 10191a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10201a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10211a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10224d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 102349b5e25fSSatish Balay break; 102449b5e25fSSatish Balay default: 102549b5e25fSSatish Balay a->row = row; 10261a3463dfSHong Zhang a->icol = col; 102749b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 102849b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 102949b5e25fSSatish Balay 103049b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 103149b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1032b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 103349b5e25fSSatish Balay 103449b5e25fSSatish Balay if (!a->solve_work) { 103587828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 103687828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 103749b5e25fSSatish Balay } 103849b5e25fSSatish Balay } 103949b5e25fSSatish Balay 10401a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 104149b5e25fSSatish Balay 104249b5e25fSSatish Balay PetscFunctionReturn(0); 104349b5e25fSSatish Balay } 1044950f1e5bSHong Zhang #endif 1045950f1e5bSHong Zhang 10464a2ae208SSatish Balay #undef __FUNCT__ 10474a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 104849b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 104949b5e25fSSatish Balay { 105049b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 105149b5e25fSSatish Balay MPI_Comm comm = A->comm; 105249b5e25fSSatish Balay int ierr; 105349b5e25fSSatish Balay 105449b5e25fSSatish Balay PetscFunctionBegin; 105549b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 105649b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 105749b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 105849b5e25fSSatish Balay PetscFunctionReturn(0); 105949b5e25fSSatish Balay } 106049b5e25fSSatish Balay 106149b5e25fSSatish Balay EXTERN_C_BEGIN 10624a2ae208SSatish Balay #undef __FUNCT__ 10634a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 106449b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 106549b5e25fSSatish Balay { 1066045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 106749b5e25fSSatish Balay int i,nz,n; 106849b5e25fSSatish Balay 106949b5e25fSSatish Balay PetscFunctionBegin; 1070045c9aa0SHong Zhang nz = baij->s_maxnz; 1071c464158bSHong Zhang n = mat->n; 107249b5e25fSSatish Balay for (i=0; i<nz; i++) { 107349b5e25fSSatish Balay baij->j[i] = indices[i]; 107449b5e25fSSatish Balay } 1075045c9aa0SHong Zhang baij->s_nz = nz; 107649b5e25fSSatish Balay for (i=0; i<n; i++) { 107749b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 107849b5e25fSSatish Balay } 107949b5e25fSSatish Balay 108049b5e25fSSatish Balay PetscFunctionReturn(0); 108149b5e25fSSatish Balay } 108249b5e25fSSatish Balay EXTERN_C_END 108349b5e25fSSatish Balay 10844a2ae208SSatish Balay #undef __FUNCT__ 10854a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 108649b5e25fSSatish Balay /*@ 108719585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 108849b5e25fSSatish Balay in the matrix. 108949b5e25fSSatish Balay 109049b5e25fSSatish Balay Input Parameters: 109119585528SSatish Balay + mat - the SeqSBAIJ matrix 109249b5e25fSSatish Balay - indices - the column indices 109349b5e25fSSatish Balay 109449b5e25fSSatish Balay Level: advanced 109549b5e25fSSatish Balay 109649b5e25fSSatish Balay Notes: 109749b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 109849b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 109949b5e25fSSatish Balay of the MatSetValues() operation. 110049b5e25fSSatish Balay 110149b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 110219585528SSatish Balay MatCreateSeqSBAIJ(). 110349b5e25fSSatish Balay 1104ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 110549b5e25fSSatish Balay 1106ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 110749b5e25fSSatish Balay @*/ 110849b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 110949b5e25fSSatish Balay { 111049b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 111149b5e25fSSatish Balay 111249b5e25fSSatish Balay PetscFunctionBegin; 111349b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1114c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 111549b5e25fSSatish Balay if (f) { 111649b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 111749b5e25fSSatish Balay } else { 111829bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 111949b5e25fSSatish Balay } 112049b5e25fSSatish Balay PetscFunctionReturn(0); 112149b5e25fSSatish Balay } 112249b5e25fSSatish Balay 11234a2ae208SSatish Balay #undef __FUNCT__ 11244a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1125273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1126273d9f13SBarry Smith { 1127273d9f13SBarry Smith int ierr; 1128273d9f13SBarry Smith 1129273d9f13SBarry Smith PetscFunctionBegin; 1130273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1131273d9f13SBarry Smith PetscFunctionReturn(0); 1132273d9f13SBarry Smith } 1133273d9f13SBarry Smith 113449b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 113549b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 113649b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 113749b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 113849b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 113949b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 114049b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 114149b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 114249b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 114349b5e25fSSatish Balay 0, 114449b5e25fSSatish Balay 0, 114549b5e25fSSatish Balay 0, 114649b5e25fSSatish Balay 0, 11475f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1148d06b337dSHong Zhang MatRelax_SeqSBAIJ, 114949b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 115049b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 115149b5e25fSSatish Balay MatEqual_SeqSBAIJ, 115249b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 115349b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 115449b5e25fSSatish Balay MatNorm_SeqSBAIJ, 115549b5e25fSSatish Balay 0, 115649b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 115749b5e25fSSatish Balay 0, 115849b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 115949b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 116049b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 116149b5e25fSSatish Balay 0, 116249b5e25fSSatish Balay 0, 11635f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 11645f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1165273d9f13SBarry Smith MatSetUpPreallocation_SeqSBAIJ, 1166c464158bSHong Zhang 0, 11674d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 116849b5e25fSSatish Balay 0, 116949b5e25fSSatish Balay 0, 117049b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 117149b5e25fSSatish Balay 0, 117249b5e25fSSatish Balay 0, 117349b5e25fSSatish Balay 0, 1174950f1e5bSHong Zhang 0, 117549b5e25fSSatish Balay 0, 117649b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 117749b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 117849b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 117949b5e25fSSatish Balay 0, 118049b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 118149b5e25fSSatish Balay MatScale_SeqSBAIJ, 118249b5e25fSSatish Balay 0, 118349b5e25fSSatish Balay 0, 118449b5e25fSSatish Balay 0, 118549b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 118649b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 118749b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 118849b5e25fSSatish Balay 0, 118949b5e25fSSatish Balay 0, 119049b5e25fSSatish Balay 0, 119149b5e25fSSatish Balay 0, 119249b5e25fSSatish Balay 0, 119349b5e25fSSatish Balay 0, 119449b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 119549b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 119649b5e25fSSatish Balay 0, 119749b5e25fSSatish Balay 0, 11988a124369SBarry Smith MatGetPetscMaps_Petsc, 1199d959ec07SHong Zhang 0, 1200d959ec07SHong Zhang 0, 1201d959ec07SHong Zhang 0, 1202d959ec07SHong Zhang 0, 1203d959ec07SHong Zhang 0, 1204d959ec07SHong Zhang 0, 120524d5174aSHong Zhang MatGetRowMax_SeqSBAIJ}; 120649b5e25fSSatish Balay 120749b5e25fSSatish Balay EXTERN_C_BEGIN 12084a2ae208SSatish Balay #undef __FUNCT__ 12094a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 121049b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 121149b5e25fSSatish Balay { 12124afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1213c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 121449b5e25fSSatish Balay int ierr; 121549b5e25fSSatish Balay 121649b5e25fSSatish Balay PetscFunctionBegin; 121749b5e25fSSatish Balay if (aij->nonew != 1) { 121829bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 121949b5e25fSSatish Balay } 122049b5e25fSSatish Balay 122149b5e25fSSatish Balay /* allocate space for values if not already there */ 122249b5e25fSSatish Balay if (!aij->saved_values) { 122387828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 122449b5e25fSSatish Balay } 122549b5e25fSSatish Balay 122649b5e25fSSatish Balay /* copy values over */ 122787828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 122849b5e25fSSatish Balay PetscFunctionReturn(0); 122949b5e25fSSatish Balay } 123049b5e25fSSatish Balay EXTERN_C_END 123149b5e25fSSatish Balay 123249b5e25fSSatish Balay EXTERN_C_BEGIN 12334a2ae208SSatish Balay #undef __FUNCT__ 12344a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 123549b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 123649b5e25fSSatish Balay { 12374afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1238c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 123949b5e25fSSatish Balay 124049b5e25fSSatish Balay PetscFunctionBegin; 124149b5e25fSSatish Balay if (aij->nonew != 1) { 124229bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 124349b5e25fSSatish Balay } 124449b5e25fSSatish Balay if (!aij->saved_values) { 124529bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 124649b5e25fSSatish Balay } 124749b5e25fSSatish Balay 124849b5e25fSSatish Balay /* copy values over */ 124987828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 125049b5e25fSSatish Balay PetscFunctionReturn(0); 125149b5e25fSSatish Balay } 125249b5e25fSSatish Balay EXTERN_C_END 125349b5e25fSSatish Balay 12548549e402SHong Zhang EXTERN_C_BEGIN 12554a2ae208SSatish Balay #undef __FUNCT__ 12564a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ" 1257c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1258c464158bSHong Zhang { 1259c464158bSHong Zhang Mat_SeqSBAIJ *b; 12600c955e93SHong Zhang int ierr,size; 1261c464158bSHong Zhang 1262c464158bSHong Zhang PetscFunctionBegin; 1263c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1264c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1265273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1266273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1267c464158bSHong Zhang 1268b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1269b0a32e0cSBarry Smith B->data = (void*)b; 1270c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1271c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1272c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1273c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1274c464158bSHong Zhang B->factor = 0; 1275c464158bSHong Zhang B->lupivotthreshold = 1.0; 1276c464158bSHong Zhang B->mapping = 0; 1277c464158bSHong Zhang b->row = 0; 1278c464158bSHong Zhang b->icol = 0; 1279c464158bSHong Zhang b->reallocs = 0; 1280c464158bSHong Zhang b->saved_values = 0; 1281c464158bSHong Zhang 12828a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 12838a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1284c464158bSHong Zhang 1285c464158bSHong Zhang b->sorted = PETSC_FALSE; 1286c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1287c464158bSHong Zhang b->nonew = 0; 1288c464158bSHong Zhang b->diag = 0; 1289c464158bSHong Zhang b->solve_work = 0; 1290c464158bSHong Zhang b->mult_work = 0; 1291*2a1b7f2aSHong Zhang B->spptr = 0; 1292c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1293c464158bSHong Zhang 1294c464158bSHong Zhang b->inew = 0; 1295c464158bSHong Zhang b->jnew = 0; 1296c464158bSHong Zhang b->anew = 0; 1297c464158bSHong Zhang b->a2anew = 0; 1298c464158bSHong Zhang b->permute = PETSC_FALSE; 1299c464158bSHong Zhang 1300c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1301c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1302c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1303c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1304c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1305c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1306c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1307c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1308c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1309c464158bSHong Zhang PetscFunctionReturn(0); 1310c464158bSHong Zhang } 13118549e402SHong Zhang EXTERN_C_END 1312c464158bSHong Zhang 13134a2ae208SSatish Balay #undef __FUNCT__ 13144a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 131549b5e25fSSatish Balay /*@C 1316273d9f13SBarry Smith MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 131749b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 131849b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 131949b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 132049b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 132149b5e25fSSatish Balay 1322c464158bSHong Zhang Collective on Mat 132349b5e25fSSatish Balay 132449b5e25fSSatish Balay Input Parameters: 1325c464158bSHong Zhang + A - the symmetric matrix 132649b5e25fSSatish Balay . bs - size of block 132749b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 132849b5e25fSSatish Balay - nnz - array containing the number of block nonzeros in the various block rows 132949b5e25fSSatish Balay (possibly different for each block row) or PETSC_NULL 133049b5e25fSSatish Balay 133149b5e25fSSatish Balay Options Database Keys: 133249b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 133349b5e25fSSatish Balay block calculations (much slower) 133449b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 133549b5e25fSSatish Balay 133649b5e25fSSatish Balay Level: intermediate 133749b5e25fSSatish Balay 133849b5e25fSSatish Balay Notes: 1339c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 134049b5e25fSSatish Balay storage. That is, the stored row and column indices can begin at 134149b5e25fSSatish Balay either one (as in Fortran) or zero. See the users' manual for details. 134249b5e25fSSatish Balay 134349b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 134449b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 134549b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 134649b5e25fSSatish Balay matrices. 134749b5e25fSSatish Balay 1348a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 134949b5e25fSSatish Balay @*/ 1350273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 135149b5e25fSSatish Balay { 1352c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13530c955e93SHong Zhang int i,len,ierr,mbs,bs2; 135449b5e25fSSatish Balay PetscTruth flg; 13554afc71dfSHong Zhang int s_nz; 135649b5e25fSSatish Balay 135749b5e25fSSatish Balay PetscFunctionBegin; 1358273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1359b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1360c464158bSHong Zhang mbs = B->m/bs; 136149b5e25fSSatish Balay bs2 = bs*bs; 136249b5e25fSSatish Balay 1363c464158bSHong Zhang if (mbs*bs != B->m) { 136429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 136549b5e25fSSatish Balay } 136649b5e25fSSatish Balay 1367435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1368435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 136949b5e25fSSatish Balay if (nnz) { 137049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 137129bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 137280fe4e49SBarry 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); 137349b5e25fSSatish Balay } 137449b5e25fSSatish Balay } 137549b5e25fSSatish Balay 1376b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 137749b5e25fSSatish Balay if (!flg) { 137849b5e25fSSatish Balay switch (bs) { 137949b5e25fSSatish Balay case 1: 13804ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 138149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 138249b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 138349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 138449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 138549b5e25fSSatish Balay break; 138649b5e25fSSatish Balay case 2: 13874ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 138849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 138949b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 139049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 139149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 139249b5e25fSSatish Balay break; 139349b5e25fSSatish Balay case 3: 13945f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 139549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 139649b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 139749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 139849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 139949b5e25fSSatish Balay break; 140049b5e25fSSatish Balay case 4: 14015f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 140249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 140349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 140449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 140549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 140649b5e25fSSatish Balay break; 140749b5e25fSSatish Balay case 5: 14085f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 140949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 141049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 141149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 141249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 141349b5e25fSSatish Balay break; 141449b5e25fSSatish Balay case 6: 14155f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 141649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 141749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 141849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 141949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 142049b5e25fSSatish Balay break; 142149b5e25fSSatish Balay case 7: 1422de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 142349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 142449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1425de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 142649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 142749b5e25fSSatish Balay break; 142849b5e25fSSatish Balay } 142949b5e25fSSatish Balay } 143049b5e25fSSatish Balay 143149b5e25fSSatish Balay b->mbs = mbs; 14324afc71dfSHong Zhang b->nbs = mbs; 1433b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 143449b5e25fSSatish Balay if (!nnz) { 1435435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 143649b5e25fSSatish Balay else if (nz <= 0) nz = 1; 143749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14388cef66ccSHong Zhang b->imax[i] = nz; 143949b5e25fSSatish Balay } 1440153ea458SHong Zhang nz = nz*mbs; /* total nz */ 144149b5e25fSSatish Balay } else { 144249b5e25fSSatish Balay nz = 0; 14438cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 144449b5e25fSSatish Balay } 14458cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 14468cef66ccSHong Zhang s_nz = nz; 144749b5e25fSSatish Balay 144849b5e25fSSatish Balay /* allocate the matrix space */ 1449c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 145082502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 145149b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 145249b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 145349b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 145449b5e25fSSatish Balay b->i = b->j + s_nz; 145549b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 145649b5e25fSSatish Balay 145749b5e25fSSatish Balay /* pointer to beginning of each row */ 145849b5e25fSSatish Balay b->i[0] = 0; 145949b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 146049b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 146149b5e25fSSatish Balay } 146249b5e25fSSatish Balay 146349b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 14645033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1465b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 146649b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 146749b5e25fSSatish Balay 146849b5e25fSSatish Balay b->bs = bs; 146949b5e25fSSatish Balay b->bs2 = bs2; 147049b5e25fSSatish Balay b->s_nz = 0; 147149b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1472153ea458SHong Zhang 147316cdd363SHong Zhang b->inew = 0; 147416cdd363SHong Zhang b->jnew = 0; 147516cdd363SHong Zhang b->anew = 0; 147616cdd363SHong Zhang b->a2anew = 0; 14771a3463dfSHong Zhang b->permute = PETSC_FALSE; 1478c464158bSHong Zhang PetscFunctionReturn(0); 1479c464158bSHong Zhang } 1480153ea458SHong Zhang 148149b5e25fSSatish Balay 14824a2ae208SSatish Balay #undef __FUNCT__ 14834a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1484c464158bSHong Zhang /*@C 1485c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1486c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1487c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1488c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1489c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 149049b5e25fSSatish Balay 1491c464158bSHong Zhang Collective on MPI_Comm 1492c464158bSHong Zhang 1493c464158bSHong Zhang Input Parameters: 1494c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1495c464158bSHong Zhang . bs - size of block 1496c464158bSHong Zhang . m - number of rows, or number of columns 1497c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1498c464158bSHong Zhang - nnz - array containing the number of block nonzeros in the various block rows 1499c464158bSHong Zhang (possibly different for each block row) or PETSC_NULL 1500c464158bSHong Zhang 1501c464158bSHong Zhang Output Parameter: 1502c464158bSHong Zhang . A - the symmetric matrix 1503c464158bSHong Zhang 1504c464158bSHong Zhang Options Database Keys: 1505c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1506c464158bSHong Zhang block calculations (much slower) 1507c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1508c464158bSHong Zhang 1509c464158bSHong Zhang Level: intermediate 1510c464158bSHong Zhang 1511c464158bSHong Zhang Notes: 1512c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 1513c464158bSHong Zhang storage. That is, the stored row and column indices can begin at 1514c464158bSHong Zhang either one (as in Fortran) or zero. See the users' manual for details. 1515c464158bSHong Zhang 1516c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1517c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1518c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1519c464158bSHong Zhang matrices. 1520c464158bSHong Zhang 1521c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1522c464158bSHong Zhang @*/ 1523c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1524c464158bSHong Zhang { 1525c464158bSHong Zhang int ierr; 1526c464158bSHong Zhang 1527c464158bSHong Zhang PetscFunctionBegin; 1528c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1529c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1530273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 153149b5e25fSSatish Balay PetscFunctionReturn(0); 153249b5e25fSSatish Balay } 153349b5e25fSSatish Balay 15344a2ae208SSatish Balay #undef __FUNCT__ 15354a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 153649b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 153749b5e25fSSatish Balay { 153849b5e25fSSatish Balay Mat C; 153949b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 154049b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 154149b5e25fSSatish Balay 154249b5e25fSSatish Balay PetscFunctionBegin; 154329bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 154449b5e25fSSatish Balay 154549b5e25fSSatish Balay *B = 0; 1546692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1547692f9cbeSHong Zhang ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 1548692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1549692f9cbeSHong Zhang 155049b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1551273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 155249b5e25fSSatish Balay C->factor = A->factor; 155349b5e25fSSatish Balay c->row = 0; 155449b5e25fSSatish Balay c->icol = 0; 155549b5e25fSSatish Balay c->saved_values = 0; 155649b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 155749b5e25fSSatish Balay C->assembled = PETSC_TRUE; 155849b5e25fSSatish Balay 155949b5e25fSSatish Balay c->bs = a->bs; 156049b5e25fSSatish Balay c->bs2 = a->bs2; 156149b5e25fSSatish Balay c->mbs = a->mbs; 156249b5e25fSSatish Balay c->nbs = a->nbs; 156349b5e25fSSatish Balay 1564b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1565b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 156649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 156749b5e25fSSatish Balay c->imax[i] = a->imax[i]; 156849b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 156949b5e25fSSatish Balay } 157049b5e25fSSatish Balay 157149b5e25fSSatish Balay /* allocate the matrix space */ 157249b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 157349b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 157482502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 157549b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 157649b5e25fSSatish Balay c->i = c->j + nz; 157749b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 157849b5e25fSSatish Balay if (mbs > 0) { 157949b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 158049b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 158149b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 158249b5e25fSSatish Balay } else { 158349b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 158449b5e25fSSatish Balay } 158549b5e25fSSatish Balay } 158649b5e25fSSatish Balay 1587b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 158849b5e25fSSatish Balay c->sorted = a->sorted; 158949b5e25fSSatish Balay c->roworiented = a->roworiented; 159049b5e25fSSatish Balay c->nonew = a->nonew; 159149b5e25fSSatish Balay 159249b5e25fSSatish Balay if (a->diag) { 1593b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1594b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 159549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 159649b5e25fSSatish Balay c->diag[i] = a->diag[i]; 159749b5e25fSSatish Balay } 159849b5e25fSSatish Balay } else c->diag = 0; 159949b5e25fSSatish Balay c->s_nz = a->s_nz; 160049b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 160149b5e25fSSatish Balay c->solve_work = 0; 1602*2a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 160349b5e25fSSatish Balay c->mult_work = 0; 160449b5e25fSSatish Balay *B = C; 1605b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 160649b5e25fSSatish Balay PetscFunctionReturn(0); 160749b5e25fSSatish Balay } 160849b5e25fSSatish Balay 16098549e402SHong Zhang EXTERN_C_BEGIN 16104a2ae208SSatish Balay #undef __FUNCT__ 16114a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1612b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A) 161349b5e25fSSatish Balay { 161449b5e25fSSatish Balay Mat_SeqSBAIJ *a; 161549b5e25fSSatish Balay Mat B; 161649b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 16173f7326a9SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 161849b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 161949b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 162087828ca2SBarry Smith PetscScalar *aa; 162149b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 162249b5e25fSSatish Balay 162349b5e25fSSatish Balay PetscFunctionBegin; 1624b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 162549b5e25fSSatish Balay bs2 = bs*bs; 162649b5e25fSSatish Balay 162749b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 162829bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1629b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 163049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1631552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 163249b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 163349b5e25fSSatish Balay 163449b5e25fSSatish Balay if (header[3] < 0) { 163529bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 163649b5e25fSSatish Balay } 163749b5e25fSSatish Balay 163829bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 163949b5e25fSSatish Balay 164049b5e25fSSatish Balay /* 164149b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 164249b5e25fSSatish Balay divisible by the blocksize 164349b5e25fSSatish Balay */ 164449b5e25fSSatish Balay mbs = M/bs; 164549b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 164649b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 164749b5e25fSSatish Balay else mbs++; 164849b5e25fSSatish Balay if (extra_rows) { 1649b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 165049b5e25fSSatish Balay } 165149b5e25fSSatish Balay 165249b5e25fSSatish Balay /* read in row lengths */ 1653b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 165449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 165549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 165649b5e25fSSatish Balay 165749b5e25fSSatish Balay /* read in column indices */ 1658b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 165949b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 166049b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 166149b5e25fSSatish Balay 166249b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 166382502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 166482502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 166549b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 166649b5e25fSSatish Balay masked = mask + mbs; 166749b5e25fSSatish Balay rowcount = 0; nzcount = 0; 166849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 166949b5e25fSSatish Balay nmask = 0; 167049b5e25fSSatish Balay for (j=0; j<bs; j++) { 167149b5e25fSSatish Balay kmax = rowlengths[rowcount]; 167249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 16732d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 167403630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 167549b5e25fSSatish Balay } 167649b5e25fSSatish Balay rowcount++; 167749b5e25fSSatish Balay } 1678574b2666SHong Zhang s_browlengths[i] += nmask; 1679574b2666SHong Zhang 168049b5e25fSSatish Balay /* zero out the mask elements we set */ 168149b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 168249b5e25fSSatish Balay } 168349b5e25fSSatish Balay 168449b5e25fSSatish Balay /* create our matrix */ 16857fe2be48SHong Zhang ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr); 168649b5e25fSSatish Balay B = *A; 168749b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 168849b5e25fSSatish Balay 168949b5e25fSSatish Balay /* set matrix "i" values */ 169049b5e25fSSatish Balay a->i[0] = 0; 169149b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1692574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1693574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 169449b5e25fSSatish Balay } 16957fe2be48SHong Zhang a->s_nz = a->i[mbs]; 169649b5e25fSSatish Balay 169749b5e25fSSatish Balay /* read in nonzero values */ 169887828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 169949b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 170049b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 170149b5e25fSSatish Balay 170249b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 170349b5e25fSSatish Balay nzcount = 0; jcount = 0; 170449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 170549b5e25fSSatish Balay nzcountb = nzcount; 170649b5e25fSSatish Balay nmask = 0; 170749b5e25fSSatish Balay for (j=0; j<bs; j++) { 170849b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 170949b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17102d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 171103630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 17122d703238SHong Zhang } 17132d703238SHong Zhang } 17142d703238SHong Zhang /* sort the masked values */ 17152d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 17162d703238SHong Zhang 17172d703238SHong Zhang /* set "j" values into matrix */ 17182d703238SHong Zhang maskcount = 1; 17192d703238SHong Zhang for (j=0; j<nmask; j++) { 172049b5e25fSSatish Balay a->j[jcount++] = masked[j]; 172149b5e25fSSatish Balay mask[masked[j]] = maskcount++; 172249b5e25fSSatish Balay } 1723574b2666SHong Zhang 172449b5e25fSSatish Balay /* set "a" values into matrix */ 172549b5e25fSSatish Balay ishift = bs2*a->i[i]; 172649b5e25fSSatish Balay for (j=0; j<bs; j++) { 172749b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 172849b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1729574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1730574b2666SHong Zhang if (tmp >= i){ 173149b5e25fSSatish Balay block = mask[tmp] - 1; 173249b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 173349b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1734574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1735574b2666SHong Zhang } 1736574b2666SHong Zhang nzcountb++; 173749b5e25fSSatish Balay } 173849b5e25fSSatish Balay } 173949b5e25fSSatish Balay /* zero out the mask elements we set */ 174049b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 174149b5e25fSSatish Balay } 174229bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 174349b5e25fSSatish Balay 174449b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1745574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 174649b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 174749b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 174849b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 174949b5e25fSSatish Balay 175049b5e25fSSatish Balay B->assembled = PETSC_TRUE; 175149b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 175249b5e25fSSatish Balay PetscFunctionReturn(0); 175349b5e25fSSatish Balay } 17548549e402SHong Zhang EXTERN_C_END 1755574b2666SHong Zhang 1756d06b337dSHong Zhang #undef __FUNCT__ 1757d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1758c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx) 1759d06b337dSHong Zhang { 1760d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1761d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1762d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 1763d06b337dSHong Zhang int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr; 1764d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1765d06b337dSHong Zhang 1766d06b337dSHong Zhang PetscFunctionBegin; 1767b965ef7fSBarry Smith its = its*lits; 176891723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 1769d06b337dSHong Zhang 1770d06b337dSHong Zhang if (bs > 1) 1771d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1772d06b337dSHong Zhang 1773d06b337dSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1774d06b337dSHong Zhang if (xx != bb) { 1775d06b337dSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1776d06b337dSHong Zhang } else { 1777d06b337dSHong Zhang b = x; 1778d06b337dSHong Zhang } 1779d06b337dSHong Zhang 1780d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1781d06b337dSHong Zhang 1782d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 17832798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1784d06b337dSHong Zhang for (i=0; i<m; i++) 1785d06b337dSHong Zhang t[i] = b[i]; 1786d06b337dSHong Zhang 1787d06b337dSHong Zhang for (i=0; i<m; i++){ 178844706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1789d06b337dSHong Zhang v = aa + ai[i] + 1; 1790d06b337dSHong Zhang vj = aj + ai[i] + 1; 1791d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1792d06b337dSHong Zhang x[i] = omega*t[i]/d; 1793d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 179444706875SHong Zhang PetscLogFlops(2*nz-1); 1795d06b337dSHong Zhang } 1796d06b337dSHong Zhang } 1797d06b337dSHong Zhang 17982798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1799d06b337dSHong Zhang for (i=0; i<m; i++) 1800d06b337dSHong Zhang t[i] = b[i]; 1801d06b337dSHong Zhang 1802d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1803d06b337dSHong Zhang v = aa + ai[i] + 1; 1804d06b337dSHong Zhang vj = aj + ai[i] + 1; 1805d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1806d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 180744706875SHong Zhang PetscLogFlops(2*nz-1); 1808d06b337dSHong Zhang } 1809d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1810d06b337dSHong Zhang d = *(aa + ai[i]); 1811d06b337dSHong Zhang v = aa + ai[i] + 1; 1812d06b337dSHong Zhang vj = aj + ai[i] + 1; 1813d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1814d06b337dSHong Zhang sum = t[i]; 1815d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 181644706875SHong Zhang PetscLogFlops(2*nz-1); 1817d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1818d06b337dSHong Zhang } 1819d06b337dSHong Zhang } 1820d06b337dSHong Zhang its--; 1821d06b337dSHong Zhang } 1822d06b337dSHong Zhang 1823d06b337dSHong Zhang while (its--) { 182444706875SHong Zhang /* 182544706875SHong Zhang forward sweep: 182644706875SHong Zhang for i=0,...,m-1: 182744706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 182844706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 182944706875SHong Zhang b = b - x[i]*U^T(i,:); 1830d06b337dSHong Zhang 183144706875SHong Zhang */ 18322798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1833d06b337dSHong Zhang for (i=0; i<m; i++) 1834d06b337dSHong Zhang t[i] = b[i]; 1835d06b337dSHong Zhang 1836d06b337dSHong Zhang for (i=0; i<m; i++){ 183744706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1838d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 1839d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 1840d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 1841d06b337dSHong Zhang sum = t[i]; 1842d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 1843d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1844d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 184544706875SHong Zhang PetscLogFlops(4*nz-2); 1846d06b337dSHong Zhang } 1847d06b337dSHong Zhang } 1848d06b337dSHong Zhang 18492798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 185044706875SHong Zhang /* 185144706875SHong Zhang backward sweep: 185244706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 185344706875SHong Zhang for i=m-1,...,0: 185444706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 185544706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 185644706875SHong Zhang */ 1857d06b337dSHong Zhang for (i=0; i<m; i++) 1858d06b337dSHong Zhang t[i] = b[i]; 1859d06b337dSHong Zhang 1860d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1861d06b337dSHong Zhang v = aa + ai[i] + 1; 1862d06b337dSHong Zhang vj = aj + ai[i] + 1; 1863d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1864d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 186544706875SHong Zhang PetscLogFlops(2*nz-1); 1866d06b337dSHong Zhang } 1867d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1868d06b337dSHong Zhang d = *(aa + ai[i]); 1869d06b337dSHong Zhang v = aa + ai[i] + 1; 1870d06b337dSHong Zhang vj = aj + ai[i] + 1; 1871d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1872d06b337dSHong Zhang sum = t[i]; 1873d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 187444706875SHong Zhang PetscLogFlops(2*nz-1); 1875d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1876d06b337dSHong Zhang } 1877d06b337dSHong Zhang } 1878d06b337dSHong Zhang } 1879d06b337dSHong Zhang 1880d06b337dSHong Zhang ierr = PetscFree(t); CHKERRQ(ierr); 1881d06b337dSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1882d06b337dSHong Zhang if (bb != xx) { 1883d06b337dSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1884d06b337dSHong Zhang } 18852798e883SHong Zhang 1886d06b337dSHong Zhang PetscFunctionReturn(0); 1887d06b337dSHong Zhang } 1888d06b337dSHong Zhang 1889d06b337dSHong Zhang 1890d06b337dSHong Zhang 1891d06b337dSHong Zhang 189249b5e25fSSatish Balay 189349b5e25fSSatish Balay 1894