1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*c3c66ee5SSatish Balay static char vcid[] = "$Id: baij.c,v 1.126 1998/03/12 23:19:14 bsmith Exp balay $"; 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" 1052d61bbb3SSatish Balay int MatDestroy_SeqBAIJ(PetscObject obj) 1062d61bbb3SSatish Balay { 1072d61bbb3SSatish Balay Mat A = (Mat) obj; 1082d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1092d61bbb3SSatish Balay 1102d61bbb3SSatish Balay #if defined(USE_PETSC_LOG) 1112d61bbb3SSatish Balay PLogObjectState(obj,"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); 1202d61bbb3SSatish Balay PetscFree(a); 1212d61bbb3SSatish Balay PLogObjectDestroy(A); 1222d61bbb3SSatish Balay PetscHeaderDestroy(A); 1232d61bbb3SSatish Balay PetscFunctionReturn(0); 1242d61bbb3SSatish Balay } 1252d61bbb3SSatish Balay 1262d61bbb3SSatish Balay #undef __FUNC__ 1272d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1282d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1292d61bbb3SSatish Balay { 1302d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1312d61bbb3SSatish Balay 1322d61bbb3SSatish Balay PetscFunctionBegin; 1332d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1342d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1352d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1362d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1372d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1382d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 1392d61bbb3SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 1402d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1412d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1422d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1432d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1442d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1452d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 146b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 147b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 148981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1492d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1502d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1512d61bbb3SSatish Balay } else { 1522d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1532d61bbb3SSatish Balay } 1542d61bbb3SSatish Balay PetscFunctionReturn(0); 1552d61bbb3SSatish Balay } 1562d61bbb3SSatish Balay 1572d61bbb3SSatish Balay 1582d61bbb3SSatish Balay #undef __FUNC__ 1592d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1602d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1612d61bbb3SSatish Balay { 1622d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1632d61bbb3SSatish Balay 1642d61bbb3SSatish Balay PetscFunctionBegin; 1652d61bbb3SSatish Balay if (m) *m = a->m; 1662d61bbb3SSatish Balay if (n) *n = a->n; 1672d61bbb3SSatish Balay PetscFunctionReturn(0); 1682d61bbb3SSatish Balay } 1692d61bbb3SSatish Balay 1702d61bbb3SSatish Balay #undef __FUNC__ 1712d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 1722d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 1732d61bbb3SSatish Balay { 1742d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1752d61bbb3SSatish Balay 1762d61bbb3SSatish Balay PetscFunctionBegin; 1772d61bbb3SSatish Balay *m = 0; *n = a->m; 1782d61bbb3SSatish Balay PetscFunctionReturn(0); 1792d61bbb3SSatish Balay } 1802d61bbb3SSatish Balay 1812d61bbb3SSatish Balay #undef __FUNC__ 1822d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 1832d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 1842d61bbb3SSatish Balay { 1852d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1862d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 1872d61bbb3SSatish Balay Scalar *aa,*v_i,*aa_i; 1882d61bbb3SSatish Balay 1892d61bbb3SSatish Balay PetscFunctionBegin; 1902d61bbb3SSatish Balay bs = a->bs; 1912d61bbb3SSatish Balay ai = a->i; 1922d61bbb3SSatish Balay aj = a->j; 1932d61bbb3SSatish Balay aa = a->a; 1942d61bbb3SSatish Balay bs2 = a->bs2; 1952d61bbb3SSatish Balay 1962d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 1972d61bbb3SSatish Balay 1982d61bbb3SSatish Balay bn = row/bs; /* Block number */ 1992d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2002d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2012d61bbb3SSatish Balay *nz = bs*M; 2022d61bbb3SSatish Balay 2032d61bbb3SSatish Balay if (v) { 2042d61bbb3SSatish Balay *v = 0; 2052d61bbb3SSatish Balay if (*nz) { 2062d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 2072d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2082d61bbb3SSatish Balay v_i = *v + i*bs; 2092d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2102d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2112d61bbb3SSatish Balay } 2122d61bbb3SSatish Balay } 2132d61bbb3SSatish Balay } 2142d61bbb3SSatish Balay 2152d61bbb3SSatish Balay if (idx) { 2162d61bbb3SSatish Balay *idx = 0; 2172d61bbb3SSatish Balay if (*nz) { 2182d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 2192d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2202d61bbb3SSatish Balay idx_i = *idx + i*bs; 2212d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2222d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2232d61bbb3SSatish Balay } 2242d61bbb3SSatish Balay } 2252d61bbb3SSatish Balay } 2262d61bbb3SSatish Balay PetscFunctionReturn(0); 2272d61bbb3SSatish Balay } 2282d61bbb3SSatish Balay 2292d61bbb3SSatish Balay #undef __FUNC__ 2302d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2312d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2322d61bbb3SSatish Balay { 2332d61bbb3SSatish Balay PetscFunctionBegin; 2342d61bbb3SSatish Balay if (idx) {if (*idx) PetscFree(*idx);} 2352d61bbb3SSatish Balay if (v) {if (*v) PetscFree(*v);} 2362d61bbb3SSatish Balay PetscFunctionReturn(0); 2372d61bbb3SSatish Balay } 2382d61bbb3SSatish Balay 2392d61bbb3SSatish Balay #undef __FUNC__ 2402d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2412d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2422d61bbb3SSatish Balay { 2432d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2442d61bbb3SSatish Balay Mat C; 2452d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2462d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2472d61bbb3SSatish Balay Scalar *array=a->a; 2482d61bbb3SSatish Balay 2492d61bbb3SSatish Balay PetscFunctionBegin; 2502d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2512d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col); 2522d61bbb3SSatish Balay PetscMemzero(col,(1+nbs)*sizeof(int)); 2532d61bbb3SSatish Balay 2542d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2552d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr); 2562d61bbb3SSatish Balay PetscFree(col); 2572d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows); 2582d61bbb3SSatish Balay cols = rows + bs; 2592d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2602d61bbb3SSatish Balay cols[0] = i*bs; 2612d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 2622d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 2632d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 2642d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 2652d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 2662d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr); 2672d61bbb3SSatish Balay array += bs2; 2682d61bbb3SSatish Balay } 2692d61bbb3SSatish Balay } 2702d61bbb3SSatish Balay PetscFree(rows); 2712d61bbb3SSatish Balay 2722d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2732d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 2742d61bbb3SSatish Balay 2752d61bbb3SSatish Balay if (B != PETSC_NULL) { 2762d61bbb3SSatish Balay *B = C; 2772d61bbb3SSatish Balay } else { 278f830108cSBarry Smith PetscOps *Abops; 279f830108cSBarry Smith struct _MatOps *Aops; 280f830108cSBarry Smith 2812d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 2822d61bbb3SSatish Balay PetscFree(a->a); 2832d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 2842d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 2852d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 2862d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 2872d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 2882d61bbb3SSatish Balay PetscFree(a); 289f830108cSBarry Smith 290f830108cSBarry Smith /* 291f830108cSBarry Smith This is horrible, horrible code. We need to keep the 292f830108cSBarry Smith A pointers for the bops and ops but copy everything 293f830108cSBarry Smith else from C. 294f830108cSBarry Smith */ 295f830108cSBarry Smith Abops = A->bops; 296f830108cSBarry Smith Aops = A->ops; 2972d61bbb3SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 298f830108cSBarry Smith A->bops = Abops; 299f830108cSBarry Smith A->ops = Aops; 300f830108cSBarry Smith 3012d61bbb3SSatish Balay PetscHeaderDestroy(C); 3022d61bbb3SSatish Balay } 3032d61bbb3SSatish Balay PetscFunctionReturn(0); 3042d61bbb3SSatish Balay } 3052d61bbb3SSatish Balay 3062d61bbb3SSatish Balay 3072d61bbb3SSatish Balay 3083b2fbd54SBarry Smith 3095615d1e5SSatish Balay #undef __FUNC__ 310d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 311b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3122593348eSBarry Smith { 313b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3149df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 315b6490206SBarry Smith Scalar *aa; 3162593348eSBarry Smith 3173a40ed3dSBarry Smith PetscFunctionBegin; 31890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 3192593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3202593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3213b2fbd54SBarry Smith 3222593348eSBarry Smith col_lens[1] = a->m; 3232593348eSBarry Smith col_lens[2] = a->n; 3247e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3252593348eSBarry Smith 3262593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 327b6490206SBarry Smith count = 0; 328b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 329b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 330b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 331b6490206SBarry Smith } 3322593348eSBarry Smith } 3330752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3342593348eSBarry Smith PetscFree(col_lens); 3352593348eSBarry Smith 3362593348eSBarry Smith /* store column indices (zero start index) */ 33766963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj); 338b6490206SBarry Smith count = 0; 339b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 340b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 341b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 342b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 343b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3442593348eSBarry Smith } 3452593348eSBarry Smith } 346b6490206SBarry Smith } 347b6490206SBarry Smith } 3480752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr); 349b6490206SBarry Smith PetscFree(jj); 3502593348eSBarry Smith 3512593348eSBarry Smith /* store nonzero values */ 35266963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa); 353b6490206SBarry Smith count = 0; 354b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 355b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 356b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 357b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 3587e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 359b6490206SBarry Smith } 360b6490206SBarry Smith } 361b6490206SBarry Smith } 362b6490206SBarry Smith } 3630752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 364b6490206SBarry Smith PetscFree(aa); 3653a40ed3dSBarry Smith PetscFunctionReturn(0); 3662593348eSBarry Smith } 3672593348eSBarry Smith 3685615d1e5SSatish Balay #undef __FUNC__ 369d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 370b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 3712593348eSBarry Smith { 372b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3739df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 3742593348eSBarry Smith FILE *fd; 3752593348eSBarry Smith char *outputname; 3762593348eSBarry Smith 3773a40ed3dSBarry Smith PetscFunctionBegin; 37890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 3792593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 38090ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 381639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 3827ddc982cSLois Curfman McInnes fprintf(fd," block size is %d\n",bs); 3832593348eSBarry Smith } 384639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 385a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 3862593348eSBarry Smith } 387639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 38844cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 38944cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 39044cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 39144cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 39244cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 3933a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 394766eeae4SLois Curfman McInnes if (imag(a->a[bs2*k + l*bs + j]) > 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 39544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 39644cd7ae7SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 397766eeae4SLois Curfman McInnes else if (imag(a->a[bs2*k + l*bs + j]) < 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 398766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 399766eeae4SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j])); 40044cd7ae7SLois Curfman McInnes else if (real(a->a[bs2*k + l*bs + j]) != 0.0) 40144cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 40244cd7ae7SLois Curfman McInnes #else 40344cd7ae7SLois Curfman McInnes if (a->a[bs2*k + l*bs + j] != 0.0) 40444cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 40544cd7ae7SLois Curfman McInnes #endif 40644cd7ae7SLois Curfman McInnes } 40744cd7ae7SLois Curfman McInnes } 40844cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 40944cd7ae7SLois Curfman McInnes } 41044cd7ae7SLois Curfman McInnes } 41144cd7ae7SLois Curfman McInnes } 4122593348eSBarry Smith else { 413b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 414b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 415b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 416b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 417b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4183a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 419766eeae4SLois Curfman McInnes if (imag(a->a[bs2*k + l*bs + j]) > 0.0) { 42088685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 4217e67e3f9SSatish Balay real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 42288685aaeSLois Curfman McInnes } 423766eeae4SLois Curfman McInnes else if (imag(a->a[bs2*k + l*bs + j]) < 0.0) { 424766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 425766eeae4SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j])); 426766eeae4SLois Curfman McInnes } 42788685aaeSLois Curfman McInnes else { 4287e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 42988685aaeSLois Curfman McInnes } 43088685aaeSLois Curfman McInnes #else 4317e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 43288685aaeSLois Curfman McInnes #endif 4332593348eSBarry Smith } 4342593348eSBarry Smith } 4352593348eSBarry Smith fprintf(fd,"\n"); 4362593348eSBarry Smith } 4372593348eSBarry Smith } 438b6490206SBarry Smith } 4392593348eSBarry Smith fflush(fd); 4403a40ed3dSBarry Smith PetscFunctionReturn(0); 4412593348eSBarry Smith } 4422593348eSBarry Smith 4435615d1e5SSatish Balay #undef __FUNC__ 444d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 4453270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 4463270192aSSatish Balay { 4473270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4483270192aSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2; 4493270192aSSatish Balay double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 4503270192aSSatish Balay Scalar *aa; 4513270192aSSatish Balay Draw draw; 4523270192aSSatish Balay DrawButton button; 4533270192aSSatish Balay PetscTruth isnull; 4543270192aSSatish Balay 4553a40ed3dSBarry Smith PetscFunctionBegin; 4563a40ed3dSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr); 4573a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 4583270192aSSatish Balay 4593270192aSSatish Balay xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 4603270192aSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 4613270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 4623270192aSSatish Balay /* loop over matrix elements drawing boxes */ 4633270192aSSatish Balay color = DRAW_BLUE; 4643270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4653270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4663270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4673270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4683270192aSSatish Balay aa = a->a + j*bs2; 4693270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4703270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4713270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 472b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4733270192aSSatish Balay } 4743270192aSSatish Balay } 4753270192aSSatish Balay } 4763270192aSSatish Balay } 4773270192aSSatish Balay color = DRAW_CYAN; 4783270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4793270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4803270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4813270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4823270192aSSatish Balay aa = a->a + j*bs2; 4833270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4843270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4853270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 486b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4873270192aSSatish Balay } 4883270192aSSatish Balay } 4893270192aSSatish Balay } 4903270192aSSatish Balay } 4913270192aSSatish Balay 4923270192aSSatish Balay color = DRAW_RED; 4933270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4943270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4953270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4963270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4973270192aSSatish Balay aa = a->a + j*bs2; 4983270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4993270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5003270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 501b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5023270192aSSatish Balay } 5033270192aSSatish Balay } 5043270192aSSatish Balay } 5053270192aSSatish Balay } 5063270192aSSatish Balay 50755843e3eSBarry Smith DrawSynchronizedFlush(draw); 5083270192aSSatish Balay DrawGetPause(draw,&pause); 5093a40ed3dSBarry Smith if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 5103270192aSSatish Balay 5113270192aSSatish Balay /* allow the matrix to zoom or shrink */ 51255843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5133270192aSSatish Balay while (button != BUTTON_RIGHT) { 51455843e3eSBarry Smith DrawSynchronizedClear(draw); 5153270192aSSatish Balay if (button == BUTTON_LEFT) scale = .5; 5163270192aSSatish Balay else if (button == BUTTON_CENTER) scale = 2.; 5173270192aSSatish Balay xl = scale*(xl + w - xc) + xc - w*scale; 5183270192aSSatish Balay xr = scale*(xr - w - xc) + xc + w*scale; 5193270192aSSatish Balay yl = scale*(yl + h - yc) + yc - h*scale; 5203270192aSSatish Balay yr = scale*(yr - h - yc) + yc + h*scale; 5213270192aSSatish Balay w *= scale; h *= scale; 5223270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 5233270192aSSatish Balay color = DRAW_BLUE; 5243270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5253270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5263270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5273270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5283270192aSSatish Balay aa = a->a + j*bs2; 5293270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5303270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5313270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 532b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5333270192aSSatish Balay } 5343270192aSSatish Balay } 5353270192aSSatish Balay } 5363270192aSSatish Balay } 5373270192aSSatish Balay color = DRAW_CYAN; 5383270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5393270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5403270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5413270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5423270192aSSatish Balay aa = a->a + j*bs2; 5433270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5443270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5453270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 546b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5473270192aSSatish Balay } 5483270192aSSatish Balay } 5493270192aSSatish Balay } 5503270192aSSatish Balay } 5513270192aSSatish Balay 5523270192aSSatish Balay color = DRAW_RED; 5533270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5543270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5553270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5563270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5573270192aSSatish Balay aa = a->a + j*bs2; 5583270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5593270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5603270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 561b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5623270192aSSatish Balay } 5633270192aSSatish Balay } 5643270192aSSatish Balay } 5653270192aSSatish Balay } 56655843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5673270192aSSatish Balay } 5683a40ed3dSBarry Smith PetscFunctionReturn(0); 5693270192aSSatish Balay } 5703270192aSSatish Balay 5715615d1e5SSatish Balay #undef __FUNC__ 572d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 5738f6be9afSLois Curfman McInnes int MatView_SeqBAIJ(PetscObject obj,Viewer viewer) 5742593348eSBarry Smith { 5752593348eSBarry Smith Mat A = (Mat) obj; 57619bcc07fSBarry Smith ViewerType vtype; 57719bcc07fSBarry Smith int ierr; 5782593348eSBarry Smith 5793a40ed3dSBarry Smith PetscFunctionBegin; 58019bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 58119bcc07fSBarry Smith if (vtype == MATLAB_VIEWER) { 582a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab viewer not supported"); 5833a40ed3dSBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 5843a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 5853a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 5863a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 5873a40ed3dSBarry Smith } else if (vtype == DRAW_VIEWER) { 5883a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 5892593348eSBarry Smith } 5903a40ed3dSBarry Smith PetscFunctionReturn(0); 5912593348eSBarry Smith } 592b6490206SBarry Smith 593cd0e1443SSatish Balay 5945615d1e5SSatish Balay #undef __FUNC__ 5952d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 5962d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 597cd0e1443SSatish Balay { 598cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 5992d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6002d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6012d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6022d61bbb3SSatish Balay Scalar *ap, *aa = a->a, zero = 0.0; 603cd0e1443SSatish Balay 6043a40ed3dSBarry Smith PetscFunctionBegin; 6052d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 606cd0e1443SSatish Balay row = im[k]; brow = row/bs; 607a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 608a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6092d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6102c3acbe9SBarry Smith nrow = ailen[brow]; 6112d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 612a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 613a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6142d61bbb3SSatish Balay col = in[l] ; 6152d61bbb3SSatish Balay bcol = col/bs; 6162d61bbb3SSatish Balay cidx = col%bs; 6172d61bbb3SSatish Balay ridx = row%bs; 6182d61bbb3SSatish Balay high = nrow; 6192d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6202d61bbb3SSatish Balay while (high-low > 5) { 621cd0e1443SSatish Balay t = (low+high)/2; 622cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 623cd0e1443SSatish Balay else low = t; 624cd0e1443SSatish Balay } 625cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 626cd0e1443SSatish Balay if (rp[i] > bcol) break; 627cd0e1443SSatish Balay if (rp[i] == bcol) { 6282d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6292d61bbb3SSatish Balay goto finished; 630cd0e1443SSatish Balay } 631cd0e1443SSatish Balay } 6322d61bbb3SSatish Balay *v++ = zero; 6332d61bbb3SSatish Balay finished:; 634cd0e1443SSatish Balay } 635cd0e1443SSatish Balay } 6363a40ed3dSBarry Smith PetscFunctionReturn(0); 637cd0e1443SSatish Balay } 638cd0e1443SSatish Balay 6392d61bbb3SSatish Balay 6405615d1e5SSatish Balay #undef __FUNC__ 64105a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 64292c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 64392c4ed94SBarry Smith { 64492c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6458a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 646dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 647dd9472c6SBarry Smith int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval; 6480e324ae4SSatish Balay Scalar *ap,*value=v,*aa=a->a,*bap; 64992c4ed94SBarry Smith 6503a40ed3dSBarry Smith PetscFunctionBegin; 6510e324ae4SSatish Balay if (roworiented) { 6520e324ae4SSatish Balay stepval = (n-1)*bs; 6530e324ae4SSatish Balay } else { 6540e324ae4SSatish Balay stepval = (m-1)*bs; 6550e324ae4SSatish Balay } 65692c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 65792c4ed94SBarry Smith row = im[k]; 6583a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 659a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 660a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 66192c4ed94SBarry Smith #endif 66292c4ed94SBarry Smith rp = aj + ai[row]; 66392c4ed94SBarry Smith ap = aa + bs2*ai[row]; 66492c4ed94SBarry Smith rmax = imax[row]; 66592c4ed94SBarry Smith nrow = ailen[row]; 66692c4ed94SBarry Smith low = 0; 66792c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6683a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 669a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 670a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 67192c4ed94SBarry Smith #endif 67292c4ed94SBarry Smith col = in[l]; 67392c4ed94SBarry Smith if (roworiented) { 6740e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6750e324ae4SSatish Balay } else { 6760e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 67792c4ed94SBarry Smith } 67892c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 67992c4ed94SBarry Smith while (high-low > 7) { 68092c4ed94SBarry Smith t = (low+high)/2; 68192c4ed94SBarry Smith if (rp[t] > col) high = t; 68292c4ed94SBarry Smith else low = t; 68392c4ed94SBarry Smith } 68492c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 68592c4ed94SBarry Smith if (rp[i] > col) break; 68692c4ed94SBarry Smith if (rp[i] == col) { 6878a84c255SSatish Balay bap = ap + bs2*i; 6880e324ae4SSatish Balay if (roworiented) { 6898a84c255SSatish Balay if (is == ADD_VALUES) { 690dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 691dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 6928a84c255SSatish Balay bap[jj] += *value++; 693dd9472c6SBarry Smith } 694dd9472c6SBarry Smith } 6950e324ae4SSatish Balay } else { 696dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 697dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 6980e324ae4SSatish Balay bap[jj] = *value++; 6998a84c255SSatish Balay } 700dd9472c6SBarry Smith } 701dd9472c6SBarry Smith } 7020e324ae4SSatish Balay } else { 7030e324ae4SSatish Balay if (is == ADD_VALUES) { 704dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 705dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7060e324ae4SSatish Balay *bap++ += *value++; 707dd9472c6SBarry Smith } 708dd9472c6SBarry Smith } 7090e324ae4SSatish Balay } else { 710dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 711dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7120e324ae4SSatish Balay *bap++ = *value++; 7130e324ae4SSatish Balay } 7148a84c255SSatish Balay } 715dd9472c6SBarry Smith } 716dd9472c6SBarry Smith } 717f1241b54SBarry Smith goto noinsert2; 71892c4ed94SBarry Smith } 71992c4ed94SBarry Smith } 72089280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 721a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 72292c4ed94SBarry Smith if (nrow >= rmax) { 72392c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 72492c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 72592c4ed94SBarry Smith Scalar *new_a; 72692c4ed94SBarry Smith 727a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 72889280ab3SLois Curfman McInnes 72992c4ed94SBarry Smith /* malloc new storage space */ 73092c4ed94SBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 73192c4ed94SBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 73292c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 73392c4ed94SBarry Smith new_i = new_j + new_nz; 73492c4ed94SBarry Smith 73592c4ed94SBarry Smith /* copy over old data into new slots */ 73692c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 73792c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 73892c4ed94SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int)); 73992c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 740dd9472c6SBarry Smith PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int)); 74192c4ed94SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(Scalar)); 74292c4ed94SBarry Smith PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 743dd9472c6SBarry Smith PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(Scalar)); 74492c4ed94SBarry Smith /* free up old matrix storage */ 74592c4ed94SBarry Smith PetscFree(a->a); 74692c4ed94SBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 74792c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 74892c4ed94SBarry Smith a->singlemalloc = 1; 74992c4ed94SBarry Smith 75092c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 75192c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 75292c4ed94SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 75392c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 75492c4ed94SBarry Smith a->reallocs++; 75592c4ed94SBarry Smith a->nz++; 75692c4ed94SBarry Smith } 75792c4ed94SBarry Smith N = nrow++ - 1; 75892c4ed94SBarry Smith /* shift up all the later entries in this row */ 75992c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 76092c4ed94SBarry Smith rp[ii+1] = rp[ii]; 76192c4ed94SBarry Smith PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 76292c4ed94SBarry Smith } 76392c4ed94SBarry Smith if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 76492c4ed94SBarry Smith rp[i] = col; 7658a84c255SSatish Balay bap = ap + bs2*i; 7660e324ae4SSatish Balay if (roworiented) { 767dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 768dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7690e324ae4SSatish Balay bap[jj] = *value++; 770dd9472c6SBarry Smith } 771dd9472c6SBarry Smith } 7720e324ae4SSatish Balay } else { 773dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 774dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7750e324ae4SSatish Balay *bap++ = *value++; 7760e324ae4SSatish Balay } 777dd9472c6SBarry Smith } 778dd9472c6SBarry Smith } 779f1241b54SBarry Smith noinsert2:; 78092c4ed94SBarry Smith low = i; 78192c4ed94SBarry Smith } 78292c4ed94SBarry Smith ailen[row] = nrow; 78392c4ed94SBarry Smith } 7843a40ed3dSBarry Smith PetscFunctionReturn(0); 78592c4ed94SBarry Smith } 78692c4ed94SBarry Smith 787f2501298SSatish Balay 7885615d1e5SSatish Balay #undef __FUNC__ 7895615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 7908f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 791584200bdSSatish Balay { 792584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 793584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 794a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 79543ee02c3SBarry Smith int mbs = a->mbs, bs2 = a->bs2,rmax = 0; 796584200bdSSatish Balay Scalar *aa = a->a, *ap; 797584200bdSSatish Balay 7983a40ed3dSBarry Smith PetscFunctionBegin; 7993a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 800584200bdSSatish Balay 80143ee02c3SBarry Smith if (m) rmax = ailen[0]; 802584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 803584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 804584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 805d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 806584200bdSSatish Balay if (fshift) { 807a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 808584200bdSSatish Balay N = ailen[i]; 809584200bdSSatish Balay for ( j=0; j<N; j++ ) { 810584200bdSSatish Balay ip[j-fshift] = ip[j]; 8117e67e3f9SSatish Balay PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar)); 812584200bdSSatish Balay } 813584200bdSSatish Balay } 814584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 815584200bdSSatish Balay } 816584200bdSSatish Balay if (mbs) { 817584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 818584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 819584200bdSSatish Balay } 820584200bdSSatish Balay /* reset ilen and imax for each row */ 821584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 822584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 823584200bdSSatish Balay } 824a7c10996SSatish Balay a->nz = ai[mbs]; 825584200bdSSatish Balay 826584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 827584200bdSSatish Balay if (fshift && a->diag) { 828584200bdSSatish Balay PetscFree(a->diag); 829584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 830584200bdSSatish Balay a->diag = 0; 831584200bdSSatish Balay } 8323d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 833219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8343d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 835584200bdSSatish Balay a->reallocs); 836d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 837e2f3b5e9SSatish Balay a->reallocs = 0; 8384e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8394e220ebcSLois Curfman McInnes 8403a40ed3dSBarry Smith PetscFunctionReturn(0); 841584200bdSSatish Balay } 842584200bdSSatish Balay 8435a838052SSatish Balay 844d9b7c43dSSatish Balay /* idx should be of length atlease bs */ 8455615d1e5SSatish Balay #undef __FUNC__ 8465615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Block" 847d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg) 848d9b7c43dSSatish Balay { 849d9b7c43dSSatish Balay int i,row; 8503a40ed3dSBarry Smith 8513a40ed3dSBarry Smith PetscFunctionBegin; 852d9b7c43dSSatish Balay row = idx[0]; 8533a40ed3dSBarry Smith if (row%bs!=0) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 854d9b7c43dSSatish Balay 855d9b7c43dSSatish Balay for ( i=1; i<bs; i++ ) { 8563a40ed3dSBarry Smith if (row+i != idx[i]) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } 857d9b7c43dSSatish Balay } 858d9b7c43dSSatish Balay *flg = PETSC_TRUE; 8593a40ed3dSBarry Smith PetscFunctionReturn(0); 860d9b7c43dSSatish Balay } 861d9b7c43dSSatish Balay 8625615d1e5SSatish Balay #undef __FUNC__ 8635615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 8648f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 865d9b7c43dSSatish Balay { 866d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 867d9b7c43dSSatish Balay IS is_local; 868d9b7c43dSSatish Balay int ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2; 869d9b7c43dSSatish Balay PetscTruth flg; 870d9b7c43dSSatish Balay Scalar zero = 0.0,*aa; 871d9b7c43dSSatish Balay 8723a40ed3dSBarry Smith PetscFunctionBegin; 873d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 874d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 875d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 876537820f0SBarry Smith ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr); 877d9b7c43dSSatish Balay ierr = ISSort(is_local); CHKERRQ(ierr); 878d9b7c43dSSatish Balay ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr); 879d9b7c43dSSatish Balay 880d9b7c43dSSatish Balay i = 0; 881d9b7c43dSSatish Balay while (i < is_n) { 882a8c6a408SBarry Smith if (rows[i]<0 || rows[i]>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 883d9b7c43dSSatish Balay flg = PETSC_FALSE; 884d9b7c43dSSatish Balay if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); } 885d9b7c43dSSatish Balay count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs; 886d9b7c43dSSatish Balay aa = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs); 8872d61bbb3SSatish Balay if (flg) { /* There exists a block of rows to be Zerowed */ 888a07cd24cSSatish Balay if (baij->ilen[rows[i]/bs] > 0) { 8892d61bbb3SSatish Balay PetscMemzero(aa,count*bs*sizeof(Scalar)); 8902d61bbb3SSatish Balay baij->ilen[rows[i]/bs] = 1; 8912d61bbb3SSatish Balay baij->j[baij->i[rows[i]/bs]] = rows[i]/bs; 892a07cd24cSSatish Balay } 8932d61bbb3SSatish Balay i += bs; 8942d61bbb3SSatish Balay } else { /* Zero out only the requested row */ 895d9b7c43dSSatish Balay for ( j=0; j<count; j++ ) { 896d9b7c43dSSatish Balay aa[0] = zero; 897d9b7c43dSSatish Balay aa+=bs; 898d9b7c43dSSatish Balay } 899d9b7c43dSSatish Balay i++; 900d9b7c43dSSatish Balay } 901d9b7c43dSSatish Balay } 902d9b7c43dSSatish Balay if (diag) { 903d9b7c43dSSatish Balay for ( j=0; j<is_n; j++ ) { 904f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 9052d61bbb3SSatish Balay /* ierr = MatSetValues(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); */ 906d9b7c43dSSatish Balay } 907d9b7c43dSSatish Balay } 908d9b7c43dSSatish Balay ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr); 909d9b7c43dSSatish Balay ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr); 910d9b7c43dSSatish Balay ierr = ISDestroy(is_local); CHKERRQ(ierr); 9119a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 9123a40ed3dSBarry Smith PetscFunctionReturn(0); 913d9b7c43dSSatish Balay } 9141c351548SSatish Balay 9155615d1e5SSatish Balay #undef __FUNC__ 9162d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9172d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9182d61bbb3SSatish Balay { 9192d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9202d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9212d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9222d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 9232d61bbb3SSatish Balay int ridx,cidx,bs2=a->bs2; 9242d61bbb3SSatish Balay Scalar *ap,value,*aa=a->a,*bap; 9252d61bbb3SSatish Balay 9262d61bbb3SSatish Balay PetscFunctionBegin; 9272d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9282d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9292d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9302d61bbb3SSatish Balay if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 9312d61bbb3SSatish Balay if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 9322d61bbb3SSatish Balay #endif 9332d61bbb3SSatish Balay rp = aj + ai[brow]; 9342d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9352d61bbb3SSatish Balay rmax = imax[brow]; 9362d61bbb3SSatish Balay nrow = ailen[brow]; 9372d61bbb3SSatish Balay low = 0; 9382d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9392d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 9402d61bbb3SSatish Balay if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 9412d61bbb3SSatish Balay if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 9422d61bbb3SSatish Balay #endif 9432d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 9442d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 9452d61bbb3SSatish Balay if (roworiented) { 9462d61bbb3SSatish Balay value = *v++; 9472d61bbb3SSatish Balay } else { 9482d61bbb3SSatish Balay value = v[k + l*m]; 9492d61bbb3SSatish Balay } 9502d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 9512d61bbb3SSatish Balay while (high-low > 7) { 9522d61bbb3SSatish Balay t = (low+high)/2; 9532d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 9542d61bbb3SSatish Balay else low = t; 9552d61bbb3SSatish Balay } 9562d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 9572d61bbb3SSatish Balay if (rp[i] > bcol) break; 9582d61bbb3SSatish Balay if (rp[i] == bcol) { 9592d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 9602d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 9612d61bbb3SSatish Balay else *bap = value; 9622d61bbb3SSatish Balay goto noinsert1; 9632d61bbb3SSatish Balay } 9642d61bbb3SSatish Balay } 9652d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 9662d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9672d61bbb3SSatish Balay if (nrow >= rmax) { 9682d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 9692d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 9702d61bbb3SSatish Balay Scalar *new_a; 9712d61bbb3SSatish Balay 9722d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 9732d61bbb3SSatish Balay 9742d61bbb3SSatish Balay /* Malloc new storage space */ 9752d61bbb3SSatish Balay len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 9762d61bbb3SSatish Balay new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 9772d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 9782d61bbb3SSatish Balay new_i = new_j + new_nz; 9792d61bbb3SSatish Balay 9802d61bbb3SSatish Balay /* copy over old data into new slots */ 9812d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 9822d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 9832d61bbb3SSatish Balay PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int)); 9842d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 9852d61bbb3SSatish Balay PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow, 9862d61bbb3SSatish Balay len*sizeof(int)); 9872d61bbb3SSatish Balay PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar)); 9882d61bbb3SSatish Balay PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 9892d61bbb3SSatish Balay PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE), 9902d61bbb3SSatish Balay aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar)); 9912d61bbb3SSatish Balay /* free up old matrix storage */ 9922d61bbb3SSatish Balay PetscFree(a->a); 9932d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 9942d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 9952d61bbb3SSatish Balay a->singlemalloc = 1; 9962d61bbb3SSatish Balay 9972d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 9982d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 9992d61bbb3SSatish Balay PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 10002d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10012d61bbb3SSatish Balay a->reallocs++; 10022d61bbb3SSatish Balay a->nz++; 10032d61bbb3SSatish Balay } 10042d61bbb3SSatish Balay N = nrow++ - 1; 10052d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10062d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10072d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 10082d61bbb3SSatish Balay PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 10092d61bbb3SSatish Balay } 10102d61bbb3SSatish Balay if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 10112d61bbb3SSatish Balay rp[i] = bcol; 10122d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10132d61bbb3SSatish Balay noinsert1:; 10142d61bbb3SSatish Balay low = i; 10152d61bbb3SSatish Balay } 10162d61bbb3SSatish Balay ailen[brow] = nrow; 10172d61bbb3SSatish Balay } 10182d61bbb3SSatish Balay PetscFunctionReturn(0); 10192d61bbb3SSatish Balay } 10202d61bbb3SSatish Balay 10212d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10222d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10232d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1024d3825aa8SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatGetSubMatrixCall,Mat*); 10252d61bbb3SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**); 10262d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10272d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10282d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10292d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10302d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10312d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10322d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10332d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10342d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10352d61bbb3SSatish Balay 10362d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10372d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 10382d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 10392d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 10402d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 10412d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 10422d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 10432d61bbb3SSatish Balay 10442d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 10452d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 10462d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 10472d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 10482d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 10492d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 10502d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4_NaturalOrdering(Mat,Vec,Vec); 10512d61bbb3SSatish Balay 10522d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 10532d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 10542d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 10552d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 10562d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 10572d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 10582d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 10592d61bbb3SSatish Balay 10602d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 10612d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 10622d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 10632d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 10642d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 10652d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 10662d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 10672d61bbb3SSatish Balay 10682d61bbb3SSatish Balay /* 10692d61bbb3SSatish Balay note: This can only work for identity for row and col. It would 10702d61bbb3SSatish Balay be good to check this and otherwise generate an error. 10712d61bbb3SSatish Balay */ 10722d61bbb3SSatish Balay #undef __FUNC__ 10732d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 10742d61bbb3SSatish Balay int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 10752d61bbb3SSatish Balay { 10762d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 10772d61bbb3SSatish Balay Mat outA; 10782d61bbb3SSatish Balay int ierr; 10792d61bbb3SSatish Balay 10802d61bbb3SSatish Balay PetscFunctionBegin; 10812d61bbb3SSatish Balay if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 10822d61bbb3SSatish Balay 10832d61bbb3SSatish Balay outA = inA; 10842d61bbb3SSatish Balay inA->factor = FACTOR_LU; 10852d61bbb3SSatish Balay a->row = row; 10862d61bbb3SSatish Balay a->col = col; 10872d61bbb3SSatish Balay 10882d61bbb3SSatish Balay if (!a->solve_work) { 10892d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 10902d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 10912d61bbb3SSatish Balay } 10922d61bbb3SSatish Balay 10932d61bbb3SSatish Balay if (!a->diag) { 10942d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 10952d61bbb3SSatish Balay } 10962d61bbb3SSatish Balay ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 10972d61bbb3SSatish Balay 10982d61bbb3SSatish Balay /* 10992d61bbb3SSatish Balay Blocksize 4 has a special faster solver for ILU(0) factorization 11002d61bbb3SSatish Balay with natural ordering 11012d61bbb3SSatish Balay */ 11022d61bbb3SSatish Balay if (a->bs == 4) { 1103f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 11042d61bbb3SSatish Balay } 11052d61bbb3SSatish Balay 11062d61bbb3SSatish Balay PetscFunctionReturn(0); 11072d61bbb3SSatish Balay } 11082d61bbb3SSatish Balay #undef __FUNC__ 1109d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1110ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1111ba4ca20aSSatish Balay { 1112ba4ca20aSSatish Balay static int called = 0; 1113ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1114ba4ca20aSSatish Balay 11153a40ed3dSBarry Smith PetscFunctionBegin; 11163a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 111776be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n"); 111876be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n"); 11193a40ed3dSBarry Smith PetscFunctionReturn(0); 1120ba4ca20aSSatish Balay } 1121d9b7c43dSSatish Balay 11222593348eSBarry Smith /* -------------------------------------------------------------------*/ 1123cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ, 11249867e207SSatish Balay MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ, 1125f44d3a6dSSatish Balay MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N, 1126faf6580fSSatish Balay MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ, 11271a6a6d98SBarry Smith MatSolve_SeqBAIJ_N,0, 1128b6490206SBarry Smith 0,0, 1129de6a44a3SBarry Smith MatLUFactor_SeqBAIJ,0, 1130b6490206SBarry Smith 0, 1131f2501298SSatish Balay MatTranspose_SeqBAIJ, 113217e48fc4SSatish Balay MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ, 11339867e207SSatish Balay MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ, 1134584200bdSSatish Balay 0,MatAssemblyEnd_SeqBAIJ, 1135b6490206SBarry Smith 0, 1136d9b7c43dSSatish Balay MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ, 11377fc0212eSBarry Smith MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 1138b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 1139de6a44a3SBarry Smith MatILUFactorSymbolic_SeqBAIJ,0, 1140d402145bSBarry Smith 0,0, 1141b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 1142b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 1143af015674SSatish Balay MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ, 11447e67e3f9SSatish Balay MatGetValues_SeqBAIJ,0, 1145ba4ca20aSSatish Balay MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 11463b2fbd54SBarry Smith 0,0,0,MatGetBlockSize_SeqBAIJ, 11473b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 114892c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 114992c4ed94SBarry Smith 0,0,0,0,0,0, 1150d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1151d3825aa8SBarry Smith MatGetSubMatrix_SeqBAIJ}; 11522593348eSBarry Smith 11535615d1e5SSatish Balay #undef __FUNC__ 11545615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 11552593348eSBarry Smith /*@C 115644cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 115744cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 115844cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 11592bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 11602bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 11612593348eSBarry Smith 11622593348eSBarry Smith Input Parameters: 1163029af93fSBarry Smith . comm - MPI communicator, set to PETSC_COMM_SELF 1164b6490206SBarry Smith . bs - size of block 11652593348eSBarry Smith . m - number of rows 11662593348eSBarry Smith . n - number of columns 1167b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 11682bd5e0b2SLois Curfman McInnes . nzz - array containing the number of block nonzeros in the various block rows 11692bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 11702593348eSBarry Smith 11712593348eSBarry Smith Output Parameter: 11722593348eSBarry Smith . A - the matrix 11732593348eSBarry Smith 11740513a670SBarry Smith Options Database Keys: 11750513a670SBarry Smith $ -mat_no_unroll - uses code that does not unroll the loops in the 11760effe34fSLois Curfman McInnes $ block calculations (much slower) 11770513a670SBarry Smith $ -mat_block_size - size of the blocks to use 11780513a670SBarry Smith 11792593348eSBarry Smith Notes: 118044cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 11812593348eSBarry Smith storage. That is, the stored row and column indices can begin at 118244cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 11832593348eSBarry Smith 11842593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 11852593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 11862593348eSBarry Smith allocation. For additional details, see the users manual chapter on 11876da5968aSLois Curfman McInnes matrices. 11882593348eSBarry Smith 118944cd7ae7SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 11902593348eSBarry Smith @*/ 1191b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 11922593348eSBarry Smith { 11932593348eSBarry Smith Mat B; 1194b6490206SBarry Smith Mat_SeqBAIJ *b; 11953b2fbd54SBarry Smith int i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size; 11963b2fbd54SBarry Smith 11973a40ed3dSBarry Smith PetscFunctionBegin; 11983b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1199a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1200b6490206SBarry Smith 12010513a670SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 12020513a670SBarry Smith 12033a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1204a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 12053a40ed3dSBarry Smith } 12062593348eSBarry Smith 12072593348eSBarry Smith *A = 0; 1208f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView); 12092593348eSBarry Smith PLogObjectCreate(B); 1210b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 121144cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqBAIJ)); 1212f830108cSBarry Smith PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps)); 12131a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr); 12141a6a6d98SBarry Smith if (!flg) { 12157fc0212eSBarry Smith switch (bs) { 12167fc0212eSBarry Smith case 1: 1217f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1218f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1219f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1220f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 12217fc0212eSBarry Smith break; 12224eeb42bcSBarry Smith case 2: 1223f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1224f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1225f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1226f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 12276d84be18SBarry Smith break; 1228f327dd97SBarry Smith case 3: 1229f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1230f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1231f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1232f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 12334eeb42bcSBarry Smith break; 12346d84be18SBarry Smith case 4: 1235f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1236f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1237f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1238f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 12396d84be18SBarry Smith break; 12406d84be18SBarry Smith case 5: 1241f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1242f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1243f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1244f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 12456d84be18SBarry Smith break; 124648e9ddb2SSatish Balay case 7: 1247f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1248f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1249f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 125048e9ddb2SSatish Balay break; 12517fc0212eSBarry Smith } 12521a6a6d98SBarry Smith } 1253b6490206SBarry Smith B->destroy = MatDestroy_SeqBAIJ; 1254b6490206SBarry Smith B->view = MatView_SeqBAIJ; 12552593348eSBarry Smith B->factor = 0; 12562593348eSBarry Smith B->lupivotthreshold = 1.0; 125790f02eecSBarry Smith B->mapping = 0; 12582593348eSBarry Smith b->row = 0; 12592593348eSBarry Smith b->col = 0; 12602593348eSBarry Smith b->reallocs = 0; 12612593348eSBarry Smith 126244cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 126344cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1264b6490206SBarry Smith b->mbs = mbs; 1265f2501298SSatish Balay b->nbs = nbs; 1266b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 12672593348eSBarry Smith if (nnz == PETSC_NULL) { 1268119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 12692593348eSBarry Smith else if (nz <= 0) nz = 1; 1270b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1271b6490206SBarry Smith nz = nz*mbs; 12723a40ed3dSBarry Smith } else { 12732593348eSBarry Smith nz = 0; 1274b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 12752593348eSBarry Smith } 12762593348eSBarry Smith 12772593348eSBarry Smith /* allocate the matrix space */ 12787e67e3f9SSatish Balay len = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int); 12792593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 12807e67e3f9SSatish Balay PetscMemzero(b->a,nz*bs2*sizeof(Scalar)); 12817e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 12822593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 12832593348eSBarry Smith b->i = b->j + nz; 12842593348eSBarry Smith b->singlemalloc = 1; 12852593348eSBarry Smith 1286de6a44a3SBarry Smith b->i[0] = 0; 1287b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 12882593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 12892593348eSBarry Smith } 12902593348eSBarry Smith 1291b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1292b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1293f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1294b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 12952593348eSBarry Smith 1296b6490206SBarry Smith b->bs = bs; 12979df24120SSatish Balay b->bs2 = bs2; 1298b6490206SBarry Smith b->mbs = mbs; 12992593348eSBarry Smith b->nz = 0; 130018e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 13012593348eSBarry Smith b->sorted = 0; 13022593348eSBarry Smith b->roworiented = 1; 13032593348eSBarry Smith b->nonew = 0; 13042593348eSBarry Smith b->diag = 0; 13052593348eSBarry Smith b->solve_work = 0; 1306de6a44a3SBarry Smith b->mult_work = 0; 13072593348eSBarry Smith b->spptr = 0; 13084e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 13094e220ebcSLois Curfman McInnes 13102593348eSBarry Smith *A = B; 13112593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 13122593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 13133a40ed3dSBarry Smith PetscFunctionReturn(0); 13142593348eSBarry Smith } 13152593348eSBarry Smith 13165615d1e5SSatish Balay #undef __FUNC__ 13175615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ" 1318b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 13192593348eSBarry Smith { 13202593348eSBarry Smith Mat C; 1321b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 13229df24120SSatish Balay int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 1323de6a44a3SBarry Smith 13243a40ed3dSBarry Smith PetscFunctionBegin; 1325a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 13262593348eSBarry Smith 13272593348eSBarry Smith *B = 0; 1328f830108cSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView); 13292593348eSBarry Smith PLogObjectCreate(C); 1330b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 1331*c3c66ee5SSatish Balay PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 1332b6490206SBarry Smith C->destroy = MatDestroy_SeqBAIJ; 1333b6490206SBarry Smith C->view = MatView_SeqBAIJ; 13342593348eSBarry Smith C->factor = A->factor; 13352593348eSBarry Smith c->row = 0; 13362593348eSBarry Smith c->col = 0; 13372593348eSBarry Smith C->assembled = PETSC_TRUE; 13382593348eSBarry Smith 133944cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 134044cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 134144cd7ae7SLois Curfman McInnes C->M = a->m; 134244cd7ae7SLois Curfman McInnes C->N = a->n; 134344cd7ae7SLois Curfman McInnes 1344b6490206SBarry Smith c->bs = a->bs; 13459df24120SSatish Balay c->bs2 = a->bs2; 1346b6490206SBarry Smith c->mbs = a->mbs; 1347b6490206SBarry Smith c->nbs = a->nbs; 13482593348eSBarry Smith 1349b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1350b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1351b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 13522593348eSBarry Smith c->imax[i] = a->imax[i]; 13532593348eSBarry Smith c->ilen[i] = a->ilen[i]; 13542593348eSBarry Smith } 13552593348eSBarry Smith 13562593348eSBarry Smith /* allocate the matrix space */ 13572593348eSBarry Smith c->singlemalloc = 1; 13587e67e3f9SSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int)); 13592593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 13607e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1361de6a44a3SBarry Smith c->i = c->j + nz; 1362b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1363b6490206SBarry Smith if (mbs > 0) { 1364de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 13652593348eSBarry Smith if (cpvalues == COPY_VALUES) { 13667e67e3f9SSatish Balay PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar)); 13672593348eSBarry Smith } 13682593348eSBarry Smith } 13692593348eSBarry Smith 1370f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 13712593348eSBarry Smith c->sorted = a->sorted; 13722593348eSBarry Smith c->roworiented = a->roworiented; 13732593348eSBarry Smith c->nonew = a->nonew; 13742593348eSBarry Smith 13752593348eSBarry Smith if (a->diag) { 1376b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1377b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1378b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 13792593348eSBarry Smith c->diag[i] = a->diag[i]; 13802593348eSBarry Smith } 13812593348eSBarry Smith } 13822593348eSBarry Smith else c->diag = 0; 13832593348eSBarry Smith c->nz = a->nz; 13842593348eSBarry Smith c->maxnz = a->maxnz; 13852593348eSBarry Smith c->solve_work = 0; 13862593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 13877fc0212eSBarry Smith c->mult_work = 0; 13882593348eSBarry Smith *B = C; 13893a40ed3dSBarry Smith PetscFunctionReturn(0); 13902593348eSBarry Smith } 13912593348eSBarry Smith 13925615d1e5SSatish Balay #undef __FUNC__ 13935615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 139419bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 13952593348eSBarry Smith { 1396b6490206SBarry Smith Mat_SeqBAIJ *a; 13972593348eSBarry Smith Mat B; 1398de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1399b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 140035aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1401a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1402b6490206SBarry Smith Scalar *aa; 140319bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 14042593348eSBarry Smith 14053a40ed3dSBarry Smith PetscFunctionBegin; 14061a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1407de6a44a3SBarry Smith bs2 = bs*bs; 1408b6490206SBarry Smith 14092593348eSBarry Smith MPI_Comm_size(comm,&size); 1410cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 141190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 14120752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 1413a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 14142593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 14152593348eSBarry Smith 1416d64ed03dSBarry Smith if (header[3] < 0) { 1417a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1418d64ed03dSBarry Smith } 1419d64ed03dSBarry Smith 1420a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 142135aab85fSBarry Smith 142235aab85fSBarry Smith /* 142335aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 142435aab85fSBarry Smith divisible by the blocksize 142535aab85fSBarry Smith */ 1426b6490206SBarry Smith mbs = M/bs; 142735aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 142835aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 142935aab85fSBarry Smith else mbs++; 143035aab85fSBarry Smith if (extra_rows) { 1431537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 143235aab85fSBarry Smith } 1433b6490206SBarry Smith 14342593348eSBarry Smith /* read in row lengths */ 143535aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 14360752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 143735aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 14382593348eSBarry Smith 1439b6490206SBarry Smith /* read in column indices */ 144035aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 14410752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr); 144235aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1443b6490206SBarry Smith 1444b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1445b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1446b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 144735aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 144835aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 144935aab85fSBarry Smith masked = mask + mbs; 1450b6490206SBarry Smith rowcount = 0; nzcount = 0; 1451b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 145235aab85fSBarry Smith nmask = 0; 1453b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1454b6490206SBarry Smith kmax = rowlengths[rowcount]; 1455b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 145635aab85fSBarry Smith tmp = jj[nzcount++]/bs; 145735aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1458b6490206SBarry Smith } 1459b6490206SBarry Smith rowcount++; 1460b6490206SBarry Smith } 146135aab85fSBarry Smith browlengths[i] += nmask; 146235aab85fSBarry Smith /* zero out the mask elements we set */ 146335aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1464b6490206SBarry Smith } 1465b6490206SBarry Smith 14662593348eSBarry Smith /* create our matrix */ 14673eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 14682593348eSBarry Smith B = *A; 1469b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 14702593348eSBarry Smith 14712593348eSBarry Smith /* set matrix "i" values */ 1472de6a44a3SBarry Smith a->i[0] = 0; 1473b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1474b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1475b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 14762593348eSBarry Smith } 1477b6490206SBarry Smith a->nz = 0; 1478b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 14792593348eSBarry Smith 1480b6490206SBarry Smith /* read in nonzero values */ 148135aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 14820752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr); 148335aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1484b6490206SBarry Smith 1485b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1486b6490206SBarry Smith nzcount = 0; jcount = 0; 1487b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1488b6490206SBarry Smith nzcountb = nzcount; 148935aab85fSBarry Smith nmask = 0; 1490b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1491b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1492b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 149335aab85fSBarry Smith tmp = jj[nzcount++]/bs; 149435aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1495b6490206SBarry Smith } 1496b6490206SBarry Smith rowcount++; 1497b6490206SBarry Smith } 1498de6a44a3SBarry Smith /* sort the masked values */ 149977c4ece6SBarry Smith PetscSortInt(nmask,masked); 1500de6a44a3SBarry Smith 1501b6490206SBarry Smith /* set "j" values into matrix */ 1502b6490206SBarry Smith maskcount = 1; 150335aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 150435aab85fSBarry Smith a->j[jcount++] = masked[j]; 1505de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1506b6490206SBarry Smith } 1507b6490206SBarry Smith /* set "a" values into matrix */ 1508de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1509b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1510b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1511b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1512de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1513de6a44a3SBarry Smith block = mask[tmp] - 1; 1514de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1515de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1516b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1517b6490206SBarry Smith } 1518b6490206SBarry Smith } 151935aab85fSBarry Smith /* zero out the mask elements we set */ 152035aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1521b6490206SBarry Smith } 1522a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1523b6490206SBarry Smith 1524b6490206SBarry Smith PetscFree(rowlengths); 1525b6490206SBarry Smith PetscFree(browlengths); 1526b6490206SBarry Smith PetscFree(aa); 1527b6490206SBarry Smith PetscFree(jj); 1528b6490206SBarry Smith PetscFree(mask); 1529b6490206SBarry Smith 1530b6490206SBarry Smith B->assembled = PETSC_TRUE; 1531de6a44a3SBarry Smith 15329c01be13SBarry Smith ierr = MatView_Private(B); CHKERRQ(ierr); 15333a40ed3dSBarry Smith PetscFunctionReturn(0); 15342593348eSBarry Smith } 15352593348eSBarry Smith 15362593348eSBarry Smith 15372593348eSBarry Smith 1538