1*8a124369SBarry Smith /*$Id: sbaij.c,v 1.59 2001/06/22 00:20:57 buschelm Exp bsmith $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 449b5e25fSSatish Balay Defines the basic matrix operations for the BAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 7ab9f2c04SSatish Balay #include "petscsys.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h" 949b5e25fSSatish Balay #include "src/vec/vecimpl.h" 1049b5e25fSSatish Balay #include "src/inline/spops.h" 11aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay #define CHUNKSIZE 10 1449b5e25fSSatish Balay 1549b5e25fSSatish Balay /* 1649b5e25fSSatish Balay Checks for missing diagonals 1749b5e25fSSatish Balay */ 184a2ae208SSatish Balay #undef __FUNCT__ 194a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 2049b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2149b5e25fSSatish Balay { 22045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2349b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2449b5e25fSSatish Balay 2549b5e25fSSatish Balay PetscFunctionBegin; 26045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2749b5e25fSSatish Balay diag = a->diag; 2849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2949b5e25fSSatish Balay if (jj[diag[i]] != i) { 3029bbc08cSBarry Smith SETERRQ1(1,"Matrix is missing diagonal number %d",i); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay } 3349b5e25fSSatish Balay PetscFunctionReturn(0); 3449b5e25fSSatish Balay } 3549b5e25fSSatish Balay 364a2ae208SSatish Balay #undef __FUNCT__ 374a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 3849b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3949b5e25fSSatish Balay { 40045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 4182502324SSatish Balay int i,j,*diag,m = a->mbs,ierr; 4249b5e25fSSatish Balay 4349b5e25fSSatish Balay PetscFunctionBegin; 4449b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4549b5e25fSSatish Balay 46b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 47b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 4849b5e25fSSatish Balay for (i=0; i<m; i++) { 4949b5e25fSSatish Balay diag[i] = a->i[i+1]; 5049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 5149b5e25fSSatish Balay if (a->j[j] == i) { 5249b5e25fSSatish Balay diag[i] = j; 5349b5e25fSSatish Balay break; 5449b5e25fSSatish Balay } 5549b5e25fSSatish Balay } 5649b5e25fSSatish Balay } 5749b5e25fSSatish Balay a->diag = diag; 5849b5e25fSSatish Balay PetscFunctionReturn(0); 5949b5e25fSSatish Balay } 6049b5e25fSSatish Balay 6149b5e25fSSatish Balay 6249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 6349b5e25fSSatish Balay 644a2ae208SSatish Balay #undef __FUNCT__ 654a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 6649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 6749b5e25fSSatish Balay { 6849b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6949b5e25fSSatish Balay 7049b5e25fSSatish Balay PetscFunctionBegin; 71045c9aa0SHong Zhang if (ia) { 7229bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 7349b5e25fSSatish Balay } 74045c9aa0SHong Zhang *nn = a->mbs; 7549b5e25fSSatish Balay PetscFunctionReturn(0); 7649b5e25fSSatish Balay } 7749b5e25fSSatish Balay 784a2ae208SSatish Balay #undef __FUNCT__ 794a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 80045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 8149b5e25fSSatish Balay { 8249b5e25fSSatish Balay PetscFunctionBegin; 8349b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 8429bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 8516cdd363SHong Zhang /* PetscFunctionReturn(0); */ 8649b5e25fSSatish Balay } 8749b5e25fSSatish Balay 884a2ae208SSatish Balay #undef __FUNCT__ 894a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 9049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9149b5e25fSSatish Balay { 9249b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 9349b5e25fSSatish Balay 9449b5e25fSSatish Balay PetscFunctionBegin; 9549b5e25fSSatish Balay *bs = sbaij->bs; 9649b5e25fSSatish Balay PetscFunctionReturn(0); 9749b5e25fSSatish Balay } 9849b5e25fSSatish Balay 994a2ae208SSatish Balay #undef __FUNCT__ 1004a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 10149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 10249b5e25fSSatish Balay { 10349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10449b5e25fSSatish Balay int ierr; 10549b5e25fSSatish Balay 10649b5e25fSSatish Balay PetscFunctionBegin; 10749b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 108b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 10949b5e25fSSatish Balay #endif 11049b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 11149b5e25fSSatish Balay if (!a->singlemalloc) { 11249b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 11363c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 11449b5e25fSSatish Balay } 11549b5e25fSSatish Balay if (a->row) { 11649b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 11749b5e25fSSatish Balay } 11849b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 11949b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 12049b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 12149b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 12249b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 12349b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 12449b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1251a3463dfSHong Zhang 1261a3463dfSHong Zhang if (a->inew){ 1271a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1281a3463dfSHong Zhang a->inew = 0; 1291a3463dfSHong Zhang } 13049b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 13149b5e25fSSatish Balay PetscFunctionReturn(0); 13249b5e25fSSatish Balay } 13349b5e25fSSatish Balay 1344a2ae208SSatish Balay #undef __FUNCT__ 1354a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 13649b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 13749b5e25fSSatish Balay { 138045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 13949b5e25fSSatish Balay 14049b5e25fSSatish Balay PetscFunctionBegin; 1414d9d31abSKris Buschelman switch (op) { 1424d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1434d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1444d9d31abSKris Buschelman break; 1454d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1464d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1474d9d31abSKris Buschelman break; 1484d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1494d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1504d9d31abSKris Buschelman break; 1514d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1524d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1554d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1584d9d31abSKris Buschelman a->nonew = 1; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1614d9d31abSKris Buschelman a->nonew = -1; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1644d9d31abSKris Buschelman a->nonew = -2; 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1674d9d31abSKris Buschelman a->nonew = 0; 1684d9d31abSKris Buschelman break; 1694d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1704d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1714d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1724d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1734d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 174d03495bdSKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 175b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1764d9d31abSKris Buschelman break; 1774d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 17829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1794d9d31abSKris Buschelman break; 1804d9d31abSKris Buschelman default: 18129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 1824d9d31abSKris Buschelman break; 18349b5e25fSSatish Balay } 18449b5e25fSSatish Balay PetscFunctionReturn(0); 18549b5e25fSSatish Balay } 18649b5e25fSSatish Balay 1874a2ae208SSatish Balay #undef __FUNCT__ 1884a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_SeqSBAIJ" 18949b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n) 19049b5e25fSSatish Balay { 19149b5e25fSSatish Balay PetscFunctionBegin; 192c464158bSHong Zhang *m = 0; 1937e8a2ce2SHong Zhang *n = A->m; 19449b5e25fSSatish Balay PetscFunctionReturn(0); 19549b5e25fSSatish Balay } 19649b5e25fSSatish Balay 1974a2ae208SSatish Balay #undef __FUNCT__ 1984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 19949b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v) 20049b5e25fSSatish Balay { 20149b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 20282502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr; 20349b5e25fSSatish Balay MatScalar *aa,*aa_i; 20449b5e25fSSatish Balay Scalar *v_i; 20549b5e25fSSatish Balay 20649b5e25fSSatish Balay PetscFunctionBegin; 20749b5e25fSSatish Balay bs = a->bs; 20849b5e25fSSatish Balay ai = a->i; 20949b5e25fSSatish Balay aj = a->j; 21049b5e25fSSatish Balay aa = a->a; 21149b5e25fSSatish Balay bs2 = a->bs2; 21249b5e25fSSatish Balay 213c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 21449b5e25fSSatish Balay 21549b5e25fSSatish Balay bn = row/bs; /* Block number */ 21649b5e25fSSatish Balay bp = row % bs; /* Block position */ 21749b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 21849b5e25fSSatish Balay *ncols = bs*M; 21949b5e25fSSatish Balay 22049b5e25fSSatish Balay if (v) { 22149b5e25fSSatish Balay *v = 0; 22249b5e25fSSatish Balay if (*ncols) { 22382502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(Scalar),v);CHKERRQ(ierr); 22449b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22549b5e25fSSatish Balay v_i = *v + i*bs; 22649b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 22749b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 22849b5e25fSSatish Balay } 22949b5e25fSSatish Balay } 23049b5e25fSSatish Balay } 23149b5e25fSSatish Balay 23249b5e25fSSatish Balay if (cols) { 23349b5e25fSSatish Balay *cols = 0; 23449b5e25fSSatish Balay if (*ncols) { 23582502324SSatish Balay ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr); 23649b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23749b5e25fSSatish Balay cols_i = *cols + i*bs; 23849b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 23949b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 24049b5e25fSSatish Balay } 24149b5e25fSSatish Balay } 24249b5e25fSSatish Balay } 24349b5e25fSSatish Balay 24449b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2455ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2465ddb2528SHong Zhang #ifdef column_search 24749b5e25fSSatish Balay v_i = *v + M*bs; 24849b5e25fSSatish Balay cols_i = *cols + M*bs; 24949b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 25049b5e25fSSatish Balay M = ai[i+1] - ai[i]; 25149b5e25fSSatish Balay for (j=0; j<M; j++){ 25249b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 25349b5e25fSSatish Balay if (itmp == bn){ 25449b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 25549b5e25fSSatish Balay for (k=0; k<bs; k++) { 25649b5e25fSSatish Balay *cols_i++ = i*bs+k; 25749b5e25fSSatish Balay *v_i++ = aa_i[k]; 25849b5e25fSSatish Balay } 25949b5e25fSSatish Balay *ncols += bs; 26049b5e25fSSatish Balay break; 26149b5e25fSSatish Balay } 26249b5e25fSSatish Balay } 26349b5e25fSSatish Balay } 2645ddb2528SHong Zhang #endif 26549b5e25fSSatish Balay 26649b5e25fSSatish Balay PetscFunctionReturn(0); 26749b5e25fSSatish Balay } 26849b5e25fSSatish Balay 2694a2ae208SSatish Balay #undef __FUNCT__ 2704a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ" 27149b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 27249b5e25fSSatish Balay { 27349b5e25fSSatish Balay int ierr; 27449b5e25fSSatish Balay 27549b5e25fSSatish Balay PetscFunctionBegin; 27649b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 27749b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 27849b5e25fSSatish Balay PetscFunctionReturn(0); 27949b5e25fSSatish Balay } 28049b5e25fSSatish Balay 2814a2ae208SSatish Balay #undef __FUNCT__ 2824a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ" 28349b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 28449b5e25fSSatish Balay { 28549b5e25fSSatish Balay PetscFunctionBegin; 28629bbc08cSBarry Smith SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called"); 28716cdd363SHong Zhang /* PetscFunctionReturn(0); */ 28849b5e25fSSatish Balay } 28949b5e25fSSatish Balay 2904a2ae208SSatish Balay #undef __FUNCT__ 2914a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary" 292b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer) 29349b5e25fSSatish Balay { 29449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 29549b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 29649b5e25fSSatish Balay Scalar *aa; 29749b5e25fSSatish Balay FILE *file; 29849b5e25fSSatish Balay 29949b5e25fSSatish Balay PetscFunctionBegin; 300b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 30182502324SSatish Balay ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 30249b5e25fSSatish Balay col_lens[0] = MAT_COOKIE; 30349b5e25fSSatish Balay 304c464158bSHong Zhang col_lens[1] = A->m; 305c464158bSHong Zhang col_lens[2] = A->m; 30649b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 30749b5e25fSSatish Balay 30849b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 30949b5e25fSSatish Balay count = 0; 31049b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31149b5e25fSSatish Balay for (j=0; j<bs; j++) { 31249b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 31349b5e25fSSatish Balay } 31449b5e25fSSatish Balay } 315c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 31649b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 31749b5e25fSSatish Balay 31849b5e25fSSatish Balay /* store column indices (zero start index) */ 31982502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 32049b5e25fSSatish Balay count = 0; 32149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 32249b5e25fSSatish Balay for (j=0; j<bs; j++) { 32349b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 32449b5e25fSSatish Balay for (l=0; l<bs; l++) { 32549b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 32649b5e25fSSatish Balay } 32749b5e25fSSatish Balay } 32849b5e25fSSatish Balay } 32949b5e25fSSatish Balay } 33049b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 33149b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 33249b5e25fSSatish Balay 33349b5e25fSSatish Balay /* store nonzero values */ 33482502324SSatish Balay ierr = PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar),&aa);CHKERRQ(ierr); 33549b5e25fSSatish Balay count = 0; 33649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33749b5e25fSSatish Balay for (j=0; j<bs; j++) { 33849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33949b5e25fSSatish Balay for (l=0; l<bs; l++) { 34049b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 34149b5e25fSSatish Balay } 34249b5e25fSSatish Balay } 34349b5e25fSSatish Balay } 34449b5e25fSSatish Balay } 34549b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 34649b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 34749b5e25fSSatish Balay 348b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 34949b5e25fSSatish Balay if (file) { 35049b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 35149b5e25fSSatish Balay } 35249b5e25fSSatish Balay PetscFunctionReturn(0); 35349b5e25fSSatish Balay } 35449b5e25fSSatish Balay 3554a2ae208SSatish Balay #undef __FUNCT__ 3564a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 357b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 35849b5e25fSSatish Balay { 35949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 360fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 361fb9695e5SSatish Balay char *name; 362f3ef73ceSBarry Smith PetscViewerFormat format; 36349b5e25fSSatish Balay 36449b5e25fSSatish Balay PetscFunctionBegin; 36580fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 366b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 367fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) { 368b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 369fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 37029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 371fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 372b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 37349b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 37449b5e25fSSatish Balay for (j=0; j<bs; j++) { 375b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 37649b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 37749b5e25fSSatish Balay for (l=0; l<bs; l++) { 37849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37949b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 380b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 38149b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38249b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 383b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 38449b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38549b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 386b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38749b5e25fSSatish Balay } 38849b5e25fSSatish Balay #else 38949b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 390b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 39149b5e25fSSatish Balay } 39249b5e25fSSatish Balay #endif 39349b5e25fSSatish Balay } 39449b5e25fSSatish Balay } 395b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 39649b5e25fSSatish Balay } 39749b5e25fSSatish Balay } 398b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 39949b5e25fSSatish Balay } else { 400b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 40149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 40249b5e25fSSatish Balay for (j=0; j<bs; j++) { 403b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 40449b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 40549b5e25fSSatish Balay for (l=0; l<bs; l++) { 40649b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 40749b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 408b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 40949b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 41049b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 411b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 41249b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 41349b5e25fSSatish Balay } else { 414b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 41549b5e25fSSatish Balay } 41649b5e25fSSatish Balay #else 417b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 41849b5e25fSSatish Balay #endif 41949b5e25fSSatish Balay } 42049b5e25fSSatish Balay } 421b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 42249b5e25fSSatish Balay } 42349b5e25fSSatish Balay } 424b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 42549b5e25fSSatish Balay } 426b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 42749b5e25fSSatish Balay PetscFunctionReturn(0); 42849b5e25fSSatish Balay } 42949b5e25fSSatish Balay 4304a2ae208SSatish Balay #undef __FUNCT__ 4314a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 432b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 43349b5e25fSSatish Balay { 43449b5e25fSSatish Balay Mat A = (Mat) Aa; 43549b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 43649b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 43749b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 43849b5e25fSSatish Balay MatScalar *aa; 43949b5e25fSSatish Balay MPI_Comm comm; 440b0a32e0cSBarry Smith PetscViewer viewer; 44149b5e25fSSatish Balay 44249b5e25fSSatish Balay PetscFunctionBegin; 44349b5e25fSSatish Balay /* 44449b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 44549b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 44649b5e25fSSatish Balay rest should return immediately. 44749b5e25fSSatish Balay */ 44849b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 44949b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 45049b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 45149b5e25fSSatish Balay 45249b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 45349b5e25fSSatish Balay 454b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 455b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 45649b5e25fSSatish Balay 45749b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 458b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 45949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 461c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 46249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 46349b5e25fSSatish Balay aa = a->a + j*bs2; 46449b5e25fSSatish Balay for (k=0; k<bs; k++) { 46549b5e25fSSatish Balay for (l=0; l<bs; l++) { 46649b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 467b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46849b5e25fSSatish Balay } 46949b5e25fSSatish Balay } 47049b5e25fSSatish Balay } 47149b5e25fSSatish Balay } 472b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 47349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 47449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 475c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 47649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47749b5e25fSSatish Balay aa = a->a + j*bs2; 47849b5e25fSSatish Balay for (k=0; k<bs; k++) { 47949b5e25fSSatish Balay for (l=0; l<bs; l++) { 48049b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 481b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 48249b5e25fSSatish Balay } 48349b5e25fSSatish Balay } 48449b5e25fSSatish Balay } 48549b5e25fSSatish Balay } 48649b5e25fSSatish Balay 487b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 48849b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 48949b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 490c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 49149b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 49249b5e25fSSatish Balay aa = a->a + j*bs2; 49349b5e25fSSatish Balay for (k=0; k<bs; k++) { 49449b5e25fSSatish Balay for (l=0; l<bs; l++) { 49549b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 496b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 49749b5e25fSSatish Balay } 49849b5e25fSSatish Balay } 49949b5e25fSSatish Balay } 50049b5e25fSSatish Balay } 50149b5e25fSSatish Balay PetscFunctionReturn(0); 50249b5e25fSSatish Balay } 50349b5e25fSSatish Balay 5044a2ae208SSatish Balay #undef __FUNCT__ 5054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 506b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 50749b5e25fSSatish Balay { 50849b5e25fSSatish Balay int ierr; 50949b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 510b0a32e0cSBarry Smith PetscDraw draw; 51149b5e25fSSatish Balay PetscTruth isnull; 51249b5e25fSSatish Balay 51349b5e25fSSatish Balay PetscFunctionBegin; 51449b5e25fSSatish Balay 515b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 516b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 51749b5e25fSSatish Balay 51849b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 519c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 52049b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 521b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 522b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 52349b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 52449b5e25fSSatish Balay PetscFunctionReturn(0); 52549b5e25fSSatish Balay } 52649b5e25fSSatish Balay 5274a2ae208SSatish Balay #undef __FUNCT__ 5284a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 529b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 53049b5e25fSSatish Balay { 53149b5e25fSSatish Balay int ierr; 53249b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 53349b5e25fSSatish Balay 53449b5e25fSSatish Balay PetscFunctionBegin; 535b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 536b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 537fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 538fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 53949b5e25fSSatish Balay if (issocket) { 54029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 54149b5e25fSSatish Balay } else if (isascii){ 54249b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 54349b5e25fSSatish Balay } else if (isbinary) { 54449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 54549b5e25fSSatish Balay } else if (isdraw) { 54649b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 54749b5e25fSSatish Balay } else { 54829bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 54949b5e25fSSatish Balay } 55049b5e25fSSatish Balay PetscFunctionReturn(0); 55149b5e25fSSatish Balay } 55249b5e25fSSatish Balay 55349b5e25fSSatish Balay 5544a2ae208SSatish Balay #undef __FUNCT__ 5554a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 55649b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 55749b5e25fSSatish Balay { 558045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 55949b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 56049b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 56149b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 56249b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 56349b5e25fSSatish Balay 56449b5e25fSSatish Balay PetscFunctionBegin; 56549b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 56649b5e25fSSatish Balay row = im[k]; brow = row/bs; 56729bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 568c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56949b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 57049b5e25fSSatish Balay nrow = ailen[brow]; 57149b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 57229bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 573c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 57449b5e25fSSatish Balay col = in[l] ; 57549b5e25fSSatish Balay bcol = col/bs; 57649b5e25fSSatish Balay cidx = col%bs; 57749b5e25fSSatish Balay ridx = row%bs; 57849b5e25fSSatish Balay high = nrow; 57949b5e25fSSatish Balay low = 0; /* assume unsorted */ 58049b5e25fSSatish Balay while (high-low > 5) { 58149b5e25fSSatish Balay t = (low+high)/2; 58249b5e25fSSatish Balay if (rp[t] > bcol) high = t; 58349b5e25fSSatish Balay else low = t; 58449b5e25fSSatish Balay } 58549b5e25fSSatish Balay for (i=low; i<high; i++) { 58649b5e25fSSatish Balay if (rp[i] > bcol) break; 58749b5e25fSSatish Balay if (rp[i] == bcol) { 58849b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 58949b5e25fSSatish Balay goto finished; 59049b5e25fSSatish Balay } 59149b5e25fSSatish Balay } 59249b5e25fSSatish Balay *v++ = zero; 59349b5e25fSSatish Balay finished:; 59449b5e25fSSatish Balay } 59549b5e25fSSatish Balay } 59649b5e25fSSatish Balay PetscFunctionReturn(0); 59749b5e25fSSatish Balay } 59849b5e25fSSatish Balay 59949b5e25fSSatish Balay 6004a2ae208SSatish Balay #undef __FUNCT__ 6014a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 60249b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 60349b5e25fSSatish Balay { 60449b5e25fSSatish Balay PetscFunctionBegin; 60529bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format"); 60616cdd363SHong Zhang /* PetscFunctionReturn(0); */ 60749b5e25fSSatish Balay } 60849b5e25fSSatish Balay 6094a2ae208SSatish Balay #undef __FUNCT__ 6104a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 61149b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 61249b5e25fSSatish Balay { 61349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 61449b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 615c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 61649b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 61749b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 61849b5e25fSSatish Balay 61949b5e25fSSatish Balay PetscFunctionBegin; 62049b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 62149b5e25fSSatish Balay 62249b5e25fSSatish Balay if (m) rmax = ailen[0]; 62349b5e25fSSatish Balay for (i=1; i<mbs; i++) { 62449b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 62549b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 62649b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 62749b5e25fSSatish Balay if (fshift) { 62849b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 62949b5e25fSSatish Balay N = ailen[i]; 63049b5e25fSSatish Balay for (j=0; j<N; j++) { 63149b5e25fSSatish Balay ip[j-fshift] = ip[j]; 63249b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 63349b5e25fSSatish Balay } 63449b5e25fSSatish Balay } 63549b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 63649b5e25fSSatish Balay } 63749b5e25fSSatish Balay if (mbs) { 63849b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 63949b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 64049b5e25fSSatish Balay } 64149b5e25fSSatish Balay /* reset ilen and imax for each row */ 64249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 64349b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 64449b5e25fSSatish Balay } 64549b5e25fSSatish Balay a->s_nz = ai[mbs]; 64649b5e25fSSatish Balay 64749b5e25fSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 64849b5e25fSSatish Balay if (fshift && a->diag) { 64949b5e25fSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 650b0a32e0cSBarry Smith PetscLogObjectMemory(A,-(m+1)*sizeof(int)); 65149b5e25fSSatish Balay a->diag = 0; 65249b5e25fSSatish Balay } 653b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 654c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 655b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 65649b5e25fSSatish Balay a->reallocs); 657b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 65849b5e25fSSatish Balay a->reallocs = 0; 65949b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 66049b5e25fSSatish Balay 66149b5e25fSSatish Balay PetscFunctionReturn(0); 66249b5e25fSSatish Balay } 66349b5e25fSSatish Balay 66449b5e25fSSatish Balay /* 66549b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 66649b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 66749b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 66849b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 66949b5e25fSSatish Balay */ 6704a2ae208SSatish Balay #undef __FUNCT__ 6714a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 67249b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 67349b5e25fSSatish Balay { 67449b5e25fSSatish Balay int i,j,k,row; 67549b5e25fSSatish Balay PetscTruth flg; 67649b5e25fSSatish Balay 67749b5e25fSSatish Balay PetscFunctionBegin; 67849b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 67949b5e25fSSatish Balay row = idx[i]; 68049b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 68149b5e25fSSatish Balay sizes[j] = 1; 68249b5e25fSSatish Balay i++; 68349b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 68449b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 68549b5e25fSSatish Balay i++; 68649b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 68749b5e25fSSatish Balay flg = PETSC_TRUE; 68849b5e25fSSatish Balay for (k=1; k<bs; k++) { 68949b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 69049b5e25fSSatish Balay flg = PETSC_FALSE; 69149b5e25fSSatish Balay break; 69249b5e25fSSatish Balay } 69349b5e25fSSatish Balay } 69449b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 69549b5e25fSSatish Balay sizes[j] = bs; 69649b5e25fSSatish Balay i+= bs; 69749b5e25fSSatish Balay } else { 69849b5e25fSSatish Balay sizes[j] = 1; 69949b5e25fSSatish Balay i++; 70049b5e25fSSatish Balay } 70149b5e25fSSatish Balay } 70249b5e25fSSatish Balay } 70349b5e25fSSatish Balay *bs_max = j; 70449b5e25fSSatish Balay PetscFunctionReturn(0); 70549b5e25fSSatish Balay } 70649b5e25fSSatish Balay 7074a2ae208SSatish Balay #undef __FUNCT__ 7084a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 70949b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag) 71049b5e25fSSatish Balay { 71149b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data; 71249b5e25fSSatish Balay int ierr,i,j,k,count,is_n,*is_idx,*rows; 71349b5e25fSSatish Balay int bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max; 71449b5e25fSSatish Balay Scalar zero = 0.0; 71549b5e25fSSatish Balay MatScalar *aa; 71649b5e25fSSatish Balay 71749b5e25fSSatish Balay PetscFunctionBegin; 71849b5e25fSSatish Balay /* Make a copy of the IS and sort it */ 71949b5e25fSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 72049b5e25fSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 72149b5e25fSSatish Balay 72249b5e25fSSatish Balay /* allocate memory for rows,sizes */ 723b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 72449b5e25fSSatish Balay sizes = rows + is_n; 72549b5e25fSSatish Balay 72649b5e25fSSatish Balay /* initialize copy IS values to rows, and sort them */ 72749b5e25fSSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 72849b5e25fSSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 72949b5e25fSSatish Balay if (sbaij->keepzeroedrows) { /* do not change nonzero structure */ 73049b5e25fSSatish Balay for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */ 73149b5e25fSSatish Balay bs_max = is_n; /* bs_max: num. of contiguous block row in the row */ 73249b5e25fSSatish Balay } else { 73349b5e25fSSatish Balay ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 73449b5e25fSSatish Balay } 73549b5e25fSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 73649b5e25fSSatish Balay 73749b5e25fSSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 73849b5e25fSSatish Balay row = rows[j]; /* row to be zeroed */ 739c464158bSHong Zhang if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 74049b5e25fSSatish Balay count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */ 74149b5e25fSSatish Balay aa = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs); 74249b5e25fSSatish Balay if (sizes[i] == bs && !sbaij->keepzeroedrows) { 74349b5e25fSSatish Balay if (diag) { 74449b5e25fSSatish Balay if (sbaij->ilen[row/bs] > 0) { 74549b5e25fSSatish Balay sbaij->ilen[row/bs] = 1; 74649b5e25fSSatish Balay sbaij->j[sbaij->i[row/bs]] = row/bs; 74749b5e25fSSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 74849b5e25fSSatish Balay } 74949b5e25fSSatish Balay /* Now insert all the diagoanl values for this bs */ 75049b5e25fSSatish Balay for (k=0; k<bs; k++) { 75149b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 75249b5e25fSSatish Balay } 75349b5e25fSSatish Balay } else { /* (!diag) */ 75449b5e25fSSatish Balay sbaij->ilen[row/bs] = 0; 75549b5e25fSSatish Balay } /* end (!diag) */ 75649b5e25fSSatish Balay } else { /* (sizes[i] != bs), broken block */ 75749b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG) 75829bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 75949b5e25fSSatish Balay #endif 76049b5e25fSSatish Balay for (k=0; k<count; k++) { 76149b5e25fSSatish Balay aa[0] = zero; 76249b5e25fSSatish Balay aa+=bs; 76349b5e25fSSatish Balay } 76449b5e25fSSatish Balay if (diag) { 76549b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 76649b5e25fSSatish Balay } 76749b5e25fSSatish Balay } 76849b5e25fSSatish Balay } 76949b5e25fSSatish Balay 77049b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 77149b5e25fSSatish Balay ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77249b5e25fSSatish Balay ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77349b5e25fSSatish Balay PetscFunctionReturn(0); 77449b5e25fSSatish Balay } 77549b5e25fSSatish Balay 77649b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 77749b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 77849b5e25fSSatish Balay */ 77949b5e25fSSatish Balay 7804a2ae208SSatish Balay #undef __FUNCT__ 7814a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 78249b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 78349b5e25fSSatish Balay { 78449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 78549b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 78649b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 78749b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 78849b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 78949b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 79049b5e25fSSatish Balay 79149b5e25fSSatish Balay PetscFunctionBegin; 79249b5e25fSSatish Balay 79349b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 79449b5e25fSSatish Balay row = im[k]; /* row number */ 79549b5e25fSSatish Balay brow = row/bs; /* block row number */ 79649b5e25fSSatish Balay if (row < 0) continue; 79749b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 798c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 79949b5e25fSSatish Balay #endif 80049b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 80149b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 80249b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 80349b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 80449b5e25fSSatish Balay low = 0; 80549b5e25fSSatish Balay 80649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 80749b5e25fSSatish Balay if (in[l] < 0) continue; 80849b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 809c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 81049b5e25fSSatish Balay #endif 81149b5e25fSSatish Balay col = in[l]; 81249b5e25fSSatish Balay bcol = col/bs; /* block col number */ 81349b5e25fSSatish Balay 81449b5e25fSSatish Balay if (brow <= bcol){ 81549b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8168549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 81749b5e25fSSatish Balay /* element value a(k,l) */ 81849b5e25fSSatish Balay if (roworiented) { 81949b5e25fSSatish Balay value = v[l + k*n]; 82049b5e25fSSatish Balay } else { 82149b5e25fSSatish Balay value = v[k + l*m]; 82249b5e25fSSatish Balay } 82349b5e25fSSatish Balay 82449b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 82549b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 82649b5e25fSSatish Balay while (high-low > 7) { 82749b5e25fSSatish Balay t = (low+high)/2; 82849b5e25fSSatish Balay if (rp[t] > bcol) high = t; 82949b5e25fSSatish Balay else low = t; 83049b5e25fSSatish Balay } 83149b5e25fSSatish Balay for (i=low; i<high; i++) { 83249b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 83349b5e25fSSatish Balay if (rp[i] > bcol) break; 83449b5e25fSSatish Balay if (rp[i] == bcol) { 83549b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 83649b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 83749b5e25fSSatish Balay else *bap = value; 8388549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8398549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8408549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8418549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8428549e402SHong Zhang else *bap = value; 8438549e402SHong Zhang } 84449b5e25fSSatish Balay goto noinsert1; 84549b5e25fSSatish Balay } 84649b5e25fSSatish Balay } 84749b5e25fSSatish Balay 84849b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 84929bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 85049b5e25fSSatish Balay if (nrow >= rmax) { 85149b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 85249b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 85349b5e25fSSatish Balay MatScalar *new_a; 85449b5e25fSSatish Balay 85529bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 85649b5e25fSSatish Balay 85749b5e25fSSatish Balay /* Malloc new storage space */ 85849b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 85982502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 86049b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 86149b5e25fSSatish Balay new_i = new_j + new_nz; 86249b5e25fSSatish Balay 86349b5e25fSSatish Balay /* copy over old data into new slots */ 86449b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 86549b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 86649b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 86749b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 86849b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 86949b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 87049b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 87149b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 87249b5e25fSSatish Balay /* free up old matrix storage */ 87349b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 87449b5e25fSSatish Balay if (!a->singlemalloc) { 87549b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 87649b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 87749b5e25fSSatish Balay } 87849b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 87949b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 88049b5e25fSSatish Balay 88149b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 88249b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 883b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 88449b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 88549b5e25fSSatish Balay a->reallocs++; 88649b5e25fSSatish Balay a->s_nz++; 88749b5e25fSSatish Balay } 88849b5e25fSSatish Balay 88949b5e25fSSatish Balay N = nrow++ - 1; 89049b5e25fSSatish Balay /* shift up all the later entries in this row */ 89149b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 89249b5e25fSSatish Balay rp[ii+1] = rp[ii]; 89349b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 89449b5e25fSSatish Balay } 89549b5e25fSSatish Balay if (N>=i) { 89649b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 89749b5e25fSSatish Balay } 89849b5e25fSSatish Balay rp[i] = bcol; 89949b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 90049b5e25fSSatish Balay noinsert1:; 90149b5e25fSSatish Balay low = i; 90249b5e25fSSatish Balay /* } */ 9038549e402SHong Zhang } 90449b5e25fSSatish Balay } /* end of if .. if.. */ 90549b5e25fSSatish Balay } /* end of loop over added columns */ 90649b5e25fSSatish Balay ailen[brow] = nrow; 90749b5e25fSSatish Balay } /* end of loop over added rows */ 90849b5e25fSSatish Balay 90949b5e25fSSatish Balay PetscFunctionReturn(0); 91049b5e25fSSatish Balay } 91149b5e25fSSatish Balay 912915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*); 913915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal); 91449b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 91549b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 91649b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 91749b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 91849b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 91949b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat); 92049b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 92149b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 92249b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 92349b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 92449b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 92549b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 92613a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 92749b5e25fSSatish Balay 92849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 92949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 93049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 93149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 93249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 93349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 93449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 93549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 93649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 93749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 93849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 93949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 94049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 94149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 94249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 94349b5e25fSSatish Balay 94407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 94507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 94607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 94707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 94807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 94907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 95007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 95107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 95207f98182SSatish Balay 9535f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9545f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9555f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9565f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9575f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9585f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9595f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 96057d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 96149b5e25fSSatish Balay 96249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 96349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 96449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 96549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 96649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 96749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 96849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 96949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 97049b5e25fSSatish Balay 97149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 97249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 97349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 97449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 97549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 97649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 97749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 97849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 97949b5e25fSSatish Balay 9804d101231SSatish Balay #ifdef HAVE_MatICCFactor 9811a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 9824a2ae208SSatish Balay #undef __FUNCT__ 9834d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 9844d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level) 98549b5e25fSSatish Balay { 9864ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 98749b5e25fSSatish Balay Mat outA; 98849b5e25fSSatish Balay int ierr; 98949b5e25fSSatish Balay PetscTruth row_identity,col_identity; 99049b5e25fSSatish Balay 99149b5e25fSSatish Balay PetscFunctionBegin; 9921a3463dfSHong Zhang /* 99329bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 99449b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 99549b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 99649b5e25fSSatish Balay if (!row_identity || !col_identity) { 99729bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 99849b5e25fSSatish Balay } 9991a3463dfSHong Zhang */ 100049b5e25fSSatish Balay 100149b5e25fSSatish Balay outA = inA; 10021260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 100349b5e25fSSatish Balay 100449b5e25fSSatish Balay if (!a->diag) { 10051a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 100649b5e25fSSatish Balay } 100749b5e25fSSatish Balay /* 100849b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 100949b5e25fSSatish Balay for ILU(0) factorization with natural ordering 101049b5e25fSSatish Balay */ 101149b5e25fSSatish Balay switch (a->bs) { 101249b5e25fSSatish Balay case 1: 10130fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 101407f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 10154d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 101649b5e25fSSatish Balay case 2: 10171a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10181a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10191a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10204d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 102149b5e25fSSatish Balay break; 102249b5e25fSSatish Balay case 3: 10231a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10241a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10251a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10264d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 102749b5e25fSSatish Balay break; 102849b5e25fSSatish Balay case 4: 10291a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10301a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10311a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10324d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 103349b5e25fSSatish Balay break; 103449b5e25fSSatish Balay case 5: 10351a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10361a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10371a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10384d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 103949b5e25fSSatish Balay break; 104049b5e25fSSatish Balay case 6: 10411a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10421a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10431a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10444d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 104549b5e25fSSatish Balay break; 104649b5e25fSSatish Balay case 7: 10471a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10481a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10491a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10504d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 105149b5e25fSSatish Balay break; 105249b5e25fSSatish Balay default: 105349b5e25fSSatish Balay a->row = row; 10541a3463dfSHong Zhang a->icol = col; 105549b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 105649b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 105749b5e25fSSatish Balay 105849b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 105949b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1060b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 106149b5e25fSSatish Balay 106249b5e25fSSatish Balay if (!a->solve_work) { 106382502324SSatish Balay ierr = PetscMalloc((A->m+a->bs)*sizeof(Scalar),&a->solve_work);CHKERRQ(ierr); 1064b0a32e0cSBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(Scalar)); 106549b5e25fSSatish Balay } 106649b5e25fSSatish Balay } 106749b5e25fSSatish Balay 10681a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 106949b5e25fSSatish Balay 107049b5e25fSSatish Balay PetscFunctionReturn(0); 107149b5e25fSSatish Balay } 1072950f1e5bSHong Zhang #endif 1073950f1e5bSHong Zhang 10744a2ae208SSatish Balay #undef __FUNCT__ 10754a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 107649b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 107749b5e25fSSatish Balay { 107849b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 107949b5e25fSSatish Balay MPI_Comm comm = A->comm; 108049b5e25fSSatish Balay int ierr; 108149b5e25fSSatish Balay 108249b5e25fSSatish Balay PetscFunctionBegin; 108349b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 108449b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 108549b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 108649b5e25fSSatish Balay PetscFunctionReturn(0); 108749b5e25fSSatish Balay } 108849b5e25fSSatish Balay 108949b5e25fSSatish Balay EXTERN_C_BEGIN 10904a2ae208SSatish Balay #undef __FUNCT__ 10914a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 109249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 109349b5e25fSSatish Balay { 1094045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 109549b5e25fSSatish Balay int i,nz,n; 109649b5e25fSSatish Balay 109749b5e25fSSatish Balay PetscFunctionBegin; 1098045c9aa0SHong Zhang nz = baij->s_maxnz; 1099c464158bSHong Zhang n = mat->n; 110049b5e25fSSatish Balay for (i=0; i<nz; i++) { 110149b5e25fSSatish Balay baij->j[i] = indices[i]; 110249b5e25fSSatish Balay } 1103045c9aa0SHong Zhang baij->s_nz = nz; 110449b5e25fSSatish Balay for (i=0; i<n; i++) { 110549b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 110649b5e25fSSatish Balay } 110749b5e25fSSatish Balay 110849b5e25fSSatish Balay PetscFunctionReturn(0); 110949b5e25fSSatish Balay } 111049b5e25fSSatish Balay EXTERN_C_END 111149b5e25fSSatish Balay 11124a2ae208SSatish Balay #undef __FUNCT__ 11134a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 111449b5e25fSSatish Balay /*@ 111519585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 111649b5e25fSSatish Balay in the matrix. 111749b5e25fSSatish Balay 111849b5e25fSSatish Balay Input Parameters: 111919585528SSatish Balay + mat - the SeqSBAIJ matrix 112049b5e25fSSatish Balay - indices - the column indices 112149b5e25fSSatish Balay 112249b5e25fSSatish Balay Level: advanced 112349b5e25fSSatish Balay 112449b5e25fSSatish Balay Notes: 112549b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 112649b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 112749b5e25fSSatish Balay of the MatSetValues() operation. 112849b5e25fSSatish Balay 112949b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 113019585528SSatish Balay MatCreateSeqSBAIJ(). 113149b5e25fSSatish Balay 1132ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 113349b5e25fSSatish Balay 1134ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 113549b5e25fSSatish Balay @*/ 113649b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 113749b5e25fSSatish Balay { 113849b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 113949b5e25fSSatish Balay 114049b5e25fSSatish Balay PetscFunctionBegin; 114149b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1142b9617806SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)())&f);CHKERRQ(ierr); 114349b5e25fSSatish Balay if (f) { 114449b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 114549b5e25fSSatish Balay } else { 114629bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 114749b5e25fSSatish Balay } 114849b5e25fSSatish Balay PetscFunctionReturn(0); 114949b5e25fSSatish Balay } 115049b5e25fSSatish Balay 11514a2ae208SSatish Balay #undef __FUNCT__ 11524a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1153273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1154273d9f13SBarry Smith { 1155273d9f13SBarry Smith int ierr; 1156273d9f13SBarry Smith 1157273d9f13SBarry Smith PetscFunctionBegin; 1158273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1159273d9f13SBarry Smith PetscFunctionReturn(0); 1160273d9f13SBarry Smith } 1161273d9f13SBarry Smith 116249b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 116349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 116449b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 116549b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 116649b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 116749b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 116849b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 116949b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 117049b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 117149b5e25fSSatish Balay 0, 117249b5e25fSSatish Balay 0, 117349b5e25fSSatish Balay 0, 117449b5e25fSSatish Balay 0, 11755f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 117649b5e25fSSatish Balay 0, 117749b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 117849b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 117949b5e25fSSatish Balay MatEqual_SeqSBAIJ, 118049b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 118149b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 118249b5e25fSSatish Balay MatNorm_SeqSBAIJ, 118349b5e25fSSatish Balay 0, 118449b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 118549b5e25fSSatish Balay 0, 118649b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 118749b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 118849b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 118949b5e25fSSatish Balay 0, 119049b5e25fSSatish Balay 0, 11915f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 11925f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1193273d9f13SBarry Smith MatSetUpPreallocation_SeqSBAIJ, 1194c464158bSHong Zhang 0, 119549b5e25fSSatish Balay MatGetOwnershipRange_SeqSBAIJ, 119649b5e25fSSatish Balay 0, 11974d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 119849b5e25fSSatish Balay 0, 119949b5e25fSSatish Balay 0, 120049b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 120149b5e25fSSatish Balay 0, 120249b5e25fSSatish Balay 0, 120349b5e25fSSatish Balay 0, 1204950f1e5bSHong Zhang 0, 120549b5e25fSSatish Balay 0, 120649b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 120749b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 120849b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 120949b5e25fSSatish Balay 0, 121049b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 121149b5e25fSSatish Balay MatScale_SeqSBAIJ, 121249b5e25fSSatish Balay 0, 121349b5e25fSSatish Balay 0, 121449b5e25fSSatish Balay 0, 121549b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 121649b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 121749b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 121849b5e25fSSatish Balay 0, 121949b5e25fSSatish Balay 0, 122049b5e25fSSatish Balay 0, 122149b5e25fSSatish Balay 0, 122249b5e25fSSatish Balay 0, 122349b5e25fSSatish Balay 0, 122449b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 122549b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 122649b5e25fSSatish Balay 0, 122749b5e25fSSatish Balay 0, 1228*8a124369SBarry Smith MatGetPetscMaps_Petsc, 1229d959ec07SHong Zhang 0, 1230d959ec07SHong Zhang 0, 1231d959ec07SHong Zhang 0, 1232d959ec07SHong Zhang 0, 1233d959ec07SHong Zhang 0, 1234d959ec07SHong Zhang 0, 123524d5174aSHong Zhang MatGetRowMax_SeqSBAIJ}; 123649b5e25fSSatish Balay 123749b5e25fSSatish Balay EXTERN_C_BEGIN 12384a2ae208SSatish Balay #undef __FUNCT__ 12394a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 124049b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 124149b5e25fSSatish Balay { 12424afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1243c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 124449b5e25fSSatish Balay int ierr; 124549b5e25fSSatish Balay 124649b5e25fSSatish Balay PetscFunctionBegin; 124749b5e25fSSatish Balay if (aij->nonew != 1) { 124829bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 124949b5e25fSSatish Balay } 125049b5e25fSSatish Balay 125149b5e25fSSatish Balay /* allocate space for values if not already there */ 125249b5e25fSSatish Balay if (!aij->saved_values) { 1253eb7adc28SSatish Balay ierr = PetscMalloc((nz+1)*sizeof(Scalar),&aij->saved_values);CHKERRQ(ierr); 125449b5e25fSSatish Balay } 125549b5e25fSSatish Balay 125649b5e25fSSatish Balay /* copy values over */ 125749b5e25fSSatish Balay ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 125849b5e25fSSatish Balay PetscFunctionReturn(0); 125949b5e25fSSatish Balay } 126049b5e25fSSatish Balay EXTERN_C_END 126149b5e25fSSatish Balay 126249b5e25fSSatish Balay EXTERN_C_BEGIN 12634a2ae208SSatish Balay #undef __FUNCT__ 12644a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 126549b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 126649b5e25fSSatish Balay { 12674afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1268c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 126949b5e25fSSatish Balay 127049b5e25fSSatish Balay PetscFunctionBegin; 127149b5e25fSSatish Balay if (aij->nonew != 1) { 127229bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 127349b5e25fSSatish Balay } 127449b5e25fSSatish Balay if (!aij->saved_values) { 127529bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 127649b5e25fSSatish Balay } 127749b5e25fSSatish Balay 127849b5e25fSSatish Balay /* copy values over */ 127949b5e25fSSatish Balay ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 128049b5e25fSSatish Balay PetscFunctionReturn(0); 128149b5e25fSSatish Balay } 128249b5e25fSSatish Balay EXTERN_C_END 128349b5e25fSSatish Balay 12848549e402SHong Zhang EXTERN_C_BEGIN 12854a2ae208SSatish Balay #undef __FUNCT__ 12864a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ" 1287c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1288c464158bSHong Zhang { 1289c464158bSHong Zhang Mat_SeqSBAIJ *b; 12900c955e93SHong Zhang int ierr,size; 1291c464158bSHong Zhang 1292c464158bSHong Zhang PetscFunctionBegin; 1293c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1294c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1295273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1296273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1297c464158bSHong Zhang 1298b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1299b0a32e0cSBarry Smith B->data = (void*)b; 1300c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1301c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1302c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1303c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1304c464158bSHong Zhang B->factor = 0; 1305c464158bSHong Zhang B->lupivotthreshold = 1.0; 1306c464158bSHong Zhang B->mapping = 0; 1307c464158bSHong Zhang b->row = 0; 1308c464158bSHong Zhang b->icol = 0; 1309c464158bSHong Zhang b->reallocs = 0; 1310c464158bSHong Zhang b->saved_values = 0; 1311c464158bSHong Zhang 1312*8a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1313*8a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1314c464158bSHong Zhang 1315c464158bSHong Zhang b->sorted = PETSC_FALSE; 1316c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1317c464158bSHong Zhang b->nonew = 0; 1318c464158bSHong Zhang b->diag = 0; 1319c464158bSHong Zhang b->solve_work = 0; 1320c464158bSHong Zhang b->mult_work = 0; 1321c464158bSHong Zhang b->spptr = 0; 1322c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1323c464158bSHong Zhang 1324c464158bSHong Zhang b->inew = 0; 1325c464158bSHong Zhang b->jnew = 0; 1326c464158bSHong Zhang b->anew = 0; 1327c464158bSHong Zhang b->a2anew = 0; 1328c464158bSHong Zhang b->permute = PETSC_FALSE; 1329c464158bSHong Zhang 1330c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1331c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1332c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1333c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1334c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1335c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1336c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1337c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1338c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1339c464158bSHong Zhang PetscFunctionReturn(0); 1340c464158bSHong Zhang } 13418549e402SHong Zhang EXTERN_C_END 1342c464158bSHong Zhang 13434a2ae208SSatish Balay #undef __FUNCT__ 13444a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 134549b5e25fSSatish Balay /*@C 1346273d9f13SBarry Smith MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 134749b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 134849b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 134949b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 135049b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 135149b5e25fSSatish Balay 1352c464158bSHong Zhang Collective on Mat 135349b5e25fSSatish Balay 135449b5e25fSSatish Balay Input Parameters: 1355c464158bSHong Zhang + A - the symmetric matrix 135649b5e25fSSatish Balay . bs - size of block 135749b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 135849b5e25fSSatish Balay - nnz - array containing the number of block nonzeros in the various block rows 135949b5e25fSSatish Balay (possibly different for each block row) or PETSC_NULL 136049b5e25fSSatish Balay 136149b5e25fSSatish Balay Options Database Keys: 136249b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 136349b5e25fSSatish Balay block calculations (much slower) 136449b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 136549b5e25fSSatish Balay 136649b5e25fSSatish Balay Level: intermediate 136749b5e25fSSatish Balay 136849b5e25fSSatish Balay Notes: 1369c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 137049b5e25fSSatish Balay storage. That is, the stored row and column indices can begin at 137149b5e25fSSatish Balay either one (as in Fortran) or zero. See the users' manual for details. 137249b5e25fSSatish Balay 137349b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 137449b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 137549b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 137649b5e25fSSatish Balay matrices. 137749b5e25fSSatish Balay 1378a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 137949b5e25fSSatish Balay @*/ 1380273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 138149b5e25fSSatish Balay { 1382c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13830c955e93SHong Zhang int i,len,ierr,mbs,bs2; 138449b5e25fSSatish Balay PetscTruth flg; 13854afc71dfSHong Zhang int s_nz; 138649b5e25fSSatish Balay 138749b5e25fSSatish Balay PetscFunctionBegin; 1388273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1389b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1390c464158bSHong Zhang mbs = B->m/bs; 139149b5e25fSSatish Balay bs2 = bs*bs; 139249b5e25fSSatish Balay 1393c464158bSHong Zhang if (mbs*bs != B->m) { 139429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 139549b5e25fSSatish Balay } 139649b5e25fSSatish Balay 1397435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1398435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 139949b5e25fSSatish Balay if (nnz) { 140049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 140129bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 140280fe4e49SBarry 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); 140349b5e25fSSatish Balay } 140449b5e25fSSatish Balay } 140549b5e25fSSatish Balay 1406b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 140749b5e25fSSatish Balay if (!flg) { 140849b5e25fSSatish Balay switch (bs) { 140949b5e25fSSatish Balay case 1: 14104ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 141149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 141249b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 141349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 141449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 141549b5e25fSSatish Balay break; 141649b5e25fSSatish Balay case 2: 14174ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 141849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 141949b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 142049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 142149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 142249b5e25fSSatish Balay break; 142349b5e25fSSatish Balay case 3: 14245f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 142549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 142649b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 142749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 142849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 142949b5e25fSSatish Balay break; 143049b5e25fSSatish Balay case 4: 14315f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 143249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 143349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 143449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 143549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 143649b5e25fSSatish Balay break; 143749b5e25fSSatish Balay case 5: 14385f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 143949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 144049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 144149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 144249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 144349b5e25fSSatish Balay break; 144449b5e25fSSatish Balay case 6: 14455f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 144649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 144749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 144849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 144949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 145049b5e25fSSatish Balay break; 145149b5e25fSSatish Balay case 7: 1452de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 145349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 145449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1455de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 145649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 145749b5e25fSSatish Balay break; 145849b5e25fSSatish Balay } 145949b5e25fSSatish Balay } 146049b5e25fSSatish Balay 146149b5e25fSSatish Balay b->mbs = mbs; 14624afc71dfSHong Zhang b->nbs = mbs; 1463b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 146449b5e25fSSatish Balay if (!nnz) { 1465435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 146649b5e25fSSatish Balay else if (nz <= 0) nz = 1; 146749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14688cef66ccSHong Zhang b->imax[i] = nz; 146949b5e25fSSatish Balay } 1470153ea458SHong Zhang nz = nz*mbs; /* total nz */ 147149b5e25fSSatish Balay } else { 147249b5e25fSSatish Balay nz = 0; 14738cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 147449b5e25fSSatish Balay } 14758cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 14768cef66ccSHong Zhang s_nz = nz; 147749b5e25fSSatish Balay 147849b5e25fSSatish Balay /* allocate the matrix space */ 1479c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 148082502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 148149b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 148249b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 148349b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 148449b5e25fSSatish Balay b->i = b->j + s_nz; 148549b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 148649b5e25fSSatish Balay 148749b5e25fSSatish Balay /* pointer to beginning of each row */ 148849b5e25fSSatish Balay b->i[0] = 0; 148949b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 149049b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 149149b5e25fSSatish Balay } 149249b5e25fSSatish Balay 149349b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 14945033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1495b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 149649b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 149749b5e25fSSatish Balay 149849b5e25fSSatish Balay b->bs = bs; 149949b5e25fSSatish Balay b->bs2 = bs2; 150049b5e25fSSatish Balay b->s_nz = 0; 150149b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1502153ea458SHong Zhang 150316cdd363SHong Zhang b->inew = 0; 150416cdd363SHong Zhang b->jnew = 0; 150516cdd363SHong Zhang b->anew = 0; 150616cdd363SHong Zhang b->a2anew = 0; 15071a3463dfSHong Zhang b->permute = PETSC_FALSE; 1508c464158bSHong Zhang PetscFunctionReturn(0); 1509c464158bSHong Zhang } 1510153ea458SHong Zhang 151149b5e25fSSatish Balay 15124a2ae208SSatish Balay #undef __FUNCT__ 15134a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1514c464158bSHong Zhang /*@C 1515c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1516c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1517c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1518c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1519c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 152049b5e25fSSatish Balay 1521c464158bSHong Zhang Collective on MPI_Comm 1522c464158bSHong Zhang 1523c464158bSHong Zhang Input Parameters: 1524c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1525c464158bSHong Zhang . bs - size of block 1526c464158bSHong Zhang . m - number of rows, or number of columns 1527c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1528c464158bSHong Zhang - nnz - array containing the number of block nonzeros in the various block rows 1529c464158bSHong Zhang (possibly different for each block row) or PETSC_NULL 1530c464158bSHong Zhang 1531c464158bSHong Zhang Output Parameter: 1532c464158bSHong Zhang . A - the symmetric matrix 1533c464158bSHong Zhang 1534c464158bSHong Zhang Options Database Keys: 1535c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1536c464158bSHong Zhang block calculations (much slower) 1537c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1538c464158bSHong Zhang 1539c464158bSHong Zhang Level: intermediate 1540c464158bSHong Zhang 1541c464158bSHong Zhang Notes: 1542c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 1543c464158bSHong Zhang storage. That is, the stored row and column indices can begin at 1544c464158bSHong Zhang either one (as in Fortran) or zero. See the users' manual for details. 1545c464158bSHong Zhang 1546c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1547c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1548c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1549c464158bSHong Zhang matrices. 1550c464158bSHong Zhang 1551c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1552c464158bSHong Zhang @*/ 1553c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1554c464158bSHong Zhang { 1555c464158bSHong Zhang int ierr; 1556c464158bSHong Zhang 1557c464158bSHong Zhang PetscFunctionBegin; 1558c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1559c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1560273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 156149b5e25fSSatish Balay PetscFunctionReturn(0); 156249b5e25fSSatish Balay } 156349b5e25fSSatish Balay 15644a2ae208SSatish Balay #undef __FUNCT__ 15654a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 156649b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 156749b5e25fSSatish Balay { 156849b5e25fSSatish Balay Mat C; 156949b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 157049b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 157149b5e25fSSatish Balay 157249b5e25fSSatish Balay PetscFunctionBegin; 157329bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 157449b5e25fSSatish Balay 157549b5e25fSSatish Balay *B = 0; 1576692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1577692f9cbeSHong Zhang ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr); 1578692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1579692f9cbeSHong Zhang 158049b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1581273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 158249b5e25fSSatish Balay C->factor = A->factor; 158349b5e25fSSatish Balay c->row = 0; 158449b5e25fSSatish Balay c->icol = 0; 158549b5e25fSSatish Balay c->saved_values = 0; 158649b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 158749b5e25fSSatish Balay C->assembled = PETSC_TRUE; 158849b5e25fSSatish Balay 158949b5e25fSSatish Balay c->bs = a->bs; 159049b5e25fSSatish Balay c->bs2 = a->bs2; 159149b5e25fSSatish Balay c->mbs = a->mbs; 159249b5e25fSSatish Balay c->nbs = a->nbs; 159349b5e25fSSatish Balay 1594b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1595b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 159649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 159749b5e25fSSatish Balay c->imax[i] = a->imax[i]; 159849b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 159949b5e25fSSatish Balay } 160049b5e25fSSatish Balay 160149b5e25fSSatish Balay /* allocate the matrix space */ 160249b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 160349b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 160482502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 160549b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 160649b5e25fSSatish Balay c->i = c->j + nz; 160749b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 160849b5e25fSSatish Balay if (mbs > 0) { 160949b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 161049b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 161149b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 161249b5e25fSSatish Balay } else { 161349b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 161449b5e25fSSatish Balay } 161549b5e25fSSatish Balay } 161649b5e25fSSatish Balay 1617b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 161849b5e25fSSatish Balay c->sorted = a->sorted; 161949b5e25fSSatish Balay c->roworiented = a->roworiented; 162049b5e25fSSatish Balay c->nonew = a->nonew; 162149b5e25fSSatish Balay 162249b5e25fSSatish Balay if (a->diag) { 1623b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1624b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 162549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 162649b5e25fSSatish Balay c->diag[i] = a->diag[i]; 162749b5e25fSSatish Balay } 162849b5e25fSSatish Balay } else c->diag = 0; 162949b5e25fSSatish Balay c->s_nz = a->s_nz; 163049b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 163149b5e25fSSatish Balay c->solve_work = 0; 163249b5e25fSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 163349b5e25fSSatish Balay c->mult_work = 0; 163449b5e25fSSatish Balay *B = C; 1635b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 163649b5e25fSSatish Balay PetscFunctionReturn(0); 163749b5e25fSSatish Balay } 163849b5e25fSSatish Balay 16398549e402SHong Zhang EXTERN_C_BEGIN 16404a2ae208SSatish Balay #undef __FUNCT__ 16414a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 1642b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A) 164349b5e25fSSatish Balay { 164449b5e25fSSatish Balay Mat_SeqSBAIJ *a; 164549b5e25fSSatish Balay Mat B; 164649b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1647574b2666SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount; 164849b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 164949b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 165049b5e25fSSatish Balay Scalar *aa; 165149b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 165249b5e25fSSatish Balay 165349b5e25fSSatish Balay PetscFunctionBegin; 1654b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 165549b5e25fSSatish Balay bs2 = bs*bs; 165649b5e25fSSatish Balay 165749b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 165829bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1659b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 166049b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 166129bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 166249b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 166349b5e25fSSatish Balay 166449b5e25fSSatish Balay if (header[3] < 0) { 166529bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 166649b5e25fSSatish Balay } 166749b5e25fSSatish Balay 166829bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 166949b5e25fSSatish Balay 167049b5e25fSSatish Balay /* 167149b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 167249b5e25fSSatish Balay divisible by the blocksize 167349b5e25fSSatish Balay */ 167449b5e25fSSatish Balay mbs = M/bs; 167549b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 167649b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 167749b5e25fSSatish Balay else mbs++; 167849b5e25fSSatish Balay if (extra_rows) { 1679b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 168049b5e25fSSatish Balay } 168149b5e25fSSatish Balay 168249b5e25fSSatish Balay /* read in row lengths */ 1683b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 168449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 168549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 168649b5e25fSSatish Balay 168749b5e25fSSatish Balay /* read in column indices */ 1688b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 168949b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 169049b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 169149b5e25fSSatish Balay 169249b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 169382502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 169482502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 1695574b2666SHong Zhang ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 169682502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 169749b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 169849b5e25fSSatish Balay masked = mask + mbs; 169949b5e25fSSatish Balay rowcount = 0; nzcount = 0; 170049b5e25fSSatish Balay for (i=0; i<mbs; i++) { 170149b5e25fSSatish Balay nmask = 0; 170249b5e25fSSatish Balay for (j=0; j<bs; j++) { 170349b5e25fSSatish Balay kmax = rowlengths[rowcount]; 170449b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17052d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 170603630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 170749b5e25fSSatish Balay } 170849b5e25fSSatish Balay rowcount++; 170949b5e25fSSatish Balay } 1710574b2666SHong Zhang s_browlengths[i] += nmask; 1711574b2666SHong Zhang browlengths[i] = 2*s_browlengths[i]; 1712574b2666SHong Zhang 171349b5e25fSSatish Balay /* zero out the mask elements we set */ 171449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 171549b5e25fSSatish Balay } 171649b5e25fSSatish Balay 171749b5e25fSSatish Balay /* create our matrix */ 171849b5e25fSSatish Balay ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 171949b5e25fSSatish Balay B = *A; 172049b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 172149b5e25fSSatish Balay 172249b5e25fSSatish Balay /* set matrix "i" values */ 172349b5e25fSSatish Balay a->i[0] = 0; 172449b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1725574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1726574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 172749b5e25fSSatish Balay } 172849b5e25fSSatish Balay a->s_nz = 0; 1729574b2666SHong Zhang for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i]; 173049b5e25fSSatish Balay 173149b5e25fSSatish Balay /* read in nonzero values */ 1732b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(Scalar),&aa);CHKERRQ(ierr); 173349b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 173449b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 173549b5e25fSSatish Balay 173649b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 173749b5e25fSSatish Balay nzcount = 0; jcount = 0; 173849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 173949b5e25fSSatish Balay nzcountb = nzcount; 174049b5e25fSSatish Balay nmask = 0; 174149b5e25fSSatish Balay for (j=0; j<bs; j++) { 174249b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 174349b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17442d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 174503630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 17462d703238SHong Zhang } 17472d703238SHong Zhang } 17482d703238SHong Zhang /* sort the masked values */ 17492d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 17502d703238SHong Zhang 17512d703238SHong Zhang /* set "j" values into matrix */ 17522d703238SHong Zhang maskcount = 1; 17532d703238SHong Zhang for (j=0; j<nmask; j++) { 175449b5e25fSSatish Balay a->j[jcount++] = masked[j]; 175549b5e25fSSatish Balay mask[masked[j]] = maskcount++; 175649b5e25fSSatish Balay } 1757574b2666SHong Zhang 175849b5e25fSSatish Balay /* set "a" values into matrix */ 175949b5e25fSSatish Balay ishift = bs2*a->i[i]; 176049b5e25fSSatish Balay for (j=0; j<bs; j++) { 176149b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 176249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1763574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1764574b2666SHong Zhang if (tmp >= i){ 176549b5e25fSSatish Balay block = mask[tmp] - 1; 176649b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 176749b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1768574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1769574b2666SHong Zhang } 1770574b2666SHong Zhang nzcountb++; 177149b5e25fSSatish Balay } 177249b5e25fSSatish Balay } 177349b5e25fSSatish Balay /* zero out the mask elements we set */ 177449b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 177549b5e25fSSatish Balay } 177629bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 177749b5e25fSSatish Balay 177849b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 177949b5e25fSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1780574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 178149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 178249b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 178349b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 178449b5e25fSSatish Balay 178549b5e25fSSatish Balay B->assembled = PETSC_TRUE; 178649b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 178749b5e25fSSatish Balay PetscFunctionReturn(0); 178849b5e25fSSatish Balay } 17898549e402SHong Zhang EXTERN_C_END 1790574b2666SHong Zhang 179149b5e25fSSatish Balay 179249b5e25fSSatish Balay 1793