1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*27a8da17SBarry Smith static char vcid[] = "$Id: baij.c,v 1.135 1998/04/27 03:52:29 curfman 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 1102d61bbb3SSatish Balay #if defined(USE_PETSC_LOG) 111e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1122d61bbb3SSatish Balay #endif 1132d61bbb3SSatish Balay PetscFree(a->a); 1142d61bbb3SSatish Balay if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 1152d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 1162d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 1172d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 1182d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 1192d61bbb3SSatish Balay if (a->mult_work) PetscFree(a->mult_work); 120e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 1212d61bbb3SSatish Balay PetscFree(a); 1222d61bbb3SSatish Balay PLogObjectDestroy(A); 1232d61bbb3SSatish Balay PetscHeaderDestroy(A); 1242d61bbb3SSatish Balay PetscFunctionReturn(0); 1252d61bbb3SSatish Balay } 1262d61bbb3SSatish Balay 1272d61bbb3SSatish Balay #undef __FUNC__ 1282d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1292d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1302d61bbb3SSatish Balay { 1312d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1322d61bbb3SSatish Balay 1332d61bbb3SSatish Balay PetscFunctionBegin; 1342d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1352d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1362d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1372d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1382d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1392d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 1402d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 1412d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1422d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1432d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1442d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1452d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1462d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 147b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 148b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 149981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1502d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1512d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1522d61bbb3SSatish Balay } else { 1532d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1542d61bbb3SSatish Balay } 1552d61bbb3SSatish Balay PetscFunctionReturn(0); 1562d61bbb3SSatish Balay } 1572d61bbb3SSatish Balay 1582d61bbb3SSatish Balay 1592d61bbb3SSatish Balay #undef __FUNC__ 1602d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1612d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1622d61bbb3SSatish Balay { 1632d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1642d61bbb3SSatish Balay 1652d61bbb3SSatish Balay PetscFunctionBegin; 1662d61bbb3SSatish Balay if (m) *m = a->m; 1672d61bbb3SSatish Balay if (n) *n = a->n; 1682d61bbb3SSatish Balay PetscFunctionReturn(0); 1692d61bbb3SSatish Balay } 1702d61bbb3SSatish Balay 1712d61bbb3SSatish Balay #undef __FUNC__ 1722d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 1732d61bbb3SSatish Balay int MatGetOwnershipRange_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 *m = 0; *n = a->m; 1792d61bbb3SSatish Balay PetscFunctionReturn(0); 1802d61bbb3SSatish Balay } 1812d61bbb3SSatish Balay 1822d61bbb3SSatish Balay #undef __FUNC__ 1832d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 1842d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 1852d61bbb3SSatish Balay { 1862d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1872d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 1882d61bbb3SSatish Balay Scalar *aa,*v_i,*aa_i; 1892d61bbb3SSatish Balay 1902d61bbb3SSatish Balay PetscFunctionBegin; 1912d61bbb3SSatish Balay bs = a->bs; 1922d61bbb3SSatish Balay ai = a->i; 1932d61bbb3SSatish Balay aj = a->j; 1942d61bbb3SSatish Balay aa = a->a; 1952d61bbb3SSatish Balay bs2 = a->bs2; 1962d61bbb3SSatish Balay 1972d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 1982d61bbb3SSatish Balay 1992d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2002d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2012d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2022d61bbb3SSatish Balay *nz = bs*M; 2032d61bbb3SSatish Balay 2042d61bbb3SSatish Balay if (v) { 2052d61bbb3SSatish Balay *v = 0; 2062d61bbb3SSatish Balay if (*nz) { 2072d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 2082d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2092d61bbb3SSatish Balay v_i = *v + i*bs; 2102d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2112d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2122d61bbb3SSatish Balay } 2132d61bbb3SSatish Balay } 2142d61bbb3SSatish Balay } 2152d61bbb3SSatish Balay 2162d61bbb3SSatish Balay if (idx) { 2172d61bbb3SSatish Balay *idx = 0; 2182d61bbb3SSatish Balay if (*nz) { 2192d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 2202d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2212d61bbb3SSatish Balay idx_i = *idx + i*bs; 2222d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2232d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2242d61bbb3SSatish Balay } 2252d61bbb3SSatish Balay } 2262d61bbb3SSatish Balay } 2272d61bbb3SSatish Balay PetscFunctionReturn(0); 2282d61bbb3SSatish Balay } 2292d61bbb3SSatish Balay 2302d61bbb3SSatish Balay #undef __FUNC__ 2312d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2322d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2332d61bbb3SSatish Balay { 2342d61bbb3SSatish Balay PetscFunctionBegin; 2352d61bbb3SSatish Balay if (idx) {if (*idx) PetscFree(*idx);} 2362d61bbb3SSatish Balay if (v) {if (*v) PetscFree(*v);} 2372d61bbb3SSatish Balay PetscFunctionReturn(0); 2382d61bbb3SSatish Balay } 2392d61bbb3SSatish Balay 2402d61bbb3SSatish Balay #undef __FUNC__ 2412d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2422d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2432d61bbb3SSatish Balay { 2442d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2452d61bbb3SSatish Balay Mat C; 2462d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2472d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2482d61bbb3SSatish Balay Scalar *array=a->a; 2492d61bbb3SSatish Balay 2502d61bbb3SSatish Balay PetscFunctionBegin; 2512d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2522d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col); 2532d61bbb3SSatish Balay PetscMemzero(col,(1+nbs)*sizeof(int)); 2542d61bbb3SSatish Balay 2552d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2562d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr); 2572d61bbb3SSatish Balay PetscFree(col); 2582d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows); 2592d61bbb3SSatish Balay cols = rows + bs; 2602d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2612d61bbb3SSatish Balay cols[0] = i*bs; 2622d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 2632d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 2642d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 2652d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 2662d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 2672d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr); 2682d61bbb3SSatish Balay array += bs2; 2692d61bbb3SSatish Balay } 2702d61bbb3SSatish Balay } 2712d61bbb3SSatish Balay PetscFree(rows); 2722d61bbb3SSatish Balay 2732d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2742d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2752d61bbb3SSatish Balay 2762d61bbb3SSatish Balay if (B != PETSC_NULL) { 2772d61bbb3SSatish Balay *B = C; 2782d61bbb3SSatish Balay } else { 279f830108cSBarry Smith PetscOps *Abops; 280f830108cSBarry Smith struct _MatOps *Aops; 281f830108cSBarry Smith 2822d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 2832d61bbb3SSatish Balay PetscFree(a->a); 2842d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 2852d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 2862d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 2872d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 2882d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 2892d61bbb3SSatish Balay PetscFree(a); 290f830108cSBarry Smith 291f830108cSBarry Smith /* 292f830108cSBarry Smith This is horrible, horrible code. We need to keep the 293f830108cSBarry Smith A pointers for the bops and ops but copy everything 294f830108cSBarry Smith else from C. 295f830108cSBarry Smith */ 296f830108cSBarry Smith Abops = A->bops; 297f830108cSBarry Smith Aops = A->ops; 2982d61bbb3SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 299f830108cSBarry Smith A->bops = Abops; 300f830108cSBarry Smith A->ops = Aops; 301*27a8da17SBarry Smith A->qlist = 0; 302f830108cSBarry Smith 3032d61bbb3SSatish Balay PetscHeaderDestroy(C); 3042d61bbb3SSatish Balay } 3052d61bbb3SSatish Balay PetscFunctionReturn(0); 3062d61bbb3SSatish Balay } 3072d61bbb3SSatish Balay 3082d61bbb3SSatish Balay 3092d61bbb3SSatish Balay 3103b2fbd54SBarry Smith 3115615d1e5SSatish Balay #undef __FUNC__ 312d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 313b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3142593348eSBarry Smith { 315b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3169df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 317b6490206SBarry Smith Scalar *aa; 3182593348eSBarry Smith 3193a40ed3dSBarry Smith PetscFunctionBegin; 32090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 3212593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3222593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3233b2fbd54SBarry Smith 3242593348eSBarry Smith col_lens[1] = a->m; 3252593348eSBarry Smith col_lens[2] = a->n; 3267e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3272593348eSBarry Smith 3282593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 329b6490206SBarry Smith count = 0; 330b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 331b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 332b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 333b6490206SBarry Smith } 3342593348eSBarry Smith } 3350752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3362593348eSBarry Smith PetscFree(col_lens); 3372593348eSBarry Smith 3382593348eSBarry Smith /* store column indices (zero start index) */ 33966963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj); 340b6490206SBarry Smith count = 0; 341b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 342b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 343b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 344b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 345b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3462593348eSBarry Smith } 3472593348eSBarry Smith } 348b6490206SBarry Smith } 349b6490206SBarry Smith } 3500752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr); 351b6490206SBarry Smith PetscFree(jj); 3522593348eSBarry Smith 3532593348eSBarry Smith /* store nonzero values */ 35466963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa); 355b6490206SBarry Smith count = 0; 356b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 357b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 358b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 359b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 3607e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 361b6490206SBarry Smith } 362b6490206SBarry Smith } 363b6490206SBarry Smith } 364b6490206SBarry Smith } 3650752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 366b6490206SBarry Smith PetscFree(aa); 3673a40ed3dSBarry Smith PetscFunctionReturn(0); 3682593348eSBarry Smith } 3692593348eSBarry Smith 3705615d1e5SSatish Balay #undef __FUNC__ 371d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 372b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 3732593348eSBarry Smith { 374b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3759df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 3762593348eSBarry Smith FILE *fd; 3772593348eSBarry Smith char *outputname; 3782593348eSBarry Smith 3793a40ed3dSBarry Smith PetscFunctionBegin; 38090ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 3812593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 38290ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 383639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 3847ddc982cSLois Curfman McInnes fprintf(fd," block size is %d\n",bs); 3852593348eSBarry Smith } 386639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 387a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 3882593348eSBarry Smith } 389639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 39044cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 39144cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 39244cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 39344cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 39444cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 3953a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 396766eeae4SLois Curfman McInnes if (imag(a->a[bs2*k + l*bs + j]) > 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 39744cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 39844cd7ae7SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 399766eeae4SLois Curfman McInnes else if (imag(a->a[bs2*k + l*bs + j]) < 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 400766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 401766eeae4SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j])); 40244cd7ae7SLois Curfman McInnes else if (real(a->a[bs2*k + l*bs + j]) != 0.0) 40344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 40444cd7ae7SLois Curfman McInnes #else 40544cd7ae7SLois Curfman McInnes if (a->a[bs2*k + l*bs + j] != 0.0) 40644cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 40744cd7ae7SLois Curfman McInnes #endif 40844cd7ae7SLois Curfman McInnes } 40944cd7ae7SLois Curfman McInnes } 41044cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 41144cd7ae7SLois Curfman McInnes } 41244cd7ae7SLois Curfman McInnes } 41344cd7ae7SLois Curfman McInnes } 4142593348eSBarry Smith else { 415b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 416b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 417b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 418b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 419b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4203a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 421766eeae4SLois Curfman McInnes if (imag(a->a[bs2*k + l*bs + j]) > 0.0) { 42288685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 4237e67e3f9SSatish Balay real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 42488685aaeSLois Curfman McInnes } 425766eeae4SLois Curfman McInnes else if (imag(a->a[bs2*k + l*bs + j]) < 0.0) { 426766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 427766eeae4SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j])); 428766eeae4SLois Curfman McInnes } 42988685aaeSLois Curfman McInnes else { 4307e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 43188685aaeSLois Curfman McInnes } 43288685aaeSLois Curfman McInnes #else 4337e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 43488685aaeSLois Curfman McInnes #endif 4352593348eSBarry Smith } 4362593348eSBarry Smith } 4372593348eSBarry Smith fprintf(fd,"\n"); 4382593348eSBarry Smith } 4392593348eSBarry Smith } 440b6490206SBarry Smith } 4412593348eSBarry Smith fflush(fd); 4423a40ed3dSBarry Smith PetscFunctionReturn(0); 4432593348eSBarry Smith } 4442593348eSBarry Smith 4455615d1e5SSatish Balay #undef __FUNC__ 446d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 4473270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 4483270192aSSatish Balay { 4493270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4503270192aSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2; 4513270192aSSatish Balay double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 4523270192aSSatish Balay Scalar *aa; 4533270192aSSatish Balay Draw draw; 4543270192aSSatish Balay DrawButton button; 4553270192aSSatish Balay PetscTruth isnull; 4563270192aSSatish Balay 4573a40ed3dSBarry Smith PetscFunctionBegin; 4583a40ed3dSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr); 4593a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 4603270192aSSatish Balay 4613270192aSSatish Balay xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 4623270192aSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 4633270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 4643270192aSSatish Balay /* loop over matrix elements drawing boxes */ 4653270192aSSatish Balay color = DRAW_BLUE; 4663270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4673270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4683270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4693270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4703270192aSSatish Balay aa = a->a + j*bs2; 4713270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4723270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4733270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 474b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4753270192aSSatish Balay } 4763270192aSSatish Balay } 4773270192aSSatish Balay } 4783270192aSSatish Balay } 4793270192aSSatish Balay color = DRAW_CYAN; 4803270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4813270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4823270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4833270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4843270192aSSatish Balay aa = a->a + j*bs2; 4853270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4863270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4873270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 488b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4893270192aSSatish Balay } 4903270192aSSatish Balay } 4913270192aSSatish Balay } 4923270192aSSatish Balay } 4933270192aSSatish Balay 4943270192aSSatish Balay color = DRAW_RED; 4953270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4963270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4973270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4983270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4993270192aSSatish Balay aa = a->a + j*bs2; 5003270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5013270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5023270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 503b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5043270192aSSatish Balay } 5053270192aSSatish Balay } 5063270192aSSatish Balay } 5073270192aSSatish Balay } 5083270192aSSatish Balay 50955843e3eSBarry Smith DrawSynchronizedFlush(draw); 5103270192aSSatish Balay DrawGetPause(draw,&pause); 5113a40ed3dSBarry Smith if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 5123270192aSSatish Balay 5133270192aSSatish Balay /* allow the matrix to zoom or shrink */ 51455843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5153270192aSSatish Balay while (button != BUTTON_RIGHT) { 51655843e3eSBarry Smith DrawSynchronizedClear(draw); 5173270192aSSatish Balay if (button == BUTTON_LEFT) scale = .5; 5183270192aSSatish Balay else if (button == BUTTON_CENTER) scale = 2.; 5193270192aSSatish Balay xl = scale*(xl + w - xc) + xc - w*scale; 5203270192aSSatish Balay xr = scale*(xr - w - xc) + xc + w*scale; 5213270192aSSatish Balay yl = scale*(yl + h - yc) + yc - h*scale; 5223270192aSSatish Balay yr = scale*(yr - h - yc) + yc + h*scale; 5233270192aSSatish Balay w *= scale; h *= scale; 5243270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 5253270192aSSatish Balay color = DRAW_BLUE; 5263270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5273270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5283270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5293270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5303270192aSSatish Balay aa = a->a + j*bs2; 5313270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5323270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5333270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 534b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5353270192aSSatish Balay } 5363270192aSSatish Balay } 5373270192aSSatish Balay } 5383270192aSSatish Balay } 5393270192aSSatish Balay color = DRAW_CYAN; 5403270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5413270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5423270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5433270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5443270192aSSatish Balay aa = a->a + j*bs2; 5453270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5463270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5473270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 548b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5493270192aSSatish Balay } 5503270192aSSatish Balay } 5513270192aSSatish Balay } 5523270192aSSatish Balay } 5533270192aSSatish Balay 5543270192aSSatish Balay color = DRAW_RED; 5553270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5563270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5573270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5583270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5593270192aSSatish Balay aa = a->a + j*bs2; 5603270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5613270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5623270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 563b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5643270192aSSatish Balay } 5653270192aSSatish Balay } 5663270192aSSatish Balay } 5673270192aSSatish Balay } 56855843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5693270192aSSatish Balay } 5703a40ed3dSBarry Smith PetscFunctionReturn(0); 5713270192aSSatish Balay } 5723270192aSSatish Balay 5735615d1e5SSatish Balay #undef __FUNC__ 574d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 575e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5762593348eSBarry Smith { 57719bcc07fSBarry Smith ViewerType vtype; 57819bcc07fSBarry Smith int ierr; 5792593348eSBarry Smith 5803a40ed3dSBarry Smith PetscFunctionBegin; 58119bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 58219bcc07fSBarry Smith if (vtype == MATLAB_VIEWER) { 583a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported"); 5843a40ed3dSBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 5853a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 5863a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 5873a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 5883a40ed3dSBarry Smith } else if (vtype == DRAW_VIEWER) { 5893a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 5905cd90555SBarry Smith } else { 5915cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5922593348eSBarry Smith } 5933a40ed3dSBarry Smith PetscFunctionReturn(0); 5942593348eSBarry Smith } 595b6490206SBarry Smith 596cd0e1443SSatish Balay 5975615d1e5SSatish Balay #undef __FUNC__ 5982d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 5992d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 600cd0e1443SSatish Balay { 601cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6022d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6032d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6042d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6052d61bbb3SSatish Balay Scalar *ap, *aa = a->a, zero = 0.0; 606cd0e1443SSatish Balay 6073a40ed3dSBarry Smith PetscFunctionBegin; 6082d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 609cd0e1443SSatish Balay row = im[k]; brow = row/bs; 610a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 611a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6122d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6132c3acbe9SBarry Smith nrow = ailen[brow]; 6142d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 615a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 616a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6172d61bbb3SSatish Balay col = in[l] ; 6182d61bbb3SSatish Balay bcol = col/bs; 6192d61bbb3SSatish Balay cidx = col%bs; 6202d61bbb3SSatish Balay ridx = row%bs; 6212d61bbb3SSatish Balay high = nrow; 6222d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6232d61bbb3SSatish Balay while (high-low > 5) { 624cd0e1443SSatish Balay t = (low+high)/2; 625cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 626cd0e1443SSatish Balay else low = t; 627cd0e1443SSatish Balay } 628cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 629cd0e1443SSatish Balay if (rp[i] > bcol) break; 630cd0e1443SSatish Balay if (rp[i] == bcol) { 6312d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6322d61bbb3SSatish Balay goto finished; 633cd0e1443SSatish Balay } 634cd0e1443SSatish Balay } 6352d61bbb3SSatish Balay *v++ = zero; 6362d61bbb3SSatish Balay finished:; 637cd0e1443SSatish Balay } 638cd0e1443SSatish Balay } 6393a40ed3dSBarry Smith PetscFunctionReturn(0); 640cd0e1443SSatish Balay } 641cd0e1443SSatish Balay 6422d61bbb3SSatish Balay 6435615d1e5SSatish Balay #undef __FUNC__ 64405a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 64592c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 64692c4ed94SBarry Smith { 64792c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6488a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 649dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 650dd9472c6SBarry Smith int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval; 6510e324ae4SSatish Balay Scalar *ap,*value=v,*aa=a->a,*bap; 65292c4ed94SBarry Smith 6533a40ed3dSBarry Smith PetscFunctionBegin; 6540e324ae4SSatish Balay if (roworiented) { 6550e324ae4SSatish Balay stepval = (n-1)*bs; 6560e324ae4SSatish Balay } else { 6570e324ae4SSatish Balay stepval = (m-1)*bs; 6580e324ae4SSatish Balay } 65992c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 66092c4ed94SBarry Smith row = im[k]; 6613a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 662a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 663a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 66492c4ed94SBarry Smith #endif 66592c4ed94SBarry Smith rp = aj + ai[row]; 66692c4ed94SBarry Smith ap = aa + bs2*ai[row]; 66792c4ed94SBarry Smith rmax = imax[row]; 66892c4ed94SBarry Smith nrow = ailen[row]; 66992c4ed94SBarry Smith low = 0; 67092c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6713a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 672a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 673a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 67492c4ed94SBarry Smith #endif 67592c4ed94SBarry Smith col = in[l]; 67692c4ed94SBarry Smith if (roworiented) { 6770e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6780e324ae4SSatish Balay } else { 6790e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 68092c4ed94SBarry Smith } 68192c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 68292c4ed94SBarry Smith while (high-low > 7) { 68392c4ed94SBarry Smith t = (low+high)/2; 68492c4ed94SBarry Smith if (rp[t] > col) high = t; 68592c4ed94SBarry Smith else low = t; 68692c4ed94SBarry Smith } 68792c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 68892c4ed94SBarry Smith if (rp[i] > col) break; 68992c4ed94SBarry Smith if (rp[i] == col) { 6908a84c255SSatish Balay bap = ap + bs2*i; 6910e324ae4SSatish Balay if (roworiented) { 6928a84c255SSatish Balay if (is == ADD_VALUES) { 693dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 694dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 6958a84c255SSatish Balay bap[jj] += *value++; 696dd9472c6SBarry Smith } 697dd9472c6SBarry Smith } 6980e324ae4SSatish Balay } else { 699dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 700dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7010e324ae4SSatish Balay bap[jj] = *value++; 7028a84c255SSatish Balay } 703dd9472c6SBarry Smith } 704dd9472c6SBarry Smith } 7050e324ae4SSatish Balay } else { 7060e324ae4SSatish Balay if (is == ADD_VALUES) { 707dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 708dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7090e324ae4SSatish Balay *bap++ += *value++; 710dd9472c6SBarry Smith } 711dd9472c6SBarry Smith } 7120e324ae4SSatish Balay } else { 713dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 714dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7150e324ae4SSatish Balay *bap++ = *value++; 7160e324ae4SSatish Balay } 7178a84c255SSatish Balay } 718dd9472c6SBarry Smith } 719dd9472c6SBarry Smith } 720f1241b54SBarry Smith goto noinsert2; 72192c4ed94SBarry Smith } 72292c4ed94SBarry Smith } 72389280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 724a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 72592c4ed94SBarry Smith if (nrow >= rmax) { 72692c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 72792c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 72892c4ed94SBarry Smith Scalar *new_a; 72992c4ed94SBarry Smith 730a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 73189280ab3SLois Curfman McInnes 73292c4ed94SBarry Smith /* malloc new storage space */ 73392c4ed94SBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 73492c4ed94SBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 73592c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 73692c4ed94SBarry Smith new_i = new_j + new_nz; 73792c4ed94SBarry Smith 73892c4ed94SBarry Smith /* copy over old data into new slots */ 73992c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 74092c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 74192c4ed94SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int)); 74292c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 743dd9472c6SBarry Smith PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int)); 74492c4ed94SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar)); 74592c4ed94SBarry Smith PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 746dd9472c6SBarry Smith PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar)); 74792c4ed94SBarry Smith /* free up old matrix storage */ 74892c4ed94SBarry Smith PetscFree(a->a); 74992c4ed94SBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 75092c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 75192c4ed94SBarry Smith a->singlemalloc = 1; 75292c4ed94SBarry Smith 75392c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 75492c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 75592c4ed94SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 75692c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 75792c4ed94SBarry Smith a->reallocs++; 75892c4ed94SBarry Smith a->nz++; 75992c4ed94SBarry Smith } 76092c4ed94SBarry Smith N = nrow++ - 1; 76192c4ed94SBarry Smith /* shift up all the later entries in this row */ 76292c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 76392c4ed94SBarry Smith rp[ii+1] = rp[ii]; 76492c4ed94SBarry Smith PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 76592c4ed94SBarry Smith } 76692c4ed94SBarry Smith if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 76792c4ed94SBarry Smith rp[i] = col; 7688a84c255SSatish Balay bap = ap + bs2*i; 7690e324ae4SSatish Balay if (roworiented) { 770dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 771dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7720e324ae4SSatish Balay bap[jj] = *value++; 773dd9472c6SBarry Smith } 774dd9472c6SBarry Smith } 7750e324ae4SSatish Balay } else { 776dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 777dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7780e324ae4SSatish Balay *bap++ = *value++; 7790e324ae4SSatish Balay } 780dd9472c6SBarry Smith } 781dd9472c6SBarry Smith } 782f1241b54SBarry Smith noinsert2:; 78392c4ed94SBarry Smith low = i; 78492c4ed94SBarry Smith } 78592c4ed94SBarry Smith ailen[row] = nrow; 78692c4ed94SBarry Smith } 7873a40ed3dSBarry Smith PetscFunctionReturn(0); 78892c4ed94SBarry Smith } 78992c4ed94SBarry Smith 790f2501298SSatish Balay 7915615d1e5SSatish Balay #undef __FUNC__ 7925615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 7938f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 794584200bdSSatish Balay { 795584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 796584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 797a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 79843ee02c3SBarry Smith int mbs = a->mbs, bs2 = a->bs2,rmax = 0; 799584200bdSSatish Balay Scalar *aa = a->a, *ap; 800584200bdSSatish Balay 8013a40ed3dSBarry Smith PetscFunctionBegin; 8023a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 803584200bdSSatish Balay 80443ee02c3SBarry Smith if (m) rmax = ailen[0]; 805584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 806584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 807584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 808d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 809584200bdSSatish Balay if (fshift) { 810a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 811584200bdSSatish Balay N = ailen[i]; 812584200bdSSatish Balay for ( j=0; j<N; j++ ) { 813584200bdSSatish Balay ip[j-fshift] = ip[j]; 8147e67e3f9SSatish Balay PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar)); 815584200bdSSatish Balay } 816584200bdSSatish Balay } 817584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 818584200bdSSatish Balay } 819584200bdSSatish Balay if (mbs) { 820584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 821584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 822584200bdSSatish Balay } 823584200bdSSatish Balay /* reset ilen and imax for each row */ 824584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 825584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 826584200bdSSatish Balay } 827a7c10996SSatish Balay a->nz = ai[mbs]; 828584200bdSSatish Balay 829584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 830584200bdSSatish Balay if (fshift && a->diag) { 831584200bdSSatish Balay PetscFree(a->diag); 832584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 833584200bdSSatish Balay a->diag = 0; 834584200bdSSatish Balay } 8353d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 836219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8373d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 838584200bdSSatish Balay a->reallocs); 839d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 840e2f3b5e9SSatish Balay a->reallocs = 0; 8414e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8424e220ebcSLois Curfman McInnes 8433a40ed3dSBarry Smith PetscFunctionReturn(0); 844584200bdSSatish Balay } 845584200bdSSatish Balay 8465a838052SSatish Balay 847d9b7c43dSSatish Balay /* idx should be of length atlease bs */ 8485615d1e5SSatish Balay #undef __FUNC__ 8495615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block" 850d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg) 851d9b7c43dSSatish Balay { 852d9b7c43dSSatish Balay int i,row; 8533a40ed3dSBarry Smith 8543a40ed3dSBarry Smith PetscFunctionBegin; 855d9b7c43dSSatish Balay row = idx[0]; 8563a40ed3dSBarry Smith if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 857d9b7c43dSSatish Balay 858d9b7c43dSSatish Balay for ( i=1; i<bs; i++ ) { 8593a40ed3dSBarry Smith if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 860d9b7c43dSSatish Balay } 861d9b7c43dSSatish Balay *flg = PETSC_TRUE; 8623a40ed3dSBarry Smith PetscFunctionReturn(0); 863d9b7c43dSSatish Balay } 864d9b7c43dSSatish Balay 8655615d1e5SSatish Balay #undef __FUNC__ 8665615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 8678f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 868d9b7c43dSSatish Balay { 869d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 870d9b7c43dSSatish Balay IS is_local; 871d9b7c43dSSatish Balay int ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2; 872d9b7c43dSSatish Balay PetscTruth flg; 873d9b7c43dSSatish Balay Scalar zero = 0.0,*aa; 874d9b7c43dSSatish Balay 8753a40ed3dSBarry Smith PetscFunctionBegin; 876d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 877d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 878d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 879537820f0SBarry Smith ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr); 880d9b7c43dSSatish Balay ierr = ISSort(is_local); CHKERRQ(ierr); 881d9b7c43dSSatish Balay ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr); 882d9b7c43dSSatish Balay 883d9b7c43dSSatish Balay i = 0; 884d9b7c43dSSatish Balay while (i < is_n) { 885a8c6a408SBarry Smith if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 886d9b7c43dSSatish Balay flg = PETSC_FALSE; 887d9b7c43dSSatish Balay if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); } 888d9b7c43dSSatish Balay count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs; 889d9b7c43dSSatish Balay aa = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs); 8902d61bbb3SSatish Balay if (flg) { /* There exists a block of rows to be Zerowed */ 891a07cd24cSSatish Balay if (baij->ilen[rows[i]/bs] > 0) { 8922d61bbb3SSatish Balay PetscMemzero(aa,count*bs*sizeof(Scalar)); 8932d61bbb3SSatish Balay baij->ilen[rows[i]/bs] = 1; 8942d61bbb3SSatish Balay baij->j[baij->i[rows[i]/bs]] = rows[i]/bs; 895a07cd24cSSatish Balay } 8962d61bbb3SSatish Balay i += bs; 8972d61bbb3SSatish Balay } else { /* Zero out only the requested row */ 898d9b7c43dSSatish Balay for ( j=0; j<count; j++ ) { 899d9b7c43dSSatish Balay aa[0] = zero; 900d9b7c43dSSatish Balay aa+=bs; 901d9b7c43dSSatish Balay } 902d9b7c43dSSatish Balay i++; 903d9b7c43dSSatish Balay } 904d9b7c43dSSatish Balay } 905d9b7c43dSSatish Balay if (diag) { 906d9b7c43dSSatish Balay for ( j=0; j<is_n; j++ ) { 907f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 9082d61bbb3SSatish Balay /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */ 909d9b7c43dSSatish Balay } 910d9b7c43dSSatish Balay } 911d9b7c43dSSatish Balay ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr); 912d9b7c43dSSatish Balay ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr); 913d9b7c43dSSatish Balay ierr = ISDestroy(is_local); CHKERRQ(ierr); 9149a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 9153a40ed3dSBarry Smith PetscFunctionReturn(0); 916d9b7c43dSSatish Balay } 9171c351548SSatish Balay 9185615d1e5SSatish Balay #undef __FUNC__ 9192d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9202d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9212d61bbb3SSatish Balay { 9222d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9232d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9242d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9252d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 9262d61bbb3SSatish Balay int ridx,cidx,bs2=a->bs2; 9272d61bbb3SSatish Balay Scalar *ap,value,*aa=a->a,*bap; 9282d61bbb3SSatish Balay 9292d61bbb3SSatish Balay PetscFunctionBegin; 9302d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9312d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9322d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9332d61bbb3SSatish Balay if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 9342d61bbb3SSatish Balay if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 9352d61bbb3SSatish Balay #endif 9362d61bbb3SSatish Balay rp = aj + ai[brow]; 9372d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9382d61bbb3SSatish Balay rmax = imax[brow]; 9392d61bbb3SSatish Balay nrow = ailen[brow]; 9402d61bbb3SSatish Balay low = 0; 9412d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9422d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9432d61bbb3SSatish Balay if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 9442d61bbb3SSatish Balay if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 9452d61bbb3SSatish Balay #endif 9462d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 9472d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 9482d61bbb3SSatish Balay if (roworiented) { 9492d61bbb3SSatish Balay value = *v++; 9502d61bbb3SSatish Balay } else { 9512d61bbb3SSatish Balay value = v[k + l*m]; 9522d61bbb3SSatish Balay } 9532d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 9542d61bbb3SSatish Balay while (high-low > 7) { 9552d61bbb3SSatish Balay t = (low+high)/2; 9562d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 9572d61bbb3SSatish Balay else low = t; 9582d61bbb3SSatish Balay } 9592d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 9602d61bbb3SSatish Balay if (rp[i] > bcol) break; 9612d61bbb3SSatish Balay if (rp[i] == bcol) { 9622d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 9632d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 9642d61bbb3SSatish Balay else *bap = value; 9652d61bbb3SSatish Balay goto noinsert1; 9662d61bbb3SSatish Balay } 9672d61bbb3SSatish Balay } 9682d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 9692d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9702d61bbb3SSatish Balay if (nrow >= rmax) { 9712d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 9722d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 9732d61bbb3SSatish Balay Scalar *new_a; 9742d61bbb3SSatish Balay 9752d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9762d61bbb3SSatish Balay 9772d61bbb3SSatish Balay /* Malloc new storage space */ 9782d61bbb3SSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 9792d61bbb3SSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 9802d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 9812d61bbb3SSatish Balay new_i = new_j + new_nz; 9822d61bbb3SSatish Balay 9832d61bbb3SSatish Balay /* copy over old data into new slots */ 9842d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 9852d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 9862d61bbb3SSatish Balay PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int)); 9872d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 9882d61bbb3SSatish Balay PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow, 9892d61bbb3SSatish Balay len*sizeof(int)); 9902d61bbb3SSatish Balay PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar)); 9912d61bbb3SSatish Balay PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 9922d61bbb3SSatish Balay PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE), 9932d61bbb3SSatish Balay aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar)); 9942d61bbb3SSatish Balay /* free up old matrix storage */ 9952d61bbb3SSatish Balay PetscFree(a->a); 9962d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 9972d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 9982d61bbb3SSatish Balay a->singlemalloc = 1; 9992d61bbb3SSatish Balay 10002d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10012d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10022d61bbb3SSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 10032d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10042d61bbb3SSatish Balay a->reallocs++; 10052d61bbb3SSatish Balay a->nz++; 10062d61bbb3SSatish Balay } 10072d61bbb3SSatish Balay N = nrow++ - 1; 10082d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10092d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10102d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 10112d61bbb3SSatish Balay PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 10122d61bbb3SSatish Balay } 10132d61bbb3SSatish Balay if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 10142d61bbb3SSatish Balay rp[i] = bcol; 10152d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10162d61bbb3SSatish Balay noinsert1:; 10172d61bbb3SSatish Balay low = i; 10182d61bbb3SSatish Balay } 10192d61bbb3SSatish Balay ailen[brow] = nrow; 10202d61bbb3SSatish Balay } 10212d61bbb3SSatish Balay PetscFunctionReturn(0); 10222d61bbb3SSatish Balay } 10232d61bbb3SSatish Balay 10242d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10252d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10262d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1027d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*); 10282d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**); 10292d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10302d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10312d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10322d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10332d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10342d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10352d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10362d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10372d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10382d61bbb3SSatish Balay 10392d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10402d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 10412d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 10422d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 10432d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 10442d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 10452d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 10462d61bbb3SSatish Balay 10472d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 10482d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 10492d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 10502d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 10512d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 10522d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 10532d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 10542d61bbb3SSatish Balay 10552d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 10562d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 10572d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 10582d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 10592d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 10602d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 10612d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 10622d61bbb3SSatish Balay 10632d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 10642d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 10652d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 10662d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 10672d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 10682d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 10692d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 10702d61bbb3SSatish Balay 10712d61bbb3SSatish Balay /* 10722d61bbb3SSatish Balay note: This can only work for identity for row and col. It would 10732d61bbb3SSatish Balay be good to check this and otherwise generate an error. 10742d61bbb3SSatish Balay */ 10752d61bbb3SSatish Balay #undef __FUNC__ 10762d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 10772d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 10782d61bbb3SSatish Balay { 10792d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 10802d61bbb3SSatish Balay Mat outA; 10812d61bbb3SSatish Balay int ierr; 10822d61bbb3SSatish Balay 10832d61bbb3SSatish Balay PetscFunctionBegin; 10842d61bbb3SSatish Balay if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 10852d61bbb3SSatish Balay 10862d61bbb3SSatish Balay outA = inA; 10872d61bbb3SSatish Balay inA->factor = FACTOR_LU; 10882d61bbb3SSatish Balay a->row = row; 10892d61bbb3SSatish Balay a->col = col; 10902d61bbb3SSatish Balay 1091e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1092e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 10931e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1094e51c0b9cSSatish Balay 10952d61bbb3SSatish Balay if (!a->solve_work) { 10962d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 10972d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 10982d61bbb3SSatish Balay } 10992d61bbb3SSatish Balay 11002d61bbb3SSatish Balay if (!a->diag) { 11012d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 11022d61bbb3SSatish Balay } 11032d61bbb3SSatish Balay ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 11042d61bbb3SSatish Balay 11052d61bbb3SSatish Balay /* 11062d61bbb3SSatish Balay Blocksize 4 has a special faster solver for ILU(0) factorization 11072d61bbb3SSatish Balay with natural ordering 11082d61bbb3SSatish Balay */ 11092d61bbb3SSatish Balay if (a->bs == 4) { 1110f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 11112d61bbb3SSatish Balay } 11122d61bbb3SSatish Balay 11132d61bbb3SSatish Balay PetscFunctionReturn(0); 11142d61bbb3SSatish Balay } 11152d61bbb3SSatish Balay #undef __FUNC__ 1116d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1117ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1118ba4ca20aSSatish Balay { 1119ba4ca20aSSatish Balay static int called = 0; 1120ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1121ba4ca20aSSatish Balay 11223a40ed3dSBarry Smith PetscFunctionBegin; 11233a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 112476be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n"); 112576be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n"); 11263a40ed3dSBarry Smith PetscFunctionReturn(0); 1127ba4ca20aSSatish Balay } 1128d9b7c43dSSatish Balay 1129*27a8da17SBarry Smith #undef __FUNC__ 1130*27a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 1131*27a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 1132*27a8da17SBarry Smith { 1133*27a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 1134*27a8da17SBarry Smith int i,nz,n; 1135*27a8da17SBarry Smith 1136*27a8da17SBarry Smith PetscFunctionBegin; 1137*27a8da17SBarry Smith nz = baij->maxnz; 1138*27a8da17SBarry Smith n = baij->n; 1139*27a8da17SBarry Smith for (i=0; i<nz; i++) { 1140*27a8da17SBarry Smith baij->j[i] = indices[i]; 1141*27a8da17SBarry Smith } 1142*27a8da17SBarry Smith baij->nz = nz; 1143*27a8da17SBarry Smith for ( i=0; i<n; i++ ) { 1144*27a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 1145*27a8da17SBarry Smith } 1146*27a8da17SBarry Smith 1147*27a8da17SBarry Smith PetscFunctionReturn(0); 1148*27a8da17SBarry Smith } 1149*27a8da17SBarry Smith 1150*27a8da17SBarry Smith #undef __FUNC__ 1151*27a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 1152*27a8da17SBarry Smith /*@ 1153*27a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 1154*27a8da17SBarry Smith in the matrix. 1155*27a8da17SBarry Smith 1156*27a8da17SBarry Smith Input Parameters: 1157*27a8da17SBarry Smith + mat - the SeqBAIJ matrix 1158*27a8da17SBarry Smith - indices - the column indices 1159*27a8da17SBarry Smith 1160*27a8da17SBarry Smith Notes: 1161*27a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 1162*27a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 1163*27a8da17SBarry Smith of the MatSetValues() operation. 1164*27a8da17SBarry Smith 1165*27a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 1166*27a8da17SBarry Smith MatCreateSeqBAIJ(). 1167*27a8da17SBarry Smith 1168*27a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 1169*27a8da17SBarry Smith 1170*27a8da17SBarry Smith @*/ 1171*27a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 1172*27a8da17SBarry Smith { 1173*27a8da17SBarry Smith int ierr,(*f)(Mat,int *); 1174*27a8da17SBarry Smith 1175*27a8da17SBarry Smith PetscFunctionBegin; 1176*27a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1177*27a8da17SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f); 1178*27a8da17SBarry Smith CHKERRQ(ierr); 1179*27a8da17SBarry Smith if (f) { 1180*27a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 1181*27a8da17SBarry Smith } else { 1182*27a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 1183*27a8da17SBarry Smith } 1184*27a8da17SBarry Smith PetscFunctionReturn(0); 1185*27a8da17SBarry Smith } 1186*27a8da17SBarry Smith 11872593348eSBarry Smith /* -------------------------------------------------------------------*/ 1188cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ, 11899867e207SSatish Balay MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ, 1190f44d3a6dSSatish Balay MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N, 1191faf6580fSSatish Balay MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ, 11921a6a6d98SBarry Smith MatSolve_SeqBAIJ_N,0, 1193b6490206SBarry Smith 0,0, 1194de6a44a3SBarry Smith MatLUFactor_SeqBAIJ,0, 1195b6490206SBarry Smith 0, 1196f2501298SSatish Balay MatTranspose_SeqBAIJ, 119717e48fc4SSatish Balay MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ, 11989867e207SSatish Balay MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ, 1199584200bdSSatish Balay 0,MatAssemblyEnd_SeqBAIJ, 1200b6490206SBarry Smith 0, 1201d9b7c43dSSatish Balay MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ, 12027fc0212eSBarry Smith MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 1203b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 1204de6a44a3SBarry Smith MatILUFactorSymbolic_SeqBAIJ,0, 1205d402145bSBarry Smith 0,0, 1206b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 1207b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 1208af015674SSatish Balay MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ, 12097e67e3f9SSatish Balay MatGetValues_SeqBAIJ,0, 1210ba4ca20aSSatish Balay MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 12113b2fbd54SBarry Smith 0,0,0,MatGetBlockSize_SeqBAIJ, 12123b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 121392c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 121492c4ed94SBarry Smith 0,0,0,0,0,0, 1215d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1216d3825aa8SBarry Smith MatGetSubMatrix_SeqBAIJ}; 12172593348eSBarry Smith 12185615d1e5SSatish Balay #undef __FUNC__ 12195615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 12202593348eSBarry Smith /*@C 122144cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 122244cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 122344cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 12242bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 12252bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 12262593348eSBarry Smith 1227db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1228db81eaa0SLois Curfman McInnes 12292593348eSBarry Smith Input Parameters: 1230db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1231b6490206SBarry Smith . bs - size of block 12322593348eSBarry Smith . m - number of rows 12332593348eSBarry Smith . n - number of columns 1234b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 1235db81eaa0SLois Curfman McInnes - nzz - array containing the number of block nonzeros in the various block rows 12362bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 12372593348eSBarry Smith 12382593348eSBarry Smith Output Parameter: 12392593348eSBarry Smith . A - the matrix 12402593348eSBarry Smith 12410513a670SBarry Smith Options Database Keys: 1242db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1243db81eaa0SLois Curfman McInnes block calculations (much slower) 1244db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 12450513a670SBarry Smith 12462593348eSBarry Smith Notes: 124744cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 12482593348eSBarry Smith storage. That is, the stored row and column indices can begin at 124944cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 12502593348eSBarry Smith 12512593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 12522593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 12532593348eSBarry Smith allocation. For additional details, see the users manual chapter on 12546da5968aSLois Curfman McInnes matrices. 12552593348eSBarry Smith 1256db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 12572593348eSBarry Smith @*/ 1258b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 12592593348eSBarry Smith { 12602593348eSBarry Smith Mat B; 1261b6490206SBarry Smith Mat_SeqBAIJ *b; 12623b2fbd54SBarry Smith int i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size; 12633b2fbd54SBarry Smith 12643a40ed3dSBarry Smith PetscFunctionBegin; 12653b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1266a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1267b6490206SBarry Smith 12680513a670SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 12690513a670SBarry Smith 12703a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1271a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 12723a40ed3dSBarry Smith } 12732593348eSBarry Smith 12742593348eSBarry Smith *A = 0; 1275f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView); 12762593348eSBarry Smith PLogObjectCreate(B); 1277b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 127844cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqBAIJ)); 1279f830108cSBarry Smith PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps)); 12801a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr); 12811a6a6d98SBarry Smith if (!flg) { 12827fc0212eSBarry Smith switch (bs) { 12837fc0212eSBarry Smith case 1: 1284f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1285f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1286f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1287f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 12887fc0212eSBarry Smith break; 12894eeb42bcSBarry Smith case 2: 1290f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1291f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1292f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1293f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 12946d84be18SBarry Smith break; 1295f327dd97SBarry Smith case 3: 1296f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1297f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1298f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1299f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 13004eeb42bcSBarry Smith break; 13016d84be18SBarry Smith case 4: 1302f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1303f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1304f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1305f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 13066d84be18SBarry Smith break; 13076d84be18SBarry Smith case 5: 1308f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1309f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1310f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1311f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 13126d84be18SBarry Smith break; 131348e9ddb2SSatish Balay case 7: 1314f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1315f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1316f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 131748e9ddb2SSatish Balay break; 13187fc0212eSBarry Smith } 13191a6a6d98SBarry Smith } 1320e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1321e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 13222593348eSBarry Smith B->factor = 0; 13232593348eSBarry Smith B->lupivotthreshold = 1.0; 132490f02eecSBarry Smith B->mapping = 0; 13252593348eSBarry Smith b->row = 0; 13262593348eSBarry Smith b->col = 0; 1327e51c0b9cSSatish Balay b->icol = 0; 13282593348eSBarry Smith b->reallocs = 0; 13292593348eSBarry Smith 133044cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 133144cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1332b6490206SBarry Smith b->mbs = mbs; 1333f2501298SSatish Balay b->nbs = nbs; 1334b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 13352593348eSBarry Smith if (nnz == PETSC_NULL) { 1336119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 13372593348eSBarry Smith else if (nz <= 0) nz = 1; 1338b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1339b6490206SBarry Smith nz = nz*mbs; 13403a40ed3dSBarry Smith } else { 13412593348eSBarry Smith nz = 0; 1342b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 13432593348eSBarry Smith } 13442593348eSBarry Smith 13452593348eSBarry Smith /* allocate the matrix space */ 13467e67e3f9SSatish Balay len = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int); 13472593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 13487e67e3f9SSatish Balay PetscMemzero(b->a,nz*bs2*sizeof(Scalar)); 13497e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 13502593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 13512593348eSBarry Smith b->i = b->j + nz; 13522593348eSBarry Smith b->singlemalloc = 1; 13532593348eSBarry Smith 1354de6a44a3SBarry Smith b->i[0] = 0; 1355b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 13562593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 13572593348eSBarry Smith } 13582593348eSBarry Smith 1359b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1360b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1361f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1362b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 13632593348eSBarry Smith 1364b6490206SBarry Smith b->bs = bs; 13659df24120SSatish Balay b->bs2 = bs2; 1366b6490206SBarry Smith b->mbs = mbs; 13672593348eSBarry Smith b->nz = 0; 136818e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 13692593348eSBarry Smith b->sorted = 0; 13702593348eSBarry Smith b->roworiented = 1; 13712593348eSBarry Smith b->nonew = 0; 13722593348eSBarry Smith b->diag = 0; 13732593348eSBarry Smith b->solve_work = 0; 1374de6a44a3SBarry Smith b->mult_work = 0; 13752593348eSBarry Smith b->spptr = 0; 13764e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 13774e220ebcSLois Curfman McInnes 13782593348eSBarry Smith *A = B; 13792593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 13802593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 1381*27a8da17SBarry Smith 1382*27a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 1383*27a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1384*27a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 1385*27a8da17SBarry Smith 13863a40ed3dSBarry Smith PetscFunctionReturn(0); 13872593348eSBarry Smith } 13882593348eSBarry Smith 13895615d1e5SSatish Balay #undef __FUNC__ 13905615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ" 1391b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 13922593348eSBarry Smith { 13932593348eSBarry Smith Mat C; 1394b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 1395*27a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1396de6a44a3SBarry Smith 13973a40ed3dSBarry Smith PetscFunctionBegin; 1398a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 13992593348eSBarry Smith 14002593348eSBarry Smith *B = 0; 1401f830108cSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView); 14022593348eSBarry Smith PLogObjectCreate(C); 1403b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 1404c3c66ee5SSatish Balay PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 1405e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1406e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 14072593348eSBarry Smith C->factor = A->factor; 14082593348eSBarry Smith c->row = 0; 14092593348eSBarry Smith c->col = 0; 14102593348eSBarry Smith C->assembled = PETSC_TRUE; 14112593348eSBarry Smith 141244cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 141344cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 141444cd7ae7SLois Curfman McInnes C->M = a->m; 141544cd7ae7SLois Curfman McInnes C->N = a->n; 141644cd7ae7SLois Curfman McInnes 1417b6490206SBarry Smith c->bs = a->bs; 14189df24120SSatish Balay c->bs2 = a->bs2; 1419b6490206SBarry Smith c->mbs = a->mbs; 1420b6490206SBarry Smith c->nbs = a->nbs; 14212593348eSBarry Smith 1422b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1423b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1424b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 14252593348eSBarry Smith c->imax[i] = a->imax[i]; 14262593348eSBarry Smith c->ilen[i] = a->ilen[i]; 14272593348eSBarry Smith } 14282593348eSBarry Smith 14292593348eSBarry Smith /* allocate the matrix space */ 14302593348eSBarry Smith c->singlemalloc = 1; 14317e67e3f9SSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int)); 14322593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 14337e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1434de6a44a3SBarry Smith c->i = c->j + nz; 1435b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1436b6490206SBarry Smith if (mbs > 0) { 1437de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 14382593348eSBarry Smith if (cpvalues == COPY_VALUES) { 14397e67e3f9SSatish Balay PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar)); 14402593348eSBarry Smith } 14412593348eSBarry Smith } 14422593348eSBarry Smith 1443f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 14442593348eSBarry Smith c->sorted = a->sorted; 14452593348eSBarry Smith c->roworiented = a->roworiented; 14462593348eSBarry Smith c->nonew = a->nonew; 14472593348eSBarry Smith 14482593348eSBarry Smith if (a->diag) { 1449b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1450b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1451b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 14522593348eSBarry Smith c->diag[i] = a->diag[i]; 14532593348eSBarry Smith } 14542593348eSBarry Smith } 14552593348eSBarry Smith else c->diag = 0; 14562593348eSBarry Smith c->nz = a->nz; 14572593348eSBarry Smith c->maxnz = a->maxnz; 14582593348eSBarry Smith c->solve_work = 0; 14592593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 14607fc0212eSBarry Smith c->mult_work = 0; 14612593348eSBarry Smith *B = C; 1462*27a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqBAIJSetColumnIndices_C", 1463*27a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1464*27a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 14653a40ed3dSBarry Smith PetscFunctionReturn(0); 14662593348eSBarry Smith } 14672593348eSBarry Smith 14685615d1e5SSatish Balay #undef __FUNC__ 14695615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 147019bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 14712593348eSBarry Smith { 1472b6490206SBarry Smith Mat_SeqBAIJ *a; 14732593348eSBarry Smith Mat B; 1474de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1475b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 147635aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1477a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1478b6490206SBarry Smith Scalar *aa; 147919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 14802593348eSBarry Smith 14813a40ed3dSBarry Smith PetscFunctionBegin; 14821a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1483de6a44a3SBarry Smith bs2 = bs*bs; 1484b6490206SBarry Smith 14852593348eSBarry Smith MPI_Comm_size(comm,&size); 1486cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 148790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 14880752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 1489a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 14902593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 14912593348eSBarry Smith 1492d64ed03dSBarry Smith if (header[3] < 0) { 1493a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1494d64ed03dSBarry Smith } 1495d64ed03dSBarry Smith 1496a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 149735aab85fSBarry Smith 149835aab85fSBarry Smith /* 149935aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 150035aab85fSBarry Smith divisible by the blocksize 150135aab85fSBarry Smith */ 1502b6490206SBarry Smith mbs = M/bs; 150335aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 150435aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 150535aab85fSBarry Smith else mbs++; 150635aab85fSBarry Smith if (extra_rows) { 1507537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 150835aab85fSBarry Smith } 1509b6490206SBarry Smith 15102593348eSBarry Smith /* read in row lengths */ 151135aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 15120752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 151335aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 15142593348eSBarry Smith 1515b6490206SBarry Smith /* read in column indices */ 151635aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 15170752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr); 151835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1519b6490206SBarry Smith 1520b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1521b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1522b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 152335aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 152435aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 152535aab85fSBarry Smith masked = mask + mbs; 1526b6490206SBarry Smith rowcount = 0; nzcount = 0; 1527b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 152835aab85fSBarry Smith nmask = 0; 1529b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1530b6490206SBarry Smith kmax = rowlengths[rowcount]; 1531b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 153235aab85fSBarry Smith tmp = jj[nzcount++]/bs; 153335aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1534b6490206SBarry Smith } 1535b6490206SBarry Smith rowcount++; 1536b6490206SBarry Smith } 153735aab85fSBarry Smith browlengths[i] += nmask; 153835aab85fSBarry Smith /* zero out the mask elements we set */ 153935aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1540b6490206SBarry Smith } 1541b6490206SBarry Smith 15422593348eSBarry Smith /* create our matrix */ 15433eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 15442593348eSBarry Smith B = *A; 1545b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 15462593348eSBarry Smith 15472593348eSBarry Smith /* set matrix "i" values */ 1548de6a44a3SBarry Smith a->i[0] = 0; 1549b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1550b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1551b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 15522593348eSBarry Smith } 1553b6490206SBarry Smith a->nz = 0; 1554b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 15552593348eSBarry Smith 1556b6490206SBarry Smith /* read in nonzero values */ 155735aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 15580752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr); 155935aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1560b6490206SBarry Smith 1561b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1562b6490206SBarry Smith nzcount = 0; jcount = 0; 1563b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1564b6490206SBarry Smith nzcountb = nzcount; 156535aab85fSBarry Smith nmask = 0; 1566b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1567b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1568b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 156935aab85fSBarry Smith tmp = jj[nzcount++]/bs; 157035aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1571b6490206SBarry Smith } 1572b6490206SBarry Smith rowcount++; 1573b6490206SBarry Smith } 1574de6a44a3SBarry Smith /* sort the masked values */ 157577c4ece6SBarry Smith PetscSortInt(nmask,masked); 1576de6a44a3SBarry Smith 1577b6490206SBarry Smith /* set "j" values into matrix */ 1578b6490206SBarry Smith maskcount = 1; 157935aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 158035aab85fSBarry Smith a->j[jcount++] = masked[j]; 1581de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1582b6490206SBarry Smith } 1583b6490206SBarry Smith /* set "a" values into matrix */ 1584de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1585b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1586b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1587b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1588de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1589de6a44a3SBarry Smith block = mask[tmp] - 1; 1590de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1591de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1592b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1593b6490206SBarry Smith } 1594b6490206SBarry Smith } 159535aab85fSBarry Smith /* zero out the mask elements we set */ 159635aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1597b6490206SBarry Smith } 1598a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1599b6490206SBarry Smith 1600b6490206SBarry Smith PetscFree(rowlengths); 1601b6490206SBarry Smith PetscFree(browlengths); 1602b6490206SBarry Smith PetscFree(aa); 1603b6490206SBarry Smith PetscFree(jj); 1604b6490206SBarry Smith PetscFree(mask); 1605b6490206SBarry Smith 1606b6490206SBarry Smith B->assembled = PETSC_TRUE; 1607de6a44a3SBarry Smith 16089c01be13SBarry Smith ierr = MatView_Private(B); CHKERRQ(ierr); 16093a40ed3dSBarry Smith PetscFunctionReturn(0); 16102593348eSBarry Smith } 16112593348eSBarry Smith 16122593348eSBarry Smith 16132593348eSBarry Smith 1614