1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3f1db9ecSBarry Smith static char vcid[] = "$Id: baij.c,v 1.150 1998/12/06 16:38:10 balay 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 #include "sys.h" 1070f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 111a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 121a6a6d98SBarry Smith #include "src/inline/spops.h" 133b547af2SSatish Balay 142d61bbb3SSatish Balay #define CHUNKSIZE 10 15de6a44a3SBarry Smith 165615d1e5SSatish Balay #undef __FUNC__ 175615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ" 18de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 19de6a44a3SBarry Smith { 20de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 217fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 22de6a44a3SBarry Smith 233a40ed3dSBarry Smith PetscFunctionBegin; 24de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 25de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 267fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 27ecc615b2SBarry Smith diag[i] = a->i[i+1]; 28de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 29de6a44a3SBarry Smith if (a->j[j] == i) { 30de6a44a3SBarry Smith diag[i] = j; 31de6a44a3SBarry Smith break; 32de6a44a3SBarry Smith } 33de6a44a3SBarry Smith } 34de6a44a3SBarry Smith } 35de6a44a3SBarry Smith a->diag = diag; 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37de6a44a3SBarry Smith } 382593348eSBarry Smith 392593348eSBarry Smith 403b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 413b2fbd54SBarry Smith 425615d1e5SSatish Balay #undef __FUNC__ 435615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 443b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 453b2fbd54SBarry Smith PetscTruth *done) 463b2fbd54SBarry Smith { 473b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 483b2fbd54SBarry Smith int ierr, n = a->mbs,i; 493b2fbd54SBarry Smith 503a40ed3dSBarry Smith PetscFunctionBegin; 513b2fbd54SBarry Smith *nn = n; 523a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 533b2fbd54SBarry Smith if (symmetric) { 543b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 553b2fbd54SBarry Smith } else if (oshift == 1) { 563b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 573b2fbd54SBarry Smith int nz = a->i[n] + 1; 583b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 593b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 603b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 613b2fbd54SBarry Smith } else { 623b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 633b2fbd54SBarry Smith } 643b2fbd54SBarry Smith 653a40ed3dSBarry Smith PetscFunctionReturn(0); 663b2fbd54SBarry Smith } 673b2fbd54SBarry Smith 685615d1e5SSatish Balay #undef __FUNC__ 69d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 703b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 713b2fbd54SBarry Smith PetscTruth *done) 723b2fbd54SBarry Smith { 733b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 743b2fbd54SBarry Smith int i,n = a->mbs; 753b2fbd54SBarry Smith 763a40ed3dSBarry Smith PetscFunctionBegin; 773a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 783b2fbd54SBarry Smith if (symmetric) { 793b2fbd54SBarry Smith PetscFree(*ia); 803b2fbd54SBarry Smith PetscFree(*ja); 81af5da2bfSBarry Smith } else if (oshift == 1) { 823b2fbd54SBarry Smith int nz = a->i[n]; 833b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 843b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 853b2fbd54SBarry Smith } 863a40ed3dSBarry Smith PetscFunctionReturn(0); 873b2fbd54SBarry Smith } 883b2fbd54SBarry Smith 892d61bbb3SSatish Balay #undef __FUNC__ 902d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 912d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 922d61bbb3SSatish Balay { 932d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 942d61bbb3SSatish Balay 952d61bbb3SSatish Balay PetscFunctionBegin; 962d61bbb3SSatish Balay *bs = baij->bs; 972d61bbb3SSatish Balay PetscFunctionReturn(0); 982d61bbb3SSatish Balay } 992d61bbb3SSatish Balay 1002d61bbb3SSatish Balay 1012d61bbb3SSatish Balay #undef __FUNC__ 1022d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 103e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1042d61bbb3SSatish Balay { 1052d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 106e51c0b9cSSatish Balay int ierr; 1072d61bbb3SSatish Balay 10894d884c6SBarry Smith if (--A->refct > 0) PetscFunctionReturn(0); 10994d884c6SBarry Smith 11094d884c6SBarry Smith if (A->mapping) { 11194d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 11294d884c6SBarry Smith } 11394d884c6SBarry Smith if (A->bmapping) { 11494d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr); 11594d884c6SBarry Smith } 11661b13de0SBarry Smith if (A->rmap) { 11761b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 11861b13de0SBarry Smith } 11961b13de0SBarry Smith if (A->cmap) { 12061b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 12161b13de0SBarry Smith } 1222d61bbb3SSatish Balay #if defined(USE_PETSC_LOG) 123e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1242d61bbb3SSatish Balay #endif 1252d61bbb3SSatish Balay PetscFree(a->a); 1262d61bbb3SSatish Balay if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 1272d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 1282d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 1292d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 1302d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 1312d61bbb3SSatish Balay if (a->mult_work) PetscFree(a->mult_work); 132e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 1332d61bbb3SSatish Balay PetscFree(a); 1342d61bbb3SSatish Balay PLogObjectDestroy(A); 1352d61bbb3SSatish Balay PetscHeaderDestroy(A); 1362d61bbb3SSatish Balay PetscFunctionReturn(0); 1372d61bbb3SSatish Balay } 1382d61bbb3SSatish Balay 1392d61bbb3SSatish Balay #undef __FUNC__ 1402d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1412d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1422d61bbb3SSatish Balay { 1432d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1442d61bbb3SSatish Balay 1452d61bbb3SSatish Balay PetscFunctionBegin; 1462d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1472d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1482d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1492d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1502d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1512d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 1522d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 1532d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1542d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1552d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1562d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1572d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1582d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 159b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 160b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 161981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1622d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1632d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1642d61bbb3SSatish Balay } else { 1652d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1662d61bbb3SSatish Balay } 1672d61bbb3SSatish Balay PetscFunctionReturn(0); 1682d61bbb3SSatish Balay } 1692d61bbb3SSatish Balay 1702d61bbb3SSatish Balay 1712d61bbb3SSatish Balay #undef __FUNC__ 1722d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1732d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1742d61bbb3SSatish Balay { 1752d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1762d61bbb3SSatish Balay 1772d61bbb3SSatish Balay PetscFunctionBegin; 1782d61bbb3SSatish Balay if (m) *m = a->m; 1792d61bbb3SSatish Balay if (n) *n = a->n; 1802d61bbb3SSatish Balay PetscFunctionReturn(0); 1812d61bbb3SSatish Balay } 1822d61bbb3SSatish Balay 1832d61bbb3SSatish Balay #undef __FUNC__ 1842d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 1852d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 1862d61bbb3SSatish Balay { 1872d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1882d61bbb3SSatish Balay 1892d61bbb3SSatish Balay PetscFunctionBegin; 1902d61bbb3SSatish Balay *m = 0; *n = a->m; 1912d61bbb3SSatish Balay PetscFunctionReturn(0); 1922d61bbb3SSatish Balay } 1932d61bbb3SSatish Balay 1942d61bbb3SSatish Balay #undef __FUNC__ 1952d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 1962d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 1972d61bbb3SSatish Balay { 1982d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1992d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 200*3f1db9ecSBarry Smith MatScalar *aa,*aa_i; 201*3f1db9ecSBarry Smith Scalar *v_i; 2022d61bbb3SSatish Balay 2032d61bbb3SSatish Balay PetscFunctionBegin; 2042d61bbb3SSatish Balay bs = a->bs; 2052d61bbb3SSatish Balay ai = a->i; 2062d61bbb3SSatish Balay aj = a->j; 2072d61bbb3SSatish Balay aa = a->a; 2082d61bbb3SSatish Balay bs2 = a->bs2; 2092d61bbb3SSatish Balay 2102d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2112d61bbb3SSatish Balay 2122d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2132d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2142d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2152d61bbb3SSatish Balay *nz = bs*M; 2162d61bbb3SSatish Balay 2172d61bbb3SSatish Balay if (v) { 2182d61bbb3SSatish Balay *v = 0; 2192d61bbb3SSatish Balay if (*nz) { 2202d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 2212d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2222d61bbb3SSatish Balay v_i = *v + i*bs; 2232d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2242d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2252d61bbb3SSatish Balay } 2262d61bbb3SSatish Balay } 2272d61bbb3SSatish Balay } 2282d61bbb3SSatish Balay 2292d61bbb3SSatish Balay if (idx) { 2302d61bbb3SSatish Balay *idx = 0; 2312d61bbb3SSatish Balay if (*nz) { 2322d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 2332d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2342d61bbb3SSatish Balay idx_i = *idx + i*bs; 2352d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2362d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2372d61bbb3SSatish Balay } 2382d61bbb3SSatish Balay } 2392d61bbb3SSatish Balay } 2402d61bbb3SSatish Balay PetscFunctionReturn(0); 2412d61bbb3SSatish Balay } 2422d61bbb3SSatish Balay 2432d61bbb3SSatish Balay #undef __FUNC__ 2442d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2452d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2462d61bbb3SSatish Balay { 2472d61bbb3SSatish Balay PetscFunctionBegin; 2482d61bbb3SSatish Balay if (idx) {if (*idx) PetscFree(*idx);} 2492d61bbb3SSatish Balay if (v) {if (*v) PetscFree(*v);} 2502d61bbb3SSatish Balay PetscFunctionReturn(0); 2512d61bbb3SSatish Balay } 2522d61bbb3SSatish Balay 2532d61bbb3SSatish Balay #undef __FUNC__ 2542d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2552d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2562d61bbb3SSatish Balay { 2572d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2582d61bbb3SSatish Balay Mat C; 2592d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2602d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 261*3f1db9ecSBarry Smith MatScalar *array = a->a; 2622d61bbb3SSatish Balay 2632d61bbb3SSatish Balay PetscFunctionBegin; 2642d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2652d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col); 2662d61bbb3SSatish Balay PetscMemzero(col,(1+nbs)*sizeof(int)); 2672d61bbb3SSatish Balay 2682d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2692d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr); 2702d61bbb3SSatish Balay PetscFree(col); 2712d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows); 2722d61bbb3SSatish Balay cols = rows + bs; 2732d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2742d61bbb3SSatish Balay cols[0] = i*bs; 2752d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 2762d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 2772d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 2782d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 2792d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 2802d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr); 2812d61bbb3SSatish Balay array += bs2; 2822d61bbb3SSatish Balay } 2832d61bbb3SSatish Balay } 2842d61bbb3SSatish Balay PetscFree(rows); 2852d61bbb3SSatish Balay 2862d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2872d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2882d61bbb3SSatish Balay 2892d61bbb3SSatish Balay if (B != PETSC_NULL) { 2902d61bbb3SSatish Balay *B = C; 2912d61bbb3SSatish Balay } else { 292f830108cSBarry Smith PetscOps *Abops; 293cc2dc46cSBarry Smith MatOps Aops; 294f830108cSBarry Smith 2952d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 2962d61bbb3SSatish Balay PetscFree(a->a); 2972d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 2982d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 2992d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 3002d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 3012d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 3022d61bbb3SSatish Balay PetscFree(a); 303f830108cSBarry Smith 3047413bad6SBarry Smith 3057413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3067413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3077413bad6SBarry Smith 308f830108cSBarry Smith /* 309f830108cSBarry Smith This is horrible, horrible code. We need to keep the 310f830108cSBarry Smith A pointers for the bops and ops but copy everything 311f830108cSBarry Smith else from C. 312f830108cSBarry Smith */ 313f830108cSBarry Smith Abops = A->bops; 314f830108cSBarry Smith Aops = A->ops; 3152d61bbb3SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 316f830108cSBarry Smith A->bops = Abops; 317f830108cSBarry Smith A->ops = Aops; 31827a8da17SBarry Smith A->qlist = 0; 319f830108cSBarry Smith 3202d61bbb3SSatish Balay PetscHeaderDestroy(C); 3212d61bbb3SSatish Balay } 3222d61bbb3SSatish Balay PetscFunctionReturn(0); 3232d61bbb3SSatish Balay } 3242d61bbb3SSatish Balay 3252d61bbb3SSatish Balay 3262d61bbb3SSatish Balay 3273b2fbd54SBarry Smith 3285615d1e5SSatish Balay #undef __FUNC__ 329d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 330b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3312593348eSBarry Smith { 332b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3339df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 334b6490206SBarry Smith Scalar *aa; 3352593348eSBarry Smith 3363a40ed3dSBarry Smith PetscFunctionBegin; 33790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 3382593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3392593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3403b2fbd54SBarry Smith 3412593348eSBarry Smith col_lens[1] = a->m; 3422593348eSBarry Smith col_lens[2] = a->n; 3437e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3442593348eSBarry Smith 3452593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 346b6490206SBarry Smith count = 0; 347b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 348b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 349b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 350b6490206SBarry Smith } 3512593348eSBarry Smith } 3520752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3532593348eSBarry Smith PetscFree(col_lens); 3542593348eSBarry Smith 3552593348eSBarry Smith /* store column indices (zero start index) */ 35666963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj); 357b6490206SBarry Smith count = 0; 358b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 359b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 360b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 361b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 362b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3632593348eSBarry Smith } 3642593348eSBarry Smith } 365b6490206SBarry Smith } 366b6490206SBarry Smith } 3670752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr); 368b6490206SBarry Smith PetscFree(jj); 3692593348eSBarry Smith 3702593348eSBarry Smith /* store nonzero values */ 37166963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa); 372b6490206SBarry Smith count = 0; 373b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 374b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 375b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 376b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 3777e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 378b6490206SBarry Smith } 379b6490206SBarry Smith } 380b6490206SBarry Smith } 381b6490206SBarry Smith } 3820752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 383b6490206SBarry Smith PetscFree(aa); 3843a40ed3dSBarry Smith PetscFunctionReturn(0); 3852593348eSBarry Smith } 3862593348eSBarry Smith 3875615d1e5SSatish Balay #undef __FUNC__ 388d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 389b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 3902593348eSBarry Smith { 391b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3929df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 3932593348eSBarry Smith FILE *fd; 3942593348eSBarry Smith char *outputname; 3952593348eSBarry Smith 3963a40ed3dSBarry Smith PetscFunctionBegin; 39790ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 39877ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 39990ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 400639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 4017ddc982cSLois Curfman McInnes fprintf(fd," block size is %d\n",bs); 4022593348eSBarry Smith } 403639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 404a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 4052593348eSBarry Smith } 406639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 40744cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 40844cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 40944cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 41044cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 41144cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 4123a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 413e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) 41444cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 415e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 416e20fef11SSatish Balay else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) 417766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 418e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 419e20fef11SSatish Balay else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) 420e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 42144cd7ae7SLois Curfman McInnes #else 42244cd7ae7SLois Curfman McInnes if (a->a[bs2*k + l*bs + j] != 0.0) 42344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 42444cd7ae7SLois Curfman McInnes #endif 42544cd7ae7SLois Curfman McInnes } 42644cd7ae7SLois Curfman McInnes } 42744cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 42844cd7ae7SLois Curfman McInnes } 42944cd7ae7SLois Curfman McInnes } 43044cd7ae7SLois Curfman McInnes } 4312593348eSBarry Smith else { 432b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 433b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 434b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 435b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 436b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4373a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 438e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 43988685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 440e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 44188685aaeSLois Curfman McInnes } 442e20fef11SSatish Balay else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 443766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 444e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 445766eeae4SLois Curfman McInnes } 44688685aaeSLois Curfman McInnes else { 447e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 44888685aaeSLois Curfman McInnes } 44988685aaeSLois Curfman McInnes #else 4507e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 45188685aaeSLois Curfman McInnes #endif 4522593348eSBarry Smith } 4532593348eSBarry Smith } 4542593348eSBarry Smith fprintf(fd,"\n"); 4552593348eSBarry Smith } 4562593348eSBarry Smith } 457b6490206SBarry Smith } 4582593348eSBarry Smith fflush(fd); 4593a40ed3dSBarry Smith PetscFunctionReturn(0); 4602593348eSBarry Smith } 4612593348eSBarry Smith 4625615d1e5SSatish Balay #undef __FUNC__ 46377ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 46477ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4653270192aSSatish Balay { 46677ed5343SBarry Smith Mat A = (Mat) Aa; 4673270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4687dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 469fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 470*3f1db9ecSBarry Smith MatScalar *aa; 47177ed5343SBarry Smith MPI_Comm comm; 47277ed5343SBarry Smith Viewer viewer; 4733270192aSSatish Balay 4743a40ed3dSBarry Smith PetscFunctionBegin; 47577ed5343SBarry Smith /* 47677ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 47777ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 47877ed5343SBarry Smith rest should return immediately. 47977ed5343SBarry Smith */ 48077ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 48177ed5343SBarry Smith MPI_Comm_rank(comm,&rank); 48277ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 4833270192aSSatish Balay 48477ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 48577ed5343SBarry Smith 48677ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr); 48777ed5343SBarry Smith 4883270192aSSatish Balay /* loop over matrix elements drawing boxes */ 4893270192aSSatish Balay color = DRAW_BLUE; 4903270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4913270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4923270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4933270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4943270192aSSatish Balay aa = a->a + j*bs2; 4953270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4963270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4973270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 498b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4993270192aSSatish Balay } 5003270192aSSatish Balay } 5013270192aSSatish Balay } 5023270192aSSatish Balay } 5033270192aSSatish Balay color = DRAW_CYAN; 5043270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5053270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5063270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5073270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5083270192aSSatish Balay aa = a->a + j*bs2; 5093270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5103270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5113270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 512b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5133270192aSSatish Balay } 5143270192aSSatish Balay } 5153270192aSSatish Balay } 5163270192aSSatish Balay } 5173270192aSSatish Balay 5183270192aSSatish Balay color = DRAW_RED; 5193270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5203270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5213270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5223270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5233270192aSSatish Balay aa = a->a + j*bs2; 5243270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5253270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5263270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 527b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5283270192aSSatish Balay } 5293270192aSSatish Balay } 5303270192aSSatish Balay } 5313270192aSSatish Balay } 53277ed5343SBarry Smith PetscFunctionReturn(0); 53377ed5343SBarry Smith } 5343270192aSSatish Balay 53577ed5343SBarry Smith #undef __FUNC__ 53677ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 53777ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 53877ed5343SBarry Smith { 53977ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5407dce120fSSatish Balay int ierr; 5417dce120fSSatish Balay double xl,yl,xr,yr,w,h; 54277ed5343SBarry Smith Draw draw; 54377ed5343SBarry Smith PetscTruth isnull; 5443270192aSSatish Balay 54577ed5343SBarry Smith PetscFunctionBegin; 54677ed5343SBarry Smith 54777ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 54877ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 54977ed5343SBarry Smith 55077ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 55177ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 55277ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5533270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 55477ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A); CHKERRQ(ierr); 55577ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5563a40ed3dSBarry Smith PetscFunctionReturn(0); 5573270192aSSatish Balay } 5583270192aSSatish Balay 5595615d1e5SSatish Balay #undef __FUNC__ 560d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 561e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5622593348eSBarry Smith { 56319bcc07fSBarry Smith ViewerType vtype; 56419bcc07fSBarry Smith int ierr; 5652593348eSBarry Smith 5663a40ed3dSBarry Smith PetscFunctionBegin; 56719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 568*3f1db9ecSBarry Smith if (PetscTypeCompare(vtype,MATLAB_VIEWER)) { 569a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported"); 570*3f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){ 5713a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 572*3f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5733a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 574*3f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 5753a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 5765cd90555SBarry Smith } else { 5775cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5782593348eSBarry Smith } 5793a40ed3dSBarry Smith PetscFunctionReturn(0); 5802593348eSBarry Smith } 581b6490206SBarry Smith 582cd0e1443SSatish Balay 5835615d1e5SSatish Balay #undef __FUNC__ 5842d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 5852d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 586cd0e1443SSatish Balay { 587cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 5882d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 5892d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 5902d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 591*3f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 592cd0e1443SSatish Balay 5933a40ed3dSBarry Smith PetscFunctionBegin; 5942d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 595cd0e1443SSatish Balay row = im[k]; brow = row/bs; 596a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 597a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 5982d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 5992c3acbe9SBarry Smith nrow = ailen[brow]; 6002d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 601a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 602a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6032d61bbb3SSatish Balay col = in[l] ; 6042d61bbb3SSatish Balay bcol = col/bs; 6052d61bbb3SSatish Balay cidx = col%bs; 6062d61bbb3SSatish Balay ridx = row%bs; 6072d61bbb3SSatish Balay high = nrow; 6082d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6092d61bbb3SSatish Balay while (high-low > 5) { 610cd0e1443SSatish Balay t = (low+high)/2; 611cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 612cd0e1443SSatish Balay else low = t; 613cd0e1443SSatish Balay } 614cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 615cd0e1443SSatish Balay if (rp[i] > bcol) break; 616cd0e1443SSatish Balay if (rp[i] == bcol) { 6172d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6182d61bbb3SSatish Balay goto finished; 619cd0e1443SSatish Balay } 620cd0e1443SSatish Balay } 6212d61bbb3SSatish Balay *v++ = zero; 6222d61bbb3SSatish Balay finished:; 623cd0e1443SSatish Balay } 624cd0e1443SSatish Balay } 6253a40ed3dSBarry Smith PetscFunctionReturn(0); 626cd0e1443SSatish Balay } 627cd0e1443SSatish Balay 6282d61bbb3SSatish Balay 6295615d1e5SSatish Balay #undef __FUNC__ 63005a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 63192c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 63292c4ed94SBarry Smith { 63392c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6348a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 635dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 636dd9472c6SBarry Smith int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval; 637*3f1db9ecSBarry Smith Scalar *value = v; 638*3f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 63992c4ed94SBarry Smith 6403a40ed3dSBarry Smith PetscFunctionBegin; 6410e324ae4SSatish Balay if (roworiented) { 6420e324ae4SSatish Balay stepval = (n-1)*bs; 6430e324ae4SSatish Balay } else { 6440e324ae4SSatish Balay stepval = (m-1)*bs; 6450e324ae4SSatish Balay } 64692c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 64792c4ed94SBarry Smith row = im[k]; 6483a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 649a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 650a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 65192c4ed94SBarry Smith #endif 65292c4ed94SBarry Smith rp = aj + ai[row]; 65392c4ed94SBarry Smith ap = aa + bs2*ai[row]; 65492c4ed94SBarry Smith rmax = imax[row]; 65592c4ed94SBarry Smith nrow = ailen[row]; 65692c4ed94SBarry Smith low = 0; 65792c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6583a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 659a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 660a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 66192c4ed94SBarry Smith #endif 66292c4ed94SBarry Smith col = in[l]; 66392c4ed94SBarry Smith if (roworiented) { 6640e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6650e324ae4SSatish Balay } else { 6660e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 66792c4ed94SBarry Smith } 66892c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 66992c4ed94SBarry Smith while (high-low > 7) { 67092c4ed94SBarry Smith t = (low+high)/2; 67192c4ed94SBarry Smith if (rp[t] > col) high = t; 67292c4ed94SBarry Smith else low = t; 67392c4ed94SBarry Smith } 67492c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 67592c4ed94SBarry Smith if (rp[i] > col) break; 67692c4ed94SBarry Smith if (rp[i] == col) { 6778a84c255SSatish Balay bap = ap + bs2*i; 6780e324ae4SSatish Balay if (roworiented) { 6798a84c255SSatish Balay if (is == ADD_VALUES) { 680dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 681dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 6828a84c255SSatish Balay bap[jj] += *value++; 683dd9472c6SBarry Smith } 684dd9472c6SBarry Smith } 6850e324ae4SSatish Balay } else { 686dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 687dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 6880e324ae4SSatish Balay bap[jj] = *value++; 6898a84c255SSatish Balay } 690dd9472c6SBarry Smith } 691dd9472c6SBarry Smith } 6920e324ae4SSatish Balay } else { 6930e324ae4SSatish Balay if (is == ADD_VALUES) { 694dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 695dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 6960e324ae4SSatish Balay *bap++ += *value++; 697dd9472c6SBarry Smith } 698dd9472c6SBarry Smith } 6990e324ae4SSatish Balay } else { 700dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 701dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7020e324ae4SSatish Balay *bap++ = *value++; 7030e324ae4SSatish Balay } 7048a84c255SSatish Balay } 705dd9472c6SBarry Smith } 706dd9472c6SBarry Smith } 707f1241b54SBarry Smith goto noinsert2; 70892c4ed94SBarry Smith } 70992c4ed94SBarry Smith } 71089280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 711a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 71292c4ed94SBarry Smith if (nrow >= rmax) { 71392c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 71492c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 715*3f1db9ecSBarry Smith MatScalar *new_a; 71692c4ed94SBarry Smith 717a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 71889280ab3SLois Curfman McInnes 71992c4ed94SBarry Smith /* malloc new storage space */ 720*3f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 721*3f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len ); CHKPTRQ(new_a); 72292c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 72392c4ed94SBarry Smith new_i = new_j + new_nz; 72492c4ed94SBarry Smith 72592c4ed94SBarry Smith /* copy over old data into new slots */ 72692c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 72792c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 72892c4ed94SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int)); 72992c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 730dd9472c6SBarry Smith PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int)); 731*3f1db9ecSBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar)); 732*3f1db9ecSBarry Smith PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar)); 733*3f1db9ecSBarry Smith PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar)); 73492c4ed94SBarry Smith /* free up old matrix storage */ 73592c4ed94SBarry Smith PetscFree(a->a); 73692c4ed94SBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 73792c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 73892c4ed94SBarry Smith a->singlemalloc = 1; 73992c4ed94SBarry Smith 74092c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 74192c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 742*3f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 74392c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 74492c4ed94SBarry Smith a->reallocs++; 74592c4ed94SBarry Smith a->nz++; 74692c4ed94SBarry Smith } 74792c4ed94SBarry Smith N = nrow++ - 1; 74892c4ed94SBarry Smith /* shift up all the later entries in this row */ 74992c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 75092c4ed94SBarry Smith rp[ii+1] = rp[ii]; 751*3f1db9ecSBarry Smith PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); 75292c4ed94SBarry Smith } 753*3f1db9ecSBarry Smith if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar)); 75492c4ed94SBarry Smith rp[i] = col; 7558a84c255SSatish Balay bap = ap + bs2*i; 7560e324ae4SSatish Balay if (roworiented) { 757dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 758dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7590e324ae4SSatish Balay bap[jj] = *value++; 760dd9472c6SBarry Smith } 761dd9472c6SBarry Smith } 7620e324ae4SSatish Balay } else { 763dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 764dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7650e324ae4SSatish Balay *bap++ = *value++; 7660e324ae4SSatish Balay } 767dd9472c6SBarry Smith } 768dd9472c6SBarry Smith } 769f1241b54SBarry Smith noinsert2:; 77092c4ed94SBarry Smith low = i; 77192c4ed94SBarry Smith } 77292c4ed94SBarry Smith ailen[row] = nrow; 77392c4ed94SBarry Smith } 7743a40ed3dSBarry Smith PetscFunctionReturn(0); 77592c4ed94SBarry Smith } 77692c4ed94SBarry Smith 777f2501298SSatish Balay 7785615d1e5SSatish Balay #undef __FUNC__ 7795615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 7808f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 781584200bdSSatish Balay { 782584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 783584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 784a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 78543ee02c3SBarry Smith int mbs = a->mbs, bs2 = a->bs2,rmax = 0; 786*3f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 787584200bdSSatish Balay 7883a40ed3dSBarry Smith PetscFunctionBegin; 7893a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 790584200bdSSatish Balay 79143ee02c3SBarry Smith if (m) rmax = ailen[0]; 792584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 793584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 794584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 795d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 796584200bdSSatish Balay if (fshift) { 797a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 798584200bdSSatish Balay N = ailen[i]; 799584200bdSSatish Balay for ( j=0; j<N; j++ ) { 800584200bdSSatish Balay ip[j-fshift] = ip[j]; 801*3f1db9ecSBarry Smith PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar)); 802584200bdSSatish Balay } 803584200bdSSatish Balay } 804584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 805584200bdSSatish Balay } 806584200bdSSatish Balay if (mbs) { 807584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 808584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 809584200bdSSatish Balay } 810584200bdSSatish Balay /* reset ilen and imax for each row */ 811584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 812584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 813584200bdSSatish Balay } 814a7c10996SSatish Balay a->nz = ai[mbs]; 815584200bdSSatish Balay 816584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 817584200bdSSatish Balay if (fshift && a->diag) { 818584200bdSSatish Balay PetscFree(a->diag); 819584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 820584200bdSSatish Balay a->diag = 0; 821584200bdSSatish Balay } 8223d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 823219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8243d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 825584200bdSSatish Balay a->reallocs); 826d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 827e2f3b5e9SSatish Balay a->reallocs = 0; 8284e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8294e220ebcSLois Curfman McInnes 8303a40ed3dSBarry Smith PetscFunctionReturn(0); 831584200bdSSatish Balay } 832584200bdSSatish Balay 8335a838052SSatish Balay 834d9b7c43dSSatish Balay /* idx should be of length atlease bs */ 8355615d1e5SSatish Balay #undef __FUNC__ 8365615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block" 837d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg) 838d9b7c43dSSatish Balay { 839d9b7c43dSSatish Balay int i,row; 8403a40ed3dSBarry Smith 8413a40ed3dSBarry Smith PetscFunctionBegin; 842d9b7c43dSSatish Balay row = idx[0]; 8433a40ed3dSBarry Smith if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 844d9b7c43dSSatish Balay 845d9b7c43dSSatish Balay for ( i=1; i<bs; i++ ) { 8463a40ed3dSBarry Smith if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 847d9b7c43dSSatish Balay } 848d9b7c43dSSatish Balay *flg = PETSC_TRUE; 8493a40ed3dSBarry Smith PetscFunctionReturn(0); 850d9b7c43dSSatish Balay } 851d9b7c43dSSatish Balay 8525615d1e5SSatish Balay #undef __FUNC__ 8535615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 8548f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 855d9b7c43dSSatish Balay { 856d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 857d9b7c43dSSatish Balay IS is_local; 858d9b7c43dSSatish Balay int ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2; 859d9b7c43dSSatish Balay PetscTruth flg; 860*3f1db9ecSBarry Smith Scalar zero = 0.0; 861*3f1db9ecSBarry Smith MatScalar *aa; 862d9b7c43dSSatish Balay 8633a40ed3dSBarry Smith PetscFunctionBegin; 864d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 865d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 866d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 867537820f0SBarry Smith ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr); 868d9b7c43dSSatish Balay ierr = ISSort(is_local); CHKERRQ(ierr); 869d9b7c43dSSatish Balay ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr); 870d9b7c43dSSatish Balay 871d9b7c43dSSatish Balay i = 0; 872d9b7c43dSSatish Balay while (i < is_n) { 873a8c6a408SBarry Smith if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 874d9b7c43dSSatish Balay flg = PETSC_FALSE; 875d9b7c43dSSatish Balay if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); } 876d9b7c43dSSatish Balay count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs; 877d9b7c43dSSatish Balay aa = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs); 8782d61bbb3SSatish Balay if (flg) { /* There exists a block of rows to be Zerowed */ 879a07cd24cSSatish Balay if (baij->ilen[rows[i]/bs] > 0) { 880*3f1db9ecSBarry Smith PetscMemzero(aa,count*bs*sizeof(MatScalar)); 8812d61bbb3SSatish Balay baij->ilen[rows[i]/bs] = 1; 8822d61bbb3SSatish Balay baij->j[baij->i[rows[i]/bs]] = rows[i]/bs; 883a07cd24cSSatish Balay } 8842d61bbb3SSatish Balay i += bs; 8852d61bbb3SSatish Balay } else { /* Zero out only the requested row */ 886d9b7c43dSSatish Balay for ( j=0; j<count; j++ ) { 887d9b7c43dSSatish Balay aa[0] = zero; 888d9b7c43dSSatish Balay aa+=bs; 889d9b7c43dSSatish Balay } 890d9b7c43dSSatish Balay i++; 891d9b7c43dSSatish Balay } 892d9b7c43dSSatish Balay } 893d9b7c43dSSatish Balay if (diag) { 894d9b7c43dSSatish Balay for ( j=0; j<is_n; j++ ) { 895f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 896d9b7c43dSSatish Balay } 897d9b7c43dSSatish Balay } 898d9b7c43dSSatish Balay ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr); 899d9b7c43dSSatish Balay ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr); 900d9b7c43dSSatish Balay ierr = ISDestroy(is_local); CHKERRQ(ierr); 9019a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 9023a40ed3dSBarry Smith PetscFunctionReturn(0); 903d9b7c43dSSatish Balay } 9041c351548SSatish Balay 9055615d1e5SSatish Balay #undef __FUNC__ 9062d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9072d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9082d61bbb3SSatish Balay { 9092d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9102d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9112d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9122d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 9132d61bbb3SSatish Balay int ridx,cidx,bs2=a->bs2; 914*3f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9152d61bbb3SSatish Balay 9162d61bbb3SSatish Balay PetscFunctionBegin; 9172d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9182d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9192d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9202d61bbb3SSatish Balay if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 9212d61bbb3SSatish Balay if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 9222d61bbb3SSatish Balay #endif 9232d61bbb3SSatish Balay rp = aj + ai[brow]; 9242d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9252d61bbb3SSatish Balay rmax = imax[brow]; 9262d61bbb3SSatish Balay nrow = ailen[brow]; 9272d61bbb3SSatish Balay low = 0; 9282d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9292d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9302d61bbb3SSatish Balay if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 9312d61bbb3SSatish Balay if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 9322d61bbb3SSatish Balay #endif 9332d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 9342d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 9352d61bbb3SSatish Balay if (roworiented) { 9362d61bbb3SSatish Balay value = *v++; 9372d61bbb3SSatish Balay } else { 9382d61bbb3SSatish Balay value = v[k + l*m]; 9392d61bbb3SSatish Balay } 9402d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 9412d61bbb3SSatish Balay while (high-low > 7) { 9422d61bbb3SSatish Balay t = (low+high)/2; 9432d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 9442d61bbb3SSatish Balay else low = t; 9452d61bbb3SSatish Balay } 9462d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 9472d61bbb3SSatish Balay if (rp[i] > bcol) break; 9482d61bbb3SSatish Balay if (rp[i] == bcol) { 9492d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 9502d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 9512d61bbb3SSatish Balay else *bap = value; 9522d61bbb3SSatish Balay goto noinsert1; 9532d61bbb3SSatish Balay } 9542d61bbb3SSatish Balay } 9552d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 9562d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9572d61bbb3SSatish Balay if (nrow >= rmax) { 9582d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 9592d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 960*3f1db9ecSBarry Smith MatScalar *new_a; 9612d61bbb3SSatish Balay 9622d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9632d61bbb3SSatish Balay 9642d61bbb3SSatish Balay /* Malloc new storage space */ 965*3f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 966*3f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len ); CHKPTRQ(new_a); 9672d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 9682d61bbb3SSatish Balay new_i = new_j + new_nz; 9692d61bbb3SSatish Balay 9702d61bbb3SSatish Balay /* copy over old data into new slots */ 9712d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 9722d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 9732d61bbb3SSatish Balay PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int)); 9742d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 9752d61bbb3SSatish Balay PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow, 9762d61bbb3SSatish Balay len*sizeof(int)); 977*3f1db9ecSBarry Smith PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar)); 978*3f1db9ecSBarry Smith PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar)); 9792d61bbb3SSatish Balay PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE), 980*3f1db9ecSBarry Smith aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar)); 9812d61bbb3SSatish Balay /* free up old matrix storage */ 9822d61bbb3SSatish Balay PetscFree(a->a); 9832d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 9842d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 9852d61bbb3SSatish Balay a->singlemalloc = 1; 9862d61bbb3SSatish Balay 9872d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 9882d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 989*3f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 9902d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 9912d61bbb3SSatish Balay a->reallocs++; 9922d61bbb3SSatish Balay a->nz++; 9932d61bbb3SSatish Balay } 9942d61bbb3SSatish Balay N = nrow++ - 1; 9952d61bbb3SSatish Balay /* shift up all the later entries in this row */ 9962d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 9972d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 998*3f1db9ecSBarry Smith PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); 9992d61bbb3SSatish Balay } 1000*3f1db9ecSBarry Smith if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar)); 10012d61bbb3SSatish Balay rp[i] = bcol; 10022d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10032d61bbb3SSatish Balay noinsert1:; 10042d61bbb3SSatish Balay low = i; 10052d61bbb3SSatish Balay } 10062d61bbb3SSatish Balay ailen[brow] = nrow; 10072d61bbb3SSatish Balay } 10082d61bbb3SSatish Balay PetscFunctionReturn(0); 10092d61bbb3SSatish Balay } 10102d61bbb3SSatish Balay 10112d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10122d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10132d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1014d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*); 10152d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**); 10162d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10172d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10182d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10192d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10202d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10212d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10222d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10232d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10242d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10252d61bbb3SSatish Balay 10262d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10272d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 10282d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 10292d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 10302d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 10312d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 10322d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 10332d61bbb3SSatish Balay 10342d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 10352d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 10362d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 10372d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 10382d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 10392d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 10402d61bbb3SSatish Balay 10412d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 10422d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 10432d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 10442d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 10452d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 10462d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 10472d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 10482d61bbb3SSatish Balay 10492d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 10502d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 10512d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 10522d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 10532d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 10542d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 10552d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 10562d61bbb3SSatish Balay 10572d61bbb3SSatish Balay #undef __FUNC__ 10582d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 10592d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 10602d61bbb3SSatish Balay { 10612d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 10622d61bbb3SSatish Balay Mat outA; 10632d61bbb3SSatish Balay int ierr; 1064667159a5SBarry Smith PetscTruth row_identity, col_identity; 10652d61bbb3SSatish Balay 10662d61bbb3SSatish Balay PetscFunctionBegin; 10672d61bbb3SSatish Balay if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 1068667159a5SBarry Smith ierr = ISIdentity(row,&row_identity); CHKERRQ(ierr); 1069667159a5SBarry Smith ierr = ISIdentity(col,&col_identity); CHKERRQ(ierr); 1070667159a5SBarry Smith if (!row_identity || !col_identity) { 1071667159a5SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity"); 1072667159a5SBarry Smith } 10732d61bbb3SSatish Balay 10742d61bbb3SSatish Balay outA = inA; 10752d61bbb3SSatish Balay inA->factor = FACTOR_LU; 10762d61bbb3SSatish Balay a->row = row; 10772d61bbb3SSatish Balay a->col = col; 10782d61bbb3SSatish Balay 1079e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1080e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 10811e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1082e51c0b9cSSatish Balay 10832d61bbb3SSatish Balay if (!a->solve_work) { 10842d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 10852d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 10862d61bbb3SSatish Balay } 10872d61bbb3SSatish Balay 10882d61bbb3SSatish Balay if (!a->diag) { 10892d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 10902d61bbb3SSatish Balay } 10912d61bbb3SSatish Balay /* 1092667159a5SBarry Smith Blocksize 4 and 5 have a special faster factorization/solver for ILU(0) factorization 10932d61bbb3SSatish Balay with natural ordering 10942d61bbb3SSatish Balay */ 10952d61bbb3SSatish Balay if (a->bs == 4) { 1096667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1097f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 1098667159a5SBarry Smith PLogInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n"); 1099667159a5SBarry Smith } else if (a->bs == 5) { 1100667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1101667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 1102667159a5SBarry Smith PLogInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n"); 11032d61bbb3SSatish Balay } 11042d61bbb3SSatish Balay 1105667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 1106667159a5SBarry Smith 1107667159a5SBarry Smith 11082d61bbb3SSatish Balay PetscFunctionReturn(0); 11092d61bbb3SSatish Balay } 11102d61bbb3SSatish Balay #undef __FUNC__ 1111d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1112ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1113ba4ca20aSSatish Balay { 1114ba4ca20aSSatish Balay static int called = 0; 1115ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1116ba4ca20aSSatish Balay 11173a40ed3dSBarry Smith PetscFunctionBegin; 11183a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 111976be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n"); 112076be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n"); 11213a40ed3dSBarry Smith PetscFunctionReturn(0); 1122ba4ca20aSSatish Balay } 1123d9b7c43dSSatish Balay 1124fb2e594dSBarry Smith EXTERN_C_BEGIN 112527a8da17SBarry Smith #undef __FUNC__ 112627a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 112727a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 112827a8da17SBarry Smith { 112927a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 113027a8da17SBarry Smith int i,nz,n; 113127a8da17SBarry Smith 113227a8da17SBarry Smith PetscFunctionBegin; 113327a8da17SBarry Smith nz = baij->maxnz; 113427a8da17SBarry Smith n = baij->n; 113527a8da17SBarry Smith for (i=0; i<nz; i++) { 113627a8da17SBarry Smith baij->j[i] = indices[i]; 113727a8da17SBarry Smith } 113827a8da17SBarry Smith baij->nz = nz; 113927a8da17SBarry Smith for ( i=0; i<n; i++ ) { 114027a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 114127a8da17SBarry Smith } 114227a8da17SBarry Smith 114327a8da17SBarry Smith PetscFunctionReturn(0); 114427a8da17SBarry Smith } 1145fb2e594dSBarry Smith EXTERN_C_END 114627a8da17SBarry Smith 114727a8da17SBarry Smith #undef __FUNC__ 114827a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 114927a8da17SBarry Smith /*@ 115027a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 115127a8da17SBarry Smith in the matrix. 115227a8da17SBarry Smith 115327a8da17SBarry Smith Input Parameters: 115427a8da17SBarry Smith + mat - the SeqBAIJ matrix 115527a8da17SBarry Smith - indices - the column indices 115627a8da17SBarry Smith 115727a8da17SBarry Smith Notes: 115827a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 115927a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 116027a8da17SBarry Smith of the MatSetValues() operation. 116127a8da17SBarry Smith 116227a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 116327a8da17SBarry Smith MatCreateSeqBAIJ(). 116427a8da17SBarry Smith 116527a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 116627a8da17SBarry Smith 116727a8da17SBarry Smith @*/ 116827a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 116927a8da17SBarry Smith { 117027a8da17SBarry Smith int ierr,(*f)(Mat,int *); 117127a8da17SBarry Smith 117227a8da17SBarry Smith PetscFunctionBegin; 117327a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 117427a8da17SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f); 117527a8da17SBarry Smith CHKERRQ(ierr); 117627a8da17SBarry Smith if (f) { 117727a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 117827a8da17SBarry Smith } else { 117927a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 118027a8da17SBarry Smith } 118127a8da17SBarry Smith PetscFunctionReturn(0); 118227a8da17SBarry Smith } 118327a8da17SBarry Smith 11842593348eSBarry Smith /* -------------------------------------------------------------------*/ 1185cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1186cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1187cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1188cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1189cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1190cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1191cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1192cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1193cc2dc46cSBarry Smith 0, 1194cc2dc46cSBarry Smith 0, 1195cc2dc46cSBarry Smith 0, 1196cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1197cc2dc46cSBarry Smith 0, 1198b6490206SBarry Smith 0, 1199f2501298SSatish Balay MatTranspose_SeqBAIJ, 1200cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1201cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1202cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1203cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1204cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1205b6490206SBarry Smith 0, 1206cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1207cc2dc46cSBarry Smith 0, 1208cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1209cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1210cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1211cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1212cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1213cc2dc46cSBarry Smith 0, 1214cc2dc46cSBarry Smith 0, 1215cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1216cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1217cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1218cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1219cc2dc46cSBarry Smith 0, 1220cc2dc46cSBarry Smith 0, 1221cc2dc46cSBarry Smith 0, 12222e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1223cc2dc46cSBarry Smith 0, 1224cc2dc46cSBarry Smith 0, 1225cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1226cc2dc46cSBarry Smith 0, 1227cc2dc46cSBarry Smith 0, 1228cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1229cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1230cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1231cc2dc46cSBarry Smith 0, 1232cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1233cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1234cc2dc46cSBarry Smith 0, 1235cc2dc46cSBarry Smith 0, 1236cc2dc46cSBarry Smith 0, 1237cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 12383b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 123992c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1240cc2dc46cSBarry Smith 0, 1241cc2dc46cSBarry Smith 0, 1242cc2dc46cSBarry Smith 0, 1243cc2dc46cSBarry Smith 0, 1244cc2dc46cSBarry Smith 0, 1245cc2dc46cSBarry Smith 0, 1246d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1247cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1248cc2dc46cSBarry Smith 0, 1249cc2dc46cSBarry Smith 0, 1250cc2dc46cSBarry Smith MatGetMaps_Petsc}; 12512593348eSBarry Smith 12525615d1e5SSatish Balay #undef __FUNC__ 12535615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 12542593348eSBarry Smith /*@C 125544cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 125644cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 125744cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 12582bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 12592bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 12602593348eSBarry Smith 1261db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1262db81eaa0SLois Curfman McInnes 12632593348eSBarry Smith Input Parameters: 1264db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1265b6490206SBarry Smith . bs - size of block 12662593348eSBarry Smith . m - number of rows 12672593348eSBarry Smith . n - number of columns 1268b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 1269db81eaa0SLois Curfman McInnes - nzz - array containing the number of block nonzeros in the various block rows 12702bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 12712593348eSBarry Smith 12722593348eSBarry Smith Output Parameter: 12732593348eSBarry Smith . A - the matrix 12742593348eSBarry Smith 12750513a670SBarry Smith Options Database Keys: 1276db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1277db81eaa0SLois Curfman McInnes block calculations (much slower) 1278db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 12790513a670SBarry Smith 12802593348eSBarry Smith Notes: 128144cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 12822593348eSBarry Smith storage. That is, the stored row and column indices can begin at 128344cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 12842593348eSBarry Smith 12852593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 12862593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 12872593348eSBarry Smith allocation. For additional details, see the users manual chapter on 12886da5968aSLois Curfman McInnes matrices. 12892593348eSBarry Smith 1290db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 12912593348eSBarry Smith @*/ 1292b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 12932593348eSBarry Smith { 12942593348eSBarry Smith Mat B; 1295b6490206SBarry Smith Mat_SeqBAIJ *b; 12963b2fbd54SBarry Smith int i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size; 12973b2fbd54SBarry Smith 12983a40ed3dSBarry Smith PetscFunctionBegin; 12993b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1300a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1301b6490206SBarry Smith 13020513a670SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 13030513a670SBarry Smith 13043a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1305a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 13063a40ed3dSBarry Smith } 13072593348eSBarry Smith 13082593348eSBarry Smith *A = 0; 1309*3f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 13102593348eSBarry Smith PLogObjectCreate(B); 1311b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 131244cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqBAIJ)); 1313cc2dc46cSBarry Smith PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)); 13141a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr); 13151a6a6d98SBarry Smith if (!flg) { 13167fc0212eSBarry Smith switch (bs) { 13177fc0212eSBarry Smith case 1: 1318f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1319f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1320f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1321f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 13227fc0212eSBarry Smith break; 13234eeb42bcSBarry Smith case 2: 1324f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1325f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1326f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1327f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 13286d84be18SBarry Smith break; 1329f327dd97SBarry Smith case 3: 1330f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1331f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1332f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1333f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 13344eeb42bcSBarry Smith break; 13356d84be18SBarry Smith case 4: 1336f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1337f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1338f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1339f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 13406d84be18SBarry Smith break; 13416d84be18SBarry Smith case 5: 1342f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1343f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1344f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1345f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 13466d84be18SBarry Smith break; 134748e9ddb2SSatish Balay case 7: 1348f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1349f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1350f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 135148e9ddb2SSatish Balay break; 13527fc0212eSBarry Smith } 13531a6a6d98SBarry Smith } 1354e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1355e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 13562593348eSBarry Smith B->factor = 0; 13572593348eSBarry Smith B->lupivotthreshold = 1.0; 135890f02eecSBarry Smith B->mapping = 0; 13592593348eSBarry Smith b->row = 0; 13602593348eSBarry Smith b->col = 0; 1361e51c0b9cSSatish Balay b->icol = 0; 13622593348eSBarry Smith b->reallocs = 0; 13632593348eSBarry Smith 136444cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 136544cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1366a5ae1ecdSBarry Smith 13677413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 13687413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1369a5ae1ecdSBarry Smith 1370b6490206SBarry Smith b->mbs = mbs; 1371f2501298SSatish Balay b->nbs = nbs; 1372b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 13732593348eSBarry Smith if (nnz == PETSC_NULL) { 1374119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 13752593348eSBarry Smith else if (nz <= 0) nz = 1; 1376b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1377b6490206SBarry Smith nz = nz*mbs; 13783a40ed3dSBarry Smith } else { 13792593348eSBarry Smith nz = 0; 1380b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 13812593348eSBarry Smith } 13822593348eSBarry Smith 13832593348eSBarry Smith /* allocate the matrix space */ 1384*3f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 1385*3f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1386*3f1db9ecSBarry Smith PetscMemzero(b->a,nz*bs2*sizeof(MatScalar)); 13877e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 13882593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 13892593348eSBarry Smith b->i = b->j + nz; 13902593348eSBarry Smith b->singlemalloc = 1; 13912593348eSBarry Smith 1392de6a44a3SBarry Smith b->i[0] = 0; 1393b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 13942593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 13952593348eSBarry Smith } 13962593348eSBarry Smith 1397b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1398b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1399f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1400b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 14012593348eSBarry Smith 1402b6490206SBarry Smith b->bs = bs; 14039df24120SSatish Balay b->bs2 = bs2; 1404b6490206SBarry Smith b->mbs = mbs; 14052593348eSBarry Smith b->nz = 0; 140618e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 14072593348eSBarry Smith b->sorted = 0; 14082593348eSBarry Smith b->roworiented = 1; 14092593348eSBarry Smith b->nonew = 0; 14102593348eSBarry Smith b->diag = 0; 14112593348eSBarry Smith b->solve_work = 0; 1412de6a44a3SBarry Smith b->mult_work = 0; 14132593348eSBarry Smith b->spptr = 0; 14144e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 14154e220ebcSLois Curfman McInnes 14162593348eSBarry Smith *A = B; 14172593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 14182593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 141927a8da17SBarry Smith 142027a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 142127a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 142227a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 142327a8da17SBarry Smith 14243a40ed3dSBarry Smith PetscFunctionReturn(0); 14252593348eSBarry Smith } 14262593348eSBarry Smith 14275615d1e5SSatish Balay #undef __FUNC__ 14282e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 14292e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 14302593348eSBarry Smith { 14312593348eSBarry Smith Mat C; 1432b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 143327a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1434de6a44a3SBarry Smith 14353a40ed3dSBarry Smith PetscFunctionBegin; 1436a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 14372593348eSBarry Smith 14382593348eSBarry Smith *B = 0; 1439*3f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 14402593348eSBarry Smith PLogObjectCreate(C); 1441b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 1442c3c66ee5SSatish Balay PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 1443e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1444e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 14452593348eSBarry Smith C->factor = A->factor; 14462593348eSBarry Smith c->row = 0; 14472593348eSBarry Smith c->col = 0; 14482593348eSBarry Smith C->assembled = PETSC_TRUE; 14492593348eSBarry Smith 145044cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 145144cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 145244cd7ae7SLois Curfman McInnes C->M = a->m; 145344cd7ae7SLois Curfman McInnes C->N = a->n; 145444cd7ae7SLois Curfman McInnes 1455b6490206SBarry Smith c->bs = a->bs; 14569df24120SSatish Balay c->bs2 = a->bs2; 1457b6490206SBarry Smith c->mbs = a->mbs; 1458b6490206SBarry Smith c->nbs = a->nbs; 14592593348eSBarry Smith 1460b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1461b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1462b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 14632593348eSBarry Smith c->imax[i] = a->imax[i]; 14642593348eSBarry Smith c->ilen[i] = a->ilen[i]; 14652593348eSBarry Smith } 14662593348eSBarry Smith 14672593348eSBarry Smith /* allocate the matrix space */ 14682593348eSBarry Smith c->singlemalloc = 1; 1469*3f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 1470*3f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len ); CHKPTRQ(c->a); 14717e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1472de6a44a3SBarry Smith c->i = c->j + nz; 1473b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1474b6490206SBarry Smith if (mbs > 0) { 1475de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 14762e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1477*3f1db9ecSBarry Smith PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar)); 14782e8a6d31SBarry Smith } else { 1479*3f1db9ecSBarry Smith PetscMemzero(c->a,bs2*nz*sizeof(MatScalar)); 14802593348eSBarry Smith } 14812593348eSBarry Smith } 14822593348eSBarry Smith 1483f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 14842593348eSBarry Smith c->sorted = a->sorted; 14852593348eSBarry Smith c->roworiented = a->roworiented; 14862593348eSBarry Smith c->nonew = a->nonew; 14872593348eSBarry Smith 14882593348eSBarry Smith if (a->diag) { 1489b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1490b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1491b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 14922593348eSBarry Smith c->diag[i] = a->diag[i]; 14932593348eSBarry Smith } 14942593348eSBarry Smith } 14952593348eSBarry Smith else c->diag = 0; 14962593348eSBarry Smith c->nz = a->nz; 14972593348eSBarry Smith c->maxnz = a->maxnz; 14982593348eSBarry Smith c->solve_work = 0; 14992593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 15007fc0212eSBarry Smith c->mult_work = 0; 15012593348eSBarry Smith *B = C; 150227a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqBAIJSetColumnIndices_C", 150327a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 150427a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 15053a40ed3dSBarry Smith PetscFunctionReturn(0); 15062593348eSBarry Smith } 15072593348eSBarry Smith 15085615d1e5SSatish Balay #undef __FUNC__ 15095615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 151019bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 15112593348eSBarry Smith { 1512b6490206SBarry Smith Mat_SeqBAIJ *a; 15132593348eSBarry Smith Mat B; 1514de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1515b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 151635aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1517a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1518b6490206SBarry Smith Scalar *aa; 151919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 15202593348eSBarry Smith 15213a40ed3dSBarry Smith PetscFunctionBegin; 15221a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1523de6a44a3SBarry Smith bs2 = bs*bs; 1524b6490206SBarry Smith 15252593348eSBarry Smith MPI_Comm_size(comm,&size); 1526cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 152790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 15280752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 1529a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 15302593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 15312593348eSBarry Smith 1532d64ed03dSBarry Smith if (header[3] < 0) { 1533a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1534d64ed03dSBarry Smith } 1535d64ed03dSBarry Smith 1536a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 153735aab85fSBarry Smith 153835aab85fSBarry Smith /* 153935aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 154035aab85fSBarry Smith divisible by the blocksize 154135aab85fSBarry Smith */ 1542b6490206SBarry Smith mbs = M/bs; 154335aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 154435aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 154535aab85fSBarry Smith else mbs++; 154635aab85fSBarry Smith if (extra_rows) { 1547537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 154835aab85fSBarry Smith } 1549b6490206SBarry Smith 15502593348eSBarry Smith /* read in row lengths */ 155135aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 15520752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 155335aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 15542593348eSBarry Smith 1555b6490206SBarry Smith /* read in column indices */ 155635aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 15570752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr); 155835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1559b6490206SBarry Smith 1560b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1561b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1562b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 156335aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 156435aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 156535aab85fSBarry Smith masked = mask + mbs; 1566b6490206SBarry Smith rowcount = 0; nzcount = 0; 1567b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 156835aab85fSBarry Smith nmask = 0; 1569b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1570b6490206SBarry Smith kmax = rowlengths[rowcount]; 1571b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 157235aab85fSBarry Smith tmp = jj[nzcount++]/bs; 157335aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1574b6490206SBarry Smith } 1575b6490206SBarry Smith rowcount++; 1576b6490206SBarry Smith } 157735aab85fSBarry Smith browlengths[i] += nmask; 157835aab85fSBarry Smith /* zero out the mask elements we set */ 157935aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1580b6490206SBarry Smith } 1581b6490206SBarry Smith 15822593348eSBarry Smith /* create our matrix */ 15833eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 15842593348eSBarry Smith B = *A; 1585b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 15862593348eSBarry Smith 15872593348eSBarry Smith /* set matrix "i" values */ 1588de6a44a3SBarry Smith a->i[0] = 0; 1589b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1590b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1591b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 15922593348eSBarry Smith } 1593b6490206SBarry Smith a->nz = 0; 1594b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 15952593348eSBarry Smith 1596b6490206SBarry Smith /* read in nonzero values */ 159735aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 15980752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr); 159935aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1600b6490206SBarry Smith 1601b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1602b6490206SBarry Smith nzcount = 0; jcount = 0; 1603b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1604b6490206SBarry Smith nzcountb = nzcount; 160535aab85fSBarry Smith nmask = 0; 1606b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1607b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1608b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 160935aab85fSBarry Smith tmp = jj[nzcount++]/bs; 161035aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1611b6490206SBarry Smith } 1612b6490206SBarry Smith rowcount++; 1613b6490206SBarry Smith } 1614de6a44a3SBarry Smith /* sort the masked values */ 161577c4ece6SBarry Smith PetscSortInt(nmask,masked); 1616de6a44a3SBarry Smith 1617b6490206SBarry Smith /* set "j" values into matrix */ 1618b6490206SBarry Smith maskcount = 1; 161935aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 162035aab85fSBarry Smith a->j[jcount++] = masked[j]; 1621de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1622b6490206SBarry Smith } 1623b6490206SBarry Smith /* set "a" values into matrix */ 1624de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1625b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1626b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1627b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1628de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1629de6a44a3SBarry Smith block = mask[tmp] - 1; 1630de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1631de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1632b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1633b6490206SBarry Smith } 1634b6490206SBarry Smith } 163535aab85fSBarry Smith /* zero out the mask elements we set */ 163635aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1637b6490206SBarry Smith } 1638a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1639b6490206SBarry Smith 1640b6490206SBarry Smith PetscFree(rowlengths); 1641b6490206SBarry Smith PetscFree(browlengths); 1642b6490206SBarry Smith PetscFree(aa); 1643b6490206SBarry Smith PetscFree(jj); 1644b6490206SBarry Smith PetscFree(mask); 1645b6490206SBarry Smith 1646b6490206SBarry Smith B->assembled = PETSC_TRUE; 1647de6a44a3SBarry Smith 16489c01be13SBarry Smith ierr = MatView_Private(B); CHKERRQ(ierr); 16493a40ed3dSBarry Smith PetscFunctionReturn(0); 16502593348eSBarry Smith } 16512593348eSBarry Smith 16522593348eSBarry Smith 16532593348eSBarry Smith 1654