1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*e1311b90SBarry Smith static char vcid[] = "$Id: baij.c,v 1.129 1998/03/26 22:32:11 balay Exp bsmith $"; 32593348eSBarry Smith #endif 42593348eSBarry Smith 52593348eSBarry Smith /* 6b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 72593348eSBarry Smith matrix storage format. 82593348eSBarry Smith */ 93369ce9aSBarry Smith 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" 105*e1311b90SBarry 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) 111*e1311b90SBarry 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; 301f830108cSBarry Smith 3022d61bbb3SSatish Balay PetscHeaderDestroy(C); 3032d61bbb3SSatish Balay } 3042d61bbb3SSatish Balay PetscFunctionReturn(0); 3052d61bbb3SSatish Balay } 3062d61bbb3SSatish Balay 3072d61bbb3SSatish Balay 3082d61bbb3SSatish Balay 3093b2fbd54SBarry Smith 3105615d1e5SSatish Balay #undef __FUNC__ 311d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 312b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3132593348eSBarry Smith { 314b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3159df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 316b6490206SBarry Smith Scalar *aa; 3172593348eSBarry Smith 3183a40ed3dSBarry Smith PetscFunctionBegin; 31990ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 3202593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3212593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3223b2fbd54SBarry Smith 3232593348eSBarry Smith col_lens[1] = a->m; 3242593348eSBarry Smith col_lens[2] = a->n; 3257e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3262593348eSBarry Smith 3272593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 328b6490206SBarry Smith count = 0; 329b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 330b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 331b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 332b6490206SBarry Smith } 3332593348eSBarry Smith } 3340752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3352593348eSBarry Smith PetscFree(col_lens); 3362593348eSBarry Smith 3372593348eSBarry Smith /* store column indices (zero start index) */ 33866963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj); 339b6490206SBarry Smith count = 0; 340b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 341b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 342b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 343b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 344b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3452593348eSBarry Smith } 3462593348eSBarry Smith } 347b6490206SBarry Smith } 348b6490206SBarry Smith } 3490752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr); 350b6490206SBarry Smith PetscFree(jj); 3512593348eSBarry Smith 3522593348eSBarry Smith /* store nonzero values */ 35366963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa); 354b6490206SBarry Smith count = 0; 355b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 356b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 357b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 358b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 3597e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 360b6490206SBarry Smith } 361b6490206SBarry Smith } 362b6490206SBarry Smith } 363b6490206SBarry Smith } 3640752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 365b6490206SBarry Smith PetscFree(aa); 3663a40ed3dSBarry Smith PetscFunctionReturn(0); 3672593348eSBarry Smith } 3682593348eSBarry Smith 3695615d1e5SSatish Balay #undef __FUNC__ 370d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 371b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 3722593348eSBarry Smith { 373b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3749df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 3752593348eSBarry Smith FILE *fd; 3762593348eSBarry Smith char *outputname; 3772593348eSBarry Smith 3783a40ed3dSBarry Smith PetscFunctionBegin; 37990ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 3802593348eSBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 38190ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 382639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 3837ddc982cSLois Curfman McInnes fprintf(fd," block size is %d\n",bs); 3842593348eSBarry Smith } 385639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 386a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 3872593348eSBarry Smith } 388639f9d9dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 38944cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 39044cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 39144cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 39244cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 39344cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 3943a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 395766eeae4SLois Curfman McInnes if (imag(a->a[bs2*k + l*bs + j]) > 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 39644cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 39744cd7ae7SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 398766eeae4SLois Curfman McInnes else if (imag(a->a[bs2*k + l*bs + j]) < 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 399766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 400766eeae4SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j])); 40144cd7ae7SLois Curfman McInnes else if (real(a->a[bs2*k + l*bs + j]) != 0.0) 40244cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 40344cd7ae7SLois Curfman McInnes #else 40444cd7ae7SLois Curfman McInnes if (a->a[bs2*k + l*bs + j] != 0.0) 40544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 40644cd7ae7SLois Curfman McInnes #endif 40744cd7ae7SLois Curfman McInnes } 40844cd7ae7SLois Curfman McInnes } 40944cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 41044cd7ae7SLois Curfman McInnes } 41144cd7ae7SLois Curfman McInnes } 41244cd7ae7SLois Curfman McInnes } 4132593348eSBarry Smith else { 414b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 415b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 416b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 417b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 418b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4193a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 420766eeae4SLois Curfman McInnes if (imag(a->a[bs2*k + l*bs + j]) > 0.0) { 42188685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 4227e67e3f9SSatish Balay real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 42388685aaeSLois Curfman McInnes } 424766eeae4SLois Curfman McInnes else if (imag(a->a[bs2*k + l*bs + j]) < 0.0) { 425766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 426766eeae4SLois Curfman McInnes real(a->a[bs2*k + l*bs + j]),-imag(a->a[bs2*k + l*bs + j])); 427766eeae4SLois Curfman McInnes } 42888685aaeSLois Curfman McInnes else { 4297e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 43088685aaeSLois Curfman McInnes } 43188685aaeSLois Curfman McInnes #else 4327e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 43388685aaeSLois Curfman McInnes #endif 4342593348eSBarry Smith } 4352593348eSBarry Smith } 4362593348eSBarry Smith fprintf(fd,"\n"); 4372593348eSBarry Smith } 4382593348eSBarry Smith } 439b6490206SBarry Smith } 4402593348eSBarry Smith fflush(fd); 4413a40ed3dSBarry Smith PetscFunctionReturn(0); 4422593348eSBarry Smith } 4432593348eSBarry Smith 4445615d1e5SSatish Balay #undef __FUNC__ 445d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 4463270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 4473270192aSSatish Balay { 4483270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4493270192aSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2; 4503270192aSSatish Balay double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 4513270192aSSatish Balay Scalar *aa; 4523270192aSSatish Balay Draw draw; 4533270192aSSatish Balay DrawButton button; 4543270192aSSatish Balay PetscTruth isnull; 4553270192aSSatish Balay 4563a40ed3dSBarry Smith PetscFunctionBegin; 4573a40ed3dSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw);CHKERRQ(ierr); 4583a40ed3dSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 4593270192aSSatish Balay 4603270192aSSatish Balay xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 4613270192aSSatish Balay xr += w; yr += h; xl = -w; yl = -h; 4623270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 4633270192aSSatish Balay /* loop over matrix elements drawing boxes */ 4643270192aSSatish Balay color = DRAW_BLUE; 4653270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4663270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4673270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4683270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4693270192aSSatish Balay aa = a->a + j*bs2; 4703270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4713270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4723270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 473b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4743270192aSSatish Balay } 4753270192aSSatish Balay } 4763270192aSSatish Balay } 4773270192aSSatish Balay } 4783270192aSSatish Balay color = DRAW_CYAN; 4793270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4803270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4813270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4823270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4833270192aSSatish Balay aa = a->a + j*bs2; 4843270192aSSatish Balay for ( k=0; k<bs; k++ ) { 4853270192aSSatish Balay for ( l=0; l<bs; l++ ) { 4863270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 487b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 4883270192aSSatish Balay } 4893270192aSSatish Balay } 4903270192aSSatish Balay } 4913270192aSSatish Balay } 4923270192aSSatish Balay 4933270192aSSatish Balay color = DRAW_RED; 4943270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 4953270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 4963270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 4973270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4983270192aSSatish Balay aa = a->a + j*bs2; 4993270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5003270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5013270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 502b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5033270192aSSatish Balay } 5043270192aSSatish Balay } 5053270192aSSatish Balay } 5063270192aSSatish Balay } 5073270192aSSatish Balay 50855843e3eSBarry Smith DrawSynchronizedFlush(draw); 5093270192aSSatish Balay DrawGetPause(draw,&pause); 5103a40ed3dSBarry Smith if (pause >= 0) { PetscSleep(pause); PetscFunctionReturn(0);} 5113270192aSSatish Balay 5123270192aSSatish Balay /* allow the matrix to zoom or shrink */ 51355843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5143270192aSSatish Balay while (button != BUTTON_RIGHT) { 51555843e3eSBarry Smith DrawSynchronizedClear(draw); 5163270192aSSatish Balay if (button == BUTTON_LEFT) scale = .5; 5173270192aSSatish Balay else if (button == BUTTON_CENTER) scale = 2.; 5183270192aSSatish Balay xl = scale*(xl + w - xc) + xc - w*scale; 5193270192aSSatish Balay xr = scale*(xr - w - xc) + xc + w*scale; 5203270192aSSatish Balay yl = scale*(yl + h - yc) + yc - h*scale; 5213270192aSSatish Balay yr = scale*(yr - h - yc) + yc + h*scale; 5223270192aSSatish Balay w *= scale; h *= scale; 5233270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 5243270192aSSatish Balay color = DRAW_BLUE; 5253270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5263270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5273270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5283270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5293270192aSSatish Balay aa = a->a + j*bs2; 5303270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5313270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5323270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 533b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5343270192aSSatish Balay } 5353270192aSSatish Balay } 5363270192aSSatish Balay } 5373270192aSSatish Balay } 5383270192aSSatish Balay color = DRAW_CYAN; 5393270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5403270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5413270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5423270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5433270192aSSatish Balay aa = a->a + j*bs2; 5443270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5453270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5463270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 547b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5483270192aSSatish Balay } 5493270192aSSatish Balay } 5503270192aSSatish Balay } 5513270192aSSatish Balay } 5523270192aSSatish Balay 5533270192aSSatish Balay color = DRAW_RED; 5543270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5553270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5563270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5573270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5583270192aSSatish Balay aa = a->a + j*bs2; 5593270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5603270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5613270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 562b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5633270192aSSatish Balay } 5643270192aSSatish Balay } 5653270192aSSatish Balay } 5663270192aSSatish Balay } 56755843e3eSBarry Smith ierr = DrawSynchronizedGetMouseButton(draw,&button,&xc,&yc,0,0); 5683270192aSSatish Balay } 5693a40ed3dSBarry Smith PetscFunctionReturn(0); 5703270192aSSatish Balay } 5713270192aSSatish Balay 5725615d1e5SSatish Balay #undef __FUNC__ 573d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 574*e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5752593348eSBarry Smith { 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 1088e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1089e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 10901e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1091e51c0b9cSSatish Balay 10922d61bbb3SSatish Balay if (!a->solve_work) { 10932d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 10942d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 10952d61bbb3SSatish Balay } 10962d61bbb3SSatish Balay 10972d61bbb3SSatish Balay if (!a->diag) { 10982d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 10992d61bbb3SSatish Balay } 11002d61bbb3SSatish Balay ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 11012d61bbb3SSatish Balay 11022d61bbb3SSatish Balay /* 11032d61bbb3SSatish Balay Blocksize 4 has a special faster solver for ILU(0) factorization 11042d61bbb3SSatish Balay with natural ordering 11052d61bbb3SSatish Balay */ 11062d61bbb3SSatish Balay if (a->bs == 4) { 1107f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 11082d61bbb3SSatish Balay } 11092d61bbb3SSatish Balay 11102d61bbb3SSatish Balay PetscFunctionReturn(0); 11112d61bbb3SSatish Balay } 11122d61bbb3SSatish Balay #undef __FUNC__ 1113d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1114ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1115ba4ca20aSSatish Balay { 1116ba4ca20aSSatish Balay static int called = 0; 1117ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1118ba4ca20aSSatish Balay 11193a40ed3dSBarry Smith PetscFunctionBegin; 11203a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 112176be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n"); 112276be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n"); 11233a40ed3dSBarry Smith PetscFunctionReturn(0); 1124ba4ca20aSSatish Balay } 1125d9b7c43dSSatish Balay 11262593348eSBarry Smith /* -------------------------------------------------------------------*/ 1127cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ, 11289867e207SSatish Balay MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ, 1129f44d3a6dSSatish Balay MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N, 1130faf6580fSSatish Balay MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ, 11311a6a6d98SBarry Smith MatSolve_SeqBAIJ_N,0, 1132b6490206SBarry Smith 0,0, 1133de6a44a3SBarry Smith MatLUFactor_SeqBAIJ,0, 1134b6490206SBarry Smith 0, 1135f2501298SSatish Balay MatTranspose_SeqBAIJ, 113617e48fc4SSatish Balay MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ, 11379867e207SSatish Balay MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ, 1138584200bdSSatish Balay 0,MatAssemblyEnd_SeqBAIJ, 1139b6490206SBarry Smith 0, 1140d9b7c43dSSatish Balay MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ, 11417fc0212eSBarry Smith MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 1142b6490206SBarry Smith MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 1143de6a44a3SBarry Smith MatILUFactorSymbolic_SeqBAIJ,0, 1144d402145bSBarry Smith 0,0, 1145b6490206SBarry Smith MatConvertSameType_SeqBAIJ,0,0, 1146b6490206SBarry Smith MatILUFactor_SeqBAIJ,0,0, 1147af015674SSatish Balay MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ, 11487e67e3f9SSatish Balay MatGetValues_SeqBAIJ,0, 1149ba4ca20aSSatish Balay MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 11503b2fbd54SBarry Smith 0,0,0,MatGetBlockSize_SeqBAIJ, 11513b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 115292c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 115392c4ed94SBarry Smith 0,0,0,0,0,0, 1154d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1155d3825aa8SBarry Smith MatGetSubMatrix_SeqBAIJ}; 11562593348eSBarry Smith 11575615d1e5SSatish Balay #undef __FUNC__ 11585615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 11592593348eSBarry Smith /*@C 116044cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 116144cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 116244cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 11632bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 11642bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 11652593348eSBarry Smith 11662593348eSBarry Smith Input Parameters: 1167029af93fSBarry Smith . comm - MPI communicator, set to PETSC_COMM_SELF 1168b6490206SBarry Smith . bs - size of block 11692593348eSBarry Smith . m - number of rows 11702593348eSBarry Smith . n - number of columns 1171b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 11722bd5e0b2SLois Curfman McInnes . nzz - array containing the number of block nonzeros in the various block rows 11732bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 11742593348eSBarry Smith 11752593348eSBarry Smith Output Parameter: 11762593348eSBarry Smith . A - the matrix 11772593348eSBarry Smith 11780513a670SBarry Smith Options Database Keys: 11790513a670SBarry Smith $ -mat_no_unroll - uses code that does not unroll the loops in the 11800effe34fSLois Curfman McInnes $ block calculations (much slower) 11810513a670SBarry Smith $ -mat_block_size - size of the blocks to use 11820513a670SBarry Smith 11832593348eSBarry Smith Notes: 118444cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 11852593348eSBarry Smith storage. That is, the stored row and column indices can begin at 118644cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 11872593348eSBarry Smith 11882593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 11892593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 11902593348eSBarry Smith allocation. For additional details, see the users manual chapter on 11916da5968aSLois Curfman McInnes matrices. 11922593348eSBarry Smith 119344cd7ae7SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 11942593348eSBarry Smith @*/ 1195b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 11962593348eSBarry Smith { 11972593348eSBarry Smith Mat B; 1198b6490206SBarry Smith Mat_SeqBAIJ *b; 11993b2fbd54SBarry Smith int i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size; 12003b2fbd54SBarry Smith 12013a40ed3dSBarry Smith PetscFunctionBegin; 12023b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1203a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1204b6490206SBarry Smith 12050513a670SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 12060513a670SBarry Smith 12073a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1208a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 12093a40ed3dSBarry Smith } 12102593348eSBarry Smith 12112593348eSBarry Smith *A = 0; 1212f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,comm,MatDestroy,MatView); 12132593348eSBarry Smith PLogObjectCreate(B); 1214b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 121544cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqBAIJ)); 1216f830108cSBarry Smith PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps)); 12171a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr); 12181a6a6d98SBarry Smith if (!flg) { 12197fc0212eSBarry Smith switch (bs) { 12207fc0212eSBarry Smith case 1: 1221f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1222f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1223f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1224f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 12257fc0212eSBarry Smith break; 12264eeb42bcSBarry Smith case 2: 1227f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1228f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1229f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1230f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 12316d84be18SBarry Smith break; 1232f327dd97SBarry Smith case 3: 1233f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1234f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1235f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1236f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 12374eeb42bcSBarry Smith break; 12386d84be18SBarry Smith case 4: 1239f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1240f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1241f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1242f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 12436d84be18SBarry Smith break; 12446d84be18SBarry Smith case 5: 1245f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1246f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1247f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1248f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 12496d84be18SBarry Smith break; 125048e9ddb2SSatish Balay case 7: 1251f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1252f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1253f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 125448e9ddb2SSatish Balay break; 12557fc0212eSBarry Smith } 12561a6a6d98SBarry Smith } 1257*e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1258*e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 12592593348eSBarry Smith B->factor = 0; 12602593348eSBarry Smith B->lupivotthreshold = 1.0; 126190f02eecSBarry Smith B->mapping = 0; 12622593348eSBarry Smith b->row = 0; 12632593348eSBarry Smith b->col = 0; 1264e51c0b9cSSatish Balay b->icol = 0; 12652593348eSBarry Smith b->reallocs = 0; 12662593348eSBarry Smith 126744cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 126844cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1269b6490206SBarry Smith b->mbs = mbs; 1270f2501298SSatish Balay b->nbs = nbs; 1271b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 12722593348eSBarry Smith if (nnz == PETSC_NULL) { 1273119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 12742593348eSBarry Smith else if (nz <= 0) nz = 1; 1275b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1276b6490206SBarry Smith nz = nz*mbs; 12773a40ed3dSBarry Smith } else { 12782593348eSBarry Smith nz = 0; 1279b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 12802593348eSBarry Smith } 12812593348eSBarry Smith 12822593348eSBarry Smith /* allocate the matrix space */ 12837e67e3f9SSatish Balay len = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int); 12842593348eSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 12857e67e3f9SSatish Balay PetscMemzero(b->a,nz*bs2*sizeof(Scalar)); 12867e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 12872593348eSBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 12882593348eSBarry Smith b->i = b->j + nz; 12892593348eSBarry Smith b->singlemalloc = 1; 12902593348eSBarry Smith 1291de6a44a3SBarry Smith b->i[0] = 0; 1292b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 12932593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 12942593348eSBarry Smith } 12952593348eSBarry Smith 1296b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1297b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1298f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1299b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 13002593348eSBarry Smith 1301b6490206SBarry Smith b->bs = bs; 13029df24120SSatish Balay b->bs2 = bs2; 1303b6490206SBarry Smith b->mbs = mbs; 13042593348eSBarry Smith b->nz = 0; 130518e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 13062593348eSBarry Smith b->sorted = 0; 13072593348eSBarry Smith b->roworiented = 1; 13082593348eSBarry Smith b->nonew = 0; 13092593348eSBarry Smith b->diag = 0; 13102593348eSBarry Smith b->solve_work = 0; 1311de6a44a3SBarry Smith b->mult_work = 0; 13122593348eSBarry Smith b->spptr = 0; 13134e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 13144e220ebcSLois Curfman McInnes 13152593348eSBarry Smith *A = B; 13162593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 13172593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 13183a40ed3dSBarry Smith PetscFunctionReturn(0); 13192593348eSBarry Smith } 13202593348eSBarry Smith 13215615d1e5SSatish Balay #undef __FUNC__ 13225615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqBAIJ" 1323b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 13242593348eSBarry Smith { 13252593348eSBarry Smith Mat C; 1326b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 13279df24120SSatish Balay int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 1328de6a44a3SBarry Smith 13293a40ed3dSBarry Smith PetscFunctionBegin; 1330a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 13312593348eSBarry Smith 13322593348eSBarry Smith *B = 0; 1333f830108cSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,A->comm,MatDestroy,MatView); 13342593348eSBarry Smith PLogObjectCreate(C); 1335b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 1336c3c66ee5SSatish Balay PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 1337*e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1338*e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 13392593348eSBarry Smith C->factor = A->factor; 13402593348eSBarry Smith c->row = 0; 13412593348eSBarry Smith c->col = 0; 13422593348eSBarry Smith C->assembled = PETSC_TRUE; 13432593348eSBarry Smith 134444cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 134544cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 134644cd7ae7SLois Curfman McInnes C->M = a->m; 134744cd7ae7SLois Curfman McInnes C->N = a->n; 134844cd7ae7SLois Curfman McInnes 1349b6490206SBarry Smith c->bs = a->bs; 13509df24120SSatish Balay c->bs2 = a->bs2; 1351b6490206SBarry Smith c->mbs = a->mbs; 1352b6490206SBarry Smith c->nbs = a->nbs; 13532593348eSBarry Smith 1354b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1355b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1356b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 13572593348eSBarry Smith c->imax[i] = a->imax[i]; 13582593348eSBarry Smith c->ilen[i] = a->ilen[i]; 13592593348eSBarry Smith } 13602593348eSBarry Smith 13612593348eSBarry Smith /* allocate the matrix space */ 13622593348eSBarry Smith c->singlemalloc = 1; 13637e67e3f9SSatish Balay len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int)); 13642593348eSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 13657e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1366de6a44a3SBarry Smith c->i = c->j + nz; 1367b6490206SBarry Smith PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1368b6490206SBarry Smith if (mbs > 0) { 1369de6a44a3SBarry Smith PetscMemcpy(c->j,a->j,nz*sizeof(int)); 13702593348eSBarry Smith if (cpvalues == COPY_VALUES) { 13717e67e3f9SSatish Balay PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar)); 13722593348eSBarry Smith } 13732593348eSBarry Smith } 13742593348eSBarry Smith 1375f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 13762593348eSBarry Smith c->sorted = a->sorted; 13772593348eSBarry Smith c->roworiented = a->roworiented; 13782593348eSBarry Smith c->nonew = a->nonew; 13792593348eSBarry Smith 13802593348eSBarry Smith if (a->diag) { 1381b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1382b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1383b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 13842593348eSBarry Smith c->diag[i] = a->diag[i]; 13852593348eSBarry Smith } 13862593348eSBarry Smith } 13872593348eSBarry Smith else c->diag = 0; 13882593348eSBarry Smith c->nz = a->nz; 13892593348eSBarry Smith c->maxnz = a->maxnz; 13902593348eSBarry Smith c->solve_work = 0; 13912593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 13927fc0212eSBarry Smith c->mult_work = 0; 13932593348eSBarry Smith *B = C; 13943a40ed3dSBarry Smith PetscFunctionReturn(0); 13952593348eSBarry Smith } 13962593348eSBarry Smith 13975615d1e5SSatish Balay #undef __FUNC__ 13985615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 139919bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 14002593348eSBarry Smith { 1401b6490206SBarry Smith Mat_SeqBAIJ *a; 14022593348eSBarry Smith Mat B; 1403de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1404b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 140535aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1406a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1407b6490206SBarry Smith Scalar *aa; 140819bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 14092593348eSBarry Smith 14103a40ed3dSBarry Smith PetscFunctionBegin; 14111a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1412de6a44a3SBarry Smith bs2 = bs*bs; 1413b6490206SBarry Smith 14142593348eSBarry Smith MPI_Comm_size(comm,&size); 1415cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 141690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 14170752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 1418a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 14192593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 14202593348eSBarry Smith 1421d64ed03dSBarry Smith if (header[3] < 0) { 1422a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1423d64ed03dSBarry Smith } 1424d64ed03dSBarry Smith 1425a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 142635aab85fSBarry Smith 142735aab85fSBarry Smith /* 142835aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 142935aab85fSBarry Smith divisible by the blocksize 143035aab85fSBarry Smith */ 1431b6490206SBarry Smith mbs = M/bs; 143235aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 143335aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 143435aab85fSBarry Smith else mbs++; 143535aab85fSBarry Smith if (extra_rows) { 1436537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 143735aab85fSBarry Smith } 1438b6490206SBarry Smith 14392593348eSBarry Smith /* read in row lengths */ 144035aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 14410752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 144235aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 14432593348eSBarry Smith 1444b6490206SBarry Smith /* read in column indices */ 144535aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 14460752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr); 144735aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1448b6490206SBarry Smith 1449b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1450b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1451b6490206SBarry Smith PetscMemzero(browlengths,mbs*sizeof(int)); 145235aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 145335aab85fSBarry Smith PetscMemzero(mask,mbs*sizeof(int)); 145435aab85fSBarry Smith masked = mask + mbs; 1455b6490206SBarry Smith rowcount = 0; nzcount = 0; 1456b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 145735aab85fSBarry Smith nmask = 0; 1458b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1459b6490206SBarry Smith kmax = rowlengths[rowcount]; 1460b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 146135aab85fSBarry Smith tmp = jj[nzcount++]/bs; 146235aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1463b6490206SBarry Smith } 1464b6490206SBarry Smith rowcount++; 1465b6490206SBarry Smith } 146635aab85fSBarry Smith browlengths[i] += nmask; 146735aab85fSBarry Smith /* zero out the mask elements we set */ 146835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1469b6490206SBarry Smith } 1470b6490206SBarry Smith 14712593348eSBarry Smith /* create our matrix */ 14723eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 14732593348eSBarry Smith B = *A; 1474b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 14752593348eSBarry Smith 14762593348eSBarry Smith /* set matrix "i" values */ 1477de6a44a3SBarry Smith a->i[0] = 0; 1478b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1479b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1480b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 14812593348eSBarry Smith } 1482b6490206SBarry Smith a->nz = 0; 1483b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 14842593348eSBarry Smith 1485b6490206SBarry Smith /* read in nonzero values */ 148635aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 14870752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr); 148835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1489b6490206SBarry Smith 1490b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1491b6490206SBarry Smith nzcount = 0; jcount = 0; 1492b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1493b6490206SBarry Smith nzcountb = nzcount; 149435aab85fSBarry Smith nmask = 0; 1495b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1496b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1497b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 149835aab85fSBarry Smith tmp = jj[nzcount++]/bs; 149935aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1500b6490206SBarry Smith } 1501b6490206SBarry Smith rowcount++; 1502b6490206SBarry Smith } 1503de6a44a3SBarry Smith /* sort the masked values */ 150477c4ece6SBarry Smith PetscSortInt(nmask,masked); 1505de6a44a3SBarry Smith 1506b6490206SBarry Smith /* set "j" values into matrix */ 1507b6490206SBarry Smith maskcount = 1; 150835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 150935aab85fSBarry Smith a->j[jcount++] = masked[j]; 1510de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1511b6490206SBarry Smith } 1512b6490206SBarry Smith /* set "a" values into matrix */ 1513de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1514b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1515b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1516b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1517de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1518de6a44a3SBarry Smith block = mask[tmp] - 1; 1519de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1520de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1521b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1522b6490206SBarry Smith } 1523b6490206SBarry Smith } 152435aab85fSBarry Smith /* zero out the mask elements we set */ 152535aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1526b6490206SBarry Smith } 1527a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1528b6490206SBarry Smith 1529b6490206SBarry Smith PetscFree(rowlengths); 1530b6490206SBarry Smith PetscFree(browlengths); 1531b6490206SBarry Smith PetscFree(aa); 1532b6490206SBarry Smith PetscFree(jj); 1533b6490206SBarry Smith PetscFree(mask); 1534b6490206SBarry Smith 1535b6490206SBarry Smith B->assembled = PETSC_TRUE; 1536de6a44a3SBarry Smith 15379c01be13SBarry Smith ierr = MatView_Private(B); CHKERRQ(ierr); 15383a40ed3dSBarry Smith PetscFunctionReturn(0); 15392593348eSBarry Smith } 15402593348eSBarry Smith 15412593348eSBarry Smith 15422593348eSBarry Smith 1543