173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 4a1373b80SHong Zhang Defines the basic matrix operations for the SBAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/inline/spops.h" 9aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1049b5e25fSSatish Balay 1149b5e25fSSatish Balay #define CHUNKSIZE 10 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay /* 1449b5e25fSSatish Balay Checks for missing diagonals 1549b5e25fSSatish Balay */ 164a2ae208SSatish Balay #undef __FUNCT__ 174a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ" 1849b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 1949b5e25fSSatish Balay { 20045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2149b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2249b5e25fSSatish Balay 2349b5e25fSSatish Balay PetscFunctionBegin; 24045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2549b5e25fSSatish Balay diag = a->diag; 2649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 27e11b0e14SHong Zhang if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i); 2849b5e25fSSatish Balay } 2949b5e25fSSatish Balay PetscFunctionReturn(0); 3049b5e25fSSatish Balay } 3149b5e25fSSatish Balay 324a2ae208SSatish Balay #undef __FUNCT__ 334a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ" 3449b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3549b5e25fSSatish Balay { 36045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 37b424e231SHong Zhang int i,mbs = a->mbs,ierr; 3849b5e25fSSatish Balay 3949b5e25fSSatish Balay PetscFunctionBegin; 4049b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4149b5e25fSSatish Balay 42b424e231SHong Zhang ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr); 43b424e231SHong Zhang PetscLogObjectMemory(A,(mbs+1)*sizeof(int)); 44b424e231SHong Zhang for (i=0; i<mbs; i++) a->diag[i] = a->i[i]; 4549b5e25fSSatish Balay PetscFunctionReturn(0); 4649b5e25fSSatish Balay } 4749b5e25fSSatish Balay 484a2ae208SSatish Balay #undef __FUNCT__ 494a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ" 504e7234bfSBarry Smith static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done) 5149b5e25fSSatish Balay { 52a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 539e37e433SSatish Balay int n = a->mbs,i; 5449b5e25fSSatish Balay 5549b5e25fSSatish Balay PetscFunctionBegin; 56d3e5a4abSHong Zhang *nn = n; 57a1373b80SHong Zhang if (!ia) PetscFunctionReturn(0); 58a1373b80SHong Zhang 59a6ece127SHong Zhang if (oshift == 1) { 60a6ece127SHong Zhang /* temporarily add 1 to i and j indices */ 616c6c5352SBarry Smith int nz = a->i[n]; 626c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 63a1373b80SHong Zhang for (i=0; i<n+1; i++) a->i[i]++; 64a1373b80SHong Zhang *ia = a->i; *ja = a->j; 65a1373b80SHong Zhang } else { 66a1373b80SHong Zhang *ia = a->i; *ja = a->j; 67a6ece127SHong Zhang } 6849b5e25fSSatish Balay PetscFunctionReturn(0); 6949b5e25fSSatish Balay } 7049b5e25fSSatish Balay 714a2ae208SSatish Balay #undef __FUNCT__ 724a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ" 734e7234bfSBarry Smith static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done) 7449b5e25fSSatish Balay { 75b7aaefc3SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 76a15ac5e4SHong Zhang int i,n = a->mbs; 77a6ece127SHong Zhang 7849b5e25fSSatish Balay PetscFunctionBegin; 7949b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 80a6ece127SHong Zhang 81a6ece127SHong Zhang if (oshift == 1) { 826c6c5352SBarry Smith int nz = a->i[n]-1; 836c6c5352SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 84a6ece127SHong Zhang for (i=0; i<n+1; i++) a->i[i]--; 85a6ece127SHong Zhang } 86a6ece127SHong Zhang PetscFunctionReturn(0); 8749b5e25fSSatish Balay } 8849b5e25fSSatish Balay 894a2ae208SSatish Balay #undef __FUNCT__ 904a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ" 9149b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9249b5e25fSSatish Balay { 9349b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 9449b5e25fSSatish Balay 9549b5e25fSSatish Balay PetscFunctionBegin; 9649b5e25fSSatish Balay *bs = sbaij->bs; 9749b5e25fSSatish Balay PetscFunctionReturn(0); 9849b5e25fSSatish Balay } 9949b5e25fSSatish Balay 1004a2ae208SSatish Balay #undef __FUNCT__ 1014a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ" 10249b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 10349b5e25fSSatish Balay { 10449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10549b5e25fSSatish Balay int ierr; 10649b5e25fSSatish Balay 10749b5e25fSSatish Balay PetscFunctionBegin; 10849b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 1096c6c5352SBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, NZ=%d",A->m,a->nz); 11049b5e25fSSatish Balay #endif 11149b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 11249b5e25fSSatish Balay if (!a->singlemalloc) { 11349b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 11463c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 11549b5e25fSSatish Balay } 11649b5e25fSSatish Balay if (a->row) { 11749b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 11849b5e25fSSatish Balay } 11949b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 12049b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 12149b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 12249b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 123d59c15a7SBarry Smith if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 124d59c15a7SBarry Smith if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);} 125d59c15a7SBarry Smith if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 12649b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 127c4319e64SHong Zhang if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 1281a3463dfSHong Zhang 1291a3463dfSHong Zhang if (a->inew){ 1301a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1311a3463dfSHong Zhang a->inew = 0; 1321a3463dfSHong Zhang } 13349b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 13449b5e25fSSatish Balay PetscFunctionReturn(0); 13549b5e25fSSatish Balay } 13649b5e25fSSatish Balay 1374a2ae208SSatish Balay #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ" 13949b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 14049b5e25fSSatish Balay { 141045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 14249b5e25fSSatish Balay 14349b5e25fSSatish Balay PetscFunctionBegin; 1444d9d31abSKris Buschelman switch (op) { 1454d9d31abSKris Buschelman case MAT_ROW_ORIENTED: 1464d9d31abSKris Buschelman a->roworiented = PETSC_TRUE; 1474d9d31abSKris Buschelman break; 1484d9d31abSKris Buschelman case MAT_COLUMN_ORIENTED: 1494d9d31abSKris Buschelman a->roworiented = PETSC_FALSE; 1504d9d31abSKris Buschelman break; 1514d9d31abSKris Buschelman case MAT_COLUMNS_SORTED: 1524d9d31abSKris Buschelman a->sorted = PETSC_TRUE; 1534d9d31abSKris Buschelman break; 1544d9d31abSKris Buschelman case MAT_COLUMNS_UNSORTED: 1554d9d31abSKris Buschelman a->sorted = PETSC_FALSE; 1564d9d31abSKris Buschelman break; 1574d9d31abSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 1584d9d31abSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 1594d9d31abSKris Buschelman break; 1604d9d31abSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 1614d9d31abSKris Buschelman a->nonew = 1; 1624d9d31abSKris Buschelman break; 1634d9d31abSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1644d9d31abSKris Buschelman a->nonew = -1; 1654d9d31abSKris Buschelman break; 1664d9d31abSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 1674d9d31abSKris Buschelman a->nonew = -2; 1684d9d31abSKris Buschelman break; 1694d9d31abSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 1704d9d31abSKris Buschelman a->nonew = 0; 1714d9d31abSKris Buschelman break; 1724d9d31abSKris Buschelman case MAT_ROWS_SORTED: 1734d9d31abSKris Buschelman case MAT_ROWS_UNSORTED: 1744d9d31abSKris Buschelman case MAT_YES_NEW_DIAGONALS: 1754d9d31abSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1764d9d31abSKris Buschelman case MAT_USE_HASH_TABLE: 177b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 1784d9d31abSKris Buschelman break; 1794d9d31abSKris Buschelman case MAT_NO_NEW_DIAGONALS: 18029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1819a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 1829a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1839a4540c5SBarry Smith case MAT_HERMITIAN: 1849a4540c5SBarry Smith SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric"); 18577e54ba9SKris Buschelman case MAT_SYMMETRIC: 18677e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 1879a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 1889a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1899a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 19077e54ba9SKris Buschelman break; 1914d9d31abSKris Buschelman default: 19229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 19349b5e25fSSatish Balay } 19449b5e25fSSatish Balay PetscFunctionReturn(0); 19549b5e25fSSatish Balay } 19649b5e25fSSatish Balay 1974a2ae208SSatish Balay #undef __FUNCT__ 1984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ" 19987828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **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; 20487828ca2SBarry Smith PetscScalar *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 213a45adfd6SMatthew Knepley if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %d out of range", row); 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) { 22387828ca2SBarry Smith ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),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" 27187828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **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 { 2858115998fSBarry Smith int ierr; 28649b5e25fSSatish Balay PetscFunctionBegin; 287999d9058SBarry Smith ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr); 2888115998fSBarry Smith PetscFunctionReturn(0); 28949b5e25fSSatish Balay } 29049b5e25fSSatish Balay 2914a2ae208SSatish Balay #undef __FUNCT__ 2924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII" 293b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer) 29449b5e25fSSatish Balay { 29549b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 296fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 297fb9695e5SSatish Balay char *name; 298f3ef73ceSBarry Smith PetscViewerFormat format; 29949b5e25fSSatish Balay 30049b5e25fSSatish Balay PetscFunctionBegin; 30180fe4e49SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 302b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 303456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 304b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 305fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 30629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 307fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 308b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 30949b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31049b5e25fSSatish Balay for (j=0; j<bs; j++) { 311b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 31249b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31349b5e25fSSatish Balay for (l=0; l<bs; l++) { 31449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 31549b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 31749b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 31849b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 31962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 32049b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32149b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 32262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 32349b5e25fSSatish Balay } 32449b5e25fSSatish Balay #else 32549b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 32662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 32749b5e25fSSatish Balay } 32849b5e25fSSatish Balay #endif 32949b5e25fSSatish Balay } 33049b5e25fSSatish Balay } 331b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 33249b5e25fSSatish Balay } 33349b5e25fSSatish Balay } 334b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 33549b5e25fSSatish Balay } else { 336b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 33749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33849b5e25fSSatish Balay for (j=0; j<bs; j++) { 339b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 34049b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 34149b5e25fSSatish Balay for (l=0; l<bs; l++) { 34249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 34349b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 34462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 34549b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34649b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 34762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 34849b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 34949b5e25fSSatish Balay } else { 35062b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 35149b5e25fSSatish Balay } 35249b5e25fSSatish Balay #else 353b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 35449b5e25fSSatish Balay #endif 35549b5e25fSSatish Balay } 35649b5e25fSSatish Balay } 357b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 35849b5e25fSSatish Balay } 35949b5e25fSSatish Balay } 360b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 36149b5e25fSSatish Balay } 362b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 36349b5e25fSSatish Balay PetscFunctionReturn(0); 36449b5e25fSSatish Balay } 36549b5e25fSSatish Balay 3664a2ae208SSatish Balay #undef __FUNCT__ 3674a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom" 368b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 36949b5e25fSSatish Balay { 37049b5e25fSSatish Balay Mat A = (Mat) Aa; 37149b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 37249b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 37349b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 37449b5e25fSSatish Balay MatScalar *aa; 37549b5e25fSSatish Balay MPI_Comm comm; 376b0a32e0cSBarry Smith PetscViewer viewer; 37749b5e25fSSatish Balay 37849b5e25fSSatish Balay PetscFunctionBegin; 37949b5e25fSSatish Balay /* 38049b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 38149b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 38249b5e25fSSatish Balay rest should return immediately. 38349b5e25fSSatish Balay */ 38449b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 38549b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 38649b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 38749b5e25fSSatish Balay 38849b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 38949b5e25fSSatish Balay 390b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 391b0a32e0cSBarry Smith PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric"); 39249b5e25fSSatish Balay 39349b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 394b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 39549b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 39649b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 397c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 39849b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 39949b5e25fSSatish Balay aa = a->a + j*bs2; 40049b5e25fSSatish Balay for (k=0; k<bs; k++) { 40149b5e25fSSatish Balay for (l=0; l<bs; l++) { 40249b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 403b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 40449b5e25fSSatish Balay } 40549b5e25fSSatish Balay } 40649b5e25fSSatish Balay } 40749b5e25fSSatish Balay } 408b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 40949b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 41049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 411c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 41249b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 41349b5e25fSSatish Balay aa = a->a + j*bs2; 41449b5e25fSSatish Balay for (k=0; k<bs; k++) { 41549b5e25fSSatish Balay for (l=0; l<bs; l++) { 41649b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 417b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 41849b5e25fSSatish Balay } 41949b5e25fSSatish Balay } 42049b5e25fSSatish Balay } 42149b5e25fSSatish Balay } 42249b5e25fSSatish Balay 423b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 42449b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 42549b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 426c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 42749b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 42849b5e25fSSatish Balay aa = a->a + j*bs2; 42949b5e25fSSatish Balay for (k=0; k<bs; k++) { 43049b5e25fSSatish Balay for (l=0; l<bs; l++) { 43149b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 432b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 43349b5e25fSSatish Balay } 43449b5e25fSSatish Balay } 43549b5e25fSSatish Balay } 43649b5e25fSSatish Balay } 43749b5e25fSSatish Balay PetscFunctionReturn(0); 43849b5e25fSSatish Balay } 43949b5e25fSSatish Balay 4404a2ae208SSatish Balay #undef __FUNCT__ 4414a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw" 442b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer) 44349b5e25fSSatish Balay { 44449b5e25fSSatish Balay int ierr; 44549b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 446b0a32e0cSBarry Smith PetscDraw draw; 44749b5e25fSSatish Balay PetscTruth isnull; 44849b5e25fSSatish Balay 44949b5e25fSSatish Balay PetscFunctionBegin; 45049b5e25fSSatish Balay 451b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 452b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 45349b5e25fSSatish Balay 45449b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 455c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 45649b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 457b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 458b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 45949b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 46049b5e25fSSatish Balay PetscFunctionReturn(0); 46149b5e25fSSatish Balay } 46249b5e25fSSatish Balay 4634a2ae208SSatish Balay #undef __FUNCT__ 4644a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ" 465b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer) 46649b5e25fSSatish Balay { 46749b5e25fSSatish Balay int ierr; 468a5e6ed63SBarry Smith PetscTruth isascii,isdraw; 46949b5e25fSSatish Balay 47049b5e25fSSatish Balay PetscFunctionBegin; 471b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 472fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 473a5e6ed63SBarry Smith if (isascii){ 47449b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 47549b5e25fSSatish Balay } else if (isdraw) { 47649b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 47749b5e25fSSatish Balay } else { 478a5e6ed63SBarry Smith Mat B; 479a5e6ed63SBarry Smith ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr); 480a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 481a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 48249b5e25fSSatish Balay } 48349b5e25fSSatish Balay PetscFunctionReturn(0); 48449b5e25fSSatish Balay } 48549b5e25fSSatish Balay 48649b5e25fSSatish Balay 4874a2ae208SSatish Balay #undef __FUNCT__ 4884a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ" 489f15d580aSBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[]) 49049b5e25fSSatish Balay { 491045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 49249b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 49349b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 49449b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 49549b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 49649b5e25fSSatish Balay 49749b5e25fSSatish Balay PetscFunctionBegin; 49849b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 49949b5e25fSSatish Balay row = im[k]; brow = row/bs; 500590ac198SBarry Smith if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row); 501590ac198SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1); 50249b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 50349b5e25fSSatish Balay nrow = ailen[brow]; 50449b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 505590ac198SBarry Smith if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]); 506590ac198SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1); 50749b5e25fSSatish Balay col = in[l] ; 50849b5e25fSSatish Balay bcol = col/bs; 50949b5e25fSSatish Balay cidx = col%bs; 51049b5e25fSSatish Balay ridx = row%bs; 51149b5e25fSSatish Balay high = nrow; 51249b5e25fSSatish Balay low = 0; /* assume unsorted */ 51349b5e25fSSatish Balay while (high-low > 5) { 51449b5e25fSSatish Balay t = (low+high)/2; 51549b5e25fSSatish Balay if (rp[t] > bcol) high = t; 51649b5e25fSSatish Balay else low = t; 51749b5e25fSSatish Balay } 51849b5e25fSSatish Balay for (i=low; i<high; i++) { 51949b5e25fSSatish Balay if (rp[i] > bcol) break; 52049b5e25fSSatish Balay if (rp[i] == bcol) { 52149b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 52249b5e25fSSatish Balay goto finished; 52349b5e25fSSatish Balay } 52449b5e25fSSatish Balay } 52549b5e25fSSatish Balay *v++ = zero; 52649b5e25fSSatish Balay finished:; 52749b5e25fSSatish Balay } 52849b5e25fSSatish Balay } 52949b5e25fSSatish Balay PetscFunctionReturn(0); 53049b5e25fSSatish Balay } 53149b5e25fSSatish Balay 53249b5e25fSSatish Balay 5334a2ae208SSatish Balay #undef __FUNCT__ 5344a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ" 535f15d580aSBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is) 53649b5e25fSSatish Balay { 5370880e062SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 5380880e062SHong Zhang int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 5390880e062SHong Zhang int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 5400880e062SHong Zhang int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 5410880e062SHong Zhang PetscTruth roworiented=a->roworiented; 542f15d580aSBarry Smith const MatScalar *value = v; 543f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 5440880e062SHong Zhang 54549b5e25fSSatish Balay PetscFunctionBegin; 5460880e062SHong Zhang if (roworiented) { 5470880e062SHong Zhang stepval = (n-1)*bs; 5480880e062SHong Zhang } else { 5490880e062SHong Zhang stepval = (m-1)*bs; 5500880e062SHong Zhang } 5510880e062SHong Zhang for (k=0; k<m; k++) { /* loop over added rows */ 5520880e062SHong Zhang row = im[k]; 5530880e062SHong Zhang if (row < 0) continue; 5540880e062SHong Zhang #if defined(PETSC_USE_BOPT_g) 555590ac198SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1); 5560880e062SHong Zhang #endif 5570880e062SHong Zhang rp = aj + ai[row]; 5580880e062SHong Zhang ap = aa + bs2*ai[row]; 5590880e062SHong Zhang rmax = imax[row]; 5600880e062SHong Zhang nrow = ailen[row]; 5610880e062SHong Zhang low = 0; 5620880e062SHong Zhang for (l=0; l<n; l++) { /* loop over added columns */ 5630880e062SHong Zhang if (in[l] < 0) continue; 5640880e062SHong Zhang col = in[l]; 565b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g) 566b1823623SSatish Balay if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",col,a->nbs-1); 567b1823623SSatish Balay #endif 5680880e062SHong Zhang if (col < row) continue; /* ignore lower triangular block */ 5690880e062SHong Zhang if (roworiented) { 5700880e062SHong Zhang value = v + k*(stepval+bs)*bs + l*bs; 5710880e062SHong Zhang } else { 5720880e062SHong Zhang value = v + l*(stepval+bs)*bs + k*bs; 5730880e062SHong Zhang } 5740880e062SHong Zhang if (!sorted) low = 0; high = nrow; 5750880e062SHong Zhang while (high-low > 7) { 5760880e062SHong Zhang t = (low+high)/2; 5770880e062SHong Zhang if (rp[t] > col) high = t; 5780880e062SHong Zhang else low = t; 5790880e062SHong Zhang } 5800880e062SHong Zhang for (i=low; i<high; i++) { 5810880e062SHong Zhang if (rp[i] > col) break; 5820880e062SHong Zhang if (rp[i] == col) { 5830880e062SHong Zhang bap = ap + bs2*i; 5840880e062SHong Zhang if (roworiented) { 5850880e062SHong Zhang if (is == ADD_VALUES) { 5860880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5870880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5880880e062SHong Zhang bap[jj] += *value++; 5890880e062SHong Zhang } 5900880e062SHong Zhang } 5910880e062SHong Zhang } else { 5920880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 5930880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 5940880e062SHong Zhang bap[jj] = *value++; 5950880e062SHong Zhang } 5960880e062SHong Zhang } 5970880e062SHong Zhang } 5980880e062SHong Zhang } else { 5990880e062SHong Zhang if (is == ADD_VALUES) { 6000880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6010880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6020880e062SHong Zhang *bap++ += *value++; 6030880e062SHong Zhang } 6040880e062SHong Zhang } 6050880e062SHong Zhang } else { 6060880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6070880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6080880e062SHong Zhang *bap++ = *value++; 6090880e062SHong Zhang } 6100880e062SHong Zhang } 6110880e062SHong Zhang } 6120880e062SHong Zhang } 6130880e062SHong Zhang goto noinsert2; 6140880e062SHong Zhang } 6150880e062SHong Zhang } 6160880e062SHong Zhang if (nonew == 1) goto noinsert2; 617a45adfd6SMatthew Knepley else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 6180880e062SHong Zhang if (nrow >= rmax) { 6190880e062SHong Zhang /* there is no extra room in row, therefore enlarge */ 6200880e062SHong Zhang int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 6210880e062SHong Zhang MatScalar *new_a; 6220880e062SHong Zhang 623a45adfd6SMatthew Knepley if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 6240880e062SHong Zhang 6250880e062SHong Zhang /* malloc new storage space */ 6260880e062SHong Zhang len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 6270880e062SHong Zhang ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 6280880e062SHong Zhang new_j = (int*)(new_a + bs2*new_nz); 6290880e062SHong Zhang new_i = new_j + new_nz; 6300880e062SHong Zhang 6310880e062SHong Zhang /* copy over old data into new slots */ 6320880e062SHong Zhang for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 6330880e062SHong Zhang for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 6340880e062SHong Zhang ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 6350880e062SHong Zhang len = (new_nz - CHUNKSIZE - ai[row] - nrow); 6360880e062SHong Zhang ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 6370880e062SHong Zhang ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 6380880e062SHong Zhang ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 6390880e062SHong Zhang ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 6400880e062SHong Zhang /* free up old matrix storage */ 6410880e062SHong Zhang ierr = PetscFree(a->a);CHKERRQ(ierr); 6420880e062SHong Zhang if (!a->singlemalloc) { 6430880e062SHong Zhang ierr = PetscFree(a->i);CHKERRQ(ierr); 6440880e062SHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 6450880e062SHong Zhang } 6460880e062SHong Zhang aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 6470880e062SHong Zhang a->singlemalloc = PETSC_TRUE; 6480880e062SHong Zhang 6490880e062SHong Zhang rp = aj + ai[row]; ap = aa + bs2*ai[row]; 6500880e062SHong Zhang rmax = imax[row] = imax[row] + CHUNKSIZE; 6510880e062SHong Zhang PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 6526c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 6530880e062SHong Zhang a->reallocs++; 6546c6c5352SBarry Smith a->nz++; 6550880e062SHong Zhang } 6560880e062SHong Zhang N = nrow++ - 1; 6570880e062SHong Zhang /* shift up all the later entries in this row */ 6580880e062SHong Zhang for (ii=N; ii>=i; ii--) { 6590880e062SHong Zhang rp[ii+1] = rp[ii]; 6600880e062SHong Zhang ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 6610880e062SHong Zhang } 6620880e062SHong Zhang if (N >= i) { 6630880e062SHong Zhang ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 6640880e062SHong Zhang } 6650880e062SHong Zhang rp[i] = col; 6660880e062SHong Zhang bap = ap + bs2*i; 6670880e062SHong Zhang if (roworiented) { 6680880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6690880e062SHong Zhang for (jj=ii; jj<bs2; jj+=bs) { 6700880e062SHong Zhang bap[jj] = *value++; 6710880e062SHong Zhang } 6720880e062SHong Zhang } 6730880e062SHong Zhang } else { 6740880e062SHong Zhang for (ii=0; ii<bs; ii++,value+=stepval) { 6750880e062SHong Zhang for (jj=0; jj<bs; jj++) { 6760880e062SHong Zhang *bap++ = *value++; 6770880e062SHong Zhang } 6780880e062SHong Zhang } 6790880e062SHong Zhang } 6800880e062SHong Zhang noinsert2:; 6810880e062SHong Zhang low = i; 6820880e062SHong Zhang } 6830880e062SHong Zhang ailen[row] = nrow; 6840880e062SHong Zhang } 6850880e062SHong Zhang PetscFunctionReturn(0); 68649b5e25fSSatish Balay } 68749b5e25fSSatish Balay 6884a2ae208SSatish Balay #undef __FUNCT__ 6894a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ" 69049b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 69149b5e25fSSatish Balay { 69249b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 69349b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 694c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 69549b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 69649b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 69749b5e25fSSatish Balay 69849b5e25fSSatish Balay PetscFunctionBegin; 69949b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 70049b5e25fSSatish Balay 70149b5e25fSSatish Balay if (m) rmax = ailen[0]; 70249b5e25fSSatish Balay for (i=1; i<mbs; i++) { 70349b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 70449b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 70549b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 70649b5e25fSSatish Balay if (fshift) { 70749b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 70849b5e25fSSatish Balay N = ailen[i]; 70949b5e25fSSatish Balay for (j=0; j<N; j++) { 71049b5e25fSSatish Balay ip[j-fshift] = ip[j]; 71149b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 71249b5e25fSSatish Balay } 71349b5e25fSSatish Balay } 71449b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 71549b5e25fSSatish Balay } 71649b5e25fSSatish Balay if (mbs) { 71749b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 71849b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 71949b5e25fSSatish Balay } 72049b5e25fSSatish Balay /* reset ilen and imax for each row */ 72149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 72249b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 72349b5e25fSSatish Balay } 7246c6c5352SBarry Smith a->nz = ai[mbs]; 72549b5e25fSSatish Balay 726b424e231SHong Zhang /* diagonals may have moved, reset it */ 727b424e231SHong Zhang if (a->diag) { 728b424e231SHong Zhang ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr); 72949b5e25fSSatish Balay } 730b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 7316c6c5352SBarry Smith m,A->m,a->bs,fshift*bs2,a->nz*bs2); 732b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 73349b5e25fSSatish Balay a->reallocs); 734b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 73549b5e25fSSatish Balay a->reallocs = 0; 73649b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 73749b5e25fSSatish Balay 73849b5e25fSSatish Balay PetscFunctionReturn(0); 73949b5e25fSSatish Balay } 74049b5e25fSSatish Balay 74149b5e25fSSatish Balay /* 74249b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 74349b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 74449b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 74549b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 74649b5e25fSSatish Balay */ 7474a2ae208SSatish Balay #undef __FUNCT__ 7484a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 749db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 75049b5e25fSSatish Balay { 75149b5e25fSSatish Balay int i,j,k,row; 75249b5e25fSSatish Balay PetscTruth flg; 75349b5e25fSSatish Balay 75449b5e25fSSatish Balay PetscFunctionBegin; 75549b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 75649b5e25fSSatish Balay row = idx[i]; 75749b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 75849b5e25fSSatish Balay sizes[j] = 1; 75949b5e25fSSatish Balay i++; 76049b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 76149b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 76249b5e25fSSatish Balay i++; 76349b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 76449b5e25fSSatish Balay flg = PETSC_TRUE; 76549b5e25fSSatish Balay for (k=1; k<bs; k++) { 76649b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 76749b5e25fSSatish Balay flg = PETSC_FALSE; 76849b5e25fSSatish Balay break; 76949b5e25fSSatish Balay } 77049b5e25fSSatish Balay } 77149b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 77249b5e25fSSatish Balay sizes[j] = bs; 77349b5e25fSSatish Balay i+= bs; 77449b5e25fSSatish Balay } else { 77549b5e25fSSatish Balay sizes[j] = 1; 77649b5e25fSSatish Balay i++; 77749b5e25fSSatish Balay } 77849b5e25fSSatish Balay } 77949b5e25fSSatish Balay } 78049b5e25fSSatish Balay *bs_max = j; 78149b5e25fSSatish Balay PetscFunctionReturn(0); 78249b5e25fSSatish Balay } 78349b5e25fSSatish Balay 7844a2ae208SSatish Balay #undef __FUNCT__ 7854a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ" 786268466fbSBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag) 78749b5e25fSSatish Balay { 78849b5e25fSSatish Balay PetscFunctionBegin; 789c0f24835SHong Zhang SETERRQ(PETSC_ERR_SUP,"No support for this function yet"); 79049b5e25fSSatish Balay } 79149b5e25fSSatish Balay 79249b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 79349b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 79449b5e25fSSatish Balay */ 79549b5e25fSSatish Balay 7964a2ae208SSatish Balay #undef __FUNCT__ 7974a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ" 798f15d580aSBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is) 79949b5e25fSSatish Balay { 80049b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 80149b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 80249b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 80349b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 80449b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 80549b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 80649b5e25fSSatish Balay 80749b5e25fSSatish Balay PetscFunctionBegin; 80849b5e25fSSatish Balay 80949b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 81049b5e25fSSatish Balay row = im[k]; /* row number */ 81149b5e25fSSatish Balay brow = row/bs; /* block row number */ 81249b5e25fSSatish Balay if (row < 0) continue; 81349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 814590ac198SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1); 81549b5e25fSSatish Balay #endif 81649b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 81749b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 81849b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 81949b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 82049b5e25fSSatish Balay low = 0; 82149b5e25fSSatish Balay 82249b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 82349b5e25fSSatish Balay if (in[l] < 0) continue; 82449b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 825590ac198SBarry Smith if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m-1); 82649b5e25fSSatish Balay #endif 82749b5e25fSSatish Balay col = in[l]; 82849b5e25fSSatish Balay bcol = col/bs; /* block col number */ 82949b5e25fSSatish Balay 83049b5e25fSSatish Balay if (brow <= bcol){ 83149b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 8328549e402SHong Zhang if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ 83349b5e25fSSatish Balay /* element value a(k,l) */ 83449b5e25fSSatish Balay if (roworiented) { 83549b5e25fSSatish Balay value = v[l + k*n]; 83649b5e25fSSatish Balay } else { 83749b5e25fSSatish Balay value = v[k + l*m]; 83849b5e25fSSatish Balay } 83949b5e25fSSatish Balay 84049b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 84149b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 84249b5e25fSSatish Balay while (high-low > 7) { 84349b5e25fSSatish Balay t = (low+high)/2; 84449b5e25fSSatish Balay if (rp[t] > bcol) high = t; 84549b5e25fSSatish Balay else low = t; 84649b5e25fSSatish Balay } 84749b5e25fSSatish Balay for (i=low; i<high; i++) { 84849b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 84949b5e25fSSatish Balay if (rp[i] > bcol) break; 85049b5e25fSSatish Balay if (rp[i] == bcol) { 85149b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 85249b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 85349b5e25fSSatish Balay else *bap = value; 8548549e402SHong Zhang /* for diag block, add/insert its symmetric element a(cidx,ridx) */ 8558549e402SHong Zhang if (brow == bcol && ridx < cidx){ 8568549e402SHong Zhang bap = ap + bs2*i + bs*ridx + cidx; 8578549e402SHong Zhang if (is == ADD_VALUES) *bap += value; 8588549e402SHong Zhang else *bap = value; 8598549e402SHong Zhang } 86049b5e25fSSatish Balay goto noinsert1; 86149b5e25fSSatish Balay } 86249b5e25fSSatish Balay } 86349b5e25fSSatish Balay 86449b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 865a45adfd6SMatthew Knepley else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 86649b5e25fSSatish Balay if (nrow >= rmax) { 86749b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 86849b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 86949b5e25fSSatish Balay MatScalar *new_a; 87049b5e25fSSatish Balay 871a45adfd6SMatthew Knepley if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 87249b5e25fSSatish Balay 87349b5e25fSSatish Balay /* Malloc new storage space */ 87449b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 87582502324SSatish Balay ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 87649b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 87749b5e25fSSatish Balay new_i = new_j + new_nz; 87849b5e25fSSatish Balay 87949b5e25fSSatish Balay /* copy over old data into new slots */ 88049b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 88149b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 88249b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 88349b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 88449b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 88549b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 88649b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 88749b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 88849b5e25fSSatish Balay /* free up old matrix storage */ 88949b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 89049b5e25fSSatish Balay if (!a->singlemalloc) { 89149b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 89249b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 89349b5e25fSSatish Balay } 89449b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 89549b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 89649b5e25fSSatish Balay 89749b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 89849b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 899b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 9006c6c5352SBarry Smith a->maxnz += bs2*CHUNKSIZE; 90149b5e25fSSatish Balay a->reallocs++; 9026c6c5352SBarry Smith a->nz++; 90349b5e25fSSatish Balay } 90449b5e25fSSatish Balay 90549b5e25fSSatish Balay N = nrow++ - 1; 90649b5e25fSSatish Balay /* shift up all the later entries in this row */ 90749b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 90849b5e25fSSatish Balay rp[ii+1] = rp[ii]; 90949b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 91049b5e25fSSatish Balay } 91149b5e25fSSatish Balay if (N>=i) { 91249b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 91349b5e25fSSatish Balay } 91449b5e25fSSatish Balay rp[i] = bcol; 91549b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 91649b5e25fSSatish Balay noinsert1:; 91749b5e25fSSatish Balay low = i; 91849b5e25fSSatish Balay /* } */ 9198549e402SHong Zhang } 92049b5e25fSSatish Balay } /* end of if .. if.. */ 92149b5e25fSSatish Balay } /* end of loop over added columns */ 92249b5e25fSSatish Balay ailen[brow] = nrow; 92349b5e25fSSatish Balay } /* end of loop over added rows */ 92449b5e25fSSatish Balay 92549b5e25fSSatish Balay PetscFunctionReturn(0); 92649b5e25fSSatish Balay } 92749b5e25fSSatish Balay 92815e8a5b3SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*); 92915e8a5b3SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*); 93052ebccd6SSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int); 93149b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 932268466fbSBarry Smith extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]); 93349b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 93449b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 935268466fbSBarry Smith extern int MatScale_SeqSBAIJ(const PetscScalar*,Mat); 93649b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 93749b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 93849b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 93949b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 94049b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 94149b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 94213a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 94349b5e25fSSatish Balay 94449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 94549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 94649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 94749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 94849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 94949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 95049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 95149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 95249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 95349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 95449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 95549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 95649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 95749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 95849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 95949b5e25fSSatish Balay 960d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs); 961d59c15a7SBarry Smith 96207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec); 96307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec); 96407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec); 96507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 96607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec); 96707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec); 96807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec); 96907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec); 97007f98182SSatish Balay 9715f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9725f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9735f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9745f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9755f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9765f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9775f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 97857d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*); 9793e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*); 98049b5e25fSSatish Balay 98149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 98249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 98349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 98449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 98549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 98649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 98749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 98849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 98949b5e25fSSatish Balay 99049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 99149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 99249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 99349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 99449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 99549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 99649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 99749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 99849b5e25fSSatish Balay 9994d101231SSatish Balay #ifdef HAVE_MatICCFactor 10001a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 10014a2ae208SSatish Balay #undef __FUNCT__ 10024d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ" 100315e8a5b3SHong Zhang int MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info) 100449b5e25fSSatish Balay { 10054ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 100649b5e25fSSatish Balay Mat outA; 100749b5e25fSSatish Balay int ierr; 100849b5e25fSSatish Balay PetscTruth row_identity,col_identity; 100949b5e25fSSatish Balay 101049b5e25fSSatish Balay PetscFunctionBegin; 10111a3463dfSHong Zhang /* 101229bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 101349b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 101449b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 101549b5e25fSSatish Balay if (!row_identity || !col_identity) { 101629bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 101749b5e25fSSatish Balay } 10181a3463dfSHong Zhang */ 101949b5e25fSSatish Balay 102049b5e25fSSatish Balay outA = inA; 10211260e1f8SHong Zhang inA->factor = FACTOR_CHOLESKY; 102249b5e25fSSatish Balay 102349b5e25fSSatish Balay if (!a->diag) { 10241a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 102549b5e25fSSatish Balay } 102649b5e25fSSatish Balay /* 102749b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 102849b5e25fSSatish Balay for ILU(0) factorization with natural ordering 102949b5e25fSSatish Balay */ 103049b5e25fSSatish Balay switch (a->bs) { 103149b5e25fSSatish Balay case 1: 10320fe9e5beSBarry Smith inA->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 103307f98182SSatish Balay inA->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1034d59c15a7SBarry Smith inA->ops->solves = MatSolves_SeqSBAIJ_1; 10354d101231SSatish Balay PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 103649b5e25fSSatish Balay case 2: 10371a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 10381a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10391a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 10404d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 104149b5e25fSSatish Balay break; 104249b5e25fSSatish Balay case 3: 10431a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10441a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10451a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10464d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 104749b5e25fSSatish Balay break; 104849b5e25fSSatish Balay case 4: 10491a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10501a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10511a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10524d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 105349b5e25fSSatish Balay break; 105449b5e25fSSatish Balay case 5: 10551a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10561a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10571a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10584d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 105949b5e25fSSatish Balay break; 106049b5e25fSSatish Balay case 6: 10611a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10621a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10631a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10644d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 106549b5e25fSSatish Balay break; 106649b5e25fSSatish Balay case 7: 10671a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10681a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10691a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10704d101231SSatish Balay PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 107149b5e25fSSatish Balay break; 107249b5e25fSSatish Balay default: 107349b5e25fSSatish Balay a->row = row; 10741a3463dfSHong Zhang a->icol = col; 107549b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 107649b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 107749b5e25fSSatish Balay 107849b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 107949b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 1080b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 108149b5e25fSSatish Balay 108249b5e25fSSatish Balay if (!a->solve_work) { 108387828ca2SBarry Smith ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 108487828ca2SBarry Smith PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar)); 108549b5e25fSSatish Balay } 108649b5e25fSSatish Balay } 108749b5e25fSSatish Balay 10881a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 108949b5e25fSSatish Balay 109049b5e25fSSatish Balay PetscFunctionReturn(0); 109149b5e25fSSatish Balay } 1092950f1e5bSHong Zhang #endif 1093950f1e5bSHong Zhang 10944a2ae208SSatish Balay #undef __FUNCT__ 10954a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ" 109649b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 109749b5e25fSSatish Balay { 109849b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 109949b5e25fSSatish Balay MPI_Comm comm = A->comm; 110049b5e25fSSatish Balay int ierr; 110149b5e25fSSatish Balay 110249b5e25fSSatish Balay PetscFunctionBegin; 110349b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 110449b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 110549b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 110649b5e25fSSatish Balay PetscFunctionReturn(0); 110749b5e25fSSatish Balay } 110849b5e25fSSatish Balay 110949b5e25fSSatish Balay EXTERN_C_BEGIN 11104a2ae208SSatish Balay #undef __FUNCT__ 11114a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 111249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 111349b5e25fSSatish Balay { 1114045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 111549b5e25fSSatish Balay int i,nz,n; 111649b5e25fSSatish Balay 111749b5e25fSSatish Balay PetscFunctionBegin; 11186c6c5352SBarry Smith nz = baij->maxnz; 1119c464158bSHong Zhang n = mat->n; 112049b5e25fSSatish Balay for (i=0; i<nz; i++) { 112149b5e25fSSatish Balay baij->j[i] = indices[i]; 112249b5e25fSSatish Balay } 11236c6c5352SBarry Smith baij->nz = nz; 112449b5e25fSSatish Balay for (i=0; i<n; i++) { 112549b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 112649b5e25fSSatish Balay } 112749b5e25fSSatish Balay 112849b5e25fSSatish Balay PetscFunctionReturn(0); 112949b5e25fSSatish Balay } 113049b5e25fSSatish Balay EXTERN_C_END 113149b5e25fSSatish Balay 11324a2ae208SSatish Balay #undef __FUNCT__ 11334a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices" 113449b5e25fSSatish Balay /*@ 113519585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 113649b5e25fSSatish Balay in the matrix. 113749b5e25fSSatish Balay 113849b5e25fSSatish Balay Input Parameters: 113919585528SSatish Balay + mat - the SeqSBAIJ matrix 114049b5e25fSSatish Balay - indices - the column indices 114149b5e25fSSatish Balay 114249b5e25fSSatish Balay Level: advanced 114349b5e25fSSatish Balay 114449b5e25fSSatish Balay Notes: 114549b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 114649b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 114749b5e25fSSatish Balay of the MatSetValues() operation. 114849b5e25fSSatish Balay 114949b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 115019585528SSatish Balay MatCreateSeqSBAIJ(). 115149b5e25fSSatish Balay 1152ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 115349b5e25fSSatish Balay 1154ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 115549b5e25fSSatish Balay @*/ 115649b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 115749b5e25fSSatish Balay { 115849b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 115949b5e25fSSatish Balay 116049b5e25fSSatish Balay PetscFunctionBegin; 116149b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 1162c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 116349b5e25fSSatish Balay if (f) { 116449b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 116549b5e25fSSatish Balay } else { 116629bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 116749b5e25fSSatish Balay } 116849b5e25fSSatish Balay PetscFunctionReturn(0); 116949b5e25fSSatish Balay } 117049b5e25fSSatish Balay 11714a2ae208SSatish Balay #undef __FUNCT__ 11724a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ" 1173273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A) 1174273d9f13SBarry Smith { 1175273d9f13SBarry Smith int ierr; 1176273d9f13SBarry Smith 1177273d9f13SBarry Smith PetscFunctionBegin; 1178273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1179273d9f13SBarry Smith PetscFunctionReturn(0); 1180273d9f13SBarry Smith } 1181273d9f13SBarry Smith 1182a6ece127SHong Zhang #undef __FUNCT__ 1183a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ" 11844e7234bfSBarry Smith int MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1185a6ece127SHong Zhang { 1186a6ece127SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1187a6ece127SHong Zhang PetscFunctionBegin; 1188a6ece127SHong Zhang *array = a->a; 1189a6ece127SHong Zhang PetscFunctionReturn(0); 1190a6ece127SHong Zhang } 1191a6ece127SHong Zhang 1192a6ece127SHong Zhang #undef __FUNCT__ 1193a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ" 11944e7234bfSBarry Smith int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[]) 1195a6ece127SHong Zhang { 1196a6ece127SHong Zhang PetscFunctionBegin; 1197a6ece127SHong Zhang PetscFunctionReturn(0); 1198a6ece127SHong Zhang } 1199a6ece127SHong Zhang 120042ee4b1aSHong Zhang #include "petscblaslapack.h" 120142ee4b1aSHong Zhang #undef __FUNCT__ 120242ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ" 1203268466fbSBarry Smith int MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 120442ee4b1aSHong Zhang { 120542ee4b1aSHong Zhang Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data; 1206c4319e64SHong Zhang int ierr,one=1,i,bs=y->bs,bs2,j; 120742ee4b1aSHong Zhang 120842ee4b1aSHong Zhang PetscFunctionBegin; 120942ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 12106c6c5352SBarry Smith BLaxpy_(&x->nz,(PetscScalar*)a,x->a,&one,y->a,&one); 1211c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1212c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1213c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1214c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1215c537a176SHong Zhang } 1216c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1217c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1218c4319e64SHong Zhang y->XtoY = X; 1219c537a176SHong Zhang } 1220c4319e64SHong Zhang bs2 = bs*bs; 12216c6c5352SBarry Smith for (i=0; i<x->nz; i++) { 1222c4319e64SHong Zhang j = 0; 1223c4319e64SHong Zhang while (j < bs2){ 12246550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1225c4319e64SHong Zhang j++; 1226c537a176SHong Zhang } 1227c4319e64SHong Zhang } 12286c6c5352SBarry Smith PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %d/%d = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz)); 122942ee4b1aSHong Zhang } else { 123042ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 123142ee4b1aSHong Zhang } 123242ee4b1aSHong Zhang PetscFunctionReturn(0); 123342ee4b1aSHong Zhang } 123442ee4b1aSHong Zhang 1235efcf0fc3SBarry Smith #undef __FUNCT__ 1236efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ" 1237efcf0fc3SBarry Smith int MatIsSymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1238efcf0fc3SBarry Smith { 1239efcf0fc3SBarry Smith PetscFunctionBegin; 1240efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1241efcf0fc3SBarry Smith PetscFunctionReturn(0); 1242efcf0fc3SBarry Smith } 1243efcf0fc3SBarry Smith 1244efcf0fc3SBarry Smith #undef __FUNCT__ 1245efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ" 1246efcf0fc3SBarry Smith int MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg) 1247efcf0fc3SBarry Smith { 1248efcf0fc3SBarry Smith PetscFunctionBegin; 1249efcf0fc3SBarry Smith *flg = PETSC_TRUE; 1250efcf0fc3SBarry Smith PetscFunctionReturn(0); 1251efcf0fc3SBarry Smith } 1252efcf0fc3SBarry Smith 1253efcf0fc3SBarry Smith #undef __FUNCT__ 1254efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ" 1255efcf0fc3SBarry Smith int MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg) 1256efcf0fc3SBarry Smith { 1257efcf0fc3SBarry Smith PetscFunctionBegin; 1258efcf0fc3SBarry Smith *flg = PETSC_FALSE; 1259efcf0fc3SBarry Smith PetscFunctionReturn(0); 1260efcf0fc3SBarry Smith } 1261efcf0fc3SBarry Smith 126249b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 126349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 126449b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 126549b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 126649b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 126797304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N, 126849b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 126949b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 127049b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 127149b5e25fSSatish Balay 0, 127249b5e25fSSatish Balay 0, 127397304618SKris Buschelman /*10*/ 0, 127449b5e25fSSatish Balay 0, 12755f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 1276d06b337dSHong Zhang MatRelax_SeqSBAIJ, 127749b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 127897304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ, 127949b5e25fSSatish Balay MatEqual_SeqSBAIJ, 128049b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 128149b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 128249b5e25fSSatish Balay MatNorm_SeqSBAIJ, 128397304618SKris Buschelman /*20*/ 0, 128449b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 128549b5e25fSSatish Balay 0, 128649b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 128749b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 128897304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ, 128949b5e25fSSatish Balay 0, 129049b5e25fSSatish Balay 0, 12915f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 12925f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 129397304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ, 1294c464158bSHong Zhang 0, 12954d101231SSatish Balay MatICCFactorSymbolic_SeqSBAIJ, 1296a6ece127SHong Zhang MatGetArray_SeqSBAIJ, 1297a6ece127SHong Zhang MatRestoreArray_SeqSBAIJ, 129897304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ, 129949b5e25fSSatish Balay 0, 130049b5e25fSSatish Balay 0, 130149b5e25fSSatish Balay 0, 1302950f1e5bSHong Zhang 0, 130397304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ, 130449b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 130549b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 130649b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 130749b5e25fSSatish Balay 0, 130897304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ, 130949b5e25fSSatish Balay MatScale_SeqSBAIJ, 131049b5e25fSSatish Balay 0, 131149b5e25fSSatish Balay 0, 131249b5e25fSSatish Balay 0, 131397304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ, 131449b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 131549b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 131649b5e25fSSatish Balay 0, 131749b5e25fSSatish Balay 0, 131897304618SKris Buschelman /*55*/ 0, 131949b5e25fSSatish Balay 0, 132049b5e25fSSatish Balay 0, 132149b5e25fSSatish Balay 0, 132249b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 132397304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ, 132449b5e25fSSatish Balay 0, 132549b5e25fSSatish Balay 0, 13268a124369SBarry Smith MatGetPetscMaps_Petsc, 1327d959ec07SHong Zhang 0, 132897304618SKris Buschelman /*65*/ 0, 1329d959ec07SHong Zhang 0, 1330d959ec07SHong Zhang 0, 1331d959ec07SHong Zhang 0, 1332d959ec07SHong Zhang 0, 133397304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ, 13343e0d88b5SBarry Smith 0, 13353e0d88b5SBarry Smith 0, 13363e0d88b5SBarry Smith 0, 13373e0d88b5SBarry Smith 0, 133897304618SKris Buschelman /*75*/ 0, 13393e0d88b5SBarry Smith 0, 13403e0d88b5SBarry Smith 0, 13413e0d88b5SBarry Smith 0, 13423e0d88b5SBarry Smith 0, 134397304618SKris Buschelman /*80*/ 0, 13443e0d88b5SBarry Smith 0, 13453e0d88b5SBarry Smith 0, 13463e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX) 134797304618SKris Buschelman MatGetInertia_SeqSBAIJ, 13483e0d88b5SBarry Smith #else 134997304618SKris Buschelman 0, 13503e0d88b5SBarry Smith #endif 1351efcf0fc3SBarry Smith /*85*/ MatLoad_SeqSBAIJ, 1352efcf0fc3SBarry Smith MatIsSymmetric_SeqSBAIJ, 1353efcf0fc3SBarry Smith MatIsStructurallySymmetric_SeqSBAIJ, 1354efcf0fc3SBarry Smith MatIsHermitian_SeqSBAIJ 13553e0d88b5SBarry Smith }; 135649b5e25fSSatish Balay 135749b5e25fSSatish Balay EXTERN_C_BEGIN 13584a2ae208SSatish Balay #undef __FUNCT__ 13594a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ" 136049b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 136149b5e25fSSatish Balay { 13624afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1363c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 136449b5e25fSSatish Balay int ierr; 136549b5e25fSSatish Balay 136649b5e25fSSatish Balay PetscFunctionBegin; 136749b5e25fSSatish Balay if (aij->nonew != 1) { 136829bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 136949b5e25fSSatish Balay } 137049b5e25fSSatish Balay 137149b5e25fSSatish Balay /* allocate space for values if not already there */ 137249b5e25fSSatish Balay if (!aij->saved_values) { 137387828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 137449b5e25fSSatish Balay } 137549b5e25fSSatish Balay 137649b5e25fSSatish Balay /* copy values over */ 137787828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 137849b5e25fSSatish Balay PetscFunctionReturn(0); 137949b5e25fSSatish Balay } 138049b5e25fSSatish Balay EXTERN_C_END 138149b5e25fSSatish Balay 138249b5e25fSSatish Balay EXTERN_C_BEGIN 13834a2ae208SSatish Balay #undef __FUNCT__ 13844a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ" 138549b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 138649b5e25fSSatish Balay { 13874afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1388c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 138949b5e25fSSatish Balay 139049b5e25fSSatish Balay PetscFunctionBegin; 139149b5e25fSSatish Balay if (aij->nonew != 1) { 139229bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 139349b5e25fSSatish Balay } 139449b5e25fSSatish Balay if (!aij->saved_values) { 139529bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 139649b5e25fSSatish Balay } 139749b5e25fSSatish Balay 139849b5e25fSSatish Balay /* copy values over */ 139987828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 140049b5e25fSSatish Balay PetscFunctionReturn(0); 140149b5e25fSSatish Balay } 140249b5e25fSSatish Balay EXTERN_C_END 140349b5e25fSSatish Balay 14048549e402SHong Zhang EXTERN_C_BEGIN 14054a2ae208SSatish Balay #undef __FUNCT__ 1406a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ" 1407a23d5eceSKris Buschelman int MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz) 140849b5e25fSSatish Balay { 1409c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 14100c955e93SHong Zhang int i,len,ierr,mbs,bs2; 141149b5e25fSSatish Balay PetscTruth flg; 141249b5e25fSSatish Balay 141349b5e25fSSatish Balay PetscFunctionBegin; 1414273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1415e82a3eeeSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1416c464158bSHong Zhang mbs = B->m/bs; 141749b5e25fSSatish Balay bs2 = bs*bs; 141849b5e25fSSatish Balay 1419c464158bSHong Zhang if (mbs*bs != B->m) { 142029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 142149b5e25fSSatish Balay } 142249b5e25fSSatish Balay 1423435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3; 1424435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 142549b5e25fSSatish Balay if (nnz) { 142649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 142729bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 142880fe4e49SBarry 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); 142949b5e25fSSatish Balay } 143049b5e25fSSatish Balay } 143149b5e25fSSatish Balay 1432e82a3eeeSBarry Smith ierr = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr); 143349b5e25fSSatish Balay if (!flg) { 143449b5e25fSSatish Balay switch (bs) { 143549b5e25fSSatish Balay case 1: 14364ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 143749b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 1438d59c15a7SBarry Smith B->ops->solves = MatSolves_SeqSBAIJ_1; 143949b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 144049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 144149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 144249b5e25fSSatish Balay break; 144349b5e25fSSatish Balay case 2: 14444ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 144549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 144649b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 144749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 144849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 144949b5e25fSSatish Balay break; 145049b5e25fSSatish Balay case 3: 14515f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 145249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 145349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 145449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 145549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 145649b5e25fSSatish Balay break; 145749b5e25fSSatish Balay case 4: 14585f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 145949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 146049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 146149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 146249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 146349b5e25fSSatish Balay break; 146449b5e25fSSatish Balay case 5: 14655f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 146649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 146749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 146849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 146949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 147049b5e25fSSatish Balay break; 147149b5e25fSSatish Balay case 6: 14725f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 147349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 147449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 147549b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 147649b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 147749b5e25fSSatish Balay break; 147849b5e25fSSatish Balay case 7: 1479de53e5efSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7; 148049b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 148149b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 1482de53e5efSHong Zhang B->ops->mult = MatMult_SeqSBAIJ_7; 148349b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 148449b5e25fSSatish Balay break; 148549b5e25fSSatish Balay } 148649b5e25fSSatish Balay } 148749b5e25fSSatish Balay 148849b5e25fSSatish Balay b->mbs = mbs; 14894afc71dfSHong Zhang b->nbs = mbs; 1490b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 149149b5e25fSSatish Balay if (!nnz) { 1492435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 149349b5e25fSSatish Balay else if (nz <= 0) nz = 1; 149449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14958cef66ccSHong Zhang b->imax[i] = nz; 149649b5e25fSSatish Balay } 1497153ea458SHong Zhang nz = nz*mbs; /* total nz */ 149849b5e25fSSatish Balay } else { 149949b5e25fSSatish Balay nz = 0; 15008cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 150149b5e25fSSatish Balay } 15026c6c5352SBarry Smith /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 150349b5e25fSSatish Balay 150449b5e25fSSatish Balay /* allocate the matrix space */ 15056c6c5352SBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 150682502324SSatish Balay ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 15076c6c5352SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15086c6c5352SBarry Smith b->j = (int*)(b->a + nz*bs2); 15096c6c5352SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15106c6c5352SBarry Smith b->i = b->j + nz; 151149b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 151249b5e25fSSatish Balay 151349b5e25fSSatish Balay /* pointer to beginning of each row */ 151449b5e25fSSatish Balay b->i[0] = 0; 151549b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 151649b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 151749b5e25fSSatish Balay } 151849b5e25fSSatish Balay 151949b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 15205033bc48SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1521b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 152249b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 152349b5e25fSSatish Balay 152449b5e25fSSatish Balay b->bs = bs; 152549b5e25fSSatish Balay b->bs2 = bs2; 15266c6c5352SBarry Smith b->nz = 0; 15276c6c5352SBarry Smith b->maxnz = nz*bs2; 1528153ea458SHong Zhang 152916cdd363SHong Zhang b->inew = 0; 153016cdd363SHong Zhang b->jnew = 0; 153116cdd363SHong Zhang b->anew = 0; 153216cdd363SHong Zhang b->a2anew = 0; 15331a3463dfSHong Zhang b->permute = PETSC_FALSE; 1534c464158bSHong Zhang PetscFunctionReturn(0); 1535c464158bSHong Zhang } 1536a23d5eceSKris Buschelman EXTERN_C_END 1537153ea458SHong Zhang 1538d769727bSBarry Smith EXTERN_C_BEGIN 1539d769727bSBarry Smith EXTERN int MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,Mat*); 1540a0e1a404SHong Zhang EXTERN int MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,Mat*); 1541d769727bSBarry Smith EXTERN_C_END 1542d769727bSBarry Smith 15430bad9183SKris Buschelman /*MC 1544fafad747SKris Buschelman MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices, 15450bad9183SKris Buschelman based on block compressed sparse row format. Only the upper triangular portion of the matrix is stored. 15460bad9183SKris Buschelman 15470bad9183SKris Buschelman Options Database Keys: 15480bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions() 15490bad9183SKris Buschelman 15500bad9183SKris Buschelman Level: beginner 15510bad9183SKris Buschelman 15520bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ 15530bad9183SKris Buschelman M*/ 15540bad9183SKris Buschelman 1555a23d5eceSKris Buschelman EXTERN_C_BEGIN 1556a23d5eceSKris Buschelman #undef __FUNCT__ 1557a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ" 1558a23d5eceSKris Buschelman int MatCreate_SeqSBAIJ(Mat B) 1559a23d5eceSKris Buschelman { 1560a23d5eceSKris Buschelman Mat_SeqSBAIJ *b; 1561a23d5eceSKris Buschelman int ierr,size; 1562a23d5eceSKris Buschelman 1563a23d5eceSKris Buschelman PetscFunctionBegin; 1564a23d5eceSKris Buschelman ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1565a23d5eceSKris Buschelman if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1566a23d5eceSKris Buschelman B->m = B->M = PetscMax(B->m,B->M); 1567a23d5eceSKris Buschelman B->n = B->N = PetscMax(B->n,B->N); 1568a23d5eceSKris Buschelman 1569a23d5eceSKris Buschelman ierr = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr); 1570a23d5eceSKris Buschelman B->data = (void*)b; 1571a23d5eceSKris Buschelman ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1572a23d5eceSKris Buschelman ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1573a23d5eceSKris Buschelman B->ops->destroy = MatDestroy_SeqSBAIJ; 1574a23d5eceSKris Buschelman B->ops->view = MatView_SeqSBAIJ; 1575a23d5eceSKris Buschelman B->factor = 0; 1576a23d5eceSKris Buschelman B->lupivotthreshold = 1.0; 1577a23d5eceSKris Buschelman B->mapping = 0; 1578a23d5eceSKris Buschelman b->row = 0; 1579a23d5eceSKris Buschelman b->icol = 0; 1580a23d5eceSKris Buschelman b->reallocs = 0; 1581a23d5eceSKris Buschelman b->saved_values = 0; 1582a23d5eceSKris Buschelman 1583a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1584a23d5eceSKris Buschelman ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1585a23d5eceSKris Buschelman 1586a23d5eceSKris Buschelman b->sorted = PETSC_FALSE; 1587a23d5eceSKris Buschelman b->roworiented = PETSC_TRUE; 1588a23d5eceSKris Buschelman b->nonew = 0; 1589a23d5eceSKris Buschelman b->diag = 0; 1590a23d5eceSKris Buschelman b->solve_work = 0; 1591a23d5eceSKris Buschelman b->mult_work = 0; 1592a23d5eceSKris Buschelman B->spptr = 0; 1593a23d5eceSKris Buschelman b->keepzeroedrows = PETSC_FALSE; 1594a23d5eceSKris Buschelman b->xtoy = 0; 1595a23d5eceSKris Buschelman b->XtoY = 0; 1596a23d5eceSKris Buschelman 1597a23d5eceSKris Buschelman b->inew = 0; 1598a23d5eceSKris Buschelman b->jnew = 0; 1599a23d5eceSKris Buschelman b->anew = 0; 1600a23d5eceSKris Buschelman b->a2anew = 0; 1601a23d5eceSKris Buschelman b->permute = PETSC_FALSE; 1602a23d5eceSKris Buschelman 1603a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1604a23d5eceSKris Buschelman "MatStoreValues_SeqSBAIJ", 1605a23d5eceSKris Buschelman MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1606a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1607a23d5eceSKris Buschelman "MatRetrieveValues_SeqSBAIJ", 1608a23d5eceSKris Buschelman (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1609a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1610a23d5eceSKris Buschelman "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1611a23d5eceSKris Buschelman MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 16124e5e7fe4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C", 16134e5e7fe4SHong Zhang "MatConvert_SeqSBAIJ_SeqAIJ", 16144e5e7fe4SHong Zhang MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr); 1615a0e1a404SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C", 1616a0e1a404SHong Zhang "MatConvert_SeqSBAIJ_SeqBAIJ", 1617a0e1a404SHong Zhang MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr); 1618a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C", 1619a23d5eceSKris Buschelman "MatSeqSBAIJSetPreallocation_SeqSBAIJ", 1620a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr); 1621a23d5eceSKris Buschelman PetscFunctionReturn(0); 1622a23d5eceSKris Buschelman } 1623a23d5eceSKris Buschelman EXTERN_C_END 1624a23d5eceSKris Buschelman 1625a23d5eceSKris Buschelman #undef __FUNCT__ 1626a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation" 1627a23d5eceSKris Buschelman /*@C 1628a23d5eceSKris Buschelman MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block 1629a23d5eceSKris Buschelman compressed row) format. For good matrix assembly performance the 1630a23d5eceSKris Buschelman user should preallocate the matrix storage by setting the parameter nz 1631a23d5eceSKris Buschelman (or the array nnz). By setting these parameters accurately, performance 1632a23d5eceSKris Buschelman during matrix assembly can be increased by more than a factor of 50. 1633a23d5eceSKris Buschelman 1634a23d5eceSKris Buschelman Collective on Mat 1635a23d5eceSKris Buschelman 1636a23d5eceSKris Buschelman Input Parameters: 1637a23d5eceSKris Buschelman + A - the symmetric matrix 1638a23d5eceSKris Buschelman . bs - size of block 1639a23d5eceSKris Buschelman . nz - number of block nonzeros per block row (same for all rows) 1640a23d5eceSKris Buschelman - nnz - array containing the number of block nonzeros in the upper triangular plus 1641a23d5eceSKris Buschelman diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1642a23d5eceSKris Buschelman 1643a23d5eceSKris Buschelman Options Database Keys: 1644a23d5eceSKris Buschelman . -mat_no_unroll - uses code that does not unroll the loops in the 1645a23d5eceSKris Buschelman block calculations (much slower) 1646a23d5eceSKris Buschelman . -mat_block_size - size of the blocks to use 1647a23d5eceSKris Buschelman 1648a23d5eceSKris Buschelman Level: intermediate 1649a23d5eceSKris Buschelman 1650a23d5eceSKris Buschelman Notes: 1651a23d5eceSKris Buschelman Specify the preallocated storage with either nz or nnz (not both). 1652a23d5eceSKris Buschelman Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1653a23d5eceSKris Buschelman allocation. For additional details, see the users manual chapter on 1654a23d5eceSKris Buschelman matrices. 1655a23d5eceSKris Buschelman 1656a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1657a23d5eceSKris Buschelman @*/ 1658ca01db9bSBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) { 1659ca01db9bSBarry Smith int ierr,(*f)(Mat,int,int,const int[]); 1660a23d5eceSKris Buschelman 1661a23d5eceSKris Buschelman PetscFunctionBegin; 1662a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1663a23d5eceSKris Buschelman if (f) { 1664a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 1665a23d5eceSKris Buschelman } 1666a23d5eceSKris Buschelman PetscFunctionReturn(0); 1667a23d5eceSKris Buschelman } 166849b5e25fSSatish Balay 16694a2ae208SSatish Balay #undef __FUNCT__ 16704a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ" 1671c464158bSHong Zhang /*@C 1672c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1673c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1674c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1675c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1676c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 167749b5e25fSSatish Balay 1678c464158bSHong Zhang Collective on MPI_Comm 1679c464158bSHong Zhang 1680c464158bSHong Zhang Input Parameters: 1681c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1682c464158bSHong Zhang . bs - size of block 1683c464158bSHong Zhang . m - number of rows, or number of columns 1684c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1685744e8345SSatish Balay - nnz - array containing the number of block nonzeros in the upper triangular plus 1686744e8345SSatish Balay diagonal portion of each block (possibly different for each block row) or PETSC_NULL 1687c464158bSHong Zhang 1688c464158bSHong Zhang Output Parameter: 1689c464158bSHong Zhang . A - the symmetric matrix 1690c464158bSHong Zhang 1691c464158bSHong Zhang Options Database Keys: 1692c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1693c464158bSHong Zhang block calculations (much slower) 1694c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1695c464158bSHong Zhang 1696c464158bSHong Zhang Level: intermediate 1697c464158bSHong Zhang 1698c464158bSHong Zhang Notes: 1699c464158bSHong Zhang 1700c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1701c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1702c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1703c464158bSHong Zhang matrices. 1704c464158bSHong Zhang 1705c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1706c464158bSHong Zhang @*/ 1707ca01db9bSBarry Smith int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A) 1708c464158bSHong Zhang { 1709c464158bSHong Zhang int ierr; 1710c464158bSHong Zhang 1711c464158bSHong Zhang PetscFunctionBegin; 1712c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1713c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1714273d9f13SBarry Smith ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 171549b5e25fSSatish Balay PetscFunctionReturn(0); 171649b5e25fSSatish Balay } 171749b5e25fSSatish Balay 17184a2ae208SSatish Balay #undef __FUNCT__ 17194a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ" 172049b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 172149b5e25fSSatish Balay { 172249b5e25fSSatish Balay Mat C; 172349b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 17246c6c5352SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 172549b5e25fSSatish Balay 172649b5e25fSSatish Balay PetscFunctionBegin; 172729bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 172849b5e25fSSatish Balay 172949b5e25fSSatish Balay *B = 0; 1730692f9cbeSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1731*be5d1d56SKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 1732692f9cbeSHong Zhang c = (Mat_SeqSBAIJ*)C->data; 1733692f9cbeSHong Zhang 173449b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1735273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 173649b5e25fSSatish Balay C->factor = A->factor; 173749b5e25fSSatish Balay c->row = 0; 173849b5e25fSSatish Balay c->icol = 0; 173949b5e25fSSatish Balay c->saved_values = 0; 174049b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 174149b5e25fSSatish Balay C->assembled = PETSC_TRUE; 174249b5e25fSSatish Balay 174382327fa8SHong Zhang C->M = A->M; 174482327fa8SHong Zhang C->N = A->N; 174549b5e25fSSatish Balay c->bs = a->bs; 174649b5e25fSSatish Balay c->bs2 = a->bs2; 174749b5e25fSSatish Balay c->mbs = a->mbs; 174849b5e25fSSatish Balay c->nbs = a->nbs; 174949b5e25fSSatish Balay 1750b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1751b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 175249b5e25fSSatish Balay for (i=0; i<mbs; i++) { 175349b5e25fSSatish Balay c->imax[i] = a->imax[i]; 175449b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 175549b5e25fSSatish Balay } 175649b5e25fSSatish Balay 175749b5e25fSSatish Balay /* allocate the matrix space */ 175849b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 175949b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 176082502324SSatish Balay ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 176149b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 176249b5e25fSSatish Balay c->i = c->j + nz; 176349b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 176449b5e25fSSatish Balay if (mbs > 0) { 176549b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 176649b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 176749b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 176849b5e25fSSatish Balay } else { 176949b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 177049b5e25fSSatish Balay } 177149b5e25fSSatish Balay } 177249b5e25fSSatish Balay 1773b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 177449b5e25fSSatish Balay c->sorted = a->sorted; 177549b5e25fSSatish Balay c->roworiented = a->roworiented; 177649b5e25fSSatish Balay c->nonew = a->nonew; 177749b5e25fSSatish Balay 177849b5e25fSSatish Balay if (a->diag) { 1779b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1780b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 178149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 178249b5e25fSSatish Balay c->diag[i] = a->diag[i]; 178349b5e25fSSatish Balay } 178449b5e25fSSatish Balay } else c->diag = 0; 17856c6c5352SBarry Smith c->nz = a->nz; 17866c6c5352SBarry Smith c->maxnz = a->maxnz; 178749b5e25fSSatish Balay c->solve_work = 0; 17882a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 178949b5e25fSSatish Balay c->mult_work = 0; 179049b5e25fSSatish Balay *B = C; 1791b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 179249b5e25fSSatish Balay PetscFunctionReturn(0); 179349b5e25fSSatish Balay } 179449b5e25fSSatish Balay 17954a2ae208SSatish Balay #undef __FUNCT__ 17964a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ" 17978e9aea5cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A) 179849b5e25fSSatish Balay { 179949b5e25fSSatish Balay Mat_SeqSBAIJ *a; 180049b5e25fSSatish Balay Mat B; 180149b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 18023f7326a9SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount; 180349b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 180449b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 180587828ca2SBarry Smith PetscScalar *aa; 180649b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 180749b5e25fSSatish Balay 180849b5e25fSSatish Balay PetscFunctionBegin; 1809b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 181049b5e25fSSatish Balay bs2 = bs*bs; 181149b5e25fSSatish Balay 181249b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 181329bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1814b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 181549b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1816552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 181749b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 181849b5e25fSSatish Balay 181949b5e25fSSatish Balay if (header[3] < 0) { 182029bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 182149b5e25fSSatish Balay } 182249b5e25fSSatish Balay 182329bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 182449b5e25fSSatish Balay 182549b5e25fSSatish Balay /* 182649b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 182749b5e25fSSatish Balay divisible by the blocksize 182849b5e25fSSatish Balay */ 182949b5e25fSSatish Balay mbs = M/bs; 183049b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 183149b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 183249b5e25fSSatish Balay else mbs++; 183349b5e25fSSatish Balay if (extra_rows) { 1834b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 183549b5e25fSSatish Balay } 183649b5e25fSSatish Balay 183749b5e25fSSatish Balay /* read in row lengths */ 1838b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 183949b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 184049b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 184149b5e25fSSatish Balay 184249b5e25fSSatish Balay /* read in column indices */ 1843b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 184449b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 184549b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 184649b5e25fSSatish Balay 184749b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 184882502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr); 1849a91a614fSBarry Smith ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 185082502324SSatish Balay ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 185149b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 185249b5e25fSSatish Balay masked = mask + mbs; 185349b5e25fSSatish Balay rowcount = 0; nzcount = 0; 185449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 185549b5e25fSSatish Balay nmask = 0; 185649b5e25fSSatish Balay for (j=0; j<bs; j++) { 185749b5e25fSSatish Balay kmax = rowlengths[rowcount]; 185849b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18592d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 186003630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 186149b5e25fSSatish Balay } 186249b5e25fSSatish Balay rowcount++; 186349b5e25fSSatish Balay } 1864574b2666SHong Zhang s_browlengths[i] += nmask; 1865574b2666SHong Zhang 186649b5e25fSSatish Balay /* zero out the mask elements we set */ 186749b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 186849b5e25fSSatish Balay } 186949b5e25fSSatish Balay 187049b5e25fSSatish Balay /* create our matrix */ 18719abb65ffSKris Buschelman ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr); 18729abb65ffSKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 187378473fd7SKris Buschelman ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr); 187449b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 187549b5e25fSSatish Balay 187649b5e25fSSatish Balay /* set matrix "i" values */ 187749b5e25fSSatish Balay a->i[0] = 0; 187849b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1879574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1880574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 188149b5e25fSSatish Balay } 18826c6c5352SBarry Smith a->nz = a->i[mbs]; 188349b5e25fSSatish Balay 188449b5e25fSSatish Balay /* read in nonzero values */ 188587828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 188649b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 188749b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 188849b5e25fSSatish Balay 188949b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 189049b5e25fSSatish Balay nzcount = 0; jcount = 0; 189149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 189249b5e25fSSatish Balay nzcountb = nzcount; 189349b5e25fSSatish Balay nmask = 0; 189449b5e25fSSatish Balay for (j=0; j<bs; j++) { 189549b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 189649b5e25fSSatish Balay for (k=0; k<kmax; k++) { 18972d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 189803630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 18992d703238SHong Zhang } 19002d703238SHong Zhang } 19012d703238SHong Zhang /* sort the masked values */ 19022d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 19032d703238SHong Zhang 19042d703238SHong Zhang /* set "j" values into matrix */ 19052d703238SHong Zhang maskcount = 1; 19062d703238SHong Zhang for (j=0; j<nmask; j++) { 190749b5e25fSSatish Balay a->j[jcount++] = masked[j]; 190849b5e25fSSatish Balay mask[masked[j]] = maskcount++; 190949b5e25fSSatish Balay } 1910574b2666SHong Zhang 191149b5e25fSSatish Balay /* set "a" values into matrix */ 191249b5e25fSSatish Balay ishift = bs2*a->i[i]; 191349b5e25fSSatish Balay for (j=0; j<bs; j++) { 191449b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 191549b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1916574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1917574b2666SHong Zhang if (tmp >= i){ 191849b5e25fSSatish Balay block = mask[tmp] - 1; 191949b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 192049b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1921574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1922574b2666SHong Zhang } 1923574b2666SHong Zhang nzcountb++; 192449b5e25fSSatish Balay } 192549b5e25fSSatish Balay } 192649b5e25fSSatish Balay /* zero out the mask elements we set */ 192749b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 192849b5e25fSSatish Balay } 19296c6c5352SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 193049b5e25fSSatish Balay 193149b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1932574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 193349b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 193449b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 193549b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 193649b5e25fSSatish Balay 19379abb65ffSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19389abb65ffSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 193949b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 19409abb65ffSKris Buschelman *A = B; 194149b5e25fSSatish Balay PetscFunctionReturn(0); 194249b5e25fSSatish Balay } 1943574b2666SHong Zhang 1944d06b337dSHong Zhang #undef __FUNCT__ 1945d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ" 1946c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx) 1947d06b337dSHong Zhang { 1948d06b337dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 1949d06b337dSHong Zhang MatScalar *aa=a->a,*v,*v1; 1950d06b337dSHong Zhang PetscScalar *x,*b,*t,sum,d; 1951d06b337dSHong Zhang int m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr; 1952d06b337dSHong Zhang int nz,nz1,*vj,*vj1,i; 1953d06b337dSHong Zhang 1954d06b337dSHong Zhang PetscFunctionBegin; 1955b965ef7fSBarry Smith its = its*lits; 195691723122SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 1957d06b337dSHong Zhang 1958d06b337dSHong Zhang if (bs > 1) 1959d06b337dSHong Zhang SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 1960d06b337dSHong Zhang 1961b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 1962d06b337dSHong Zhang if (xx != bb) { 1963b1d4fb26SBarry Smith ierr = VecGetArrayFast(bb,&b);CHKERRQ(ierr); 1964d06b337dSHong Zhang } else { 1965d06b337dSHong Zhang b = x; 1966d06b337dSHong Zhang } 1967d06b337dSHong Zhang 1968d06b337dSHong Zhang ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr); 1969d06b337dSHong Zhang 1970d06b337dSHong Zhang if (flag & SOR_ZERO_INITIAL_GUESS) { 19712798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 1972d06b337dSHong Zhang for (i=0; i<m; i++) 1973d06b337dSHong Zhang t[i] = b[i]; 1974d06b337dSHong Zhang 1975d06b337dSHong Zhang for (i=0; i<m; i++){ 197644706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 1977d06b337dSHong Zhang v = aa + ai[i] + 1; 1978d06b337dSHong Zhang vj = aj + ai[i] + 1; 1979d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1980d06b337dSHong Zhang x[i] = omega*t[i]/d; 1981d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */ 198244706875SHong Zhang PetscLogFlops(2*nz-1); 1983d06b337dSHong Zhang } 1984d06b337dSHong Zhang } 1985d06b337dSHong Zhang 19862798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 1987d06b337dSHong Zhang for (i=0; i<m; i++) 1988d06b337dSHong Zhang t[i] = b[i]; 1989d06b337dSHong Zhang 1990d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 1991d06b337dSHong Zhang v = aa + ai[i] + 1; 1992d06b337dSHong Zhang vj = aj + ai[i] + 1; 1993d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 1994d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 199544706875SHong Zhang PetscLogFlops(2*nz-1); 1996d06b337dSHong Zhang } 1997d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 1998d06b337dSHong Zhang d = *(aa + ai[i]); 1999d06b337dSHong Zhang v = aa + ai[i] + 1; 2000d06b337dSHong Zhang vj = aj + ai[i] + 1; 2001d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2002d06b337dSHong Zhang sum = t[i]; 2003d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 200444706875SHong Zhang PetscLogFlops(2*nz-1); 2005d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2006d06b337dSHong Zhang } 2007d06b337dSHong Zhang } 2008d06b337dSHong Zhang its--; 2009d06b337dSHong Zhang } 2010d06b337dSHong Zhang 2011d06b337dSHong Zhang while (its--) { 201244706875SHong Zhang /* 201344706875SHong Zhang forward sweep: 201444706875SHong Zhang for i=0,...,m-1: 201544706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 201644706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 201744706875SHong Zhang b = b - x[i]*U^T(i,:); 2018d06b337dSHong Zhang 201944706875SHong Zhang */ 20202798e883SHong Zhang if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 2021d06b337dSHong Zhang for (i=0; i<m; i++) 2022d06b337dSHong Zhang t[i] = b[i]; 2023d06b337dSHong Zhang 2024d06b337dSHong Zhang for (i=0; i<m; i++){ 202544706875SHong Zhang d = *(aa + ai[i]); /* diag[i] */ 2026d06b337dSHong Zhang v = aa + ai[i] + 1; v1=v; 2027d06b337dSHong Zhang vj = aj + ai[i] + 1; vj1=vj; 2028d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; nz1=nz; 2029d06b337dSHong Zhang sum = t[i]; 2030d06b337dSHong Zhang while (nz1--) sum -= (*v1++)*x[*vj1++]; 2031d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2032d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 203344706875SHong Zhang PetscLogFlops(4*nz-2); 2034d06b337dSHong Zhang } 2035d06b337dSHong Zhang } 2036d06b337dSHong Zhang 20372798e883SHong Zhang if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 203844706875SHong Zhang /* 203944706875SHong Zhang backward sweep: 204044706875SHong Zhang b = b - x[i]*U^T(i,:), i=0,...,n-2 204144706875SHong Zhang for i=m-1,...,0: 204244706875SHong Zhang sum[i] = (b[i] - U(i,:)x )/d[i]; 204344706875SHong Zhang x[i] = (1-omega)x[i] + omega*sum[i]; 204444706875SHong Zhang */ 2045d06b337dSHong Zhang for (i=0; i<m; i++) 2046d06b337dSHong Zhang t[i] = b[i]; 2047d06b337dSHong Zhang 2048d06b337dSHong Zhang for (i=0; i<m-1; i++){ /* update rhs */ 2049d06b337dSHong Zhang v = aa + ai[i] + 1; 2050d06b337dSHong Zhang vj = aj + ai[i] + 1; 2051d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2052d06b337dSHong Zhang while (nz--) t[*vj++] -= x[i]*(*v++); 205344706875SHong Zhang PetscLogFlops(2*nz-1); 2054d06b337dSHong Zhang } 2055d06b337dSHong Zhang for (i=m-1; i>=0; i--){ 2056d06b337dSHong Zhang d = *(aa + ai[i]); 2057d06b337dSHong Zhang v = aa + ai[i] + 1; 2058d06b337dSHong Zhang vj = aj + ai[i] + 1; 2059d06b337dSHong Zhang nz = ai[i+1] - ai[i] - 1; 2060d06b337dSHong Zhang sum = t[i]; 2061d06b337dSHong Zhang while (nz--) sum -= x[*vj++]*(*v++); 206244706875SHong Zhang PetscLogFlops(2*nz-1); 2063d06b337dSHong Zhang x[i] = (1-omega)*x[i] + omega*sum/d; 2064d06b337dSHong Zhang } 2065d06b337dSHong Zhang } 2066d06b337dSHong Zhang } 2067d06b337dSHong Zhang 2068d06b337dSHong Zhang ierr = PetscFree(t); CHKERRQ(ierr); 2069b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 2070d06b337dSHong Zhang if (bb != xx) { 2071b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(bb,&b);CHKERRQ(ierr); 2072d06b337dSHong Zhang } 20732798e883SHong Zhang 2074d06b337dSHong Zhang PetscFunctionReturn(0); 2075d06b337dSHong Zhang } 2076d06b337dSHong Zhang 2077d06b337dSHong Zhang 2078d06b337dSHong Zhang 2079d06b337dSHong Zhang 208049b5e25fSSatish Balay 208149b5e25fSSatish Balay 2082