1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*7413bad6SBarry Smith static char vcid[] = "$Id: baij.c,v 1.142 1998/07/14 03:05:01 bsmith Exp bsmith $"; 32593348eSBarry Smith #endif 42593348eSBarry Smith 52593348eSBarry Smith /* 6b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 72593348eSBarry Smith matrix storage format. 82593348eSBarry Smith */ 93369ce9aSBarry Smith 103369ce9aSBarry Smith #include "pinclude/pviewer.h" 113369ce9aSBarry Smith #include "sys.h" 1270f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 131a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 141a6a6d98SBarry Smith #include "src/inline/spops.h" 152593348eSBarry Smith #include "petsc.h" 163b547af2SSatish Balay 172d61bbb3SSatish Balay #define CHUNKSIZE 10 18de6a44a3SBarry Smith 195615d1e5SSatish Balay #undef __FUNC__ 205615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ" 21de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 22de6a44a3SBarry Smith { 23de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 247fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 25de6a44a3SBarry Smith 263a40ed3dSBarry Smith PetscFunctionBegin; 27de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 28de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 297fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 30de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 31de6a44a3SBarry Smith if (a->j[j] == i) { 32de6a44a3SBarry Smith diag[i] = j; 33de6a44a3SBarry Smith break; 34de6a44a3SBarry Smith } 35de6a44a3SBarry Smith } 36de6a44a3SBarry Smith } 37de6a44a3SBarry Smith a->diag = diag; 383a40ed3dSBarry Smith PetscFunctionReturn(0); 39de6a44a3SBarry Smith } 402593348eSBarry Smith 412593348eSBarry Smith 423b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 433b2fbd54SBarry Smith 445615d1e5SSatish Balay #undef __FUNC__ 455615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 463b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 473b2fbd54SBarry Smith PetscTruth *done) 483b2fbd54SBarry Smith { 493b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 503b2fbd54SBarry Smith int ierr, n = a->mbs,i; 513b2fbd54SBarry Smith 523a40ed3dSBarry Smith PetscFunctionBegin; 533b2fbd54SBarry Smith *nn = n; 543a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 553b2fbd54SBarry Smith if (symmetric) { 563b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 573b2fbd54SBarry Smith } else if (oshift == 1) { 583b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 593b2fbd54SBarry Smith int nz = a->i[n] + 1; 603b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 613b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 623b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 633b2fbd54SBarry Smith } else { 643b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 653b2fbd54SBarry Smith } 663b2fbd54SBarry Smith 673a40ed3dSBarry Smith PetscFunctionReturn(0); 683b2fbd54SBarry Smith } 693b2fbd54SBarry Smith 705615d1e5SSatish Balay #undef __FUNC__ 71d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 723b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 733b2fbd54SBarry Smith PetscTruth *done) 743b2fbd54SBarry Smith { 753b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 763b2fbd54SBarry Smith int i,n = a->mbs; 773b2fbd54SBarry Smith 783a40ed3dSBarry Smith PetscFunctionBegin; 793a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 803b2fbd54SBarry Smith if (symmetric) { 813b2fbd54SBarry Smith PetscFree(*ia); 823b2fbd54SBarry Smith PetscFree(*ja); 83af5da2bfSBarry Smith } else if (oshift == 1) { 843b2fbd54SBarry Smith int nz = a->i[n]; 853b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 863b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 873b2fbd54SBarry Smith } 883a40ed3dSBarry Smith PetscFunctionReturn(0); 893b2fbd54SBarry Smith } 903b2fbd54SBarry Smith 912d61bbb3SSatish Balay #undef __FUNC__ 922d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 932d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 942d61bbb3SSatish Balay { 952d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 962d61bbb3SSatish Balay 972d61bbb3SSatish Balay PetscFunctionBegin; 982d61bbb3SSatish Balay *bs = baij->bs; 992d61bbb3SSatish Balay PetscFunctionReturn(0); 1002d61bbb3SSatish Balay } 1012d61bbb3SSatish Balay 1022d61bbb3SSatish Balay 1032d61bbb3SSatish Balay #undef __FUNC__ 1042d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 105e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1062d61bbb3SSatish Balay { 1072d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 108e51c0b9cSSatish Balay int ierr; 1092d61bbb3SSatish Balay 11094d884c6SBarry Smith if (--A->refct > 0) PetscFunctionReturn(0); 11194d884c6SBarry Smith 11294d884c6SBarry Smith if (A->mapping) { 11394d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 11494d884c6SBarry Smith } 11594d884c6SBarry Smith if (A->bmapping) { 11694d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr); 11794d884c6SBarry Smith } 11861b13de0SBarry Smith if (A->rmap) { 11961b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 12061b13de0SBarry Smith } 12161b13de0SBarry Smith if (A->cmap) { 12261b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 12361b13de0SBarry Smith } 1242d61bbb3SSatish Balay #if defined(USE_PETSC_LOG) 125e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1262d61bbb3SSatish Balay #endif 1272d61bbb3SSatish Balay PetscFree(a->a); 1282d61bbb3SSatish Balay if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 1292d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 1302d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 1312d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 1322d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 1332d61bbb3SSatish Balay if (a->mult_work) PetscFree(a->mult_work); 134e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 1352d61bbb3SSatish Balay PetscFree(a); 1362d61bbb3SSatish Balay PLogObjectDestroy(A); 1372d61bbb3SSatish Balay PetscHeaderDestroy(A); 1382d61bbb3SSatish Balay PetscFunctionReturn(0); 1392d61bbb3SSatish Balay } 1402d61bbb3SSatish Balay 1412d61bbb3SSatish Balay #undef __FUNC__ 1422d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1432d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1442d61bbb3SSatish Balay { 1452d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1462d61bbb3SSatish Balay 1472d61bbb3SSatish Balay PetscFunctionBegin; 1482d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1492d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1502d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1512d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1522d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1532d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 1542d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 1552d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1562d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1572d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1582d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1592d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1602d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 161b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 162b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 163981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1642d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1652d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1662d61bbb3SSatish Balay } else { 1672d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1682d61bbb3SSatish Balay } 1692d61bbb3SSatish Balay PetscFunctionReturn(0); 1702d61bbb3SSatish Balay } 1712d61bbb3SSatish Balay 1722d61bbb3SSatish Balay 1732d61bbb3SSatish Balay #undef __FUNC__ 1742d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1752d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1762d61bbb3SSatish Balay { 1772d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1782d61bbb3SSatish Balay 1792d61bbb3SSatish Balay PetscFunctionBegin; 1802d61bbb3SSatish Balay if (m) *m = a->m; 1812d61bbb3SSatish Balay if (n) *n = a->n; 1822d61bbb3SSatish Balay PetscFunctionReturn(0); 1832d61bbb3SSatish Balay } 1842d61bbb3SSatish Balay 1852d61bbb3SSatish Balay #undef __FUNC__ 1862d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 1872d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 1882d61bbb3SSatish Balay { 1892d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1902d61bbb3SSatish Balay 1912d61bbb3SSatish Balay PetscFunctionBegin; 1922d61bbb3SSatish Balay *m = 0; *n = a->m; 1932d61bbb3SSatish Balay PetscFunctionReturn(0); 1942d61bbb3SSatish Balay } 1952d61bbb3SSatish Balay 1962d61bbb3SSatish Balay #undef __FUNC__ 1972d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 1982d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 1992d61bbb3SSatish Balay { 2002d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2012d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2022d61bbb3SSatish Balay Scalar *aa,*v_i,*aa_i; 2032d61bbb3SSatish Balay 2042d61bbb3SSatish Balay PetscFunctionBegin; 2052d61bbb3SSatish Balay bs = a->bs; 2062d61bbb3SSatish Balay ai = a->i; 2072d61bbb3SSatish Balay aj = a->j; 2082d61bbb3SSatish Balay aa = a->a; 2092d61bbb3SSatish Balay bs2 = a->bs2; 2102d61bbb3SSatish Balay 2112d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2122d61bbb3SSatish Balay 2132d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2142d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2152d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2162d61bbb3SSatish Balay *nz = bs*M; 2172d61bbb3SSatish Balay 2182d61bbb3SSatish Balay if (v) { 2192d61bbb3SSatish Balay *v = 0; 2202d61bbb3SSatish Balay if (*nz) { 2212d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 2222d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2232d61bbb3SSatish Balay v_i = *v + i*bs; 2242d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2252d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2262d61bbb3SSatish Balay } 2272d61bbb3SSatish Balay } 2282d61bbb3SSatish Balay } 2292d61bbb3SSatish Balay 2302d61bbb3SSatish Balay if (idx) { 2312d61bbb3SSatish Balay *idx = 0; 2322d61bbb3SSatish Balay if (*nz) { 2332d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 2342d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2352d61bbb3SSatish Balay idx_i = *idx + i*bs; 2362d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2372d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2382d61bbb3SSatish Balay } 2392d61bbb3SSatish Balay } 2402d61bbb3SSatish Balay } 2412d61bbb3SSatish Balay PetscFunctionReturn(0); 2422d61bbb3SSatish Balay } 2432d61bbb3SSatish Balay 2442d61bbb3SSatish Balay #undef __FUNC__ 2452d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2462d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2472d61bbb3SSatish Balay { 2482d61bbb3SSatish Balay PetscFunctionBegin; 2492d61bbb3SSatish Balay if (idx) {if (*idx) PetscFree(*idx);} 2502d61bbb3SSatish Balay if (v) {if (*v) PetscFree(*v);} 2512d61bbb3SSatish Balay PetscFunctionReturn(0); 2522d61bbb3SSatish Balay } 2532d61bbb3SSatish Balay 2542d61bbb3SSatish Balay #undef __FUNC__ 2552d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2562d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2572d61bbb3SSatish Balay { 2582d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2592d61bbb3SSatish Balay Mat C; 2602d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2612d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2622d61bbb3SSatish Balay Scalar *array=a->a; 2632d61bbb3SSatish Balay 2642d61bbb3SSatish Balay PetscFunctionBegin; 2652d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2662d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col); 2672d61bbb3SSatish Balay PetscMemzero(col,(1+nbs)*sizeof(int)); 2682d61bbb3SSatish Balay 2692d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2702d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr); 2712d61bbb3SSatish Balay PetscFree(col); 2722d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows); 2732d61bbb3SSatish Balay cols = rows + bs; 2742d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2752d61bbb3SSatish Balay cols[0] = i*bs; 2762d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 2772d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 2782d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 2792d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 2802d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 2812d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr); 2822d61bbb3SSatish Balay array += bs2; 2832d61bbb3SSatish Balay } 2842d61bbb3SSatish Balay } 2852d61bbb3SSatish Balay PetscFree(rows); 2862d61bbb3SSatish Balay 2872d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2882d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2892d61bbb3SSatish Balay 2902d61bbb3SSatish Balay if (B != PETSC_NULL) { 2912d61bbb3SSatish Balay *B = C; 2922d61bbb3SSatish Balay } else { 293f830108cSBarry Smith PetscOps *Abops; 294cc2dc46cSBarry Smith MatOps Aops; 295f830108cSBarry Smith 2962d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 2972d61bbb3SSatish Balay PetscFree(a->a); 2982d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 2992d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 3002d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 3012d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 3022d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 3032d61bbb3SSatish Balay PetscFree(a); 304f830108cSBarry Smith 305*7413bad6SBarry Smith 306*7413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 307*7413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 308*7413bad6SBarry Smith 309f830108cSBarry Smith /* 310f830108cSBarry Smith This is horrible, horrible code. We need to keep the 311f830108cSBarry Smith A pointers for the bops and ops but copy everything 312f830108cSBarry Smith else from C. 313f830108cSBarry Smith */ 314f830108cSBarry Smith Abops = A->bops; 315f830108cSBarry Smith Aops = A->ops; 3162d61bbb3SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 317f830108cSBarry Smith A->bops = Abops; 318f830108cSBarry Smith A->ops = Aops; 31927a8da17SBarry Smith A->qlist = 0; 320f830108cSBarry Smith 3212d61bbb3SSatish Balay PetscHeaderDestroy(C); 3222d61bbb3SSatish Balay } 3232d61bbb3SSatish Balay PetscFunctionReturn(0); 3242d61bbb3SSatish Balay } 3252d61bbb3SSatish Balay 3262d61bbb3SSatish Balay 3272d61bbb3SSatish Balay 3283b2fbd54SBarry Smith 3295615d1e5SSatish Balay #undef __FUNC__ 330d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 331b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3322593348eSBarry Smith { 333b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3349df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 335b6490206SBarry Smith Scalar *aa; 3362593348eSBarry Smith 3373a40ed3dSBarry Smith PetscFunctionBegin; 33890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 3392593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3402593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3413b2fbd54SBarry Smith 3422593348eSBarry Smith col_lens[1] = a->m; 3432593348eSBarry Smith col_lens[2] = a->n; 3447e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3452593348eSBarry Smith 3462593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 347b6490206SBarry Smith count = 0; 348b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 349b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 350b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 351b6490206SBarry Smith } 3522593348eSBarry Smith } 3530752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3542593348eSBarry Smith PetscFree(col_lens); 3552593348eSBarry Smith 3562593348eSBarry Smith /* store column indices (zero start index) */ 35766963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj); 358b6490206SBarry Smith count = 0; 359b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 360b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 361b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 362b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 363b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3642593348eSBarry Smith } 3652593348eSBarry Smith } 366b6490206SBarry Smith } 367b6490206SBarry Smith } 3680752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr); 369b6490206SBarry Smith PetscFree(jj); 3702593348eSBarry Smith 3712593348eSBarry Smith /* store nonzero values */ 37266963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa); 373b6490206SBarry Smith count = 0; 374b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 375b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 376b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 377b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 3787e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 379b6490206SBarry Smith } 380b6490206SBarry Smith } 381b6490206SBarry Smith } 382b6490206SBarry Smith } 3830752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 384b6490206SBarry Smith PetscFree(aa); 3853a40ed3dSBarry Smith PetscFunctionReturn(0); 3862593348eSBarry Smith } 3872593348eSBarry Smith 3885615d1e5SSatish Balay #undef __FUNC__ 389d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 390b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 3912593348eSBarry Smith { 392b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3939df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 3942593348eSBarry Smith FILE *fd; 3952593348eSBarry Smith char *outputname; 3962593348eSBarry Smith 3973a40ed3dSBarry Smith PetscFunctionBegin; 39890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 3992593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 40090ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 401639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 4027ddc982cSLois Curfman McInnes fprintf(fd," block size is %d\n",bs); 4032593348eSBarry Smith } 404639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 405a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 4062593348eSBarry Smith } 407639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 40844cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 40944cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 41044cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 41144cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 41244cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 4133a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 414e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) 41544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 416e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 417e20fef11SSatish Balay else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) 418766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 419e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 420e20fef11SSatish Balay else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) 421e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 42244cd7ae7SLois Curfman McInnes #else 42344cd7ae7SLois Curfman McInnes if (a->a[bs2*k + l*bs + j] != 0.0) 42444cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 42544cd7ae7SLois Curfman McInnes #endif 42644cd7ae7SLois Curfman McInnes } 42744cd7ae7SLois Curfman McInnes } 42844cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 42944cd7ae7SLois Curfman McInnes } 43044cd7ae7SLois Curfman McInnes } 43144cd7ae7SLois Curfman McInnes } 4322593348eSBarry Smith else { 433b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 434b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 435b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 436b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 437b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4383a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 439e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 44088685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 441e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 44288685aaeSLois Curfman McInnes } 443e20fef11SSatish Balay else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 444766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 445e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 446766eeae4SLois Curfman McInnes } 44788685aaeSLois Curfman McInnes else { 448e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 44988685aaeSLois Curfman McInnes } 45088685aaeSLois Curfman McInnes #else 4517e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 45288685aaeSLois Curfman McInnes #endif 4532593348eSBarry Smith } 4542593348eSBarry Smith } 4552593348eSBarry Smith fprintf(fd,"\n"); 4562593348eSBarry Smith } 4572593348eSBarry Smith } 458b6490206SBarry Smith } 4592593348eSBarry Smith fflush(fd); 4603a40ed3dSBarry Smith PetscFunctionReturn(0); 4612593348eSBarry Smith } 4622593348eSBarry Smith 4635615d1e5SSatish Balay #undef __FUNC__ 464d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 4653270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 4663270192aSSatish Balay { 4673270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4683270192aSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2; 4693270192aSSatish Balay double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 4703270192aSSatish Balay Scalar *aa; 4713270192aSSatish Balay Draw draw; 4723270192aSSatish Balay DrawButton button; 4733270192aSSatish Balay PetscTruth isnull; 4743270192aSSatish Balay 4753a40ed3dSBarry Smith PetscFunctionBegin; 4763a40ed3dSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr); 4773a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 4783270192aSSatish Balay 4793270192aSSatish Balay xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 4803270192aSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 4813270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 4823270192aSSatish Balay /* loop over matrix elements drawing boxes */ 4833270192aSSatish Balay color = DRAW_BLUE; 4843270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4853270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4863270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4873270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4883270192aSSatish Balay aa = a->a + j*bs2; 4893270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4903270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4913270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 492b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4933270192aSSatish Balay } 4943270192aSSatish Balay } 4953270192aSSatish Balay } 4963270192aSSatish Balay } 4973270192aSSatish Balay color = DRAW_CYAN; 4983270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4993270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5003270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5013270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5023270192aSSatish Balay aa = a->a + j*bs2; 5033270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5043270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5053270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 506b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5073270192aSSatish Balay } 5083270192aSSatish Balay } 5093270192aSSatish Balay } 5103270192aSSatish Balay } 5113270192aSSatish Balay 5123270192aSSatish Balay color = DRAW_RED; 5133270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5143270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5153270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5163270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5173270192aSSatish Balay aa = a->a + j*bs2; 5183270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5193270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5203270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 521b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5223270192aSSatish Balay } 5233270192aSSatish Balay } 5243270192aSSatish Balay } 5253270192aSSatish Balay } 5263270192aSSatish Balay 52755843e3eSBarry Smith DrawSynchronizedFlush(draw); 5283270192aSSatish Balay DrawGetPause(draw,&pause); 5293a40ed3dSBarry Smith if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 5303270192aSSatish Balay 5313270192aSSatish Balay /* allow the matrix to zoom or shrink */ 53255843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5333270192aSSatish Balay while (button != BUTTON_RIGHT) { 53455843e3eSBarry Smith DrawSynchronizedClear(draw); 5353270192aSSatish Balay if (button == BUTTON_LEFT) scale = .5; 5363270192aSSatish Balay else if (button == BUTTON_CENTER) scale = 2.; 5373270192aSSatish Balay xl = scale*(xl + w - xc) + xc - w*scale; 5383270192aSSatish Balay xr = scale*(xr - w - xc) + xc + w*scale; 5393270192aSSatish Balay yl = scale*(yl + h - yc) + yc - h*scale; 5403270192aSSatish Balay yr = scale*(yr - h - yc) + yc + h*scale; 5413270192aSSatish Balay w *= scale; h *= scale; 5423270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 5433270192aSSatish Balay color = DRAW_BLUE; 5443270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5453270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5463270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5473270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5483270192aSSatish Balay aa = a->a + j*bs2; 5493270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5503270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5513270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 552b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5533270192aSSatish Balay } 5543270192aSSatish Balay } 5553270192aSSatish Balay } 5563270192aSSatish Balay } 5573270192aSSatish Balay color = DRAW_CYAN; 5583270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5593270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5603270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5613270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5623270192aSSatish Balay aa = a->a + j*bs2; 5633270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5643270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5653270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 566b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5673270192aSSatish Balay } 5683270192aSSatish Balay } 5693270192aSSatish Balay } 5703270192aSSatish Balay } 5713270192aSSatish Balay 5723270192aSSatish Balay color = DRAW_RED; 5733270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5743270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5753270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5763270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5773270192aSSatish Balay aa = a->a + j*bs2; 5783270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5793270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5803270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 581b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5823270192aSSatish Balay } 5833270192aSSatish Balay } 5843270192aSSatish Balay } 5853270192aSSatish Balay } 58655843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5873270192aSSatish Balay } 5883a40ed3dSBarry Smith PetscFunctionReturn(0); 5893270192aSSatish Balay } 5903270192aSSatish Balay 5915615d1e5SSatish Balay #undef __FUNC__ 592d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 593e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5942593348eSBarry Smith { 59519bcc07fSBarry Smith ViewerType vtype; 59619bcc07fSBarry Smith int ierr; 5972593348eSBarry Smith 5983a40ed3dSBarry Smith PetscFunctionBegin; 59919bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 60019bcc07fSBarry Smith if (vtype == MATLAB_VIEWER) { 601a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported"); 6023a40ed3dSBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 6033a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6043a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 6053a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6063a40ed3dSBarry Smith } else if (vtype == DRAW_VIEWER) { 6073a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6085cd90555SBarry Smith } else { 6095cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 6102593348eSBarry Smith } 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6122593348eSBarry Smith } 613b6490206SBarry Smith 614cd0e1443SSatish Balay 6155615d1e5SSatish Balay #undef __FUNC__ 6162d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6172d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 618cd0e1443SSatish Balay { 619cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6202d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6212d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6222d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6232d61bbb3SSatish Balay Scalar *ap, *aa = a->a, zero = 0.0; 624cd0e1443SSatish Balay 6253a40ed3dSBarry Smith PetscFunctionBegin; 6262d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 627cd0e1443SSatish Balay row = im[k]; brow = row/bs; 628a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 629a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6302d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6312c3acbe9SBarry Smith nrow = ailen[brow]; 6322d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 633a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 634a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6352d61bbb3SSatish Balay col = in[l] ; 6362d61bbb3SSatish Balay bcol = col/bs; 6372d61bbb3SSatish Balay cidx = col%bs; 6382d61bbb3SSatish Balay ridx = row%bs; 6392d61bbb3SSatish Balay high = nrow; 6402d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6412d61bbb3SSatish Balay while (high-low > 5) { 642cd0e1443SSatish Balay t = (low+high)/2; 643cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 644cd0e1443SSatish Balay else low = t; 645cd0e1443SSatish Balay } 646cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 647cd0e1443SSatish Balay if (rp[i] > bcol) break; 648cd0e1443SSatish Balay if (rp[i] == bcol) { 6492d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6502d61bbb3SSatish Balay goto finished; 651cd0e1443SSatish Balay } 652cd0e1443SSatish Balay } 6532d61bbb3SSatish Balay *v++ = zero; 6542d61bbb3SSatish Balay finished:; 655cd0e1443SSatish Balay } 656cd0e1443SSatish Balay } 6573a40ed3dSBarry Smith PetscFunctionReturn(0); 658cd0e1443SSatish Balay } 659cd0e1443SSatish Balay 6602d61bbb3SSatish Balay 6615615d1e5SSatish Balay #undef __FUNC__ 66205a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 66392c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 66492c4ed94SBarry Smith { 66592c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6668a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 667dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 668dd9472c6SBarry Smith int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval; 6690e324ae4SSatish Balay Scalar *ap,*value=v,*aa=a->a,*bap; 67092c4ed94SBarry Smith 6713a40ed3dSBarry Smith PetscFunctionBegin; 6720e324ae4SSatish Balay if (roworiented) { 6730e324ae4SSatish Balay stepval = (n-1)*bs; 6740e324ae4SSatish Balay } else { 6750e324ae4SSatish Balay stepval = (m-1)*bs; 6760e324ae4SSatish Balay } 67792c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 67892c4ed94SBarry Smith row = im[k]; 6793a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 680a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 681a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 68292c4ed94SBarry Smith #endif 68392c4ed94SBarry Smith rp = aj + ai[row]; 68492c4ed94SBarry Smith ap = aa + bs2*ai[row]; 68592c4ed94SBarry Smith rmax = imax[row]; 68692c4ed94SBarry Smith nrow = ailen[row]; 68792c4ed94SBarry Smith low = 0; 68892c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6893a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 690a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 691a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 69292c4ed94SBarry Smith #endif 69392c4ed94SBarry Smith col = in[l]; 69492c4ed94SBarry Smith if (roworiented) { 6950e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6960e324ae4SSatish Balay } else { 6970e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 69892c4ed94SBarry Smith } 69992c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 70092c4ed94SBarry Smith while (high-low > 7) { 70192c4ed94SBarry Smith t = (low+high)/2; 70292c4ed94SBarry Smith if (rp[t] > col) high = t; 70392c4ed94SBarry Smith else low = t; 70492c4ed94SBarry Smith } 70592c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 70692c4ed94SBarry Smith if (rp[i] > col) break; 70792c4ed94SBarry Smith if (rp[i] == col) { 7088a84c255SSatish Balay bap = ap + bs2*i; 7090e324ae4SSatish Balay if (roworiented) { 7108a84c255SSatish Balay if (is == ADD_VALUES) { 711dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 712dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7138a84c255SSatish Balay bap[jj] += *value++; 714dd9472c6SBarry Smith } 715dd9472c6SBarry Smith } 7160e324ae4SSatish Balay } else { 717dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 718dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7190e324ae4SSatish Balay bap[jj] = *value++; 7208a84c255SSatish Balay } 721dd9472c6SBarry Smith } 722dd9472c6SBarry Smith } 7230e324ae4SSatish Balay } else { 7240e324ae4SSatish Balay if (is == ADD_VALUES) { 725dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 726dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7270e324ae4SSatish Balay *bap++ += *value++; 728dd9472c6SBarry Smith } 729dd9472c6SBarry Smith } 7300e324ae4SSatish Balay } else { 731dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 732dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7330e324ae4SSatish Balay *bap++ = *value++; 7340e324ae4SSatish Balay } 7358a84c255SSatish Balay } 736dd9472c6SBarry Smith } 737dd9472c6SBarry Smith } 738f1241b54SBarry Smith goto noinsert2; 73992c4ed94SBarry Smith } 74092c4ed94SBarry Smith } 74189280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 742a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74392c4ed94SBarry Smith if (nrow >= rmax) { 74492c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 74592c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 74692c4ed94SBarry Smith Scalar *new_a; 74792c4ed94SBarry Smith 748a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74989280ab3SLois Curfman McInnes 75092c4ed94SBarry Smith /* malloc new storage space */ 75192c4ed94SBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 75292c4ed94SBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 75392c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 75492c4ed94SBarry Smith new_i = new_j + new_nz; 75592c4ed94SBarry Smith 75692c4ed94SBarry Smith /* copy over old data into new slots */ 75792c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 75892c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 75992c4ed94SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int)); 76092c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 761dd9472c6SBarry Smith PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int)); 76292c4ed94SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar)); 76392c4ed94SBarry Smith PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 764dd9472c6SBarry Smith PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar)); 76592c4ed94SBarry Smith /* free up old matrix storage */ 76692c4ed94SBarry Smith PetscFree(a->a); 76792c4ed94SBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 76892c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 76992c4ed94SBarry Smith a->singlemalloc = 1; 77092c4ed94SBarry Smith 77192c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 77292c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 77392c4ed94SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 77492c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 77592c4ed94SBarry Smith a->reallocs++; 77692c4ed94SBarry Smith a->nz++; 77792c4ed94SBarry Smith } 77892c4ed94SBarry Smith N = nrow++ - 1; 77992c4ed94SBarry Smith /* shift up all the later entries in this row */ 78092c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 78192c4ed94SBarry Smith rp[ii+1] = rp[ii]; 78292c4ed94SBarry Smith PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 78392c4ed94SBarry Smith } 78492c4ed94SBarry Smith if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 78592c4ed94SBarry Smith rp[i] = col; 7868a84c255SSatish Balay bap = ap + bs2*i; 7870e324ae4SSatish Balay if (roworiented) { 788dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 789dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7900e324ae4SSatish Balay bap[jj] = *value++; 791dd9472c6SBarry Smith } 792dd9472c6SBarry Smith } 7930e324ae4SSatish Balay } else { 794dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 795dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7960e324ae4SSatish Balay *bap++ = *value++; 7970e324ae4SSatish Balay } 798dd9472c6SBarry Smith } 799dd9472c6SBarry Smith } 800f1241b54SBarry Smith noinsert2:; 80192c4ed94SBarry Smith low = i; 80292c4ed94SBarry Smith } 80392c4ed94SBarry Smith ailen[row] = nrow; 80492c4ed94SBarry Smith } 8053a40ed3dSBarry Smith PetscFunctionReturn(0); 80692c4ed94SBarry Smith } 80792c4ed94SBarry Smith 808f2501298SSatish Balay 8095615d1e5SSatish Balay #undef __FUNC__ 8105615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8118f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 812584200bdSSatish Balay { 813584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 814584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 815a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 81643ee02c3SBarry Smith int mbs = a->mbs, bs2 = a->bs2,rmax = 0; 817584200bdSSatish Balay Scalar *aa = a->a, *ap; 818584200bdSSatish Balay 8193a40ed3dSBarry Smith PetscFunctionBegin; 8203a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 821584200bdSSatish Balay 82243ee02c3SBarry Smith if (m) rmax = ailen[0]; 823584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 824584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 825584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 826d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 827584200bdSSatish Balay if (fshift) { 828a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 829584200bdSSatish Balay N = ailen[i]; 830584200bdSSatish Balay for ( j=0; j<N; j++ ) { 831584200bdSSatish Balay ip[j-fshift] = ip[j]; 8327e67e3f9SSatish Balay PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar)); 833584200bdSSatish Balay } 834584200bdSSatish Balay } 835584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 836584200bdSSatish Balay } 837584200bdSSatish Balay if (mbs) { 838584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 839584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 840584200bdSSatish Balay } 841584200bdSSatish Balay /* reset ilen and imax for each row */ 842584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 843584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 844584200bdSSatish Balay } 845a7c10996SSatish Balay a->nz = ai[mbs]; 846584200bdSSatish Balay 847584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 848584200bdSSatish Balay if (fshift && a->diag) { 849584200bdSSatish Balay PetscFree(a->diag); 850584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 851584200bdSSatish Balay a->diag = 0; 852584200bdSSatish Balay } 8533d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 854219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8553d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 856584200bdSSatish Balay a->reallocs); 857d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 858e2f3b5e9SSatish Balay a->reallocs = 0; 8594e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8604e220ebcSLois Curfman McInnes 8613a40ed3dSBarry Smith PetscFunctionReturn(0); 862584200bdSSatish Balay } 863584200bdSSatish Balay 8645a838052SSatish Balay 865d9b7c43dSSatish Balay /* idx should be of length atlease bs */ 8665615d1e5SSatish Balay #undef __FUNC__ 8675615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block" 868d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg) 869d9b7c43dSSatish Balay { 870d9b7c43dSSatish Balay int i,row; 8713a40ed3dSBarry Smith 8723a40ed3dSBarry Smith PetscFunctionBegin; 873d9b7c43dSSatish Balay row = idx[0]; 8743a40ed3dSBarry Smith if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 875d9b7c43dSSatish Balay 876d9b7c43dSSatish Balay for ( i=1; i<bs; i++ ) { 8773a40ed3dSBarry Smith if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 878d9b7c43dSSatish Balay } 879d9b7c43dSSatish Balay *flg = PETSC_TRUE; 8803a40ed3dSBarry Smith PetscFunctionReturn(0); 881d9b7c43dSSatish Balay } 882d9b7c43dSSatish Balay 8835615d1e5SSatish Balay #undef __FUNC__ 8845615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 8858f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 886d9b7c43dSSatish Balay { 887d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 888d9b7c43dSSatish Balay IS is_local; 889d9b7c43dSSatish Balay int ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2; 890d9b7c43dSSatish Balay PetscTruth flg; 891d9b7c43dSSatish Balay Scalar zero = 0.0,*aa; 892d9b7c43dSSatish Balay 8933a40ed3dSBarry Smith PetscFunctionBegin; 894d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 895d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 896d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 897537820f0SBarry Smith ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr); 898d9b7c43dSSatish Balay ierr = ISSort(is_local); CHKERRQ(ierr); 899d9b7c43dSSatish Balay ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr); 900d9b7c43dSSatish Balay 901d9b7c43dSSatish Balay i = 0; 902d9b7c43dSSatish Balay while (i < is_n) { 903a8c6a408SBarry Smith if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 904d9b7c43dSSatish Balay flg = PETSC_FALSE; 905d9b7c43dSSatish Balay if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); } 906d9b7c43dSSatish Balay count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs; 907d9b7c43dSSatish Balay aa = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs); 9082d61bbb3SSatish Balay if (flg) { /* There exists a block of rows to be Zerowed */ 909a07cd24cSSatish Balay if (baij->ilen[rows[i]/bs] > 0) { 9102d61bbb3SSatish Balay PetscMemzero(aa,count*bs*sizeof(Scalar)); 9112d61bbb3SSatish Balay baij->ilen[rows[i]/bs] = 1; 9122d61bbb3SSatish Balay baij->j[baij->i[rows[i]/bs]] = rows[i]/bs; 913a07cd24cSSatish Balay } 9142d61bbb3SSatish Balay i += bs; 9152d61bbb3SSatish Balay } else { /* Zero out only the requested row */ 916d9b7c43dSSatish Balay for ( j=0; j<count; j++ ) { 917d9b7c43dSSatish Balay aa[0] = zero; 918d9b7c43dSSatish Balay aa+=bs; 919d9b7c43dSSatish Balay } 920d9b7c43dSSatish Balay i++; 921d9b7c43dSSatish Balay } 922d9b7c43dSSatish Balay } 923d9b7c43dSSatish Balay if (diag) { 924d9b7c43dSSatish Balay for ( j=0; j<is_n; j++ ) { 925f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 9262d61bbb3SSatish Balay /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */ 927d9b7c43dSSatish Balay } 928d9b7c43dSSatish Balay } 929d9b7c43dSSatish Balay ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr); 930d9b7c43dSSatish Balay ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr); 931d9b7c43dSSatish Balay ierr = ISDestroy(is_local); CHKERRQ(ierr); 9329a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 9333a40ed3dSBarry Smith PetscFunctionReturn(0); 934d9b7c43dSSatish Balay } 9351c351548SSatish Balay 9365615d1e5SSatish Balay #undef __FUNC__ 9372d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9382d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9392d61bbb3SSatish Balay { 9402d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9412d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9422d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9432d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 9442d61bbb3SSatish Balay int ridx,cidx,bs2=a->bs2; 9452d61bbb3SSatish Balay Scalar *ap,value,*aa=a->a,*bap; 9462d61bbb3SSatish Balay 9472d61bbb3SSatish Balay PetscFunctionBegin; 9482d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9492d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9502d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9512d61bbb3SSatish Balay if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 9522d61bbb3SSatish Balay if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 9532d61bbb3SSatish Balay #endif 9542d61bbb3SSatish Balay rp = aj + ai[brow]; 9552d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9562d61bbb3SSatish Balay rmax = imax[brow]; 9572d61bbb3SSatish Balay nrow = ailen[brow]; 9582d61bbb3SSatish Balay low = 0; 9592d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9602d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9612d61bbb3SSatish Balay if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 9622d61bbb3SSatish Balay if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 9632d61bbb3SSatish Balay #endif 9642d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 9652d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 9662d61bbb3SSatish Balay if (roworiented) { 9672d61bbb3SSatish Balay value = *v++; 9682d61bbb3SSatish Balay } else { 9692d61bbb3SSatish Balay value = v[k + l*m]; 9702d61bbb3SSatish Balay } 9712d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 9722d61bbb3SSatish Balay while (high-low > 7) { 9732d61bbb3SSatish Balay t = (low+high)/2; 9742d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 9752d61bbb3SSatish Balay else low = t; 9762d61bbb3SSatish Balay } 9772d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 9782d61bbb3SSatish Balay if (rp[i] > bcol) break; 9792d61bbb3SSatish Balay if (rp[i] == bcol) { 9802d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 9812d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 9822d61bbb3SSatish Balay else *bap = value; 9832d61bbb3SSatish Balay goto noinsert1; 9842d61bbb3SSatish Balay } 9852d61bbb3SSatish Balay } 9862d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 9872d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9882d61bbb3SSatish Balay if (nrow >= rmax) { 9892d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 9902d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 9912d61bbb3SSatish Balay Scalar *new_a; 9922d61bbb3SSatish Balay 9932d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9942d61bbb3SSatish Balay 9952d61bbb3SSatish Balay /* Malloc new storage space */ 9962d61bbb3SSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 9972d61bbb3SSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 9982d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 9992d61bbb3SSatish Balay new_i = new_j + new_nz; 10002d61bbb3SSatish Balay 10012d61bbb3SSatish Balay /* copy over old data into new slots */ 10022d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10032d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 10042d61bbb3SSatish Balay PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int)); 10052d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 10062d61bbb3SSatish Balay PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow, 10072d61bbb3SSatish Balay len*sizeof(int)); 10082d61bbb3SSatish Balay PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar)); 10092d61bbb3SSatish Balay PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 10102d61bbb3SSatish Balay PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE), 10112d61bbb3SSatish Balay aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar)); 10122d61bbb3SSatish Balay /* free up old matrix storage */ 10132d61bbb3SSatish Balay PetscFree(a->a); 10142d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 10152d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 10162d61bbb3SSatish Balay a->singlemalloc = 1; 10172d61bbb3SSatish Balay 10182d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10192d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10202d61bbb3SSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 10212d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10222d61bbb3SSatish Balay a->reallocs++; 10232d61bbb3SSatish Balay a->nz++; 10242d61bbb3SSatish Balay } 10252d61bbb3SSatish Balay N = nrow++ - 1; 10262d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10272d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10282d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 10292d61bbb3SSatish Balay PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 10302d61bbb3SSatish Balay } 10312d61bbb3SSatish Balay if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 10322d61bbb3SSatish Balay rp[i] = bcol; 10332d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10342d61bbb3SSatish Balay noinsert1:; 10352d61bbb3SSatish Balay low = i; 10362d61bbb3SSatish Balay } 10372d61bbb3SSatish Balay ailen[brow] = nrow; 10382d61bbb3SSatish Balay } 10392d61bbb3SSatish Balay PetscFunctionReturn(0); 10402d61bbb3SSatish Balay } 10412d61bbb3SSatish Balay 10422d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10432d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10442d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1045d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*); 10462d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**); 10472d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10482d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10492d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10502d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10512d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10522d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10532d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10542d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10552d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10562d61bbb3SSatish Balay 10572d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10582d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 10592d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 10602d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 10612d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 10622d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 10632d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 10642d61bbb3SSatish Balay 10652d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 10662d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 10672d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 10682d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 10692d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 10702d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 10712d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 10722d61bbb3SSatish Balay 10732d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 10742d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 10752d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 10762d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 10772d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 10782d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 10792d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 10802d61bbb3SSatish Balay 10812d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 10822d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 10832d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 10842d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 10852d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 10862d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 10872d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 10882d61bbb3SSatish Balay 10892d61bbb3SSatish Balay /* 10902d61bbb3SSatish Balay note: This can only work for identity for row and col. It would 10912d61bbb3SSatish Balay be good to check this and otherwise generate an error. 10922d61bbb3SSatish Balay */ 10932d61bbb3SSatish Balay #undef __FUNC__ 10942d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 10952d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 10962d61bbb3SSatish Balay { 10972d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 10982d61bbb3SSatish Balay Mat outA; 10992d61bbb3SSatish Balay int ierr; 11002d61bbb3SSatish Balay 11012d61bbb3SSatish Balay PetscFunctionBegin; 11022d61bbb3SSatish Balay if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 11032d61bbb3SSatish Balay 11042d61bbb3SSatish Balay outA = inA; 11052d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11062d61bbb3SSatish Balay a->row = row; 11072d61bbb3SSatish Balay a->col = col; 11082d61bbb3SSatish Balay 1109e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1110e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 11111e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1112e51c0b9cSSatish Balay 11132d61bbb3SSatish Balay if (!a->solve_work) { 11142d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 11152d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 11162d61bbb3SSatish Balay } 11172d61bbb3SSatish Balay 11182d61bbb3SSatish Balay if (!a->diag) { 11192d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 11202d61bbb3SSatish Balay } 11212d61bbb3SSatish Balay ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 11222d61bbb3SSatish Balay 11232d61bbb3SSatish Balay /* 11242d61bbb3SSatish Balay Blocksize 4 has a special faster solver for ILU(0) factorization 11252d61bbb3SSatish Balay with natural ordering 11262d61bbb3SSatish Balay */ 11272d61bbb3SSatish Balay if (a->bs == 4) { 1128f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 11292d61bbb3SSatish Balay } 11302d61bbb3SSatish Balay 11312d61bbb3SSatish Balay PetscFunctionReturn(0); 11322d61bbb3SSatish Balay } 11332d61bbb3SSatish Balay #undef __FUNC__ 1134d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1135ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1136ba4ca20aSSatish Balay { 1137ba4ca20aSSatish Balay static int called = 0; 1138ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1139ba4ca20aSSatish Balay 11403a40ed3dSBarry Smith PetscFunctionBegin; 11413a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 114276be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n"); 114376be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n"); 11443a40ed3dSBarry Smith PetscFunctionReturn(0); 1145ba4ca20aSSatish Balay } 1146d9b7c43dSSatish Balay 114727a8da17SBarry Smith #undef __FUNC__ 114827a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 114927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 115027a8da17SBarry Smith { 115127a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 115227a8da17SBarry Smith int i,nz,n; 115327a8da17SBarry Smith 115427a8da17SBarry Smith PetscFunctionBegin; 115527a8da17SBarry Smith nz = baij->maxnz; 115627a8da17SBarry Smith n = baij->n; 115727a8da17SBarry Smith for (i=0; i<nz; i++) { 115827a8da17SBarry Smith baij->j[i] = indices[i]; 115927a8da17SBarry Smith } 116027a8da17SBarry Smith baij->nz = nz; 116127a8da17SBarry Smith for ( i=0; i<n; i++ ) { 116227a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 116327a8da17SBarry Smith } 116427a8da17SBarry Smith 116527a8da17SBarry Smith PetscFunctionReturn(0); 116627a8da17SBarry Smith } 116727a8da17SBarry Smith 116827a8da17SBarry Smith #undef __FUNC__ 116927a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 117027a8da17SBarry Smith /*@ 117127a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 117227a8da17SBarry Smith in the matrix. 117327a8da17SBarry Smith 117427a8da17SBarry Smith Input Parameters: 117527a8da17SBarry Smith + mat - the SeqBAIJ matrix 117627a8da17SBarry Smith - indices - the column indices 117727a8da17SBarry Smith 117827a8da17SBarry Smith Notes: 117927a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 118027a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 118127a8da17SBarry Smith of the MatSetValues() operation. 118227a8da17SBarry Smith 118327a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 118427a8da17SBarry Smith MatCreateSeqBAIJ(). 118527a8da17SBarry Smith 118627a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 118727a8da17SBarry Smith 118827a8da17SBarry Smith @*/ 118927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 119027a8da17SBarry Smith { 119127a8da17SBarry Smith int ierr,(*f)(Mat,int *); 119227a8da17SBarry Smith 119327a8da17SBarry Smith PetscFunctionBegin; 119427a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 119527a8da17SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f); 119627a8da17SBarry Smith CHKERRQ(ierr); 119727a8da17SBarry Smith if (f) { 119827a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 119927a8da17SBarry Smith } else { 120027a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 120127a8da17SBarry Smith } 120227a8da17SBarry Smith PetscFunctionReturn(0); 120327a8da17SBarry Smith } 120427a8da17SBarry Smith 12052593348eSBarry Smith /* -------------------------------------------------------------------*/ 1206cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1207cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1208cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1209cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1210cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1211cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1212cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1213cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1214cc2dc46cSBarry Smith 0, 1215cc2dc46cSBarry Smith 0, 1216cc2dc46cSBarry Smith 0, 1217cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1218cc2dc46cSBarry Smith 0, 1219b6490206SBarry Smith 0, 1220f2501298SSatish Balay MatTranspose_SeqBAIJ, 1221cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1222cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1223cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1224cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1225cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1226b6490206SBarry Smith 0, 1227cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1228cc2dc46cSBarry Smith 0, 1229cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1230cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1231cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1232cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1233cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1234cc2dc46cSBarry Smith 0, 1235cc2dc46cSBarry Smith 0, 1236cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1237cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1238cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1239cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1240cc2dc46cSBarry Smith 0, 1241cc2dc46cSBarry Smith 0, 1242cc2dc46cSBarry Smith 0, 1243cc2dc46cSBarry Smith MatConvertSameType_SeqBAIJ, 1244cc2dc46cSBarry Smith 0, 1245cc2dc46cSBarry Smith 0, 1246cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1247cc2dc46cSBarry Smith 0, 1248cc2dc46cSBarry Smith 0, 1249cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1250cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1251cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1252cc2dc46cSBarry Smith 0, 1253cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1254cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1255cc2dc46cSBarry Smith 0, 1256cc2dc46cSBarry Smith 0, 1257cc2dc46cSBarry Smith 0, 1258cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 12593b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 126092c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1261cc2dc46cSBarry Smith 0, 1262cc2dc46cSBarry Smith 0, 1263cc2dc46cSBarry Smith 0, 1264cc2dc46cSBarry Smith 0, 1265cc2dc46cSBarry Smith 0, 1266cc2dc46cSBarry Smith 0, 1267d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1268cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1269cc2dc46cSBarry Smith 0, 1270cc2dc46cSBarry Smith 0, 1271cc2dc46cSBarry Smith MatGetMaps_Petsc}; 12722593348eSBarry Smith 12735615d1e5SSatish Balay #undef __FUNC__ 12745615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 12752593348eSBarry Smith /*@C 127644cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 127744cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 127844cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 12792bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 12802bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 12812593348eSBarry Smith 1282db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1283db81eaa0SLois Curfman McInnes 12842593348eSBarry Smith Input Parameters: 1285db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1286b6490206SBarry Smith . bs - size of block 12872593348eSBarry Smith . m - number of rows 12882593348eSBarry Smith . n - number of columns 1289b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 1290db81eaa0SLois Curfman McInnes - nzz - array containing the number of block nonzeros in the various block rows 12912bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 12922593348eSBarry Smith 12932593348eSBarry Smith Output Parameter: 12942593348eSBarry Smith . A - the matrix 12952593348eSBarry Smith 12960513a670SBarry Smith Options Database Keys: 1297db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1298db81eaa0SLois Curfman McInnes block calculations (much slower) 1299db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 13000513a670SBarry Smith 13012593348eSBarry Smith Notes: 130244cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 13032593348eSBarry Smith storage. That is, the stored row and column indices can begin at 130444cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 13052593348eSBarry Smith 13062593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 13072593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 13082593348eSBarry Smith allocation. For additional details, see the users manual chapter on 13096da5968aSLois Curfman McInnes matrices. 13102593348eSBarry Smith 1311db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 13122593348eSBarry Smith @*/ 1313b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 13142593348eSBarry Smith { 13152593348eSBarry Smith Mat B; 1316b6490206SBarry Smith Mat_SeqBAIJ *b; 13173b2fbd54SBarry Smith int i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size; 13183b2fbd54SBarry Smith 13193a40ed3dSBarry Smith PetscFunctionBegin; 13203b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1321a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1322b6490206SBarry Smith 13230513a670SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 13240513a670SBarry Smith 13253a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1326a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 13273a40ed3dSBarry Smith } 13282593348eSBarry Smith 13292593348eSBarry Smith *A = 0; 1330f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView); 13312593348eSBarry Smith PLogObjectCreate(B); 1332b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 133344cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqBAIJ)); 1334cc2dc46cSBarry Smith PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)); 13351a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr); 13361a6a6d98SBarry Smith if (!flg) { 13377fc0212eSBarry Smith switch (bs) { 13387fc0212eSBarry Smith case 1: 1339f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1340f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1341f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1342f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 13437fc0212eSBarry Smith break; 13444eeb42bcSBarry Smith case 2: 1345f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1346f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1347f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1348f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 13496d84be18SBarry Smith break; 1350f327dd97SBarry Smith case 3: 1351f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1352f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1353f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1354f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 13554eeb42bcSBarry Smith break; 13566d84be18SBarry Smith case 4: 1357f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1358f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1359f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1360f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 13616d84be18SBarry Smith break; 13626d84be18SBarry Smith case 5: 1363f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1364f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1365f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1366f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 13676d84be18SBarry Smith break; 136848e9ddb2SSatish Balay case 7: 1369f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1370f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1371f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 137248e9ddb2SSatish Balay break; 13737fc0212eSBarry Smith } 13741a6a6d98SBarry Smith } 1375e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1376e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 13772593348eSBarry Smith B->factor = 0; 13782593348eSBarry Smith B->lupivotthreshold = 1.0; 137990f02eecSBarry Smith B->mapping = 0; 13802593348eSBarry Smith b->row = 0; 13812593348eSBarry Smith b->col = 0; 1382e51c0b9cSSatish Balay b->icol = 0; 13832593348eSBarry Smith b->reallocs = 0; 13842593348eSBarry Smith 138544cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 138644cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1387a5ae1ecdSBarry Smith 1388*7413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 1389*7413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1390a5ae1ecdSBarry Smith 1391b6490206SBarry Smith b->mbs = mbs; 1392f2501298SSatish Balay b->nbs = nbs; 1393b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 13942593348eSBarry Smith if (nnz == PETSC_NULL) { 1395119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 13962593348eSBarry Smith else if (nz <= 0) nz = 1; 1397b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1398b6490206SBarry Smith nz = nz*mbs; 13993a40ed3dSBarry Smith } else { 14002593348eSBarry Smith nz = 0; 1401b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 14022593348eSBarry Smith } 14032593348eSBarry Smith 14042593348eSBarry Smith /* allocate the matrix space */ 14057e67e3f9SSatish Balay len = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int); 14062593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 14077e67e3f9SSatish Balay PetscMemzero(b->a,nz*bs2*sizeof(Scalar)); 14087e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 14092593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 14102593348eSBarry Smith b->i = b->j + nz; 14112593348eSBarry Smith b->singlemalloc = 1; 14122593348eSBarry Smith 1413de6a44a3SBarry Smith b->i[0] = 0; 1414b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 14152593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 14162593348eSBarry Smith } 14172593348eSBarry Smith 1418b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1419b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1420f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1421b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 14222593348eSBarry Smith 1423b6490206SBarry Smith b->bs = bs; 14249df24120SSatish Balay b->bs2 = bs2; 1425b6490206SBarry Smith b->mbs = mbs; 14262593348eSBarry Smith b->nz = 0; 142718e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 14282593348eSBarry Smith b->sorted = 0; 14292593348eSBarry Smith b->roworiented = 1; 14302593348eSBarry Smith b->nonew = 0; 14312593348eSBarry Smith b->diag = 0; 14322593348eSBarry Smith b->solve_work = 0; 1433de6a44a3SBarry Smith b->mult_work = 0; 14342593348eSBarry Smith b->spptr = 0; 14354e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 14364e220ebcSLois Curfman McInnes 14372593348eSBarry Smith *A = B; 14382593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 14392593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 144027a8da17SBarry Smith 144127a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 144227a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 144327a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 144427a8da17SBarry Smith 14453a40ed3dSBarry Smith PetscFunctionReturn(0); 14462593348eSBarry Smith } 14472593348eSBarry Smith 14485615d1e5SSatish Balay #undef __FUNC__ 14495615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ" 1450b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 14512593348eSBarry Smith { 14522593348eSBarry Smith Mat C; 1453b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 145427a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1455de6a44a3SBarry Smith 14563a40ed3dSBarry Smith PetscFunctionBegin; 1457a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 14582593348eSBarry Smith 14592593348eSBarry Smith *B = 0; 1460f830108cSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView); 14612593348eSBarry Smith PLogObjectCreate(C); 1462b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 1463c3c66ee5SSatish Balay PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 1464e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1465e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 14662593348eSBarry Smith C->factor = A->factor; 14672593348eSBarry Smith c->row = 0; 14682593348eSBarry Smith c->col = 0; 14692593348eSBarry Smith C->assembled = PETSC_TRUE; 14702593348eSBarry Smith 147144cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 147244cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 147344cd7ae7SLois Curfman McInnes C->M = a->m; 147444cd7ae7SLois Curfman McInnes C->N = a->n; 147544cd7ae7SLois Curfman McInnes 1476b6490206SBarry Smith c->bs = a->bs; 14779df24120SSatish Balay c->bs2 = a->bs2; 1478b6490206SBarry Smith c->mbs = a->mbs; 1479b6490206SBarry Smith c->nbs = a->nbs; 14802593348eSBarry Smith 1481b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1482b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1483b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 14842593348eSBarry Smith c->imax[i] = a->imax[i]; 14852593348eSBarry Smith c->ilen[i] = a->ilen[i]; 14862593348eSBarry Smith } 14872593348eSBarry Smith 14882593348eSBarry Smith /* allocate the matrix space */ 14892593348eSBarry Smith c->singlemalloc = 1; 14907e67e3f9SSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int)); 14912593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 14927e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1493de6a44a3SBarry Smith c->i = c->j + nz; 1494b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1495b6490206SBarry Smith if (mbs > 0) { 1496de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 14972593348eSBarry Smith if (cpvalues == COPY_VALUES) { 14987e67e3f9SSatish Balay PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar)); 14992593348eSBarry Smith } 15002593348eSBarry Smith } 15012593348eSBarry Smith 1502f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 15032593348eSBarry Smith c->sorted = a->sorted; 15042593348eSBarry Smith c->roworiented = a->roworiented; 15052593348eSBarry Smith c->nonew = a->nonew; 15062593348eSBarry Smith 15072593348eSBarry Smith if (a->diag) { 1508b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1509b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1510b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 15112593348eSBarry Smith c->diag[i] = a->diag[i]; 15122593348eSBarry Smith } 15132593348eSBarry Smith } 15142593348eSBarry Smith else c->diag = 0; 15152593348eSBarry Smith c->nz = a->nz; 15162593348eSBarry Smith c->maxnz = a->maxnz; 15172593348eSBarry Smith c->solve_work = 0; 15182593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 15197fc0212eSBarry Smith c->mult_work = 0; 15202593348eSBarry Smith *B = C; 152127a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqBAIJSetColumnIndices_C", 152227a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 152327a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 15243a40ed3dSBarry Smith PetscFunctionReturn(0); 15252593348eSBarry Smith } 15262593348eSBarry Smith 15275615d1e5SSatish Balay #undef __FUNC__ 15285615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 152919bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 15302593348eSBarry Smith { 1531b6490206SBarry Smith Mat_SeqBAIJ *a; 15322593348eSBarry Smith Mat B; 1533de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1534b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 153535aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1536a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1537b6490206SBarry Smith Scalar *aa; 153819bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 15392593348eSBarry Smith 15403a40ed3dSBarry Smith PetscFunctionBegin; 15411a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1542de6a44a3SBarry Smith bs2 = bs*bs; 1543b6490206SBarry Smith 15442593348eSBarry Smith MPI_Comm_size(comm,&size); 1545cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 154690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 15470752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 1548a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 15492593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 15502593348eSBarry Smith 1551d64ed03dSBarry Smith if (header[3] < 0) { 1552a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1553d64ed03dSBarry Smith } 1554d64ed03dSBarry Smith 1555a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 155635aab85fSBarry Smith 155735aab85fSBarry Smith /* 155835aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 155935aab85fSBarry Smith divisible by the blocksize 156035aab85fSBarry Smith */ 1561b6490206SBarry Smith mbs = M/bs; 156235aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 156335aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 156435aab85fSBarry Smith else mbs++; 156535aab85fSBarry Smith if (extra_rows) { 1566537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 156735aab85fSBarry Smith } 1568b6490206SBarry Smith 15692593348eSBarry Smith /* read in row lengths */ 157035aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 15710752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 157235aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 15732593348eSBarry Smith 1574b6490206SBarry Smith /* read in column indices */ 157535aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 15760752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr); 157735aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1578b6490206SBarry Smith 1579b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1580b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1581b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 158235aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 158335aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 158435aab85fSBarry Smith masked = mask + mbs; 1585b6490206SBarry Smith rowcount = 0; nzcount = 0; 1586b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 158735aab85fSBarry Smith nmask = 0; 1588b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1589b6490206SBarry Smith kmax = rowlengths[rowcount]; 1590b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 159135aab85fSBarry Smith tmp = jj[nzcount++]/bs; 159235aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1593b6490206SBarry Smith } 1594b6490206SBarry Smith rowcount++; 1595b6490206SBarry Smith } 159635aab85fSBarry Smith browlengths[i] += nmask; 159735aab85fSBarry Smith /* zero out the mask elements we set */ 159835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1599b6490206SBarry Smith } 1600b6490206SBarry Smith 16012593348eSBarry Smith /* create our matrix */ 16023eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 16032593348eSBarry Smith B = *A; 1604b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 16052593348eSBarry Smith 16062593348eSBarry Smith /* set matrix "i" values */ 1607de6a44a3SBarry Smith a->i[0] = 0; 1608b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1609b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1610b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 16112593348eSBarry Smith } 1612b6490206SBarry Smith a->nz = 0; 1613b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 16142593348eSBarry Smith 1615b6490206SBarry Smith /* read in nonzero values */ 161635aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 16170752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr); 161835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1619b6490206SBarry Smith 1620b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1621b6490206SBarry Smith nzcount = 0; jcount = 0; 1622b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1623b6490206SBarry Smith nzcountb = nzcount; 162435aab85fSBarry Smith nmask = 0; 1625b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1626b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1627b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 162835aab85fSBarry Smith tmp = jj[nzcount++]/bs; 162935aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1630b6490206SBarry Smith } 1631b6490206SBarry Smith rowcount++; 1632b6490206SBarry Smith } 1633de6a44a3SBarry Smith /* sort the masked values */ 163477c4ece6SBarry Smith PetscSortInt(nmask,masked); 1635de6a44a3SBarry Smith 1636b6490206SBarry Smith /* set "j" values into matrix */ 1637b6490206SBarry Smith maskcount = 1; 163835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 163935aab85fSBarry Smith a->j[jcount++] = masked[j]; 1640de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1641b6490206SBarry Smith } 1642b6490206SBarry Smith /* set "a" values into matrix */ 1643de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1644b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1645b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1646b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1647de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1648de6a44a3SBarry Smith block = mask[tmp] - 1; 1649de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1650de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1651b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1652b6490206SBarry Smith } 1653b6490206SBarry Smith } 165435aab85fSBarry Smith /* zero out the mask elements we set */ 165535aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1656b6490206SBarry Smith } 1657a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1658b6490206SBarry Smith 1659b6490206SBarry Smith PetscFree(rowlengths); 1660b6490206SBarry Smith PetscFree(browlengths); 1661b6490206SBarry Smith PetscFree(aa); 1662b6490206SBarry Smith PetscFree(jj); 1663b6490206SBarry Smith PetscFree(mask); 1664b6490206SBarry Smith 1665b6490206SBarry Smith B->assembled = PETSC_TRUE; 1666de6a44a3SBarry Smith 16679c01be13SBarry Smith ierr = MatView_Private(B); CHKERRQ(ierr); 16683a40ed3dSBarry Smith PetscFunctionReturn(0); 16692593348eSBarry Smith } 16702593348eSBarry Smith 16712593348eSBarry Smith 16722593348eSBarry Smith 1673