173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 4a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/vec/vecimpl.h" 949b5e25fSSatish Balay #include "src/inline/spops.h" 10aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1149b5e25fSSatish Balay 1249b5e25fSSatish Balay #define CHUNKSIZE 10 1349b5e25fSSatish Balay 1449b5e25fSSatish Balay /* 1549b5e25fSSatish Balay Checks for missing diagonals 1649b5e25fSSatish Balay */ 174a2ae208SSatish Balay #undef __FUNCT__ 184a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 1949b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2049b5e25fSSatish Balay { 21045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2249b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2349b5e25fSSatish Balay 2449b5e25fSSatish Balay PetscFunctionBegin; 25045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2649b5e25fSSatish Balay diag = a->diag; 2749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 28e11b0e14SHong Zhang if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i); 2949b5e25fSSatish Balay } 3049b5e25fSSatish Balay PetscFunctionReturn(0); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay 334a2ae208SSatish Balay #undef __FUNCT__ 344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 3549b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3649b5e25fSSatish Balay { 37045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 38b424e231SHong Zhang int i,mbs = a->mbs,ierr; 3949b5e25fSSatish Balay 4049b5e25fSSatish Balay PetscFunctionBegin; 4149b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4249b5e25fSSatish Balay 43b424e231SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr); 44b424e231SHong Zhang PetscLogObjectMemory(A,(mbs+1)*sizeof(int)); 45b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4649b5e25fSSatish Balay PetscFunctionReturn(0); 4749b5e25fSSatish Balay } 4849b5e25fSSatish Balay 49a1373b80SHong Zhang /* 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; 58a1373b80SHong Zhang 59045c9aa0SHong Zhang if (ia) { 6029bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 6149b5e25fSSatish Balay } 62a1373b80SHong Zhang 63a1373b80SHong Zhang #ifdef NEW 64a1373b80SHong Zhang int ierr; 65a1373b80SHong Zhang 66a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 67a1373b80SHong Zhang 68a1373b80SHong Zhang /* 69a1373b80SHong Zhang if (symmetric) { 70a1373b80SHong Zhang ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 71a1373b80SHong Zhang } else if (oshift == 1) { 72a1373b80SHong Zhang int nz = a->i[n]; 73a1373b80SHong Zhang for (i=0; i<nz; i++) a->j[i]++; 74a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 75a1373b80SHong Zhang *ia = a->i; *ja = a->j; 76a1373b80SHong Zhang } else { 77a1373b80SHong Zhang */ 78a1373b80SHong Zhang *ia = a->i; *ja = a->j; 79a1373b80SHong Zhang /* } */ 80a1373b80SHong Zhang #endif 81045c9aa0SHong Zhang *nn = a->mbs; 8249b5e25fSSatish Balay PetscFunctionReturn(0); 8349b5e25fSSatish Balay } 8449b5e25fSSatish Balay 854a2ae208SSatish Balay #undef __FUNCT__ 864a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 87045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 8849b5e25fSSatish Balay { 8949b5e25fSSatish Balay PetscFunctionBegin; 9049b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 9129bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 9216cdd363SHong Zhang /* PetscFunctionReturn(0); */ 9349b5e25fSSatish Balay } 9449b5e25fSSatish Balay 954a2ae208SSatish Balay #undef __FUNCT__ 964a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 9749b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9849b5e25fSSatish Balay { 9949b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 10049b5e25fSSatish Balay 10149b5e25fSSatish Balay PetscFunctionBegin; 10249b5e25fSSatish Balay *bs = sbaij->bs; 10349b5e25fSSatish Balay PetscFunctionReturn(0); 10449b5e25fSSatish Balay } 10549b5e25fSSatish Balay 1064a2ae208SSatish Balay #undef __FUNCT__ 1074a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 10849b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 10949b5e25fSSatish Balay { 11049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 11149b5e25fSSatish Balay int ierr; 11249b5e25fSSatish Balay 11349b5e25fSSatish Balay PetscFunctionBegin; 11449b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 115b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 11649b5e25fSSatish Balay #endif 11749b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 11849b5e25fSSatish Balay if (!a->singlemalloc) { 11949b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 12063c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 12149b5e25fSSatish Balay } 12249b5e25fSSatish Balay if (a->row) { 12349b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 12449b5e25fSSatish Balay } 12549b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 12649b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 12749b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 12849b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 12949b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 13049b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 13149b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1321a3463dfSHong Zhang 1331a3463dfSHong Zhang if (a->inew){ 1341a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1351a3463dfSHong Zhang a->inew = 0; 1361a3463dfSHong Zhang } 13749b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 13849b5e25fSSatish Balay PetscFunctionReturn(0); 13949b5e25fSSatish Balay } 14049b5e25fSSatish Balay 1414a2ae208SSatish Balay #undef __FUNCT__ 1424a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 14349b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 14449b5e25fSSatish Balay { 145045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 14649b5e25fSSatish Balay 14749b5e25fSSatish Balay PetscFunctionBegin; 1484d9d31abSKris Buschelman switch (op) { 1494d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1504d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1514d9d31abSKris Buschelman break; 1524d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1534d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1544d9d31abSKris Buschelman break; 1554d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1564d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1574d9d31abSKris Buschelman break; 1584d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1594d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1604d9d31abSKris Buschelman break; 1614d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1624d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1634d9d31abSKris Buschelman break; 1644d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1654d9d31abSKris Buschelman a->nonew = 1; 1664d9d31abSKris Buschelman break; 1674d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1684d9d31abSKris Buschelman a->nonew = -1; 1694d9d31abSKris Buschelman break; 1704d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1714d9d31abSKris Buschelman a->nonew = -2; 1724d9d31abSKris Buschelman break; 1734d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1744d9d31abSKris Buschelman a->nonew = 0; 1754d9d31abSKris Buschelman break; 1764d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1774d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1784d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1794d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1804d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 181d03495bdSKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 182b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1834d9d31abSKris Buschelman break; 1844d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 18529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1864d9d31abSKris Buschelman default: 18729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 18849b5e25fSSatish Balay } 18949b5e25fSSatish Balay PetscFunctionReturn(0); 19049b5e25fSSatish Balay } 19149b5e25fSSatish Balay 1924a2ae208SSatish Balay #undef __FUNCT__ 1934a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 19487828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v) 19549b5e25fSSatish Balay { 19649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 19782502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr; 19849b5e25fSSatish Balay MatScalar *aa,*aa_i; 19987828ca2SBarry Smith PetscScalar *v_i; 20049b5e25fSSatish Balay 20149b5e25fSSatish Balay PetscFunctionBegin; 20249b5e25fSSatish Balay bs = a->bs; 20349b5e25fSSatish Balay ai = a->i; 20449b5e25fSSatish Balay aj = a->j; 20549b5e25fSSatish Balay aa = a->a; 20649b5e25fSSatish Balay bs2 = a->bs2; 20749b5e25fSSatish Balay 208c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 20949b5e25fSSatish Balay 21049b5e25fSSatish Balay bn = row/bs; /* Block number */ 21149b5e25fSSatish Balay bp = row % bs; /* Block position */ 21249b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 21349b5e25fSSatish Balay *ncols = bs*M; 21449b5e25fSSatish Balay 21549b5e25fSSatish Balay if (v) { 21649b5e25fSSatish Balay *v = 0; 21749b5e25fSSatish Balay if (*ncols) { 21887828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr); 21949b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22049b5e25fSSatish Balay v_i = *v + i*bs; 22149b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 22249b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 22349b5e25fSSatish Balay } 22449b5e25fSSatish Balay } 22549b5e25fSSatish Balay } 22649b5e25fSSatish Balay 22749b5e25fSSatish Balay if (cols) { 22849b5e25fSSatish Balay *cols = 0; 22949b5e25fSSatish Balay if (*ncols) { 23082502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 23149b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23249b5e25fSSatish Balay cols_i = *cols + i*bs; 23349b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 23449b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 23549b5e25fSSatish Balay } 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay } 23849b5e25fSSatish Balay 23949b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2405ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2415ddb2528SHong Zhang #ifdef column_search 24249b5e25fSSatish Balay v_i = *v + M*bs; 24349b5e25fSSatish Balay cols_i = *cols + M*bs; 24449b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 24549b5e25fSSatish Balay M = ai[i+1] - ai[i]; 24649b5e25fSSatish Balay for (j=0; j<M; j++){ 24749b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 24849b5e25fSSatish Balay if (itmp == bn){ 24949b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 25049b5e25fSSatish Balay for (k=0; k<bs; k++) { 25149b5e25fSSatish Balay *cols_i++ = i*bs+k; 25249b5e25fSSatish Balay *v_i++ = aa_i[k]; 25349b5e25fSSatish Balay } 25449b5e25fSSatish Balay *ncols += bs; 25549b5e25fSSatish Balay break; 25649b5e25fSSatish Balay } 25749b5e25fSSatish Balay } 25849b5e25fSSatish Balay } 2595ddb2528SHong Zhang #endif 26049b5e25fSSatish Balay 26149b5e25fSSatish Balay PetscFunctionReturn(0); 26249b5e25fSSatish Balay } 26349b5e25fSSatish Balay 2644a2ae208SSatish Balay #undef __FUNCT__ 2654a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 26687828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 26749b5e25fSSatish Balay { 26849b5e25fSSatish Balay int ierr; 26949b5e25fSSatish Balay 27049b5e25fSSatish Balay PetscFunctionBegin; 27149b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 27249b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 27349b5e25fSSatish Balay PetscFunctionReturn(0); 27449b5e25fSSatish Balay } 27549b5e25fSSatish Balay 2764a2ae208SSatish Balay #undef __FUNCT__ 2774a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 27849b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 27949b5e25fSSatish Balay { 2808115998fSBarry Smith int ierr; 28149b5e25fSSatish Balay PetscFunctionBegin; 282999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2838115998fSBarry Smith PetscFunctionReturn(0); 28449b5e25fSSatish Balay } 28549b5e25fSSatish Balay 2864a2ae208SSatish Balay #undef __FUNCT__ 2874a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary" 288b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer) 28949b5e25fSSatish Balay { 29049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 29149b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 29287828ca2SBarry Smith PetscScalar *aa; 29349b5e25fSSatish Balay FILE *file; 29449b5e25fSSatish Balay 29549b5e25fSSatish Balay PetscFunctionBegin; 296b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 29782502324SSatish Balay ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 298552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 29949b5e25fSSatish Balay 300c464158bSHong Zhang col_lens[1] = A->m; 301c464158bSHong Zhang col_lens[2] = A->m; 30249b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 30349b5e25fSSatish Balay 30449b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 30549b5e25fSSatish Balay count = 0; 30649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30749b5e25fSSatish Balay for (j=0; j<bs; j++) { 30849b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 30949b5e25fSSatish Balay } 31049b5e25fSSatish Balay } 311c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 31249b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 31349b5e25fSSatish Balay 31449b5e25fSSatish Balay /* store column indices (zero start index) */ 31582502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 31649b5e25fSSatish Balay count = 0; 31749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31849b5e25fSSatish Balay for (j=0; j<bs; j++) { 31949b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 32049b5e25fSSatish Balay for (l=0; l<bs; l++) { 32149b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 32249b5e25fSSatish Balay } 32349b5e25fSSatish Balay } 32449b5e25fSSatish Balay } 32549b5e25fSSatish Balay } 32649b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 32749b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 32849b5e25fSSatish Balay 32949b5e25fSSatish Balay /* store nonzero values */ 33087828ca2SBarry Smith ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 33149b5e25fSSatish Balay count = 0; 33249b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33349b5e25fSSatish Balay for (j=0; j<bs; j++) { 33449b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33549b5e25fSSatish Balay for (l=0; l<bs; l++) { 33649b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 33749b5e25fSSatish Balay } 33849b5e25fSSatish Balay } 33949b5e25fSSatish Balay } 34049b5e25fSSatish Balay } 34149b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 34249b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 34349b5e25fSSatish Balay 344b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 34549b5e25fSSatish Balay if (file) { 34649b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 34749b5e25fSSatish Balay } 34849b5e25fSSatish Balay PetscFunctionReturn(0); 34949b5e25fSSatish Balay } 35049b5e25fSSatish Balay 3514a2ae208SSatish Balay #undef __FUNCT__ 3524a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 353b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 35449b5e25fSSatish Balay { 35549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 356fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 357fb9695e5SSatish Balay char *name; 358f3ef73ceSBarry Smith PetscViewerFormat format; 35949b5e25fSSatish Balay 36049b5e25fSSatish Balay PetscFunctionBegin; 36180fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 362b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 363456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 364b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 365fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 36629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 367fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 368b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36949b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 37049b5e25fSSatish Balay for (j=0; j<bs; j++) { 371b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 37249b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 37349b5e25fSSatish Balay for (l=0; l<bs; l++) { 37449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37549b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 37749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37849b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 38049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38149b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 38262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38349b5e25fSSatish Balay } 38449b5e25fSSatish Balay #else 38549b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 38662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38749b5e25fSSatish Balay } 38849b5e25fSSatish Balay #endif 38949b5e25fSSatish Balay } 39049b5e25fSSatish Balay } 391b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 39249b5e25fSSatish Balay } 39349b5e25fSSatish Balay } 394b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 39549b5e25fSSatish Balay } else { 396b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 39749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 39849b5e25fSSatish Balay for (j=0; j<bs; j++) { 399b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 40049b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 40149b5e25fSSatish Balay for (l=0; l<bs; l++) { 40249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 40349b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 40462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 40549b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40649b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 40762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 40849b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40949b5e25fSSatish Balay } else { 41062b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 41149b5e25fSSatish Balay } 41249b5e25fSSatish Balay #else 413b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 41449b5e25fSSatish Balay #endif 41549b5e25fSSatish Balay } 41649b5e25fSSatish Balay } 417b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 41849b5e25fSSatish Balay } 41949b5e25fSSatish Balay } 420b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 42149b5e25fSSatish Balay } 422b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 42349b5e25fSSatish Balay PetscFunctionReturn(0); 42449b5e25fSSatish Balay } 42549b5e25fSSatish Balay 4264a2ae208SSatish Balay #undef __FUNCT__ 4274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 428b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 42949b5e25fSSatish Balay { 43049b5e25fSSatish Balay Mat A = (Mat) Aa; 43149b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 43249b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 43349b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 43449b5e25fSSatish Balay MatScalar *aa; 43549b5e25fSSatish Balay MPI_Comm comm; 436b0a32e0cSBarry Smith PetscViewer viewer; 43749b5e25fSSatish Balay 43849b5e25fSSatish Balay PetscFunctionBegin; 43949b5e25fSSatish Balay /* 44049b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 44149b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 44249b5e25fSSatish Balay rest should return immediately. 44349b5e25fSSatish Balay */ 44449b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 44549b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 44649b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 44749b5e25fSSatish Balay 44849b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 44949b5e25fSSatish Balay 450b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 451b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 45249b5e25fSSatish Balay 45349b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 454b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 45549b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45649b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 457c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45849b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45949b5e25fSSatish Balay aa = a->a + j*bs2; 46049b5e25fSSatish Balay for (k=0; k<bs; k++) { 46149b5e25fSSatish Balay for (l=0; l<bs; l++) { 46249b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 463b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46449b5e25fSSatish Balay } 46549b5e25fSSatish Balay } 46649b5e25fSSatish Balay } 46749b5e25fSSatish Balay } 468b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 46949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 47049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 471c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 47249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47349b5e25fSSatish Balay aa = a->a + j*bs2; 47449b5e25fSSatish Balay for (k=0; k<bs; k++) { 47549b5e25fSSatish Balay for (l=0; l<bs; l++) { 47649b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 477b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 47849b5e25fSSatish Balay } 47949b5e25fSSatish Balay } 48049b5e25fSSatish Balay } 48149b5e25fSSatish Balay } 48249b5e25fSSatish Balay 483b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 48449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 48549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 486c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 48749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 48849b5e25fSSatish Balay aa = a->a + j*bs2; 48949b5e25fSSatish Balay for (k=0; k<bs; k++) { 49049b5e25fSSatish Balay for (l=0; l<bs; l++) { 49149b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 492b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 49349b5e25fSSatish Balay } 49449b5e25fSSatish Balay } 49549b5e25fSSatish Balay } 49649b5e25fSSatish Balay } 49749b5e25fSSatish Balay PetscFunctionReturn(0); 49849b5e25fSSatish Balay } 49949b5e25fSSatish Balay 5004a2ae208SSatish Balay #undef __FUNCT__ 5014a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 502b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 50349b5e25fSSatish Balay { 50449b5e25fSSatish Balay int ierr; 50549b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 506b0a32e0cSBarry Smith PetscDraw draw; 50749b5e25fSSatish Balay PetscTruth isnull; 50849b5e25fSSatish Balay 50949b5e25fSSatish Balay PetscFunctionBegin; 51049b5e25fSSatish Balay 511b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 512b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 51349b5e25fSSatish Balay 51449b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 515c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 51649b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 517b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 518b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 51949b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 52049b5e25fSSatish Balay PetscFunctionReturn(0); 52149b5e25fSSatish Balay } 52249b5e25fSSatish Balay 5234a2ae208SSatish Balay #undef __FUNCT__ 5244a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 525b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 52649b5e25fSSatish Balay { 52749b5e25fSSatish Balay int ierr; 52849b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 52949b5e25fSSatish Balay 53049b5e25fSSatish Balay PetscFunctionBegin; 531b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 532b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 533fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 534fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 53549b5e25fSSatish Balay if (issocket) { 53629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 53749b5e25fSSatish Balay } else if (isascii){ 53849b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 53949b5e25fSSatish Balay } else if (isbinary) { 54049b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 54149b5e25fSSatish Balay } else if (isdraw) { 54249b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 54349b5e25fSSatish Balay } else { 54429bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 54549b5e25fSSatish Balay } 54649b5e25fSSatish Balay PetscFunctionReturn(0); 54749b5e25fSSatish Balay } 54849b5e25fSSatish Balay 54949b5e25fSSatish Balay 5504a2ae208SSatish Balay #undef __FUNCT__ 5514a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 55287828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v) 55349b5e25fSSatish Balay { 554045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 55549b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 55649b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 55749b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 55849b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 55949b5e25fSSatish Balay 56049b5e25fSSatish Balay PetscFunctionBegin; 56149b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 56249b5e25fSSatish Balay row = im[k]; brow = row/bs; 56329bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 564c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56549b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 56649b5e25fSSatish Balay nrow = ailen[brow]; 56749b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 56829bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 569c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 57049b5e25fSSatish Balay col = in[l] ; 57149b5e25fSSatish Balay bcol = col/bs; 57249b5e25fSSatish Balay cidx = col%bs; 57349b5e25fSSatish Balay ridx = row%bs; 57449b5e25fSSatish Balay high = nrow; 57549b5e25fSSatish Balay low = 0; /* assume unsorted */ 57649b5e25fSSatish Balay while (high-low > 5) { 57749b5e25fSSatish Balay t = (low+high)/2; 57849b5e25fSSatish Balay if (rp[t] > bcol) high = t; 57949b5e25fSSatish Balay else low = t; 58049b5e25fSSatish Balay } 58149b5e25fSSatish Balay for (i=low; i<high; i++) { 58249b5e25fSSatish Balay if (rp[i] > bcol) break; 58349b5e25fSSatish Balay if (rp[i] == bcol) { 58449b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 58549b5e25fSSatish Balay goto finished; 58649b5e25fSSatish Balay } 58749b5e25fSSatish Balay } 58849b5e25fSSatish Balay *v++ = zero; 58949b5e25fSSatish Balay finished:; 59049b5e25fSSatish Balay } 59149b5e25fSSatish Balay } 59249b5e25fSSatish Balay PetscFunctionReturn(0); 59349b5e25fSSatish Balay } 59449b5e25fSSatish Balay 59549b5e25fSSatish Balay 5964a2ae208SSatish Balay #undef __FUNCT__ 5974a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 59887828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 59949b5e25fSSatish Balay { 6000880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 6010880e062SHong Zhang int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 6020880e062SHong Zhang int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 6030880e062SHong Zhang int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6040880e062SHong Zhang PetscTruth roworiented=a->roworiented; 6050880e062SHong Zhang MatScalar *value = v,*ap,*aa = a->a,*bap; 6060880e062SHong Zhang 60749b5e25fSSatish Balay PetscFunctionBegin; 6080880e062SHong Zhang if (roworiented) { 6090880e062SHong Zhang stepval = (n-1)*bs; 6100880e062SHong Zhang } else { 6110880e062SHong Zhang stepval = (m-1)*bs; 6120880e062SHong Zhang } 6130880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 6140880e062SHong Zhang row = im[k]; 6150880e062SHong Zhang if (row < 0) continue; 6160880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 6170880e062SHong Zhang if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 6180880e062SHong Zhang #endif 6190880e062SHong Zhang rp = aj + ai[row]; 6200880e062SHong Zhang ap = aa + bs2*ai[row]; 6210880e062SHong Zhang rmax = imax[row]; 6220880e062SHong Zhang nrow = ailen[row]; 6230880e062SHong Zhang low = 0; 6240880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 6250880e062SHong Zhang if (in[l] < 0) continue; 6260880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 6270880e062SHong Zhang if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 6280880e062SHong Zhang #endif 6290880e062SHong Zhang col = in[l]; 6300880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 6310880e062SHong Zhang if (roworiented) { 6320880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 6330880e062SHong Zhang } else { 6340880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 6350880e062SHong Zhang } 6360880e062SHong Zhang if (!sorted) low = 0; high = nrow; 6370880e062SHong Zhang while (high-low > 7) { 6380880e062SHong Zhang t = (low+high)/2; 6390880e062SHong Zhang if (rp[t] > col) high = t; 6400880e062SHong Zhang else low = t; 6410880e062SHong Zhang } 6420880e062SHong Zhang for (i=low; i<high; i++) { 6430880e062SHong Zhang if (rp[i] > col) break; 6440880e062SHong Zhang if (rp[i] == col) { 6450880e062SHong Zhang bap = ap + bs2*i; 6460880e062SHong Zhang if (roworiented) { 6470880e062SHong Zhang if (is == ADD_VALUES) { 6480880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6490880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6500880e062SHong Zhang bap[jj] += *value++; 6510880e062SHong Zhang } 6520880e062SHong Zhang } 6530880e062SHong Zhang } else { 6540880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6550880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6560880e062SHong Zhang bap[jj] = *value++; 6570880e062SHong Zhang } 6580880e062SHong Zhang } 6590880e062SHong Zhang } 6600880e062SHong Zhang } else { 6610880e062SHong Zhang if (is == ADD_VALUES) { 6620880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6630880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6640880e062SHong Zhang *bap++ += *value++; 6650880e062SHong Zhang } 6660880e062SHong Zhang } 6670880e062SHong Zhang } else { 6680880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6690880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6700880e062SHong Zhang *bap++ = *value++; 6710880e062SHong Zhang } 6720880e062SHong Zhang } 6730880e062SHong Zhang } 6740880e062SHong Zhang } 6750880e062SHong Zhang goto noinsert2; 6760880e062SHong Zhang } 6770880e062SHong Zhang } 6780880e062SHong Zhang if (nonew == 1) goto noinsert2; 6790880e062SHong Zhang else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 6800880e062SHong Zhang if (nrow >= rmax) { 6810880e062SHong Zhang /* there is no extra room in row, therefore enlarge */ 6820880e062SHong Zhang int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6830880e062SHong Zhang MatScalar *new_a; 6840880e062SHong Zhang 6850880e062SHong Zhang if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 6860880e062SHong Zhang 6870880e062SHong Zhang /* malloc new storage space */ 6880880e062SHong Zhang len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 6890880e062SHong Zhang ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 6900880e062SHong Zhang new_j = (int*)(new_a + bs2*new_nz); 6910880e062SHong Zhang new_i = new_j + new_nz; 6920880e062SHong Zhang 6930880e062SHong Zhang /* copy over old data into new slots */ 6940880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6950880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 6960880e062SHong Zhang ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 6970880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 6980880e062SHong Zhang ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 6990880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 7000880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 7010880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 7020880e062SHong Zhang /* free up old matrix storage */ 7030880e062SHong Zhang ierr = PetscFree(a->a);CHKERRQ(ierr); 7040880e062SHong Zhang if (!a->singlemalloc) { 7050880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 7060880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 7070880e062SHong Zhang } 7080880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 7090880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 7100880e062SHong Zhang 7110880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 7120880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 7130880e062SHong Zhang PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 7140880e062SHong Zhang a->s_maxnz += bs2*CHUNKSIZE; 7150880e062SHong Zhang a->reallocs++; 7160880e062SHong Zhang a->s_nz++; 7170880e062SHong Zhang } 7180880e062SHong Zhang N = nrow++ - 1; 7190880e062SHong Zhang /* shift up all the later entries in this row */ 7200880e062SHong Zhang for (ii=N; ii>=i; ii--) { 7210880e062SHong Zhang rp[ii+1] = rp[ii]; 7220880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 7230880e062SHong Zhang } 7240880e062SHong Zhang if (N >= i) { 7250880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 7260880e062SHong Zhang } 7270880e062SHong Zhang rp[i] = col; 7280880e062SHong Zhang bap = ap + bs2*i; 7290880e062SHong Zhang if (roworiented) { 7300880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7310880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 7320880e062SHong Zhang bap[jj] = *value++; 7330880e062SHong Zhang } 7340880e062SHong Zhang } 7350880e062SHong Zhang } else { 7360880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 7370880e062SHong Zhang for (jj=0; jj<bs; jj++) { 7380880e062SHong Zhang *bap++ = *value++; 7390880e062SHong Zhang } 7400880e062SHong Zhang } 7410880e062SHong Zhang } 7420880e062SHong Zhang noinsert2:; 7430880e062SHong Zhang low = i; 7440880e062SHong Zhang } 7450880e062SHong Zhang ailen[row] = nrow; 7460880e062SHong Zhang } 7470880e062SHong Zhang PetscFunctionReturn(0); 74849b5e25fSSatish Balay } 74949b5e25fSSatish Balay 7504a2ae208SSatish Balay #undef __FUNCT__ 7514a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 75249b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 75349b5e25fSSatish Balay { 75449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 75549b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 756c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 75749b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 75849b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 759cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 760cf4441caSHong Zhang PetscTruth flag; 761cf4441caSHong Zhang #endif 76249b5e25fSSatish Balay 76349b5e25fSSatish Balay PetscFunctionBegin; 76449b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 76549b5e25fSSatish Balay 76649b5e25fSSatish Balay if (m) rmax = ailen[0]; 76749b5e25fSSatish Balay for (i=1; i<mbs; i++) { 76849b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 76949b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 77049b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 77149b5e25fSSatish Balay if (fshift) { 77249b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 77349b5e25fSSatish Balay N = ailen[i]; 77449b5e25fSSatish Balay for (j=0; j<N; j++) { 77549b5e25fSSatish Balay ip[j-fshift] = ip[j]; 77649b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 77749b5e25fSSatish Balay } 77849b5e25fSSatish Balay } 77949b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 78049b5e25fSSatish Balay } 78149b5e25fSSatish Balay if (mbs) { 78249b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 78349b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 78449b5e25fSSatish Balay } 78549b5e25fSSatish Balay /* reset ilen and imax for each row */ 78649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 78749b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 78849b5e25fSSatish Balay } 78949b5e25fSSatish Balay a->s_nz = ai[mbs]; 79049b5e25fSSatish Balay 791b424e231SHong Zhang /* diagonals may have moved, reset it */ 792b424e231SHong Zhang if (a->diag) { 793b424e231SHong Zhang ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr); 79449b5e25fSSatish Balay } 795b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 796c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 797b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 79849b5e25fSSatish Balay a->reallocs); 799b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 80049b5e25fSSatish Balay a->reallocs = 0; 80149b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 80249b5e25fSSatish Balay 803cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES) 8042c535e4dSHong Zhang ierr = PetscOptionsHasName(A->prefix,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr); 805cf4441caSHong Zhang if (flag) { ierr = MatUseSpooles_SeqSBAIJ(A);CHKERRQ(ierr); } 806cf4441caSHong Zhang #endif 807cf4441caSHong Zhang 80849b5e25fSSatish Balay PetscFunctionReturn(0); 80949b5e25fSSatish Balay } 81049b5e25fSSatish Balay 81149b5e25fSSatish Balay /* 81249b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 81349b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 81449b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 81549b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 81649b5e25fSSatish Balay */ 8174a2ae208SSatish Balay #undef __FUNCT__ 8184a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 819*db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 82049b5e25fSSatish Balay { 82149b5e25fSSatish Balay int i,j,k,row; 82249b5e25fSSatish Balay PetscTruth flg; 82349b5e25fSSatish Balay 82449b5e25fSSatish Balay PetscFunctionBegin; 82549b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 82649b5e25fSSatish Balay row = idx[i]; 82749b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 82849b5e25fSSatish Balay sizes[j] = 1; 82949b5e25fSSatish Balay i++; 83049b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 83149b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 83249b5e25fSSatish Balay i++; 83349b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 83449b5e25fSSatish Balay flg = PETSC_TRUE; 83549b5e25fSSatish Balay for (k=1; k<bs; k++) { 83649b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 83749b5e25fSSatish Balay flg = PETSC_FALSE; 83849b5e25fSSatish Balay break; 83949b5e25fSSatish Balay } 84049b5e25fSSatish Balay } 84149b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 84249b5e25fSSatish Balay sizes[j] = bs; 84349b5e25fSSatish Balay i+= bs; 84449b5e25fSSatish Balay } else { 84549b5e25fSSatish Balay sizes[j] = 1; 84649b5e25fSSatish Balay i++; 84749b5e25fSSatish Balay } 84849b5e25fSSatish Balay } 84949b5e25fSSatish Balay } 85049b5e25fSSatish Balay *bs_max = j; 85149b5e25fSSatish Balay PetscFunctionReturn(0); 85249b5e25fSSatish Balay } 85349b5e25fSSatish Balay 8544a2ae208SSatish Balay #undef __FUNCT__ 8554a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 85687828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag) 85749b5e25fSSatish Balay { 85849b5e25fSSatish Balay PetscFunctionBegin; 859c0f24835SHong Zhang SETERRQ(PETSC_ERR_SUP,"No support for this function yet"); 86049b5e25fSSatish Balay } 86149b5e25fSSatish Balay 86249b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 86349b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 86449b5e25fSSatish Balay */ 86549b5e25fSSatish Balay 8664a2ae208SSatish Balay #undef __FUNCT__ 8674a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 86887828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is) 86949b5e25fSSatish Balay { 87049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 87149b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 87249b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 87349b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 87449b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 87549b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 87649b5e25fSSatish Balay 87749b5e25fSSatish Balay PetscFunctionBegin; 87849b5e25fSSatish Balay 87949b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 88049b5e25fSSatish Balay row = im[k]; /* row number */ 88149b5e25fSSatish Balay brow = row/bs; /* block row number */ 88249b5e25fSSatish Balay if (row < 0) continue; 88349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 884c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 88549b5e25fSSatish Balay #endif 88649b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 88749b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 88849b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 88949b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 89049b5e25fSSatish Balay low = 0; 89149b5e25fSSatish Balay 89249b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 89349b5e25fSSatish Balay if (in[l] < 0) continue; 89449b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 895c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 89649b5e25fSSatish Balay #endif 89749b5e25fSSatish Balay col = in[l]; 89849b5e25fSSatish Balay bcol = col/bs; /* block col number */ 89949b5e25fSSatish Balay 90049b5e25fSSatish Balay if (brow <= bcol){ 90149b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 9028549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 90349b5e25fSSatish Balay /* element value a(k,l) */ 90449b5e25fSSatish Balay if (roworiented) { 90549b5e25fSSatish Balay value = v[l + k*n]; 90649b5e25fSSatish Balay } else { 90749b5e25fSSatish Balay value = v[k + l*m]; 90849b5e25fSSatish Balay } 90949b5e25fSSatish Balay 91049b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 91149b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 91249b5e25fSSatish Balay while (high-low > 7) { 91349b5e25fSSatish Balay t = (low+high)/2; 91449b5e25fSSatish Balay if (rp[t] > bcol) high = t; 91549b5e25fSSatish Balay else low = t; 91649b5e25fSSatish Balay } 91749b5e25fSSatish Balay for (i=low; i<high; i++) { 91849b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 91949b5e25fSSatish Balay if (rp[i] > bcol) break; 92049b5e25fSSatish Balay if (rp[i] == bcol) { 92149b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 92249b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 92349b5e25fSSatish Balay else *bap = value; 9248549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 9258549e402SHong Zhang if (brow == bcol && ridx < cidx){ 9268549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 9278549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 9288549e402SHong Zhang else *bap = value; 9298549e402SHong Zhang } 93049b5e25fSSatish Balay goto noinsert1; 93149b5e25fSSatish Balay } 93249b5e25fSSatish Balay } 93349b5e25fSSatish Balay 93449b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 93529bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 93649b5e25fSSatish Balay if (nrow >= rmax) { 93749b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 93849b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 93949b5e25fSSatish Balay MatScalar *new_a; 94049b5e25fSSatish Balay 94129bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 94249b5e25fSSatish Balay 94349b5e25fSSatish Balay /* Malloc new storage space */ 94449b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 94582502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 94649b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 94749b5e25fSSatish Balay new_i = new_j + new_nz; 94849b5e25fSSatish Balay 94949b5e25fSSatish Balay /* copy over old data into new slots */ 95049b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 95149b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 95249b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 95349b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 95449b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 95549b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 95649b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 95749b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 95849b5e25fSSatish Balay /* free up old matrix storage */ 95949b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 96049b5e25fSSatish Balay if (!a->singlemalloc) { 96149b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 96249b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 96349b5e25fSSatish Balay } 96449b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 96549b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 96649b5e25fSSatish Balay 96749b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 96849b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 969b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 97049b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 97149b5e25fSSatish Balay a->reallocs++; 97249b5e25fSSatish Balay a->s_nz++; 97349b5e25fSSatish Balay } 97449b5e25fSSatish Balay 97549b5e25fSSatish Balay N = nrow++ - 1; 97649b5e25fSSatish Balay /* shift up all the later entries in this row */ 97749b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 97849b5e25fSSatish Balay rp[ii+1] = rp[ii]; 97949b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 98049b5e25fSSatish Balay } 98149b5e25fSSatish Balay if (N>=i) { 98249b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 98349b5e25fSSatish Balay } 98449b5e25fSSatish Balay rp[i] = bcol; 98549b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 98649b5e25fSSatish Balay noinsert1:; 98749b5e25fSSatish Balay low = i; 98849b5e25fSSatish Balay /* } */ 9898549e402SHong Zhang } 99049b5e25fSSatish Balay } /* end of if .. if.. */ 99149b5e25fSSatish Balay } /* end of loop over added columns */ 99249b5e25fSSatish Balay ailen[brow] = nrow; 99349b5e25fSSatish Balay } /* end of loop over added rows */ 99449b5e25fSSatish Balay 99549b5e25fSSatish Balay PetscFunctionReturn(0); 99649b5e25fSSatish Balay } 99749b5e25fSSatish Balay 998915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*); 999915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal); 100049b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 100149b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 100249b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 100349b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 100449b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 100587828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat); 100649b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 100749b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 100849b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 100949b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 101049b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 101149b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 101213a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 101349b5e25fSSatish Balay 101449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 101549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 101649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 101749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 101849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 101949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 102049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 102149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 102249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 102349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 102449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 102549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 102649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 102749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 102849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 102949b5e25fSSatish Balay 103007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 103107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 103207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 103307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 103407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 103507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 103607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 103707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 103807f98182SSatish Balay 10395f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 10405f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 10415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 10425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 10435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 10445f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 10455f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 104657d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 104749b5e25fSSatish Balay 104849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 104949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 105049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 105149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 105249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 105349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 105449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 105549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 105649b5e25fSSatish Balay 105749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 105849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 105949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 106049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 106149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 106249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 106349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 106449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 106549b5e25fSSatish Balay 10664d101231SSatish Balay #ifdef HAVE_MatICCFactor 10671a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 10684a2ae208SSatish Balay #undef __FUNCT__ 10694d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 10704d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level) 107149b5e25fSSatish Balay { 10724ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 107349b5e25fSSatish Balay Mat outA; 107449b5e25fSSatish Balay int ierr; 107549b5e25fSSatish Balay PetscTruth row_identity,col_identity; 107649b5e25fSSatish Balay 107749b5e25fSSatish Balay PetscFunctionBegin; 10781a3463dfSHong Zhang /* 107929bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 108049b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 108149b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 108249b5e25fSSatish Balay if (!row_identity || !col_identity) { 108329bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 108449b5e25fSSatish Balay } 10851a3463dfSHong Zhang */ 108649b5e25fSSatish Balay 108749b5e25fSSatish Balay outA = inA; 10881260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 108949b5e25fSSatish Balay 109049b5e25fSSatish Balay if (!a->diag) { 10911a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 109249b5e25fSSatish Balay } 109349b5e25fSSatish Balay /* 109449b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 109549b5e25fSSatish Balay for ILU(0) factorization with natural ordering 109649b5e25fSSatish Balay */ 109749b5e25fSSatish Balay switch (a->bs) { 109849b5e25fSSatish Balay case 1: 10990fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 110007f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 11014d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 110249b5e25fSSatish Balay case 2: 11031a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 11041a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 11051a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 11064d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 110749b5e25fSSatish Balay break; 110849b5e25fSSatish Balay case 3: 11091a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 11101a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 11111a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 11124d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 111349b5e25fSSatish Balay break; 111449b5e25fSSatish Balay case 4: 11151a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 11161a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 11171a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 11184d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 111949b5e25fSSatish Balay break; 112049b5e25fSSatish Balay case 5: 11211a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 11221a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 11231a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 11244d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 112549b5e25fSSatish Balay break; 112649b5e25fSSatish Balay case 6: 11271a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 11281a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 11291a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 11304d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 113149b5e25fSSatish Balay break; 113249b5e25fSSatish Balay case 7: 11331a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 11341a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 11351a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 11364d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 113749b5e25fSSatish Balay break; 113849b5e25fSSatish Balay default: 113949b5e25fSSatish Balay a->row = row; 11401a3463dfSHong Zhang a->icol = col; 114149b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 114249b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 114349b5e25fSSatish Balay 114449b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 114549b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1146b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 114749b5e25fSSatish Balay 114849b5e25fSSatish Balay if (!a->solve_work) { 114987828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 115087828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 115149b5e25fSSatish Balay } 115249b5e25fSSatish Balay } 115349b5e25fSSatish Balay 11541a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 115549b5e25fSSatish Balay 115649b5e25fSSatish Balay PetscFunctionReturn(0); 115749b5e25fSSatish Balay } 1158950f1e5bSHong Zhang #endif 1159950f1e5bSHong Zhang 11604a2ae208SSatish Balay #undef __FUNCT__ 11614a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 116249b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 116349b5e25fSSatish Balay { 116449b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 116549b5e25fSSatish Balay MPI_Comm comm = A->comm; 116649b5e25fSSatish Balay int ierr; 116749b5e25fSSatish Balay 116849b5e25fSSatish Balay PetscFunctionBegin; 116949b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 117049b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 117149b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 117249b5e25fSSatish Balay PetscFunctionReturn(0); 117349b5e25fSSatish Balay } 117449b5e25fSSatish Balay 117549b5e25fSSatish Balay EXTERN_C_BEGIN 11764a2ae208SSatish Balay #undef __FUNCT__ 11774a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 117849b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 117949b5e25fSSatish Balay { 1180045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 118149b5e25fSSatish Balay int i,nz,n; 118249b5e25fSSatish Balay 118349b5e25fSSatish Balay PetscFunctionBegin; 1184045c9aa0SHong Zhang nz = baij->s_maxnz; 1185c464158bSHong Zhang n = mat->n; 118649b5e25fSSatish Balay for (i=0; i<nz; i++) { 118749b5e25fSSatish Balay baij->j[i] = indices[i]; 118849b5e25fSSatish Balay } 1189045c9aa0SHong Zhang baij->s_nz = nz; 119049b5e25fSSatish Balay for (i=0; i<n; i++) { 119149b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 119249b5e25fSSatish Balay } 119349b5e25fSSatish Balay 119449b5e25fSSatish Balay PetscFunctionReturn(0); 119549b5e25fSSatish Balay } 119649b5e25fSSatish Balay EXTERN_C_END 119749b5e25fSSatish Balay 11984a2ae208SSatish Balay #undef __FUNCT__ 11994a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 120049b5e25fSSatish Balay /*@ 120119585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 120249b5e25fSSatish Balay in the matrix. 120349b5e25fSSatish Balay 120449b5e25fSSatish Balay Input Parameters: 120519585528SSatish Balay + mat - the SeqSBAIJ matrix 120649b5e25fSSatish Balay - indices - the column indices 120749b5e25fSSatish Balay 120849b5e25fSSatish Balay Level: advanced 120949b5e25fSSatish Balay 121049b5e25fSSatish Balay Notes: 121149b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 121249b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 121349b5e25fSSatish Balay of the MatSetValues() operation. 121449b5e25fSSatish Balay 121549b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 121619585528SSatish Balay MatCreateSeqSBAIJ(). 121749b5e25fSSatish Balay 1218ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 121949b5e25fSSatish Balay 1220ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 122149b5e25fSSatish Balay @*/ 122249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 122349b5e25fSSatish Balay { 122449b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 122549b5e25fSSatish Balay 122649b5e25fSSatish Balay PetscFunctionBegin; 122749b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1228c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 122949b5e25fSSatish Balay if (f) { 123049b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 123149b5e25fSSatish Balay } else { 123229bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 123349b5e25fSSatish Balay } 123449b5e25fSSatish Balay PetscFunctionReturn(0); 123549b5e25fSSatish Balay } 123649b5e25fSSatish Balay 12374a2ae208SSatish Balay #undef __FUNCT__ 12384a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1239273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1240273d9f13SBarry Smith { 1241273d9f13SBarry Smith int ierr; 1242273d9f13SBarry Smith 1243273d9f13SBarry Smith PetscFunctionBegin; 1244273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1245273d9f13SBarry Smith PetscFunctionReturn(0); 1246273d9f13SBarry Smith } 1247273d9f13SBarry Smith 124842ee4b1aSHong Zhang #include "petscblaslapack.h" 124942ee4b1aSHong Zhang #undef __FUNCT__ 125042ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 125142ee4b1aSHong Zhang int MatAXPY_SeqSBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str) 125242ee4b1aSHong Zhang { 125342ee4b1aSHong Zhang int ierr,one=1; 125442ee4b1aSHong Zhang Mat_SeqSBAIJ *x = (Mat_SeqSBAIJ *)X->data,*y = (Mat_SeqSBAIJ *)Y->data; 125542ee4b1aSHong Zhang 125642ee4b1aSHong Zhang PetscFunctionBegin; 125742ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 125842ee4b1aSHong Zhang BLaxpy_(&x->s_nz,a,x->a,&one,y->a,&one); 125942ee4b1aSHong Zhang } else { 126042ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 126142ee4b1aSHong Zhang } 126242ee4b1aSHong Zhang PetscFunctionReturn(0); 126342ee4b1aSHong Zhang } 126442ee4b1aSHong Zhang 126549b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 126649b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 126749b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 126849b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 126949b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 127049b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 127149b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 127249b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 127349b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 127449b5e25fSSatish Balay 0, 127549b5e25fSSatish Balay 0, 127649b5e25fSSatish Balay 0, 127749b5e25fSSatish Balay 0, 12785f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1279d06b337dSHong Zhang MatRelax_SeqSBAIJ, 128049b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 128149b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 128249b5e25fSSatish Balay MatEqual_SeqSBAIJ, 128349b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 128449b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 128549b5e25fSSatish Balay MatNorm_SeqSBAIJ, 128649b5e25fSSatish Balay 0, 128749b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 128849b5e25fSSatish Balay 0, 128949b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 129049b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 129149b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 129249b5e25fSSatish Balay 0, 129349b5e25fSSatish Balay 0, 12945f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12955f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1296273d9f13SBarry Smith MatSetUpPreallocation_SeqSBAIJ, 1297c464158bSHong Zhang 0, 12984d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 129949b5e25fSSatish Balay 0, 130049b5e25fSSatish Balay 0, 130149b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 130249b5e25fSSatish Balay 0, 130349b5e25fSSatish Balay 0, 130449b5e25fSSatish Balay 0, 1305950f1e5bSHong Zhang 0, 130642ee4b1aSHong Zhang MatAXPY_SeqSBAIJ, 130749b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 130849b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 130949b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 131049b5e25fSSatish Balay 0, 131149b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 131249b5e25fSSatish Balay MatScale_SeqSBAIJ, 131349b5e25fSSatish Balay 0, 131449b5e25fSSatish Balay 0, 131549b5e25fSSatish Balay 0, 131649b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 131749b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 131849b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 131949b5e25fSSatish Balay 0, 132049b5e25fSSatish Balay 0, 132149b5e25fSSatish Balay 0, 132249b5e25fSSatish Balay 0, 132349b5e25fSSatish Balay 0, 132449b5e25fSSatish Balay 0, 132549b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 132649b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 132749b5e25fSSatish Balay 0, 132849b5e25fSSatish Balay 0, 13298a124369SBarry Smith MatGetPetscMaps_Petsc, 1330d959ec07SHong Zhang 0, 1331d959ec07SHong Zhang 0, 1332d959ec07SHong Zhang 0, 1333d959ec07SHong Zhang 0, 1334d959ec07SHong Zhang 0, 1335d959ec07SHong Zhang 0, 133624d5174aSHong Zhang MatGetRowMax_SeqSBAIJ}; 133749b5e25fSSatish Balay 133849b5e25fSSatish Balay EXTERN_C_BEGIN 13394a2ae208SSatish Balay #undef __FUNCT__ 13404a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 134149b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 134249b5e25fSSatish Balay { 13434afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1344c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 134549b5e25fSSatish Balay int ierr; 134649b5e25fSSatish Balay 134749b5e25fSSatish Balay PetscFunctionBegin; 134849b5e25fSSatish Balay if (aij->nonew != 1) { 134929bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 135049b5e25fSSatish Balay } 135149b5e25fSSatish Balay 135249b5e25fSSatish Balay /* allocate space for values if not already there */ 135349b5e25fSSatish Balay if (!aij->saved_values) { 135487828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 135549b5e25fSSatish Balay } 135649b5e25fSSatish Balay 135749b5e25fSSatish Balay /* copy values over */ 135887828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 135949b5e25fSSatish Balay PetscFunctionReturn(0); 136049b5e25fSSatish Balay } 136149b5e25fSSatish Balay EXTERN_C_END 136249b5e25fSSatish Balay 136349b5e25fSSatish Balay EXTERN_C_BEGIN 13644a2ae208SSatish Balay #undef __FUNCT__ 13654a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 136649b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 136749b5e25fSSatish Balay { 13684afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1369c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 137049b5e25fSSatish Balay 137149b5e25fSSatish Balay PetscFunctionBegin; 137249b5e25fSSatish Balay if (aij->nonew != 1) { 137329bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 137449b5e25fSSatish Balay } 137549b5e25fSSatish Balay if (!aij->saved_values) { 137629bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 137749b5e25fSSatish Balay } 137849b5e25fSSatish Balay 137949b5e25fSSatish Balay /* copy values over */ 138087828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 138149b5e25fSSatish Balay PetscFunctionReturn(0); 138249b5e25fSSatish Balay } 138349b5e25fSSatish Balay EXTERN_C_END 138449b5e25fSSatish Balay 13858549e402SHong Zhang EXTERN_C_BEGIN 13864a2ae208SSatish Balay #undef __FUNCT__ 13874a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ" 1388c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1389c464158bSHong Zhang { 1390c464158bSHong Zhang Mat_SeqSBAIJ *b; 13910c955e93SHong Zhang int ierr,size; 1392c464158bSHong Zhang 1393c464158bSHong Zhang PetscFunctionBegin; 1394c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1395c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1396273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1397273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1398c464158bSHong Zhang 1399b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1400b0a32e0cSBarry Smith B->data = (void*)b; 1401c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1402c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1403c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1404c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1405c464158bSHong Zhang B->factor = 0; 1406c464158bSHong Zhang B->lupivotthreshold = 1.0; 1407c464158bSHong Zhang B->mapping = 0; 1408c464158bSHong Zhang b->row = 0; 1409c464158bSHong Zhang b->icol = 0; 1410c464158bSHong Zhang b->reallocs = 0; 1411c464158bSHong Zhang b->saved_values = 0; 1412c464158bSHong Zhang 14138a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 14148a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1415c464158bSHong Zhang 1416c464158bSHong Zhang b->sorted = PETSC_FALSE; 1417c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1418c464158bSHong Zhang b->nonew = 0; 1419c464158bSHong Zhang b->diag = 0; 1420c464158bSHong Zhang b->solve_work = 0; 1421c464158bSHong Zhang b->mult_work = 0; 14222a1b7f2aSHong Zhang B->spptr = 0; 1423c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1424c464158bSHong Zhang 1425c464158bSHong Zhang b->inew = 0; 1426c464158bSHong Zhang b->jnew = 0; 1427c464158bSHong Zhang b->anew = 0; 1428c464158bSHong Zhang b->a2anew = 0; 1429c464158bSHong Zhang b->permute = PETSC_FALSE; 1430c464158bSHong Zhang 1431c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1432c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1433c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1434c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1435c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1436c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1437c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1438c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1439c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1440c464158bSHong Zhang PetscFunctionReturn(0); 1441c464158bSHong Zhang } 14428549e402SHong Zhang EXTERN_C_END 1443c464158bSHong Zhang 14444a2ae208SSatish Balay #undef __FUNCT__ 14454a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 144649b5e25fSSatish Balay /*@C 1447273d9f13SBarry Smith MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 144849b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 144949b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 145049b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 145149b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 145249b5e25fSSatish Balay 1453c464158bSHong Zhang Collective on Mat 145449b5e25fSSatish Balay 145549b5e25fSSatish Balay Input Parameters: 1456c464158bSHong Zhang + A - the symmetric matrix 145749b5e25fSSatish Balay . bs - size of block 145849b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 1459744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1460744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 146149b5e25fSSatish Balay 146249b5e25fSSatish Balay Options Database Keys: 146349b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 146449b5e25fSSatish Balay block calculations (much slower) 146549b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 146649b5e25fSSatish Balay 146749b5e25fSSatish Balay Level: intermediate 146849b5e25fSSatish Balay 146949b5e25fSSatish Balay Notes: 147049b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 147149b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 147249b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 147349b5e25fSSatish Balay matrices. 147449b5e25fSSatish Balay 1475a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 147649b5e25fSSatish Balay @*/ 1477273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 147849b5e25fSSatish Balay { 1479c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14800c955e93SHong Zhang int i,len,ierr,mbs,bs2; 148149b5e25fSSatish Balay PetscTruth flg; 14824afc71dfSHong Zhang int s_nz; 148349b5e25fSSatish Balay 148449b5e25fSSatish Balay PetscFunctionBegin; 1485273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1486e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1487c464158bSHong Zhang mbs = B->m/bs; 148849b5e25fSSatish Balay bs2 = bs*bs; 148949b5e25fSSatish Balay 1490c464158bSHong Zhang if (mbs*bs != B->m) { 149129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 149249b5e25fSSatish Balay } 149349b5e25fSSatish Balay 1494435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1495435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 149649b5e25fSSatish Balay if (nnz) { 149749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 149829bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 149980fe4e49SBarry 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); 150049b5e25fSSatish Balay } 150149b5e25fSSatish Balay } 150249b5e25fSSatish Balay 1503e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 150449b5e25fSSatish Balay if (!flg) { 150549b5e25fSSatish Balay switch (bs) { 150649b5e25fSSatish Balay case 1: 15074ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 150849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 150949b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 151049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 151149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 151249b5e25fSSatish Balay break; 151349b5e25fSSatish Balay case 2: 15144ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 151549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 151649b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 151749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 151849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 151949b5e25fSSatish Balay break; 152049b5e25fSSatish Balay case 3: 15215f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 152249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 152349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 152449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 152549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 152649b5e25fSSatish Balay break; 152749b5e25fSSatish Balay case 4: 15285f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 152949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 153049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 153149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 153249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 153349b5e25fSSatish Balay break; 153449b5e25fSSatish Balay case 5: 15355f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 153649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 153749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 153849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 153949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 154049b5e25fSSatish Balay break; 154149b5e25fSSatish Balay case 6: 15425f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 154349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 154449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 154549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 154649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 154749b5e25fSSatish Balay break; 154849b5e25fSSatish Balay case 7: 1549de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 155049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 155149b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1552de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 155349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 155449b5e25fSSatish Balay break; 155549b5e25fSSatish Balay } 155649b5e25fSSatish Balay } 155749b5e25fSSatish Balay 155849b5e25fSSatish Balay b->mbs = mbs; 15594afc71dfSHong Zhang b->nbs = mbs; 1560b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 156149b5e25fSSatish Balay if (!nnz) { 1562435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 156349b5e25fSSatish Balay else if (nz <= 0) nz = 1; 156449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 15658cef66ccSHong Zhang b->imax[i] = nz; 156649b5e25fSSatish Balay } 1567153ea458SHong Zhang nz = nz*mbs; /* total nz */ 156849b5e25fSSatish Balay } else { 156949b5e25fSSatish Balay nz = 0; 15708cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 157149b5e25fSSatish Balay } 15728cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 15738cef66ccSHong Zhang s_nz = nz; 157449b5e25fSSatish Balay 157549b5e25fSSatish Balay /* allocate the matrix space */ 1576c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 157782502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 157849b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 157949b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 158049b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 158149b5e25fSSatish Balay b->i = b->j + s_nz; 158249b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 158349b5e25fSSatish Balay 158449b5e25fSSatish Balay /* pointer to beginning of each row */ 158549b5e25fSSatish Balay b->i[0] = 0; 158649b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 158749b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 158849b5e25fSSatish Balay } 158949b5e25fSSatish Balay 159049b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 15915033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1592b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 159349b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 159449b5e25fSSatish Balay 159549b5e25fSSatish Balay b->bs = bs; 159649b5e25fSSatish Balay b->bs2 = bs2; 159749b5e25fSSatish Balay b->s_nz = 0; 159849b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1599153ea458SHong Zhang 160016cdd363SHong Zhang b->inew = 0; 160116cdd363SHong Zhang b->jnew = 0; 160216cdd363SHong Zhang b->anew = 0; 160316cdd363SHong Zhang b->a2anew = 0; 16041a3463dfSHong Zhang b->permute = PETSC_FALSE; 1605c464158bSHong Zhang PetscFunctionReturn(0); 1606c464158bSHong Zhang } 1607153ea458SHong Zhang 160849b5e25fSSatish Balay 16094a2ae208SSatish Balay #undef __FUNCT__ 16104a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1611c464158bSHong Zhang /*@C 1612c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1613c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1614c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1615c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1616c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 161749b5e25fSSatish Balay 1618c464158bSHong Zhang Collective on MPI_Comm 1619c464158bSHong Zhang 1620c464158bSHong Zhang Input Parameters: 1621c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1622c464158bSHong Zhang . bs - size of block 1623c464158bSHong Zhang . m - number of rows, or number of columns 1624c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1625744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1626744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1627c464158bSHong Zhang 1628c464158bSHong Zhang Output Parameter: 1629c464158bSHong Zhang . A - the symmetric matrix 1630c464158bSHong Zhang 1631c464158bSHong Zhang Options Database Keys: 1632c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1633c464158bSHong Zhang block calculations (much slower) 1634c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1635c464158bSHong Zhang 1636c464158bSHong Zhang Level: intermediate 1637c464158bSHong Zhang 1638c464158bSHong Zhang Notes: 1639c464158bSHong Zhang 1640c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1641c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1642c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1643c464158bSHong Zhang matrices. 1644c464158bSHong Zhang 1645c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1646c464158bSHong Zhang @*/ 1647c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1648c464158bSHong Zhang { 1649c464158bSHong Zhang int ierr; 1650c464158bSHong Zhang 1651c464158bSHong Zhang PetscFunctionBegin; 1652c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1653c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1654273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 165549b5e25fSSatish Balay PetscFunctionReturn(0); 165649b5e25fSSatish Balay } 165749b5e25fSSatish Balay 16584a2ae208SSatish Balay #undef __FUNCT__ 16594a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 166049b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 166149b5e25fSSatish Balay { 166249b5e25fSSatish Balay Mat C; 166349b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 166449b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 166549b5e25fSSatish Balay 166649b5e25fSSatish Balay PetscFunctionBegin; 166729bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 166849b5e25fSSatish Balay 166949b5e25fSSatish Balay *B = 0; 1670692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1671692f9cbeSHong Zhang ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 1672692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1673692f9cbeSHong Zhang 167449b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1675273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 167649b5e25fSSatish Balay C->factor = A->factor; 167749b5e25fSSatish Balay c->row = 0; 167849b5e25fSSatish Balay c->icol = 0; 167949b5e25fSSatish Balay c->saved_values = 0; 168049b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 168149b5e25fSSatish Balay C->assembled = PETSC_TRUE; 168249b5e25fSSatish Balay 168349b5e25fSSatish Balay c->bs = a->bs; 168449b5e25fSSatish Balay c->bs2 = a->bs2; 168549b5e25fSSatish Balay c->mbs = a->mbs; 168649b5e25fSSatish Balay c->nbs = a->nbs; 168749b5e25fSSatish Balay 1688b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1689b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 169049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 169149b5e25fSSatish Balay c->imax[i] = a->imax[i]; 169249b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 169349b5e25fSSatish Balay } 169449b5e25fSSatish Balay 169549b5e25fSSatish Balay /* allocate the matrix space */ 169649b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 169749b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 169882502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 169949b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 170049b5e25fSSatish Balay c->i = c->j + nz; 170149b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 170249b5e25fSSatish Balay if (mbs > 0) { 170349b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 170449b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 170549b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 170649b5e25fSSatish Balay } else { 170749b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 170849b5e25fSSatish Balay } 170949b5e25fSSatish Balay } 171049b5e25fSSatish Balay 1711b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 171249b5e25fSSatish Balay c->sorted = a->sorted; 171349b5e25fSSatish Balay c->roworiented = a->roworiented; 171449b5e25fSSatish Balay c->nonew = a->nonew; 171549b5e25fSSatish Balay 171649b5e25fSSatish Balay if (a->diag) { 1717b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1718b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 171949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 172049b5e25fSSatish Balay c->diag[i] = a->diag[i]; 172149b5e25fSSatish Balay } 172249b5e25fSSatish Balay } else c->diag = 0; 172349b5e25fSSatish Balay c->s_nz = a->s_nz; 172449b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 172549b5e25fSSatish Balay c->solve_work = 0; 17262a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 172749b5e25fSSatish Balay c->mult_work = 0; 172849b5e25fSSatish Balay *B = C; 1729b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 173049b5e25fSSatish Balay PetscFunctionReturn(0); 173149b5e25fSSatish Balay } 173249b5e25fSSatish Balay 17338549e402SHong Zhang EXTERN_C_BEGIN 17344a2ae208SSatish Balay #undef __FUNCT__ 17354a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1736b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A) 173749b5e25fSSatish Balay { 173849b5e25fSSatish Balay Mat_SeqSBAIJ *a; 173949b5e25fSSatish Balay Mat B; 174049b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 17413f7326a9SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 174249b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 174349b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 174487828ca2SBarry Smith PetscScalar *aa; 174549b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 174649b5e25fSSatish Balay 174749b5e25fSSatish Balay PetscFunctionBegin; 1748b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 174949b5e25fSSatish Balay bs2 = bs*bs; 175049b5e25fSSatish Balay 175149b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 175229bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1753b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 175449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1755552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 175649b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 175749b5e25fSSatish Balay 175849b5e25fSSatish Balay if (header[3] < 0) { 175929bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 176049b5e25fSSatish Balay } 176149b5e25fSSatish Balay 176229bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 176349b5e25fSSatish Balay 176449b5e25fSSatish Balay /* 176549b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 176649b5e25fSSatish Balay divisible by the blocksize 176749b5e25fSSatish Balay */ 176849b5e25fSSatish Balay mbs = M/bs; 176949b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 177049b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 177149b5e25fSSatish Balay else mbs++; 177249b5e25fSSatish Balay if (extra_rows) { 1773b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 177449b5e25fSSatish Balay } 177549b5e25fSSatish Balay 177649b5e25fSSatish Balay /* read in row lengths */ 1777b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 177849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 177949b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 178049b5e25fSSatish Balay 178149b5e25fSSatish Balay /* read in column indices */ 1782b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 178349b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 178449b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 178549b5e25fSSatish Balay 178649b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 178782502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 1788a91a614fSBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 178982502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 179049b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 179149b5e25fSSatish Balay masked = mask + mbs; 179249b5e25fSSatish Balay rowcount = 0; nzcount = 0; 179349b5e25fSSatish Balay for (i=0; i<mbs; i++) { 179449b5e25fSSatish Balay nmask = 0; 179549b5e25fSSatish Balay for (j=0; j<bs; j++) { 179649b5e25fSSatish Balay kmax = rowlengths[rowcount]; 179749b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17982d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 179903630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 180049b5e25fSSatish Balay } 180149b5e25fSSatish Balay rowcount++; 180249b5e25fSSatish Balay } 1803574b2666SHong Zhang s_browlengths[i] += nmask; 1804574b2666SHong Zhang 180549b5e25fSSatish Balay /* zero out the mask elements we set */ 180649b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 180749b5e25fSSatish Balay } 180849b5e25fSSatish Balay 180949b5e25fSSatish Balay /* create our matrix */ 18107fe2be48SHong Zhang ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr); 181149b5e25fSSatish Balay B = *A; 181249b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 181349b5e25fSSatish Balay 181449b5e25fSSatish Balay /* set matrix "i" values */ 181549b5e25fSSatish Balay a->i[0] = 0; 181649b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1817574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1818574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 181949b5e25fSSatish Balay } 18207fe2be48SHong Zhang a->s_nz = a->i[mbs]; 182149b5e25fSSatish Balay 182249b5e25fSSatish Balay /* read in nonzero values */ 182387828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 182449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 182549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 182649b5e25fSSatish Balay 182749b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 182849b5e25fSSatish Balay nzcount = 0; jcount = 0; 182949b5e25fSSatish Balay for (i=0; i<mbs; i++) { 183049b5e25fSSatish Balay nzcountb = nzcount; 183149b5e25fSSatish Balay nmask = 0; 183249b5e25fSSatish Balay for (j=0; j<bs; j++) { 183349b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 183449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18352d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 183603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 18372d703238SHong Zhang } 18382d703238SHong Zhang } 18392d703238SHong Zhang /* sort the masked values */ 18402d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 18412d703238SHong Zhang 18422d703238SHong Zhang /* set "j" values into matrix */ 18432d703238SHong Zhang maskcount = 1; 18442d703238SHong Zhang for (j=0; j<nmask; j++) { 184549b5e25fSSatish Balay a->j[jcount++] = masked[j]; 184649b5e25fSSatish Balay mask[masked[j]] = maskcount++; 184749b5e25fSSatish Balay } 1848574b2666SHong Zhang 184949b5e25fSSatish Balay /* set "a" values into matrix */ 185049b5e25fSSatish Balay ishift = bs2*a->i[i]; 185149b5e25fSSatish Balay for (j=0; j<bs; j++) { 185249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 185349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1854574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1855574b2666SHong Zhang if (tmp >= i){ 185649b5e25fSSatish Balay block = mask[tmp] - 1; 185749b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 185849b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1859574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1860574b2666SHong Zhang } 1861574b2666SHong Zhang nzcountb++; 186249b5e25fSSatish Balay } 186349b5e25fSSatish Balay } 186449b5e25fSSatish Balay /* zero out the mask elements we set */ 186549b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 186649b5e25fSSatish Balay } 186729bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 186849b5e25fSSatish Balay 186949b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1870574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 187149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 187249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 187349b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 187449b5e25fSSatish Balay 187549b5e25fSSatish Balay B->assembled = PETSC_TRUE; 187649b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 187749b5e25fSSatish Balay PetscFunctionReturn(0); 187849b5e25fSSatish Balay } 18798549e402SHong Zhang EXTERN_C_END 1880574b2666SHong Zhang 1881d06b337dSHong Zhang #undef __FUNCT__ 1882d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1883c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx) 1884d06b337dSHong Zhang { 1885d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1886d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1887d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 1888d06b337dSHong Zhang int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr; 1889d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1890d06b337dSHong Zhang 1891d06b337dSHong Zhang PetscFunctionBegin; 1892b965ef7fSBarry Smith its = its*lits; 189391723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 1894d06b337dSHong Zhang 1895d06b337dSHong Zhang if (bs > 1) 1896d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1897d06b337dSHong Zhang 1898d06b337dSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1899d06b337dSHong Zhang if (xx != bb) { 1900d06b337dSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1901d06b337dSHong Zhang } else { 1902d06b337dSHong Zhang b = x; 1903d06b337dSHong Zhang } 1904d06b337dSHong Zhang 1905d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1906d06b337dSHong Zhang 1907d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19082798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1909d06b337dSHong Zhang for (i=0; i<m; i++) 1910d06b337dSHong Zhang t[i] = b[i]; 1911d06b337dSHong Zhang 1912d06b337dSHong Zhang for (i=0; i<m; i++){ 191344706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1914d06b337dSHong Zhang v = aa + ai[i] + 1; 1915d06b337dSHong Zhang vj = aj + ai[i] + 1; 1916d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1917d06b337dSHong Zhang x[i] = omega*t[i]/d; 1918d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 191944706875SHong Zhang PetscLogFlops(2*nz-1); 1920d06b337dSHong Zhang } 1921d06b337dSHong Zhang } 1922d06b337dSHong Zhang 19232798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1924d06b337dSHong Zhang for (i=0; i<m; i++) 1925d06b337dSHong Zhang t[i] = b[i]; 1926d06b337dSHong Zhang 1927d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1928d06b337dSHong Zhang v = aa + ai[i] + 1; 1929d06b337dSHong Zhang vj = aj + ai[i] + 1; 1930d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1931d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 193244706875SHong Zhang PetscLogFlops(2*nz-1); 1933d06b337dSHong Zhang } 1934d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1935d06b337dSHong Zhang d = *(aa + ai[i]); 1936d06b337dSHong Zhang v = aa + ai[i] + 1; 1937d06b337dSHong Zhang vj = aj + ai[i] + 1; 1938d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1939d06b337dSHong Zhang sum = t[i]; 1940d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 194144706875SHong Zhang PetscLogFlops(2*nz-1); 1942d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1943d06b337dSHong Zhang } 1944d06b337dSHong Zhang } 1945d06b337dSHong Zhang its--; 1946d06b337dSHong Zhang } 1947d06b337dSHong Zhang 1948d06b337dSHong Zhang while (its--) { 194944706875SHong Zhang /* 195044706875SHong Zhang forward sweep: 195144706875SHong Zhang for i=0,...,m-1: 195244706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 195344706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 195444706875SHong Zhang b = b - x[i]*U^T(i,:); 1955d06b337dSHong Zhang 195644706875SHong Zhang */ 19572798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1958d06b337dSHong Zhang for (i=0; i<m; i++) 1959d06b337dSHong Zhang t[i] = b[i]; 1960d06b337dSHong Zhang 1961d06b337dSHong Zhang for (i=0; i<m; i++){ 196244706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1963d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 1964d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 1965d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 1966d06b337dSHong Zhang sum = t[i]; 1967d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 1968d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 1969d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 197044706875SHong Zhang PetscLogFlops(4*nz-2); 1971d06b337dSHong Zhang } 1972d06b337dSHong Zhang } 1973d06b337dSHong Zhang 19742798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 197544706875SHong Zhang /* 197644706875SHong Zhang backward sweep: 197744706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 197844706875SHong Zhang for i=m-1,...,0: 197944706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 198044706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 198144706875SHong Zhang */ 1982d06b337dSHong Zhang for (i=0; i<m; i++) 1983d06b337dSHong Zhang t[i] = b[i]; 1984d06b337dSHong Zhang 1985d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1986d06b337dSHong Zhang v = aa + ai[i] + 1; 1987d06b337dSHong Zhang vj = aj + ai[i] + 1; 1988d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1989d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 199044706875SHong Zhang PetscLogFlops(2*nz-1); 1991d06b337dSHong Zhang } 1992d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1993d06b337dSHong Zhang d = *(aa + ai[i]); 1994d06b337dSHong Zhang v = aa + ai[i] + 1; 1995d06b337dSHong Zhang vj = aj + ai[i] + 1; 1996d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1997d06b337dSHong Zhang sum = t[i]; 1998d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 199944706875SHong Zhang PetscLogFlops(2*nz-1); 2000d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2001d06b337dSHong Zhang } 2002d06b337dSHong Zhang } 2003d06b337dSHong Zhang } 2004d06b337dSHong Zhang 2005d06b337dSHong Zhang ierr = PetscFree(t); CHKERRQ(ierr); 2006d06b337dSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2007d06b337dSHong Zhang if (bb != xx) { 2008d06b337dSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 2009d06b337dSHong Zhang } 20102798e883SHong Zhang 2011d06b337dSHong Zhang PetscFunctionReturn(0); 2012d06b337dSHong Zhang } 2013d06b337dSHong Zhang 2014d06b337dSHong Zhang 2015d06b337dSHong Zhang 2016d06b337dSHong Zhang 201749b5e25fSSatish Balay 201849b5e25fSSatish Balay 2019