1*7e8a2ce2SHong Zhang /*$Id: sbaij.c,v 1.35 2000/10/12 14:03:03 hzhang Exp hzhang $*/ 249b5e25fSSatish Balay 349b5e25fSSatish Balay /* 449b5e25fSSatish Balay Defines the basic matrix operations for the BAIJ (compressed row) 549b5e25fSSatish Balay matrix storage format. 649b5e25fSSatish Balay */ 7ab9f2c04SSatish Balay #include "petscsys.h" /*I "petscmat.h" I*/ 849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h" 949b5e25fSSatish Balay #include "src/vec/vecimpl.h" 1049b5e25fSSatish Balay #include "src/inline/spops.h" 11aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h" 1249b5e25fSSatish Balay 1349b5e25fSSatish Balay #define CHUNKSIZE 10 1449b5e25fSSatish Balay 1549b5e25fSSatish Balay /* 1649b5e25fSSatish Balay Checks for missing diagonals 1749b5e25fSSatish Balay */ 1849b5e25fSSatish Balay #undef __FUNC__ 1949b5e25fSSatish Balay #define __FUNC__ "MatMissingDiagonal_SeqSBAIJ" 2049b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A) 2149b5e25fSSatish Balay { 22045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 2349b5e25fSSatish Balay int *diag,*jj = a->j,i,ierr; 2449b5e25fSSatish Balay 2549b5e25fSSatish Balay PetscFunctionBegin; 26045c9aa0SHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr); 2749b5e25fSSatish Balay diag = a->diag; 2849b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 2949b5e25fSSatish Balay if (jj[diag[i]] != i) { 3029bbc08cSBarry Smith SETERRQ1(1,"Matrix is missing diagonal number %d",i); 3149b5e25fSSatish Balay } 3249b5e25fSSatish Balay } 3349b5e25fSSatish Balay PetscFunctionReturn(0); 3449b5e25fSSatish Balay } 3549b5e25fSSatish Balay 3649b5e25fSSatish Balay #undef __FUNC__ 3749b5e25fSSatish Balay #define __FUNC__ "MatMarkDiagonal_SeqSBAIJ" 3849b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A) 3949b5e25fSSatish Balay { 40045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 4149b5e25fSSatish Balay int i,j,*diag,m = a->mbs; 4249b5e25fSSatish Balay 4349b5e25fSSatish Balay PetscFunctionBegin; 4449b5e25fSSatish Balay if (a->diag) PetscFunctionReturn(0); 4549b5e25fSSatish Balay 4649b5e25fSSatish Balay diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag); 4749b5e25fSSatish Balay PLogObjectMemory(A,(m+1)*sizeof(int)); 4849b5e25fSSatish Balay for (i=0; i<m; i++) { 4949b5e25fSSatish Balay diag[i] = a->i[i+1]; 5049b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 5149b5e25fSSatish Balay if (a->j[j] == i) { 5249b5e25fSSatish Balay diag[i] = j; 5349b5e25fSSatish Balay break; 5449b5e25fSSatish Balay } 5549b5e25fSSatish Balay } 5649b5e25fSSatish Balay } 5749b5e25fSSatish Balay a->diag = diag; 5849b5e25fSSatish Balay PetscFunctionReturn(0); 5949b5e25fSSatish Balay } 6049b5e25fSSatish Balay 6149b5e25fSSatish Balay 6249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 6349b5e25fSSatish Balay 6449b5e25fSSatish Balay #undef __FUNC__ 6549b5e25fSSatish Balay #define __FUNC__ "MatGetRowIJ_SeqSBAIJ" 6649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 6749b5e25fSSatish Balay { 6849b5e25fSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6949b5e25fSSatish Balay 7049b5e25fSSatish Balay PetscFunctionBegin; 71045c9aa0SHong Zhang if (ia) { 7229bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 7349b5e25fSSatish Balay } 74045c9aa0SHong Zhang *nn = a->mbs; 7549b5e25fSSatish Balay PetscFunctionReturn(0); 7649b5e25fSSatish Balay } 7749b5e25fSSatish Balay 7849b5e25fSSatish Balay #undef __FUNC__ 7949b5e25fSSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqSBAIJ" 80045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 8149b5e25fSSatish Balay { 8249b5e25fSSatish Balay PetscFunctionBegin; 8349b5e25fSSatish Balay if (!ia) PetscFunctionReturn(0); 8429bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering"); 8516cdd363SHong Zhang /* PetscFunctionReturn(0); */ 8649b5e25fSSatish Balay } 8749b5e25fSSatish Balay 8849b5e25fSSatish Balay #undef __FUNC__ 8949b5e25fSSatish Balay #define __FUNC__ "MatGetBlockSize_SeqSBAIJ" 9049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs) 9149b5e25fSSatish Balay { 9249b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data; 9349b5e25fSSatish Balay 9449b5e25fSSatish Balay PetscFunctionBegin; 9549b5e25fSSatish Balay *bs = sbaij->bs; 9649b5e25fSSatish Balay PetscFunctionReturn(0); 9749b5e25fSSatish Balay } 9849b5e25fSSatish Balay 9949b5e25fSSatish Balay #undef __FUNC__ 10049b5e25fSSatish Balay #define __FUNC__ "MatDestroy_SeqSBAIJ" 10149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A) 10249b5e25fSSatish Balay { 10349b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 10449b5e25fSSatish Balay int ierr; 10549b5e25fSSatish Balay 10649b5e25fSSatish Balay PetscFunctionBegin; 10749b5e25fSSatish Balay if (A->mapping) { 10849b5e25fSSatish Balay ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 10949b5e25fSSatish Balay } 11049b5e25fSSatish Balay if (A->bmapping) { 11149b5e25fSSatish Balay ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 11249b5e25fSSatish Balay } 11349b5e25fSSatish Balay if (A->rmap) { 11449b5e25fSSatish Balay ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 11549b5e25fSSatish Balay } 11649b5e25fSSatish Balay if (A->cmap) { 11749b5e25fSSatish Balay ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 11849b5e25fSSatish Balay } 11949b5e25fSSatish Balay #if defined(PETSC_USE_LOG) 120c464158bSHong Zhang PLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz); 12149b5e25fSSatish Balay #endif 12249b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 12349b5e25fSSatish Balay if (!a->singlemalloc) { 12449b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 12563c8bf9fSHong Zhang ierr = PetscFree(a->j);CHKERRQ(ierr); 12649b5e25fSSatish Balay } 12749b5e25fSSatish Balay if (a->row) { 12849b5e25fSSatish Balay ierr = ISDestroy(a->row);CHKERRQ(ierr); 12949b5e25fSSatish Balay } 13049b5e25fSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 13149b5e25fSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 13249b5e25fSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 13349b5e25fSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 13449b5e25fSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 13549b5e25fSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 13649b5e25fSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 1371a3463dfSHong Zhang 1381a3463dfSHong Zhang if (a->inew){ 1391a3463dfSHong Zhang ierr = PetscFree(a->inew);CHKERRQ(ierr); 1401a3463dfSHong Zhang a->inew = 0; 1411a3463dfSHong Zhang } 14249b5e25fSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 14349b5e25fSSatish Balay PLogObjectDestroy(A); 14449b5e25fSSatish Balay PetscHeaderDestroy(A); 14549b5e25fSSatish Balay PetscFunctionReturn(0); 14649b5e25fSSatish Balay } 14749b5e25fSSatish Balay 14849b5e25fSSatish Balay #undef __FUNC__ 14949b5e25fSSatish Balay #define __FUNC__ "MatSetOption_SeqSBAIJ" 15049b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op) 15149b5e25fSSatish Balay { 152045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 15349b5e25fSSatish Balay 15449b5e25fSSatish Balay PetscFunctionBegin; 15549b5e25fSSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = PETSC_TRUE; 15649b5e25fSSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = PETSC_FALSE; 15749b5e25fSSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = PETSC_TRUE; 15849b5e25fSSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = PETSC_FALSE; 15949b5e25fSSatish Balay else if (op == MAT_KEEP_ZEROED_ROWS) a->keepzeroedrows = PETSC_TRUE; 16049b5e25fSSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 16149b5e25fSSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 16249b5e25fSSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 16349b5e25fSSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 16449b5e25fSSatish Balay else if (op == MAT_ROWS_SORTED || 16549b5e25fSSatish Balay op == MAT_ROWS_UNSORTED || 16649b5e25fSSatish Balay op == MAT_SYMMETRIC || 16749b5e25fSSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 16849b5e25fSSatish Balay op == MAT_YES_NEW_DIAGONALS || 16949b5e25fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 17049b5e25fSSatish Balay op == MAT_USE_HASH_TABLE) { 171045c9aa0SHong Zhang PLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n"); 17249b5e25fSSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 17329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 17449b5e25fSSatish Balay } else { 17529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 17649b5e25fSSatish Balay } 17749b5e25fSSatish Balay PetscFunctionReturn(0); 17849b5e25fSSatish Balay } 17949b5e25fSSatish Balay 18049b5e25fSSatish Balay #undef __FUNC__ 18149b5e25fSSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqSBAIJ" 18249b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n) 18349b5e25fSSatish Balay { 18449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 18549b5e25fSSatish Balay 18649b5e25fSSatish Balay PetscFunctionBegin; 187c464158bSHong Zhang *m = 0; 188*7e8a2ce2SHong Zhang *n = A->m; 18949b5e25fSSatish Balay PetscFunctionReturn(0); 19049b5e25fSSatish Balay } 19149b5e25fSSatish Balay 19249b5e25fSSatish Balay #undef __FUNC__ 19349b5e25fSSatish Balay #define __FUNC__ "MatGetRow_SeqSBAIJ" 19449b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v) 19549b5e25fSSatish Balay { 19649b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 19749b5e25fSSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2; 19849b5e25fSSatish Balay MatScalar *aa,*aa_i; 19949b5e25fSSatish Balay Scalar *v_i; 20049b5e25fSSatish Balay 20149b5e25fSSatish Balay PetscFunctionBegin; 20249b5e25fSSatish Balay bs = a->bs; 20349b5e25fSSatish Balay ai = a->i; 20449b5e25fSSatish Balay aj = a->j; 20549b5e25fSSatish Balay aa = a->a; 20649b5e25fSSatish Balay bs2 = a->bs2; 20749b5e25fSSatish Balay 208c464158bSHong Zhang if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 20949b5e25fSSatish Balay 21049b5e25fSSatish Balay bn = row/bs; /* Block number */ 21149b5e25fSSatish Balay bp = row % bs; /* Block position */ 21249b5e25fSSatish Balay M = ai[bn+1] - ai[bn]; 21349b5e25fSSatish Balay *ncols = bs*M; 21449b5e25fSSatish Balay 21549b5e25fSSatish Balay if (v) { 21649b5e25fSSatish Balay *v = 0; 21749b5e25fSSatish Balay if (*ncols) { 21849b5e25fSSatish Balay *v = (Scalar*)PetscMalloc((*ncols+row)*sizeof(Scalar));CHKPTRQ(*v); 21949b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 22049b5e25fSSatish Balay v_i = *v + i*bs; 22149b5e25fSSatish Balay aa_i = aa + bs2*(ai[bn] + i); 22249b5e25fSSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 22349b5e25fSSatish Balay } 22449b5e25fSSatish Balay } 22549b5e25fSSatish Balay } 22649b5e25fSSatish Balay 22749b5e25fSSatish Balay if (cols) { 22849b5e25fSSatish Balay *cols = 0; 22949b5e25fSSatish Balay if (*ncols) { 23049b5e25fSSatish Balay *cols = (int*)PetscMalloc((*ncols+row)*sizeof(int));CHKPTRQ(*cols); 23149b5e25fSSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 23249b5e25fSSatish Balay cols_i = *cols + i*bs; 23349b5e25fSSatish Balay itmp = bs*aj[ai[bn] + i]; 23449b5e25fSSatish Balay for (j=0; j<bs; j++) {cols_i[j] = itmp++;} 23549b5e25fSSatish Balay } 23649b5e25fSSatish Balay } 23749b5e25fSSatish Balay } 23849b5e25fSSatish Balay 23949b5e25fSSatish Balay /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */ 2405ddb2528SHong Zhang /* this segment is currently removed, so only entries in the upper triangle are obtained */ 2415ddb2528SHong Zhang #ifdef column_search 24249b5e25fSSatish Balay v_i = *v + M*bs; 24349b5e25fSSatish Balay cols_i = *cols + M*bs; 24449b5e25fSSatish Balay for (i=0; i<bn; i++){ /* for each block row */ 24549b5e25fSSatish Balay M = ai[i+1] - ai[i]; 24649b5e25fSSatish Balay for (j=0; j<M; j++){ 24749b5e25fSSatish Balay itmp = aj[ai[i] + j]; /* block column value */ 24849b5e25fSSatish Balay if (itmp == bn){ 24949b5e25fSSatish Balay aa_i = aa + bs2*(ai[i] + j) + bs*bp; 25049b5e25fSSatish Balay for (k=0; k<bs; k++) { 25149b5e25fSSatish Balay *cols_i++ = i*bs+k; 25249b5e25fSSatish Balay *v_i++ = aa_i[k]; 25349b5e25fSSatish Balay } 25449b5e25fSSatish Balay *ncols += bs; 25549b5e25fSSatish Balay break; 25649b5e25fSSatish Balay } 25749b5e25fSSatish Balay } 25849b5e25fSSatish Balay } 2595ddb2528SHong Zhang #endif 26049b5e25fSSatish Balay 26149b5e25fSSatish Balay PetscFunctionReturn(0); 26249b5e25fSSatish Balay } 26349b5e25fSSatish Balay 26449b5e25fSSatish Balay #undef __FUNC__ 26549b5e25fSSatish Balay #define __FUNC__ "MatRestoreRow_SeqSBAIJ" 26649b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 26749b5e25fSSatish Balay { 26849b5e25fSSatish Balay int ierr; 26949b5e25fSSatish Balay 27049b5e25fSSatish Balay PetscFunctionBegin; 27149b5e25fSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 27249b5e25fSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 27349b5e25fSSatish Balay PetscFunctionReturn(0); 27449b5e25fSSatish Balay } 27549b5e25fSSatish Balay 27649b5e25fSSatish Balay #undef __FUNC__ 27749b5e25fSSatish Balay #define __FUNC__ "MatTranspose_SeqSBAIJ" 27849b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B) 27949b5e25fSSatish Balay { 28049b5e25fSSatish Balay PetscFunctionBegin; 28129bbc08cSBarry Smith SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called"); 28216cdd363SHong Zhang /* PetscFunctionReturn(0); */ 28349b5e25fSSatish Balay } 28449b5e25fSSatish Balay 28549b5e25fSSatish Balay #undef __FUNC__ 28649b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Binary" 28749b5e25fSSatish Balay static int MatView_SeqSBAIJ_Binary(Mat A,Viewer viewer) 28849b5e25fSSatish Balay { 28949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 29049b5e25fSSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 29149b5e25fSSatish Balay Scalar *aa; 29249b5e25fSSatish Balay FILE *file; 29349b5e25fSSatish Balay 29449b5e25fSSatish Balay PetscFunctionBegin; 29549b5e25fSSatish Balay ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 296c464158bSHong Zhang col_lens = (int*)PetscMalloc((4+A->m)*sizeof(int));CHKPTRQ(col_lens); 29749b5e25fSSatish Balay col_lens[0] = MAT_COOKIE; 29849b5e25fSSatish Balay 299c464158bSHong Zhang col_lens[1] = A->m; 300c464158bSHong Zhang col_lens[2] = A->m; 30149b5e25fSSatish Balay col_lens[3] = a->s_nz*bs2; 30249b5e25fSSatish Balay 30349b5e25fSSatish Balay /* store lengths of each row and write (including header) to file */ 30449b5e25fSSatish Balay count = 0; 30549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 30649b5e25fSSatish Balay for (j=0; j<bs; j++) { 30749b5e25fSSatish Balay col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 30849b5e25fSSatish Balay } 30949b5e25fSSatish Balay } 310c464158bSHong Zhang ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 31149b5e25fSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 31249b5e25fSSatish Balay 31349b5e25fSSatish Balay /* store column indices (zero start index) */ 31449b5e25fSSatish Balay jj = (int*)PetscMalloc((a->s_nz+1)*bs2*sizeof(int));CHKPTRQ(jj); 31549b5e25fSSatish Balay count = 0; 31649b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 31749b5e25fSSatish Balay for (j=0; j<bs; j++) { 31849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 31949b5e25fSSatish Balay for (l=0; l<bs; l++) { 32049b5e25fSSatish Balay jj[count++] = bs*a->j[k] + l; 32149b5e25fSSatish Balay } 32249b5e25fSSatish Balay } 32349b5e25fSSatish Balay } 32449b5e25fSSatish Balay } 32549b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr); 32649b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 32749b5e25fSSatish Balay 32849b5e25fSSatish Balay /* store nonzero values */ 32949b5e25fSSatish Balay aa = (Scalar*)PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 33049b5e25fSSatish Balay count = 0; 33149b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 33249b5e25fSSatish Balay for (j=0; j<bs; j++) { 33349b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 33449b5e25fSSatish Balay for (l=0; l<bs; l++) { 33549b5e25fSSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 33649b5e25fSSatish Balay } 33749b5e25fSSatish Balay } 33849b5e25fSSatish Balay } 33949b5e25fSSatish Balay } 34049b5e25fSSatish Balay ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr); 34149b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 34249b5e25fSSatish Balay 34349b5e25fSSatish Balay ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 34449b5e25fSSatish Balay if (file) { 34549b5e25fSSatish Balay fprintf(file,"-matload_block_size %d\n",a->bs); 34649b5e25fSSatish Balay } 34749b5e25fSSatish Balay PetscFunctionReturn(0); 34849b5e25fSSatish Balay } 34949b5e25fSSatish Balay 35049b5e25fSSatish Balay #undef __FUNC__ 35149b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_ASCII" 35249b5e25fSSatish Balay static int MatView_SeqSBAIJ_ASCII(Mat A,Viewer viewer) 35349b5e25fSSatish Balay { 35449b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 35549b5e25fSSatish Balay int ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2; 35649b5e25fSSatish Balay char *outputname; 35749b5e25fSSatish Balay 35849b5e25fSSatish Balay PetscFunctionBegin; 35949b5e25fSSatish Balay ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 36049b5e25fSSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 36149b5e25fSSatish Balay if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 36249b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 36349b5e25fSSatish Balay } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 36429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 36549b5e25fSSatish Balay } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 36649b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 36749b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 36849b5e25fSSatish Balay for (j=0; j<bs; j++) { 36949b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 37049b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 37149b5e25fSSatish Balay for (l=0; l<bs; l++) { 37249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 37349b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37449b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 37549b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37649b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 37749b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 37849b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 37949b5e25fSSatish Balay } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 38049b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 38149b5e25fSSatish Balay } 38249b5e25fSSatish Balay #else 38349b5e25fSSatish Balay if (a->a[bs2*k + l*bs + j] != 0.0) { 38449b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 38549b5e25fSSatish Balay } 38649b5e25fSSatish Balay #endif 38749b5e25fSSatish Balay } 38849b5e25fSSatish Balay } 38949b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 39049b5e25fSSatish Balay } 39149b5e25fSSatish Balay } 39249b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 39349b5e25fSSatish Balay } else { 39449b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 39549b5e25fSSatish Balay for (i=0; i<a->mbs; i++) { 39649b5e25fSSatish Balay for (j=0; j<bs; j++) { 39749b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 39849b5e25fSSatish Balay for (k=a->i[i]; k<a->i[i+1]; k++) { 39949b5e25fSSatish Balay for (l=0; l<bs; l++) { 40049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX) 40149b5e25fSSatish Balay if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 40249b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 40349b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40449b5e25fSSatish Balay } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 40549b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 40649b5e25fSSatish Balay PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40749b5e25fSSatish Balay } else { 40849b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 40949b5e25fSSatish Balay } 41049b5e25fSSatish Balay #else 41149b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 41249b5e25fSSatish Balay #endif 41349b5e25fSSatish Balay } 41449b5e25fSSatish Balay } 41549b5e25fSSatish Balay ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 41649b5e25fSSatish Balay } 41749b5e25fSSatish Balay } 41849b5e25fSSatish Balay ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 41949b5e25fSSatish Balay } 42049b5e25fSSatish Balay ierr = ViewerFlush(viewer);CHKERRQ(ierr); 42149b5e25fSSatish Balay PetscFunctionReturn(0); 42249b5e25fSSatish Balay } 42349b5e25fSSatish Balay 42449b5e25fSSatish Balay #undef __FUNC__ 42549b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw_Zoom" 42649b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw_Zoom(Draw draw,void *Aa) 42749b5e25fSSatish Balay { 42849b5e25fSSatish Balay Mat A = (Mat) Aa; 42949b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 43049b5e25fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 43149b5e25fSSatish Balay PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 43249b5e25fSSatish Balay MatScalar *aa; 43349b5e25fSSatish Balay MPI_Comm comm; 43449b5e25fSSatish Balay Viewer viewer; 43549b5e25fSSatish Balay 43649b5e25fSSatish Balay PetscFunctionBegin; 43749b5e25fSSatish Balay /* 43849b5e25fSSatish Balay This is nasty. If this is called from an originally parallel matrix 43949b5e25fSSatish Balay then all processes call this,but only the first has the matrix so the 44049b5e25fSSatish Balay rest should return immediately. 44149b5e25fSSatish Balay */ 44249b5e25fSSatish Balay ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 44349b5e25fSSatish Balay ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 44449b5e25fSSatish Balay if (rank) PetscFunctionReturn(0); 44549b5e25fSSatish Balay 44649b5e25fSSatish Balay ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 44749b5e25fSSatish Balay 44849b5e25fSSatish Balay ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 44949b5e25fSSatish Balay DrawString(draw, .3*(xl+xr), .3*(yl+yr), DRAW_BLACK, "symmetric"); 45049b5e25fSSatish Balay 45149b5e25fSSatish Balay /* loop over matrix elements drawing boxes */ 45249b5e25fSSatish Balay color = DRAW_BLUE; 45349b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 45449b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 455c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 45649b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 45749b5e25fSSatish Balay aa = a->a + j*bs2; 45849b5e25fSSatish Balay for (k=0; k<bs; k++) { 45949b5e25fSSatish Balay for (l=0; l<bs; l++) { 46049b5e25fSSatish Balay if (PetscRealPart(*aa++) >= 0.) continue; 46149b5e25fSSatish Balay ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 46249b5e25fSSatish Balay } 46349b5e25fSSatish Balay } 46449b5e25fSSatish Balay } 46549b5e25fSSatish Balay } 46649b5e25fSSatish Balay color = DRAW_CYAN; 46749b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 46849b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 469c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 47049b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 47149b5e25fSSatish Balay aa = a->a + j*bs2; 47249b5e25fSSatish Balay for (k=0; k<bs; k++) { 47349b5e25fSSatish Balay for (l=0; l<bs; l++) { 47449b5e25fSSatish Balay if (PetscRealPart(*aa++) != 0.) continue; 47549b5e25fSSatish Balay ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 47649b5e25fSSatish Balay } 47749b5e25fSSatish Balay } 47849b5e25fSSatish Balay } 47949b5e25fSSatish Balay } 48049b5e25fSSatish Balay 48149b5e25fSSatish Balay color = DRAW_RED; 48249b5e25fSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 48349b5e25fSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 484c464158bSHong Zhang y_l = A->m - row - 1.0; y_r = y_l + 1.0; 48549b5e25fSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 48649b5e25fSSatish Balay aa = a->a + j*bs2; 48749b5e25fSSatish Balay for (k=0; k<bs; k++) { 48849b5e25fSSatish Balay for (l=0; l<bs; l++) { 48949b5e25fSSatish Balay if (PetscRealPart(*aa++) <= 0.) continue; 49049b5e25fSSatish Balay ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 49149b5e25fSSatish Balay } 49249b5e25fSSatish Balay } 49349b5e25fSSatish Balay } 49449b5e25fSSatish Balay } 49549b5e25fSSatish Balay PetscFunctionReturn(0); 49649b5e25fSSatish Balay } 49749b5e25fSSatish Balay 49849b5e25fSSatish Balay #undef __FUNC__ 49949b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw" 50049b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw(Mat A,Viewer viewer) 50149b5e25fSSatish Balay { 50249b5e25fSSatish Balay Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 50349b5e25fSSatish Balay int ierr; 50449b5e25fSSatish Balay PetscReal xl,yl,xr,yr,w,h; 50549b5e25fSSatish Balay Draw draw; 50649b5e25fSSatish Balay PetscTruth isnull; 50749b5e25fSSatish Balay 50849b5e25fSSatish Balay PetscFunctionBegin; 50949b5e25fSSatish Balay 51049b5e25fSSatish Balay ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 51149b5e25fSSatish Balay ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 51249b5e25fSSatish Balay 51349b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 514c464158bSHong Zhang xr = A->m; yr = A->m; h = yr/10.0; w = xr/10.0; 51549b5e25fSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 51649b5e25fSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 51749b5e25fSSatish Balay ierr = DrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 51849b5e25fSSatish Balay ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 51949b5e25fSSatish Balay PetscFunctionReturn(0); 52049b5e25fSSatish Balay } 52149b5e25fSSatish Balay 52249b5e25fSSatish Balay #undef __FUNC__ 52349b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ" 52449b5e25fSSatish Balay int MatView_SeqSBAIJ(Mat A,Viewer viewer) 52549b5e25fSSatish Balay { 52649b5e25fSSatish Balay int ierr; 52749b5e25fSSatish Balay PetscTruth issocket,isascii,isbinary,isdraw; 52849b5e25fSSatish Balay 52949b5e25fSSatish Balay PetscFunctionBegin; 53049b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 53149b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 53249b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 53349b5e25fSSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 53449b5e25fSSatish Balay if (issocket) { 53529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 53649b5e25fSSatish Balay } else if (isascii){ 53749b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 53849b5e25fSSatish Balay } else if (isbinary) { 53949b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr); 54049b5e25fSSatish Balay } else if (isdraw) { 54149b5e25fSSatish Balay ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr); 54249b5e25fSSatish Balay } else { 54329bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 54449b5e25fSSatish Balay } 54549b5e25fSSatish Balay PetscFunctionReturn(0); 54649b5e25fSSatish Balay } 54749b5e25fSSatish Balay 54849b5e25fSSatish Balay 54949b5e25fSSatish Balay #undef __FUNC__ 55049b5e25fSSatish Balay #define __FUNC__ "MatGetValues_SeqSBAIJ" 55149b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 55249b5e25fSSatish Balay { 553045c9aa0SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 55449b5e25fSSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 55549b5e25fSSatish Balay int *ai = a->i,*ailen = a->ilen; 55649b5e25fSSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 55749b5e25fSSatish Balay MatScalar *ap,*aa = a->a,zero = 0.0; 55849b5e25fSSatish Balay 55949b5e25fSSatish Balay PetscFunctionBegin; 56049b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 56149b5e25fSSatish Balay row = im[k]; brow = row/bs; 56229bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 563c464158bSHong Zhang if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 56449b5e25fSSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 56549b5e25fSSatish Balay nrow = ailen[brow]; 56649b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 56729bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 568c464158bSHong Zhang if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 56949b5e25fSSatish Balay col = in[l] ; 57049b5e25fSSatish Balay bcol = col/bs; 57149b5e25fSSatish Balay cidx = col%bs; 57249b5e25fSSatish Balay ridx = row%bs; 57349b5e25fSSatish Balay high = nrow; 57449b5e25fSSatish Balay low = 0; /* assume unsorted */ 57549b5e25fSSatish Balay while (high-low > 5) { 57649b5e25fSSatish Balay t = (low+high)/2; 57749b5e25fSSatish Balay if (rp[t] > bcol) high = t; 57849b5e25fSSatish Balay else low = t; 57949b5e25fSSatish Balay } 58049b5e25fSSatish Balay for (i=low; i<high; i++) { 58149b5e25fSSatish Balay if (rp[i] > bcol) break; 58249b5e25fSSatish Balay if (rp[i] == bcol) { 58349b5e25fSSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 58449b5e25fSSatish Balay goto finished; 58549b5e25fSSatish Balay } 58649b5e25fSSatish Balay } 58749b5e25fSSatish Balay *v++ = zero; 58849b5e25fSSatish Balay finished:; 58949b5e25fSSatish Balay } 59049b5e25fSSatish Balay } 59149b5e25fSSatish Balay PetscFunctionReturn(0); 59249b5e25fSSatish Balay } 59349b5e25fSSatish Balay 59449b5e25fSSatish Balay 59549b5e25fSSatish Balay #undef __FUNC__ 59649b5e25fSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqSBAIJ" 59749b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 59849b5e25fSSatish Balay { 59949b5e25fSSatish Balay PetscFunctionBegin; 60029bbc08cSBarry Smith SETERRQ(1,"Function not yet written for SBAIJ format"); 60116cdd363SHong Zhang /* PetscFunctionReturn(0); */ 60249b5e25fSSatish Balay } 60349b5e25fSSatish Balay 60449b5e25fSSatish Balay #undef __FUNC__ 60549b5e25fSSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqSBAIJ" 60649b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode) 60749b5e25fSSatish Balay { 60849b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 60949b5e25fSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 610c464158bSHong Zhang int m = A->m,*ip,N,*ailen = a->ilen; 61149b5e25fSSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 61249b5e25fSSatish Balay MatScalar *aa = a->a,*ap; 61349b5e25fSSatish Balay 61449b5e25fSSatish Balay PetscFunctionBegin; 61549b5e25fSSatish Balay if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 61649b5e25fSSatish Balay 61749b5e25fSSatish Balay if (m) rmax = ailen[0]; 61849b5e25fSSatish Balay for (i=1; i<mbs; i++) { 61949b5e25fSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 62049b5e25fSSatish Balay fshift += imax[i-1] - ailen[i-1]; 62149b5e25fSSatish Balay rmax = PetscMax(rmax,ailen[i]); 62249b5e25fSSatish Balay if (fshift) { 62349b5e25fSSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 62449b5e25fSSatish Balay N = ailen[i]; 62549b5e25fSSatish Balay for (j=0; j<N; j++) { 62649b5e25fSSatish Balay ip[j-fshift] = ip[j]; 62749b5e25fSSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 62849b5e25fSSatish Balay } 62949b5e25fSSatish Balay } 63049b5e25fSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 63149b5e25fSSatish Balay } 63249b5e25fSSatish Balay if (mbs) { 63349b5e25fSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 63449b5e25fSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 63549b5e25fSSatish Balay } 63649b5e25fSSatish Balay /* reset ilen and imax for each row */ 63749b5e25fSSatish Balay for (i=0; i<mbs; i++) { 63849b5e25fSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 63949b5e25fSSatish Balay } 64049b5e25fSSatish Balay a->s_nz = ai[mbs]; 64149b5e25fSSatish Balay 64249b5e25fSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 64349b5e25fSSatish Balay if (fshift && a->diag) { 64449b5e25fSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 64549b5e25fSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 64649b5e25fSSatish Balay a->diag = 0; 64749b5e25fSSatish Balay } 64849b5e25fSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 649c464158bSHong Zhang m,A->m,a->bs,fshift*bs2,a->s_nz*bs2); 65049b5e25fSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n", 65149b5e25fSSatish Balay a->reallocs); 65249b5e25fSSatish Balay PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 65349b5e25fSSatish Balay a->reallocs = 0; 65449b5e25fSSatish Balay A->info.nz_unneeded = (PetscReal)fshift*bs2; 65549b5e25fSSatish Balay 65649b5e25fSSatish Balay PetscFunctionReturn(0); 65749b5e25fSSatish Balay } 65849b5e25fSSatish Balay 65949b5e25fSSatish Balay /* 66049b5e25fSSatish Balay This function returns an array of flags which indicate the locations of contiguous 66149b5e25fSSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 66249b5e25fSSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 66349b5e25fSSatish Balay Assume: sizes should be long enough to hold all the values. 66449b5e25fSSatish Balay */ 66549b5e25fSSatish Balay #undef __FUNC__ 66649b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ_Check_Blocks" 66749b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 66849b5e25fSSatish Balay { 66949b5e25fSSatish Balay int i,j,k,row; 67049b5e25fSSatish Balay PetscTruth flg; 67149b5e25fSSatish Balay 67249b5e25fSSatish Balay PetscFunctionBegin; 67349b5e25fSSatish Balay for (i=0,j=0; i<n; j++) { 67449b5e25fSSatish Balay row = idx[i]; 67549b5e25fSSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 67649b5e25fSSatish Balay sizes[j] = 1; 67749b5e25fSSatish Balay i++; 67849b5e25fSSatish Balay } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */ 67949b5e25fSSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 68049b5e25fSSatish Balay i++; 68149b5e25fSSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 68249b5e25fSSatish Balay flg = PETSC_TRUE; 68349b5e25fSSatish Balay for (k=1; k<bs; k++) { 68449b5e25fSSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 68549b5e25fSSatish Balay flg = PETSC_FALSE; 68649b5e25fSSatish Balay break; 68749b5e25fSSatish Balay } 68849b5e25fSSatish Balay } 68949b5e25fSSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 69049b5e25fSSatish Balay sizes[j] = bs; 69149b5e25fSSatish Balay i+= bs; 69249b5e25fSSatish Balay } else { 69349b5e25fSSatish Balay sizes[j] = 1; 69449b5e25fSSatish Balay i++; 69549b5e25fSSatish Balay } 69649b5e25fSSatish Balay } 69749b5e25fSSatish Balay } 69849b5e25fSSatish Balay *bs_max = j; 69949b5e25fSSatish Balay PetscFunctionReturn(0); 70049b5e25fSSatish Balay } 70149b5e25fSSatish Balay 70249b5e25fSSatish Balay #undef __FUNC__ 70349b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ" 70449b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag) 70549b5e25fSSatish Balay { 70649b5e25fSSatish Balay Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data; 70749b5e25fSSatish Balay int ierr,i,j,k,count,is_n,*is_idx,*rows; 70849b5e25fSSatish Balay int bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max; 70949b5e25fSSatish Balay Scalar zero = 0.0; 71049b5e25fSSatish Balay MatScalar *aa; 71149b5e25fSSatish Balay 71249b5e25fSSatish Balay PetscFunctionBegin; 71349b5e25fSSatish Balay /* Make a copy of the IS and sort it */ 71449b5e25fSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 71549b5e25fSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 71649b5e25fSSatish Balay 71749b5e25fSSatish Balay /* allocate memory for rows,sizes */ 71849b5e25fSSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 71949b5e25fSSatish Balay sizes = rows + is_n; 72049b5e25fSSatish Balay 72149b5e25fSSatish Balay /* initialize copy IS values to rows, and sort them */ 72249b5e25fSSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 72349b5e25fSSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 72449b5e25fSSatish Balay if (sbaij->keepzeroedrows) { /* do not change nonzero structure */ 72549b5e25fSSatish Balay for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */ 72649b5e25fSSatish Balay bs_max = is_n; /* bs_max: num. of contiguous block row in the row */ 72749b5e25fSSatish Balay } else { 72849b5e25fSSatish Balay ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 72949b5e25fSSatish Balay } 73049b5e25fSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 73149b5e25fSSatish Balay 73249b5e25fSSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 73349b5e25fSSatish Balay row = rows[j]; /* row to be zeroed */ 734c464158bSHong Zhang if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 73549b5e25fSSatish Balay count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */ 73649b5e25fSSatish Balay aa = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs); 73749b5e25fSSatish Balay if (sizes[i] == bs && !sbaij->keepzeroedrows) { 73849b5e25fSSatish Balay if (diag) { 73949b5e25fSSatish Balay if (sbaij->ilen[row/bs] > 0) { 74049b5e25fSSatish Balay sbaij->ilen[row/bs] = 1; 74149b5e25fSSatish Balay sbaij->j[sbaij->i[row/bs]] = row/bs; 74249b5e25fSSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 74349b5e25fSSatish Balay } 74449b5e25fSSatish Balay /* Now insert all the diagoanl values for this bs */ 74549b5e25fSSatish Balay for (k=0; k<bs; k++) { 74649b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 74749b5e25fSSatish Balay } 74849b5e25fSSatish Balay } else { /* (!diag) */ 74949b5e25fSSatish Balay sbaij->ilen[row/bs] = 0; 75049b5e25fSSatish Balay } /* end (!diag) */ 75149b5e25fSSatish Balay } else { /* (sizes[i] != bs), broken block */ 75249b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG) 75329bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 75449b5e25fSSatish Balay #endif 75549b5e25fSSatish Balay for (k=0; k<count; k++) { 75649b5e25fSSatish Balay aa[0] = zero; 75749b5e25fSSatish Balay aa+=bs; 75849b5e25fSSatish Balay } 75949b5e25fSSatish Balay if (diag) { 76049b5e25fSSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 76149b5e25fSSatish Balay } 76249b5e25fSSatish Balay } 76349b5e25fSSatish Balay } 76449b5e25fSSatish Balay 76549b5e25fSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 76649b5e25fSSatish Balay ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76749b5e25fSSatish Balay ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76849b5e25fSSatish Balay PetscFunctionReturn(0); 76949b5e25fSSatish Balay } 77049b5e25fSSatish Balay 77149b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks). 77249b5e25fSSatish Balay Any a(i,j) with i>j input by user is ingored. 77349b5e25fSSatish Balay */ 77449b5e25fSSatish Balay 77549b5e25fSSatish Balay #undef __FUNC__ 77649b5e25fSSatish Balay #define __FUNC__ "MatSetValues_SeqSBAIJ" 77749b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 77849b5e25fSSatish Balay { 77949b5e25fSSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 78049b5e25fSSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 78149b5e25fSSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 78249b5e25fSSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 78349b5e25fSSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 78449b5e25fSSatish Balay MatScalar *ap,value,*aa=a->a,*bap; 78549b5e25fSSatish Balay 78649b5e25fSSatish Balay PetscFunctionBegin; 78749b5e25fSSatish Balay 78849b5e25fSSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 78949b5e25fSSatish Balay row = im[k]; /* row number */ 79049b5e25fSSatish Balay brow = row/bs; /* block row number */ 79149b5e25fSSatish Balay if (row < 0) continue; 79249b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 793c464158bSHong Zhang if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 79449b5e25fSSatish Balay #endif 79549b5e25fSSatish Balay rp = aj + ai[brow]; /*ptr to beginning of column value of the row block*/ 79649b5e25fSSatish Balay ap = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/ 79749b5e25fSSatish Balay rmax = imax[brow]; /* maximum space allocated for this row */ 79849b5e25fSSatish Balay nrow = ailen[brow]; /* actual length of this row */ 79949b5e25fSSatish Balay low = 0; 80049b5e25fSSatish Balay 80149b5e25fSSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 80249b5e25fSSatish Balay if (in[l] < 0) continue; 80349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g) 804c464158bSHong Zhang if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m); 80549b5e25fSSatish Balay #endif 80649b5e25fSSatish Balay col = in[l]; 80749b5e25fSSatish Balay bcol = col/bs; /* block col number */ 80849b5e25fSSatish Balay 80949b5e25fSSatish Balay if (brow <= bcol){ 81049b5e25fSSatish Balay ridx = row % bs; cidx = col % bs; /*row and col index inside the block */ 81149b5e25fSSatish Balay /* if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ */ 81249b5e25fSSatish Balay /* element value a(k,l) */ 81349b5e25fSSatish Balay if (roworiented) { 81449b5e25fSSatish Balay value = v[l + k*n]; 81549b5e25fSSatish Balay } else { 81649b5e25fSSatish Balay value = v[k + l*m]; 81749b5e25fSSatish Balay } 81849b5e25fSSatish Balay 81949b5e25fSSatish Balay /* move pointer bap to a(k,l) quickly and add/insert value */ 82049b5e25fSSatish Balay if (!sorted) low = 0; high = nrow; 82149b5e25fSSatish Balay while (high-low > 7) { 82249b5e25fSSatish Balay t = (low+high)/2; 82349b5e25fSSatish Balay if (rp[t] > bcol) high = t; 82449b5e25fSSatish Balay else low = t; 82549b5e25fSSatish Balay } 82649b5e25fSSatish Balay for (i=low; i<high; i++) { 82749b5e25fSSatish Balay /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */ 82849b5e25fSSatish Balay if (rp[i] > bcol) break; 82949b5e25fSSatish Balay if (rp[i] == bcol) { 83049b5e25fSSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 83149b5e25fSSatish Balay if (is == ADD_VALUES) *bap += value; 83249b5e25fSSatish Balay else *bap = value; 83349b5e25fSSatish Balay goto noinsert1; 83449b5e25fSSatish Balay } 83549b5e25fSSatish Balay } 83649b5e25fSSatish Balay 83749b5e25fSSatish Balay if (nonew == 1) goto noinsert1; 83829bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 83949b5e25fSSatish Balay if (nrow >= rmax) { 84049b5e25fSSatish Balay /* there is no extra room in row, therefore enlarge */ 84149b5e25fSSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 84249b5e25fSSatish Balay MatScalar *new_a; 84349b5e25fSSatish Balay 84429bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 84549b5e25fSSatish Balay 84649b5e25fSSatish Balay /* Malloc new storage space */ 84749b5e25fSSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 84849b5e25fSSatish Balay new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 84949b5e25fSSatish Balay new_j = (int*)(new_a + bs2*new_nz); 85049b5e25fSSatish Balay new_i = new_j + new_nz; 85149b5e25fSSatish Balay 85249b5e25fSSatish Balay /* copy over old data into new slots */ 85349b5e25fSSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 85449b5e25fSSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 85549b5e25fSSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 85649b5e25fSSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 85749b5e25fSSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 85849b5e25fSSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 85949b5e25fSSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 86049b5e25fSSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 86149b5e25fSSatish Balay /* free up old matrix storage */ 86249b5e25fSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 86349b5e25fSSatish Balay if (!a->singlemalloc) { 86449b5e25fSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 86549b5e25fSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 86649b5e25fSSatish Balay } 86749b5e25fSSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 86849b5e25fSSatish Balay a->singlemalloc = PETSC_TRUE; 86949b5e25fSSatish Balay 87049b5e25fSSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 87149b5e25fSSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 87249b5e25fSSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 87349b5e25fSSatish Balay a->s_maxnz += bs2*CHUNKSIZE; 87449b5e25fSSatish Balay a->reallocs++; 87549b5e25fSSatish Balay a->s_nz++; 87649b5e25fSSatish Balay } 87749b5e25fSSatish Balay 87849b5e25fSSatish Balay N = nrow++ - 1; 87949b5e25fSSatish Balay /* shift up all the later entries in this row */ 88049b5e25fSSatish Balay for (ii=N; ii>=i; ii--) { 88149b5e25fSSatish Balay rp[ii+1] = rp[ii]; 88249b5e25fSSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 88349b5e25fSSatish Balay } 88449b5e25fSSatish Balay if (N>=i) { 88549b5e25fSSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 88649b5e25fSSatish Balay } 88749b5e25fSSatish Balay rp[i] = bcol; 88849b5e25fSSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 88949b5e25fSSatish Balay noinsert1:; 89049b5e25fSSatish Balay low = i; 89149b5e25fSSatish Balay /* } */ 89249b5e25fSSatish Balay } /* end of if .. if.. */ 89349b5e25fSSatish Balay } /* end of loop over added columns */ 89449b5e25fSSatish Balay ailen[brow] = nrow; 89549b5e25fSSatish Balay } /* end of loop over added rows */ 89649b5e25fSSatish Balay 89749b5e25fSSatish Balay PetscFunctionReturn(0); 89849b5e25fSSatish Balay } 89949b5e25fSSatish Balay 900915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*); 901915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal); 90249b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int); 90349b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 90449b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 90549b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec); 90649b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec); 90749b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat); 90849b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *); 90949b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*); 91049b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec); 91149b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec); 91249b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *); 91349b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat); 91413a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec); 91549b5e25fSSatish Balay 91649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec); 91749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec); 91849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec); 91949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec); 92049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec); 92149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec); 92249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec); 92349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec); 92449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec); 92549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec); 92649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec); 92749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec); 92849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec); 92949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec); 93049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec); 93149b5e25fSSatish Balay 9325f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*); 9335f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*); 9345f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*); 9355f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*); 9365f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*); 9375f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*); 9385f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*); 93949b5e25fSSatish Balay 94049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec); 94149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec); 94249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec); 94349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec); 94449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec); 94549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec); 94649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec); 94749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec); 94849b5e25fSSatish Balay 94949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec); 95049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec); 95149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec); 95249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec); 95349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec); 95449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec); 95549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec); 95649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec); 95749b5e25fSSatish Balay 9581a3463dfSHong Zhang #ifdef HAVE_MatIncompleteCholeskyFactor 9591a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work! */ 96049b5e25fSSatish Balay #undef __FUNC__ 9614ccecd49SHong Zhang #define __FUNC__ "MatIncompleteCholeskyFactor_SeqSBAIJ" 9621a3463dfSHong Zhang int MatIncompleteCholeskyFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level) 96349b5e25fSSatish Balay { 9644ccecd49SHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data; 96549b5e25fSSatish Balay Mat outA; 96649b5e25fSSatish Balay int ierr; 96749b5e25fSSatish Balay PetscTruth row_identity,col_identity; 96849b5e25fSSatish Balay 96949b5e25fSSatish Balay PetscFunctionBegin; 9701a3463dfSHong Zhang /* 97129bbc08cSBarry Smith if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 97249b5e25fSSatish Balay ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 97349b5e25fSSatish Balay ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 97449b5e25fSSatish Balay if (!row_identity || !col_identity) { 97529bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ICC"); 97649b5e25fSSatish Balay } 9771a3463dfSHong Zhang */ 97849b5e25fSSatish Balay 97949b5e25fSSatish Balay outA = inA; 98049b5e25fSSatish Balay inA->factor = FACTOR_LU; 98149b5e25fSSatish Balay 98249b5e25fSSatish Balay if (!a->diag) { 9831a3463dfSHong Zhang ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr); 98449b5e25fSSatish Balay } 98549b5e25fSSatish Balay /* 98649b5e25fSSatish Balay Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 98749b5e25fSSatish Balay for ILU(0) factorization with natural ordering 98849b5e25fSSatish Balay */ 98949b5e25fSSatish Balay switch (a->bs) { 99049b5e25fSSatish Balay case 1: 9911a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; /* 2? */ 9921a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 99349b5e25fSSatish Balay case 2: 9941a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 9951a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9961a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering; 9971a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 99849b5e25fSSatish Balay break; 99949b5e25fSSatish Balay case 3: 10001a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 10011a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10021a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering; 10031a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 100449b5e25fSSatish Balay break; 100549b5e25fSSatish Balay case 4: 10061a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 10071a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10081a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering; 10091a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 101049b5e25fSSatish Balay break; 101149b5e25fSSatish Balay case 5: 10121a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 10131a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10141a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_5_NaturalOrdering; 10151a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 101649b5e25fSSatish Balay break; 101749b5e25fSSatish Balay case 6: 10181a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 10191a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10201a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_6_NaturalOrdering; 10211a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 102249b5e25fSSatish Balay break; 102349b5e25fSSatish Balay case 7: 10241a3463dfSHong Zhang inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 10251a3463dfSHong Zhang inA->ops->solvetranspose = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10261a3463dfSHong Zhang inA->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 10271a3463dfSHong Zhang PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 102849b5e25fSSatish Balay break; 102949b5e25fSSatish Balay default: 103049b5e25fSSatish Balay a->row = row; 10311a3463dfSHong Zhang a->icol = col; 103249b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 103349b5e25fSSatish Balay ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 103449b5e25fSSatish Balay 103549b5e25fSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 103649b5e25fSSatish Balay ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr); 103749b5e25fSSatish Balay PLogObjectParent(inA,a->icol); 103849b5e25fSSatish Balay 103949b5e25fSSatish Balay if (!a->solve_work) { 1040c464158bSHong Zhang a->solve_work = (Scalar*)PetscMalloc((A->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1041c464158bSHong Zhang PLogObjectMemory(inA,(A->m+a->bs)*sizeof(Scalar)); 104249b5e25fSSatish Balay } 104349b5e25fSSatish Balay } 104449b5e25fSSatish Balay 10451a3463dfSHong Zhang ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr); 104649b5e25fSSatish Balay 104749b5e25fSSatish Balay PetscFunctionReturn(0); 104849b5e25fSSatish Balay } 1049950f1e5bSHong Zhang #endif 1050950f1e5bSHong Zhang 105149b5e25fSSatish Balay #undef __FUNC__ 105249b5e25fSSatish Balay #define __FUNC__ "MatPrintHelp_SeqSBAIJ" 105349b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A) 105449b5e25fSSatish Balay { 105549b5e25fSSatish Balay static PetscTruth called = PETSC_FALSE; 105649b5e25fSSatish Balay MPI_Comm comm = A->comm; 105749b5e25fSSatish Balay int ierr; 105849b5e25fSSatish Balay 105949b5e25fSSatish Balay PetscFunctionBegin; 106049b5e25fSSatish Balay if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 106149b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 106249b5e25fSSatish Balay ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 106349b5e25fSSatish Balay PetscFunctionReturn(0); 106449b5e25fSSatish Balay } 106549b5e25fSSatish Balay 106649b5e25fSSatish Balay EXTERN_C_BEGIN 106749b5e25fSSatish Balay #undef __FUNC__ 106849b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ" 106949b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices) 107049b5e25fSSatish Balay { 1071045c9aa0SHong Zhang Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data; 107249b5e25fSSatish Balay int i,nz,n; 107349b5e25fSSatish Balay 107449b5e25fSSatish Balay PetscFunctionBegin; 1075045c9aa0SHong Zhang nz = baij->s_maxnz; 1076c464158bSHong Zhang n = mat->n; 107749b5e25fSSatish Balay for (i=0; i<nz; i++) { 107849b5e25fSSatish Balay baij->j[i] = indices[i]; 107949b5e25fSSatish Balay } 1080045c9aa0SHong Zhang baij->s_nz = nz; 108149b5e25fSSatish Balay for (i=0; i<n; i++) { 108249b5e25fSSatish Balay baij->ilen[i] = baij->imax[i]; 108349b5e25fSSatish Balay } 108449b5e25fSSatish Balay 108549b5e25fSSatish Balay PetscFunctionReturn(0); 108649b5e25fSSatish Balay } 108749b5e25fSSatish Balay EXTERN_C_END 108849b5e25fSSatish Balay 108949b5e25fSSatish Balay #undef __FUNC__ 109049b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices" 109149b5e25fSSatish Balay /*@ 109219585528SSatish Balay MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows 109349b5e25fSSatish Balay in the matrix. 109449b5e25fSSatish Balay 109549b5e25fSSatish Balay Input Parameters: 109619585528SSatish Balay + mat - the SeqSBAIJ matrix 109749b5e25fSSatish Balay - indices - the column indices 109849b5e25fSSatish Balay 109949b5e25fSSatish Balay Level: advanced 110049b5e25fSSatish Balay 110149b5e25fSSatish Balay Notes: 110249b5e25fSSatish Balay This can be called if you have precomputed the nonzero structure of the 110349b5e25fSSatish Balay matrix and want to provide it to the matrix object to improve the performance 110449b5e25fSSatish Balay of the MatSetValues() operation. 110549b5e25fSSatish Balay 110649b5e25fSSatish Balay You MUST have set the correct numbers of nonzeros per row in the call to 110719585528SSatish Balay MatCreateSeqSBAIJ(). 110849b5e25fSSatish Balay 1109ab9f2c04SSatish Balay MUST be called before any calls to MatSetValues() 111049b5e25fSSatish Balay 1111ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ 111249b5e25fSSatish Balay @*/ 111349b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices) 111449b5e25fSSatish Balay { 111549b5e25fSSatish Balay int ierr,(*f)(Mat,int *); 111649b5e25fSSatish Balay 111749b5e25fSSatish Balay PetscFunctionBegin; 111849b5e25fSSatish Balay PetscValidHeaderSpecific(mat,MAT_COOKIE); 11194afc71dfSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 112049b5e25fSSatish Balay if (f) { 112149b5e25fSSatish Balay ierr = (*f)(mat,indices);CHKERRQ(ierr); 112249b5e25fSSatish Balay } else { 112329bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 112449b5e25fSSatish Balay } 112549b5e25fSSatish Balay PetscFunctionReturn(0); 112649b5e25fSSatish Balay } 112749b5e25fSSatish Balay 112849b5e25fSSatish Balay /* -------------------------------------------------------------------*/ 112949b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ, 113049b5e25fSSatish Balay MatGetRow_SeqSBAIJ, 113149b5e25fSSatish Balay MatRestoreRow_SeqSBAIJ, 113249b5e25fSSatish Balay MatMult_SeqSBAIJ_N, 113349b5e25fSSatish Balay MatMultAdd_SeqSBAIJ_N, 113449b5e25fSSatish Balay MatMultTranspose_SeqSBAIJ, 113549b5e25fSSatish Balay MatMultTransposeAdd_SeqSBAIJ, 113649b5e25fSSatish Balay MatSolve_SeqSBAIJ_N, 113749b5e25fSSatish Balay 0, 113849b5e25fSSatish Balay 0, 113949b5e25fSSatish Balay 0, 114049b5e25fSSatish Balay 0, 11415f4ef2deSHong Zhang MatCholeskyFactor_SeqSBAIJ, 114249b5e25fSSatish Balay 0, 114349b5e25fSSatish Balay MatTranspose_SeqSBAIJ, 114449b5e25fSSatish Balay MatGetInfo_SeqSBAIJ, 114549b5e25fSSatish Balay MatEqual_SeqSBAIJ, 114649b5e25fSSatish Balay MatGetDiagonal_SeqSBAIJ, 114749b5e25fSSatish Balay MatDiagonalScale_SeqSBAIJ, 114849b5e25fSSatish Balay MatNorm_SeqSBAIJ, 114949b5e25fSSatish Balay 0, 115049b5e25fSSatish Balay MatAssemblyEnd_SeqSBAIJ, 115149b5e25fSSatish Balay 0, 115249b5e25fSSatish Balay MatSetOption_SeqSBAIJ, 115349b5e25fSSatish Balay MatZeroEntries_SeqSBAIJ, 115449b5e25fSSatish Balay MatZeroRows_SeqSBAIJ, 115549b5e25fSSatish Balay 0, 115649b5e25fSSatish Balay 0, 11575f4ef2deSHong Zhang MatCholeskyFactorSymbolic_SeqSBAIJ, 11585f4ef2deSHong Zhang MatCholeskyFactorNumeric_SeqSBAIJ_N, 1159c464158bSHong Zhang 0, 1160c464158bSHong Zhang 0, 116149b5e25fSSatish Balay MatGetOwnershipRange_SeqSBAIJ, 116249b5e25fSSatish Balay 0, 11635f4ef2deSHong Zhang MatIncompleteCholeskyFactorSymbolic_SeqSBAIJ, 116449b5e25fSSatish Balay 0, 116549b5e25fSSatish Balay 0, 116649b5e25fSSatish Balay MatDuplicate_SeqSBAIJ, 116749b5e25fSSatish Balay 0, 116849b5e25fSSatish Balay 0, 116949b5e25fSSatish Balay 0, 1170950f1e5bSHong Zhang 0, 117149b5e25fSSatish Balay 0, 117249b5e25fSSatish Balay MatGetSubMatrices_SeqSBAIJ, 117349b5e25fSSatish Balay MatIncreaseOverlap_SeqSBAIJ, 117449b5e25fSSatish Balay MatGetValues_SeqSBAIJ, 117549b5e25fSSatish Balay 0, 117649b5e25fSSatish Balay MatPrintHelp_SeqSBAIJ, 117749b5e25fSSatish Balay MatScale_SeqSBAIJ, 117849b5e25fSSatish Balay 0, 117949b5e25fSSatish Balay 0, 118049b5e25fSSatish Balay 0, 118149b5e25fSSatish Balay MatGetBlockSize_SeqSBAIJ, 118249b5e25fSSatish Balay MatGetRowIJ_SeqSBAIJ, 118349b5e25fSSatish Balay MatRestoreRowIJ_SeqSBAIJ, 118449b5e25fSSatish Balay 0, 118549b5e25fSSatish Balay 0, 118649b5e25fSSatish Balay 0, 118749b5e25fSSatish Balay 0, 118849b5e25fSSatish Balay 0, 118949b5e25fSSatish Balay 0, 119049b5e25fSSatish Balay MatSetValuesBlocked_SeqSBAIJ, 119149b5e25fSSatish Balay MatGetSubMatrix_SeqSBAIJ, 119249b5e25fSSatish Balay 0, 119349b5e25fSSatish Balay 0, 1194d959ec07SHong Zhang MatGetMaps_Petsc, 1195d959ec07SHong Zhang 0, 1196d959ec07SHong Zhang 0, 1197d959ec07SHong Zhang 0, 1198d959ec07SHong Zhang 0, 1199d959ec07SHong Zhang 0, 1200d959ec07SHong Zhang 0, 120124d5174aSHong Zhang MatGetRowMax_SeqSBAIJ}; 120249b5e25fSSatish Balay 120349b5e25fSSatish Balay EXTERN_C_BEGIN 120449b5e25fSSatish Balay #undef __FUNC__ 120549b5e25fSSatish Balay #define __FUNC__ "MatStoreValues_SeqSBAIJ" 120649b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat) 120749b5e25fSSatish Balay { 12084afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1209c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2; 121049b5e25fSSatish Balay int ierr; 121149b5e25fSSatish Balay 121249b5e25fSSatish Balay PetscFunctionBegin; 121349b5e25fSSatish Balay if (aij->nonew != 1) { 121429bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 121549b5e25fSSatish Balay } 121649b5e25fSSatish Balay 121749b5e25fSSatish Balay /* allocate space for values if not already there */ 121849b5e25fSSatish Balay if (!aij->saved_values) { 121949b5e25fSSatish Balay aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 122049b5e25fSSatish Balay } 122149b5e25fSSatish Balay 122249b5e25fSSatish Balay /* copy values over */ 122349b5e25fSSatish Balay ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 122449b5e25fSSatish Balay PetscFunctionReturn(0); 122549b5e25fSSatish Balay } 122649b5e25fSSatish Balay EXTERN_C_END 122749b5e25fSSatish Balay 122849b5e25fSSatish Balay EXTERN_C_BEGIN 122949b5e25fSSatish Balay #undef __FUNC__ 123049b5e25fSSatish Balay #define __FUNC__ "MatRetrieveValues_SeqSBAIJ" 123149b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat) 123249b5e25fSSatish Balay { 12334afc71dfSHong Zhang Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data; 1234c464158bSHong Zhang int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 123549b5e25fSSatish Balay 123649b5e25fSSatish Balay PetscFunctionBegin; 123749b5e25fSSatish Balay if (aij->nonew != 1) { 123829bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 123949b5e25fSSatish Balay } 124049b5e25fSSatish Balay if (!aij->saved_values) { 124129bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 124249b5e25fSSatish Balay } 124349b5e25fSSatish Balay 124449b5e25fSSatish Balay /* copy values over */ 124549b5e25fSSatish Balay ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 124649b5e25fSSatish Balay PetscFunctionReturn(0); 124749b5e25fSSatish Balay } 124849b5e25fSSatish Balay EXTERN_C_END 124949b5e25fSSatish Balay 125049b5e25fSSatish Balay #undef __FUNC__ 1251c464158bSHong Zhang #define __FUNC__ "MatCreate_SeqSBAIJ" 1252c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B) 1253c464158bSHong Zhang { 1254c464158bSHong Zhang Mat_SeqSBAIJ *b; 1255c464158bSHong Zhang int i,ierr,size; 1256c464158bSHong Zhang PetscTruth flg; 1257c464158bSHong Zhang int s_nz; 1258c464158bSHong Zhang 1259c464158bSHong Zhang PetscFunctionBegin; 1260c464158bSHong Zhang ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1261c464158bSHong Zhang if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1262c464158bSHong Zhang 1263c464158bSHong Zhang B->data = (void*)(b = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(b); 1264c464158bSHong Zhang ierr = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr); 1265c464158bSHong Zhang ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1266c464158bSHong Zhang B->ops->destroy = MatDestroy_SeqSBAIJ; 1267c464158bSHong Zhang B->ops->view = MatView_SeqSBAIJ; 1268c464158bSHong Zhang B->factor = 0; 1269c464158bSHong Zhang B->lupivotthreshold = 1.0; 1270c464158bSHong Zhang B->mapping = 0; 1271c464158bSHong Zhang b->row = 0; 1272c464158bSHong Zhang b->icol = 0; 1273c464158bSHong Zhang b->reallocs = 0; 1274c464158bSHong Zhang b->saved_values = 0; 1275c464158bSHong Zhang 1276c464158bSHong Zhang ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 1277c464158bSHong Zhang ierr = MapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr); 1278c464158bSHong Zhang 1279c464158bSHong Zhang b->sorted = PETSC_FALSE; 1280c464158bSHong Zhang b->roworiented = PETSC_TRUE; 1281c464158bSHong Zhang b->nonew = 0; 1282c464158bSHong Zhang b->diag = 0; 1283c464158bSHong Zhang b->solve_work = 0; 1284c464158bSHong Zhang b->mult_work = 0; 1285c464158bSHong Zhang b->spptr = 0; 1286c464158bSHong Zhang b->keepzeroedrows = PETSC_FALSE; 1287c464158bSHong Zhang 1288c464158bSHong Zhang b->inew = 0; 1289c464158bSHong Zhang b->jnew = 0; 1290c464158bSHong Zhang b->anew = 0; 1291c464158bSHong Zhang b->a2anew = 0; 1292c464158bSHong Zhang b->permute = PETSC_FALSE; 1293c464158bSHong Zhang 1294c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 1295c464158bSHong Zhang "MatStoreValues_SeqSBAIJ", 1296c464158bSHong Zhang (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr); 1297c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 1298c464158bSHong Zhang "MatRetrieveValues_SeqSBAIJ", 1299c464158bSHong Zhang (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr); 1300c464158bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C", 1301c464158bSHong Zhang "MatSeqSBAIJSetColumnIndices_SeqSBAIJ", 1302c464158bSHong Zhang (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr); 1303c464158bSHong Zhang PetscFunctionReturn(0); 1304c464158bSHong Zhang } 1305c464158bSHong Zhang 1306c464158bSHong Zhang #undef __FUNC__ 1307c464158bSHong Zhang #define __FUNC__ "MatSeqSBAIJSetNonzeroStructure" 130849b5e25fSSatish Balay /*@C 1309c464158bSHong Zhang MatSeqSBAIJSetNonzeroStructure - Creates a sparse symmetric matrix in block AIJ (block 131049b5e25fSSatish Balay compressed row) format. For good matrix assembly performance the 131149b5e25fSSatish Balay user should preallocate the matrix storage by setting the parameter nz 131249b5e25fSSatish Balay (or the array nnz). By setting these parameters accurately, performance 131349b5e25fSSatish Balay during matrix assembly can be increased by more than a factor of 50. 131449b5e25fSSatish Balay 1315c464158bSHong Zhang Collective on Mat 131649b5e25fSSatish Balay 131749b5e25fSSatish Balay Input Parameters: 1318c464158bSHong Zhang + A - the symmetric matrix 131949b5e25fSSatish Balay . bs - size of block 132049b5e25fSSatish Balay . nz - number of block nonzeros per block row (same for all rows) 132149b5e25fSSatish Balay - nnz - array containing the number of block nonzeros in the various block rows 132249b5e25fSSatish Balay (possibly different for each block row) or PETSC_NULL 132349b5e25fSSatish Balay 132449b5e25fSSatish Balay Options Database Keys: 132549b5e25fSSatish Balay . -mat_no_unroll - uses code that does not unroll the loops in the 132649b5e25fSSatish Balay block calculations (much slower) 132749b5e25fSSatish Balay . -mat_block_size - size of the blocks to use 132849b5e25fSSatish Balay 132949b5e25fSSatish Balay Level: intermediate 133049b5e25fSSatish Balay 133149b5e25fSSatish Balay Notes: 1332c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 133349b5e25fSSatish Balay storage. That is, the stored row and column indices can begin at 133449b5e25fSSatish Balay either one (as in Fortran) or zero. See the users' manual for details. 133549b5e25fSSatish Balay 133649b5e25fSSatish Balay Specify the preallocated storage with either nz or nnz (not both). 133749b5e25fSSatish Balay Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 133849b5e25fSSatish Balay allocation. For additional details, see the users manual chapter on 133949b5e25fSSatish Balay matrices. 134049b5e25fSSatish Balay 1341a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 134249b5e25fSSatish Balay @*/ 1343c464158bSHong Zhang int MatSeqSBAIJSetNonzeroStructure(Mat B,int bs,int nz,int *nnz) 134449b5e25fSSatish Balay { 1345c464158bSHong Zhang Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data; 13464afc71dfSHong Zhang int i,len,ierr,mbs,bs2,size; 134749b5e25fSSatish Balay PetscTruth flg; 13484afc71dfSHong Zhang int s_nz; 134949b5e25fSSatish Balay 135049b5e25fSSatish Balay PetscFunctionBegin; 135149b5e25fSSatish Balay ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1352c464158bSHong Zhang mbs = B->m/bs; 135349b5e25fSSatish Balay bs2 = bs*bs; 135449b5e25fSSatish Balay 1355c464158bSHong Zhang if (mbs*bs != B->m) { 135629bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize"); 135749b5e25fSSatish Balay } 135849b5e25fSSatish Balay 135929bbc08cSBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than -2: value %d",nz); 136049b5e25fSSatish Balay if (nnz) { 136149b5e25fSSatish Balay for (i=0; i<mbs; i++) { 136229bbc08cSBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 136349b5e25fSSatish Balay } 136449b5e25fSSatish Balay } 136549b5e25fSSatish Balay 136649b5e25fSSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 136749b5e25fSSatish Balay if (!flg) { 136849b5e25fSSatish Balay switch (bs) { 136949b5e25fSSatish Balay case 1: 13704ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1; 137149b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_1; 137249b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_1; 137349b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_1; 137449b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_1; 137549b5e25fSSatish Balay break; 137649b5e25fSSatish Balay case 2: 13774ccecd49SHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2; 137849b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_2; 137949b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_2; 138049b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_2; 138149b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_2; 138249b5e25fSSatish Balay break; 138349b5e25fSSatish Balay case 3: 13845f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3; 138549b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_3; 138649b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_3; 138749b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_3; 138849b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_3; 138949b5e25fSSatish Balay break; 139049b5e25fSSatish Balay case 4: 13915f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4; 139249b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_4; 139349b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_4; 139449b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_4; 139549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_4; 139649b5e25fSSatish Balay break; 139749b5e25fSSatish Balay case 5: 13985f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5; 139949b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_5; 140049b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_5; 140149b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_5; 140249b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_5; 140349b5e25fSSatish Balay break; 140449b5e25fSSatish Balay case 6: 14055f4ef2deSHong Zhang B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6; 140649b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_6; 140749b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_6; 140849b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_6; 140949b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_6; 141049b5e25fSSatish Balay break; 141149b5e25fSSatish Balay case 7: 141249b5e25fSSatish Balay B->ops->mult = MatMult_SeqSBAIJ_7; 141349b5e25fSSatish Balay B->ops->solve = MatSolve_SeqSBAIJ_7; 141449b5e25fSSatish Balay B->ops->solvetranspose = MatSolveTranspose_SeqSBAIJ_7; 141549b5e25fSSatish Balay B->ops->multadd = MatMultAdd_SeqSBAIJ_7; 141649b5e25fSSatish Balay break; 141749b5e25fSSatish Balay } 141849b5e25fSSatish Balay } 141949b5e25fSSatish Balay 142049b5e25fSSatish Balay b->mbs = mbs; 14214afc71dfSHong Zhang b->nbs = mbs; 142249b5e25fSSatish Balay b->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax); 142349b5e25fSSatish Balay if (!nnz) { 142449b5e25fSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 142549b5e25fSSatish Balay else if (nz <= 0) nz = 1; 142649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 14278cef66ccSHong Zhang b->imax[i] = nz; 142849b5e25fSSatish Balay } 1429153ea458SHong Zhang nz = nz*mbs; /* total nz */ 143049b5e25fSSatish Balay } else { 143149b5e25fSSatish Balay nz = 0; 14328cef66ccSHong Zhang for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 143349b5e25fSSatish Balay } 14348cef66ccSHong Zhang /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */ 14358cef66ccSHong Zhang s_nz = nz; 143649b5e25fSSatish Balay 143749b5e25fSSatish Balay /* allocate the matrix space */ 1438c464158bSHong Zhang len = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 143949b5e25fSSatish Balay b->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a); 144049b5e25fSSatish Balay ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 144149b5e25fSSatish Balay b->j = (int*)(b->a + s_nz*bs2); 144249b5e25fSSatish Balay ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr); 144349b5e25fSSatish Balay b->i = b->j + s_nz; 144449b5e25fSSatish Balay b->singlemalloc = PETSC_TRUE; 144549b5e25fSSatish Balay 144649b5e25fSSatish Balay /* pointer to beginning of each row */ 144749b5e25fSSatish Balay b->i[0] = 0; 144849b5e25fSSatish Balay for (i=1; i<mbs+1; i++) { 144949b5e25fSSatish Balay b->i[i] = b->i[i-1] + b->imax[i-1]; 145049b5e25fSSatish Balay } 145149b5e25fSSatish Balay 145249b5e25fSSatish Balay /* b->ilen will count nonzeros in each block row so far. */ 145349b5e25fSSatish Balay b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int)); 145449b5e25fSSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 145549b5e25fSSatish Balay for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 145649b5e25fSSatish Balay 145749b5e25fSSatish Balay b->bs = bs; 145849b5e25fSSatish Balay b->bs2 = bs2; 145949b5e25fSSatish Balay b->s_nz = 0; 146049b5e25fSSatish Balay b->s_maxnz = s_nz*bs2; 1461153ea458SHong Zhang 146216cdd363SHong Zhang b->inew = 0; 146316cdd363SHong Zhang b->jnew = 0; 146416cdd363SHong Zhang b->anew = 0; 146516cdd363SHong Zhang b->a2anew = 0; 14661a3463dfSHong Zhang b->permute = PETSC_FALSE; 1467c464158bSHong Zhang PetscFunctionReturn(0); 1468c464158bSHong Zhang } 1469153ea458SHong Zhang 147049b5e25fSSatish Balay 1471c464158bSHong Zhang #undef __FUNC__ 1472c464158bSHong Zhang #define __FUNC__ "MatCreateSeqSBAIJ" 1473c464158bSHong Zhang /*@C 1474c464158bSHong Zhang MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block 1475c464158bSHong Zhang compressed row) format. For good matrix assembly performance the 1476c464158bSHong Zhang user should preallocate the matrix storage by setting the parameter nz 1477c464158bSHong Zhang (or the array nnz). By setting these parameters accurately, performance 1478c464158bSHong Zhang during matrix assembly can be increased by more than a factor of 50. 147949b5e25fSSatish Balay 1480c464158bSHong Zhang Collective on MPI_Comm 1481c464158bSHong Zhang 1482c464158bSHong Zhang Input Parameters: 1483c464158bSHong Zhang + comm - MPI communicator, set to PETSC_COMM_SELF 1484c464158bSHong Zhang . bs - size of block 1485c464158bSHong Zhang . m - number of rows, or number of columns 1486c464158bSHong Zhang . nz - number of block nonzeros per block row (same for all rows) 1487c464158bSHong Zhang - nnz - array containing the number of block nonzeros in the various block rows 1488c464158bSHong Zhang (possibly different for each block row) or PETSC_NULL 1489c464158bSHong Zhang 1490c464158bSHong Zhang Output Parameter: 1491c464158bSHong Zhang . A - the symmetric matrix 1492c464158bSHong Zhang 1493c464158bSHong Zhang Options Database Keys: 1494c464158bSHong Zhang . -mat_no_unroll - uses code that does not unroll the loops in the 1495c464158bSHong Zhang block calculations (much slower) 1496c464158bSHong Zhang . -mat_block_size - size of the blocks to use 1497c464158bSHong Zhang 1498c464158bSHong Zhang Level: intermediate 1499c464158bSHong Zhang 1500c464158bSHong Zhang Notes: 1501c464158bSHong Zhang The block AIJ format is compatible with standard Fortran 77 1502c464158bSHong Zhang storage. That is, the stored row and column indices can begin at 1503c464158bSHong Zhang either one (as in Fortran) or zero. See the users' manual for details. 1504c464158bSHong Zhang 1505c464158bSHong Zhang Specify the preallocated storage with either nz or nnz (not both). 1506c464158bSHong Zhang Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1507c464158bSHong Zhang allocation. For additional details, see the users manual chapter on 1508c464158bSHong Zhang matrices. 1509c464158bSHong Zhang 1510c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ() 1511c464158bSHong Zhang @*/ 1512c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1513c464158bSHong Zhang { 1514c464158bSHong Zhang int ierr; 1515c464158bSHong Zhang 1516c464158bSHong Zhang PetscFunctionBegin; 1517c464158bSHong Zhang ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1518c464158bSHong Zhang ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr); 1519c464158bSHong Zhang ierr = MatSeqSBAIJSetNonzeroStructure(*A,bs,nz,nnz);CHKERRQ(ierr); 152049b5e25fSSatish Balay PetscFunctionReturn(0); 152149b5e25fSSatish Balay } 152249b5e25fSSatish Balay 152349b5e25fSSatish Balay #undef __FUNC__ 152449b5e25fSSatish Balay #define __FUNC__ "MatDuplicate_SeqSBAIJ" 152549b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 152649b5e25fSSatish Balay { 152749b5e25fSSatish Balay Mat C; 152849b5e25fSSatish Balay Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data; 152949b5e25fSSatish Balay int i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr; 153049b5e25fSSatish Balay 153149b5e25fSSatish Balay PetscFunctionBegin; 153229bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 153349b5e25fSSatish Balay 153449b5e25fSSatish Balay *B = 0; 1535c464158bSHong Zhang ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&C);CHKERRQ(ierr); 1536c464158bSHong Zhang ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr); 153749b5e25fSSatish Balay C->data = (void*)(c = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(c); 153849b5e25fSSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 153949b5e25fSSatish Balay C->ops->destroy = MatDestroy_SeqSBAIJ; 154049b5e25fSSatish Balay C->ops->view = MatView_SeqSBAIJ; 154149b5e25fSSatish Balay C->factor = A->factor; 154249b5e25fSSatish Balay c->row = 0; 154349b5e25fSSatish Balay c->icol = 0; 154449b5e25fSSatish Balay c->saved_values = 0; 154549b5e25fSSatish Balay c->keepzeroedrows = a->keepzeroedrows; 154649b5e25fSSatish Balay C->assembled = PETSC_TRUE; 154749b5e25fSSatish Balay 154849b5e25fSSatish Balay c->bs = a->bs; 154949b5e25fSSatish Balay c->bs2 = a->bs2; 155049b5e25fSSatish Balay c->mbs = a->mbs; 155149b5e25fSSatish Balay c->nbs = a->nbs; 155249b5e25fSSatish Balay 155349b5e25fSSatish Balay c->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 155449b5e25fSSatish Balay c->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 155549b5e25fSSatish Balay for (i=0; i<mbs; i++) { 155649b5e25fSSatish Balay c->imax[i] = a->imax[i]; 155749b5e25fSSatish Balay c->ilen[i] = a->ilen[i]; 155849b5e25fSSatish Balay } 155949b5e25fSSatish Balay 156049b5e25fSSatish Balay /* allocate the matrix space */ 156149b5e25fSSatish Balay c->singlemalloc = PETSC_TRUE; 156249b5e25fSSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 156349b5e25fSSatish Balay c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a); 156449b5e25fSSatish Balay c->j = (int*)(c->a + nz*bs2); 156549b5e25fSSatish Balay c->i = c->j + nz; 156649b5e25fSSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 156749b5e25fSSatish Balay if (mbs > 0) { 156849b5e25fSSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 156949b5e25fSSatish Balay if (cpvalues == MAT_COPY_VALUES) { 157049b5e25fSSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 157149b5e25fSSatish Balay } else { 157249b5e25fSSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 157349b5e25fSSatish Balay } 157449b5e25fSSatish Balay } 157549b5e25fSSatish Balay 157649b5e25fSSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ)); 157749b5e25fSSatish Balay c->sorted = a->sorted; 157849b5e25fSSatish Balay c->roworiented = a->roworiented; 157949b5e25fSSatish Balay c->nonew = a->nonew; 158049b5e25fSSatish Balay 158149b5e25fSSatish Balay if (a->diag) { 158249b5e25fSSatish Balay c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag); 158349b5e25fSSatish Balay PLogObjectMemory(C,(mbs+1)*sizeof(int)); 158449b5e25fSSatish Balay for (i=0; i<mbs; i++) { 158549b5e25fSSatish Balay c->diag[i] = a->diag[i]; 158649b5e25fSSatish Balay } 158749b5e25fSSatish Balay } else c->diag = 0; 158849b5e25fSSatish Balay c->s_nz = a->s_nz; 158949b5e25fSSatish Balay c->s_maxnz = a->s_maxnz; 159049b5e25fSSatish Balay c->solve_work = 0; 159149b5e25fSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 159249b5e25fSSatish Balay c->mult_work = 0; 159349b5e25fSSatish Balay *B = C; 159449b5e25fSSatish Balay ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 159549b5e25fSSatish Balay PetscFunctionReturn(0); 159649b5e25fSSatish Balay } 159749b5e25fSSatish Balay 159849b5e25fSSatish Balay #undef __FUNC__ 159949b5e25fSSatish Balay #define __FUNC__ "MatLoad_SeqSBAIJ" 160049b5e25fSSatish Balay int MatLoad_SeqSBAIJ(Viewer viewer,MatType type,Mat *A) 160149b5e25fSSatish Balay { 160249b5e25fSSatish Balay Mat_SeqSBAIJ *a; 160349b5e25fSSatish Balay Mat B; 160449b5e25fSSatish Balay int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1605574b2666SHong Zhang int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount; 160649b5e25fSSatish Balay int kmax,jcount,block,idx,point,nzcountb,extra_rows; 160749b5e25fSSatish Balay int *masked,nmask,tmp,bs2,ishift; 160849b5e25fSSatish Balay Scalar *aa; 160949b5e25fSSatish Balay MPI_Comm comm = ((PetscObject)viewer)->comm; 161049b5e25fSSatish Balay 161149b5e25fSSatish Balay PetscFunctionBegin; 161249b5e25fSSatish Balay ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 161349b5e25fSSatish Balay bs2 = bs*bs; 161449b5e25fSSatish Balay 161549b5e25fSSatish Balay ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 161629bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 161749b5e25fSSatish Balay ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 161849b5e25fSSatish Balay ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 161929bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 162049b5e25fSSatish Balay M = header[1]; N = header[2]; nz = header[3]; 162149b5e25fSSatish Balay 162249b5e25fSSatish Balay if (header[3] < 0) { 162329bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ"); 162449b5e25fSSatish Balay } 162549b5e25fSSatish Balay 162629bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 162749b5e25fSSatish Balay 162849b5e25fSSatish Balay /* 162949b5e25fSSatish Balay This code adds extra rows to make sure the number of rows is 163049b5e25fSSatish Balay divisible by the blocksize 163149b5e25fSSatish Balay */ 163249b5e25fSSatish Balay mbs = M/bs; 163349b5e25fSSatish Balay extra_rows = bs - M + bs*(mbs); 163449b5e25fSSatish Balay if (extra_rows == bs) extra_rows = 0; 163549b5e25fSSatish Balay else mbs++; 163649b5e25fSSatish Balay if (extra_rows) { 163749b5e25fSSatish Balay PLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"); 163849b5e25fSSatish Balay } 163949b5e25fSSatish Balay 164049b5e25fSSatish Balay /* read in row lengths */ 164149b5e25fSSatish Balay rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 164249b5e25fSSatish Balay ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 164349b5e25fSSatish Balay for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 164449b5e25fSSatish Balay 164549b5e25fSSatish Balay /* read in column indices */ 164649b5e25fSSatish Balay jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj); 164749b5e25fSSatish Balay ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 164849b5e25fSSatish Balay for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 164949b5e25fSSatish Balay 165049b5e25fSSatish Balay /* loop over row lengths determining block row lengths */ 165149b5e25fSSatish Balay browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1652574b2666SHong Zhang s_browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(s_browlengths); 1653574b2666SHong Zhang ierr = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr); 165449b5e25fSSatish Balay mask = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask); 165549b5e25fSSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 165649b5e25fSSatish Balay masked = mask + mbs; 165749b5e25fSSatish Balay rowcount = 0; nzcount = 0; 165849b5e25fSSatish Balay for (i=0; i<mbs; i++) { 165949b5e25fSSatish Balay nmask = 0; 166049b5e25fSSatish Balay for (j=0; j<bs; j++) { 166149b5e25fSSatish Balay kmax = rowlengths[rowcount]; 166249b5e25fSSatish Balay for (k=0; k<kmax; k++) { 16632d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 166403630b6eSHong Zhang if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;} 166549b5e25fSSatish Balay } 166649b5e25fSSatish Balay rowcount++; 166749b5e25fSSatish Balay } 1668574b2666SHong Zhang s_browlengths[i] += nmask; 1669574b2666SHong Zhang browlengths[i] = 2*s_browlengths[i]; 1670574b2666SHong Zhang 167149b5e25fSSatish Balay /* zero out the mask elements we set */ 167249b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 167349b5e25fSSatish Balay } 167449b5e25fSSatish Balay 167549b5e25fSSatish Balay /* create our matrix */ 167649b5e25fSSatish Balay ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 167749b5e25fSSatish Balay B = *A; 167849b5e25fSSatish Balay a = (Mat_SeqSBAIJ*)B->data; 167949b5e25fSSatish Balay 168049b5e25fSSatish Balay /* set matrix "i" values */ 168149b5e25fSSatish Balay a->i[0] = 0; 168249b5e25fSSatish Balay for (i=1; i<= mbs; i++) { 1683574b2666SHong Zhang a->i[i] = a->i[i-1] + s_browlengths[i-1]; 1684574b2666SHong Zhang a->ilen[i-1] = s_browlengths[i-1]; 168549b5e25fSSatish Balay } 168649b5e25fSSatish Balay a->s_nz = 0; 1687574b2666SHong Zhang for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i]; 168849b5e25fSSatish Balay 168949b5e25fSSatish Balay /* read in nonzero values */ 169049b5e25fSSatish Balay aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 169149b5e25fSSatish Balay ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 169249b5e25fSSatish Balay for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 169349b5e25fSSatish Balay 169449b5e25fSSatish Balay /* set "a" and "j" values into matrix */ 169549b5e25fSSatish Balay nzcount = 0; jcount = 0; 169649b5e25fSSatish Balay for (i=0; i<mbs; i++) { 169749b5e25fSSatish Balay nzcountb = nzcount; 169849b5e25fSSatish Balay nmask = 0; 169949b5e25fSSatish Balay for (j=0; j<bs; j++) { 170049b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 170149b5e25fSSatish Balay for (k=0; k<kmax; k++) { 17022d703238SHong Zhang tmp = jj[nzcount++]/bs; /* block col. index */ 170303630b6eSHong Zhang if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;} 17042d703238SHong Zhang } 17052d703238SHong Zhang rowcount++; 17062d703238SHong Zhang } 17072d703238SHong Zhang /* sort the masked values */ 17082d703238SHong Zhang ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 17092d703238SHong Zhang 17102d703238SHong Zhang /* set "j" values into matrix */ 17112d703238SHong Zhang maskcount = 1; 17122d703238SHong Zhang for (j=0; j<nmask; j++) { 171349b5e25fSSatish Balay a->j[jcount++] = masked[j]; 171449b5e25fSSatish Balay mask[masked[j]] = maskcount++; 171549b5e25fSSatish Balay } 1716574b2666SHong Zhang 171749b5e25fSSatish Balay /* set "a" values into matrix */ 171849b5e25fSSatish Balay ishift = bs2*a->i[i]; 171949b5e25fSSatish Balay for (j=0; j<bs; j++) { 172049b5e25fSSatish Balay kmax = rowlengths[i*bs+j]; 172149b5e25fSSatish Balay for (k=0; k<kmax; k++) { 1722574b2666SHong Zhang tmp = jj[nzcountb]/bs ; /* block col. index */ 1723574b2666SHong Zhang if (tmp >= i){ 172449b5e25fSSatish Balay block = mask[tmp] - 1; 172549b5e25fSSatish Balay point = jj[nzcountb] - bs*tmp; 172649b5e25fSSatish Balay idx = ishift + bs2*block + j + bs*point; 1727574b2666SHong Zhang a->a[idx] = aa[nzcountb]; 1728574b2666SHong Zhang } 1729574b2666SHong Zhang nzcountb++; 173049b5e25fSSatish Balay } 173149b5e25fSSatish Balay } 173249b5e25fSSatish Balay /* zero out the mask elements we set */ 173349b5e25fSSatish Balay for (j=0; j<nmask; j++) mask[masked[j]] = 0; 173449b5e25fSSatish Balay } 173529bbc08cSBarry Smith if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 173649b5e25fSSatish Balay 173749b5e25fSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 173849b5e25fSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1739574b2666SHong Zhang ierr = PetscFree(s_browlengths);CHKERRQ(ierr); 174049b5e25fSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 174149b5e25fSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 174249b5e25fSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 174349b5e25fSSatish Balay 174449b5e25fSSatish Balay B->assembled = PETSC_TRUE; 174549b5e25fSSatish Balay ierr = MatView_Private(B);CHKERRQ(ierr); 174649b5e25fSSatish Balay PetscFunctionReturn(0); 174749b5e25fSSatish Balay } 1748574b2666SHong Zhang 1749574b2666SHong Zhang 175049b5e25fSSatish Balay 175149b5e25fSSatish Balay 1752