1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*302e33e3SBarry Smith static char vcid[] = "$Id: baij.c,v 1.174 1999/05/04 20:32:27 balay Exp bsmith $"; 32593348eSBarry Smith #endif 42593348eSBarry Smith 52593348eSBarry Smith /* 6b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 72593348eSBarry Smith matrix storage format. 82593348eSBarry Smith */ 93369ce9aSBarry Smith #include "sys.h" 1070f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 111a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 121a6a6d98SBarry Smith #include "src/inline/spops.h" 133b547af2SSatish Balay 142d61bbb3SSatish Balay #define CHUNKSIZE 10 15de6a44a3SBarry Smith 16be5855fcSBarry Smith /* 17be5855fcSBarry Smith Checks for missing diagonals 18be5855fcSBarry Smith */ 19be5855fcSBarry Smith #undef __FUNC__ 20be5855fcSBarry Smith #define __FUNC__ "MatMissingDiag_SeqBAIJ" 21be5855fcSBarry Smith int MatMissingDiag_SeqBAIJ(Mat A) 22be5855fcSBarry Smith { 23be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 24be5855fcSBarry Smith int *diag = a->diag, *jj = a->j,i; 25be5855fcSBarry Smith 26be5855fcSBarry Smith PetscFunctionBegin; 270e8e8aceSBarry Smith for ( i=0; i<a->mbs; i++ ) { 28be5855fcSBarry Smith if (jj[diag[i]] != i) { 29be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 30be5855fcSBarry Smith } 31be5855fcSBarry Smith } 32be5855fcSBarry Smith PetscFunctionReturn(0); 33be5855fcSBarry Smith } 34be5855fcSBarry Smith 355615d1e5SSatish Balay #undef __FUNC__ 365615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ" 37de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 38de6a44a3SBarry Smith { 39de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 407fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 41de6a44a3SBarry Smith 423a40ed3dSBarry Smith PetscFunctionBegin; 43de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 44de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 457fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 46ecc615b2SBarry Smith diag[i] = a->i[i+1]; 47de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 48de6a44a3SBarry Smith if (a->j[j] == i) { 49de6a44a3SBarry Smith diag[i] = j; 50de6a44a3SBarry Smith break; 51de6a44a3SBarry Smith } 52de6a44a3SBarry Smith } 53de6a44a3SBarry Smith } 54de6a44a3SBarry Smith a->diag = diag; 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56de6a44a3SBarry Smith } 572593348eSBarry Smith 582593348eSBarry Smith 593b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 603b2fbd54SBarry Smith 615615d1e5SSatish Balay #undef __FUNC__ 625615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 633b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 643b2fbd54SBarry Smith PetscTruth *done) 653b2fbd54SBarry Smith { 663b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 673b2fbd54SBarry Smith int ierr, n = a->mbs,i; 683b2fbd54SBarry Smith 693a40ed3dSBarry Smith PetscFunctionBegin; 703b2fbd54SBarry Smith *nn = n; 713a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 723b2fbd54SBarry Smith if (symmetric) { 733b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 743b2fbd54SBarry Smith } else if (oshift == 1) { 753b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 763b2fbd54SBarry Smith int nz = a->i[n] + 1; 773b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 783b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 793b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 803b2fbd54SBarry Smith } else { 813b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 823b2fbd54SBarry Smith } 833b2fbd54SBarry Smith 843a40ed3dSBarry Smith PetscFunctionReturn(0); 853b2fbd54SBarry Smith } 863b2fbd54SBarry Smith 875615d1e5SSatish Balay #undef __FUNC__ 88d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 893b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 903b2fbd54SBarry Smith PetscTruth *done) 913b2fbd54SBarry Smith { 923b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 933b2fbd54SBarry Smith int i,n = a->mbs; 943b2fbd54SBarry Smith 953a40ed3dSBarry Smith PetscFunctionBegin; 963a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 973b2fbd54SBarry Smith if (symmetric) { 983b2fbd54SBarry Smith PetscFree(*ia); 993b2fbd54SBarry Smith PetscFree(*ja); 100af5da2bfSBarry Smith } else if (oshift == 1) { 1013b2fbd54SBarry Smith int nz = a->i[n]; 1023b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 1033b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 1043b2fbd54SBarry Smith } 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 1063b2fbd54SBarry Smith } 1073b2fbd54SBarry Smith 1082d61bbb3SSatish Balay #undef __FUNC__ 1092d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 1102d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 1112d61bbb3SSatish Balay { 1122d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 1132d61bbb3SSatish Balay 1142d61bbb3SSatish Balay PetscFunctionBegin; 1152d61bbb3SSatish Balay *bs = baij->bs; 1162d61bbb3SSatish Balay PetscFunctionReturn(0); 1172d61bbb3SSatish Balay } 1182d61bbb3SSatish Balay 1192d61bbb3SSatish Balay 1202d61bbb3SSatish Balay #undef __FUNC__ 1212d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 122e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1232d61bbb3SSatish Balay { 1242d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 125e51c0b9cSSatish Balay int ierr; 1262d61bbb3SSatish Balay 12794d884c6SBarry Smith if (--A->refct > 0) PetscFunctionReturn(0); 12894d884c6SBarry Smith 12994d884c6SBarry Smith if (A->mapping) { 13094d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 13194d884c6SBarry Smith } 13294d884c6SBarry Smith if (A->bmapping) { 13394d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 13494d884c6SBarry Smith } 13561b13de0SBarry Smith if (A->rmap) { 13661b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 13761b13de0SBarry Smith } 13861b13de0SBarry Smith if (A->cmap) { 13961b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 14061b13de0SBarry Smith } 1412d61bbb3SSatish Balay #if defined(USE_PETSC_LOG) 142e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1432d61bbb3SSatish Balay #endif 1442d61bbb3SSatish Balay PetscFree(a->a); 1452d61bbb3SSatish Balay if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 1462d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 1472d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 1482d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 1492d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 1502d61bbb3SSatish Balay if (a->mult_work) PetscFree(a->mult_work); 151e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 1523e90b805SBarry Smith if (a->saved_values) PetscFree(a->saved_values); 1532d61bbb3SSatish Balay PetscFree(a); 1542d61bbb3SSatish Balay PLogObjectDestroy(A); 1552d61bbb3SSatish Balay PetscHeaderDestroy(A); 1562d61bbb3SSatish Balay PetscFunctionReturn(0); 1572d61bbb3SSatish Balay } 1582d61bbb3SSatish Balay 1592d61bbb3SSatish Balay #undef __FUNC__ 1602d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1612d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1622d61bbb3SSatish Balay { 1632d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1642d61bbb3SSatish Balay 1652d61bbb3SSatish Balay PetscFunctionBegin; 1662d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1672d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1682d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1692d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1702d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1714787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1724787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1732d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1742d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1752d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1762d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1772d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1782d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 179b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 180b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 181981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1822d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1832d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1842d61bbb3SSatish Balay } else { 1852d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1862d61bbb3SSatish Balay } 1872d61bbb3SSatish Balay PetscFunctionReturn(0); 1882d61bbb3SSatish Balay } 1892d61bbb3SSatish Balay 1902d61bbb3SSatish Balay 1912d61bbb3SSatish Balay #undef __FUNC__ 1922d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1932d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1942d61bbb3SSatish Balay { 1952d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1962d61bbb3SSatish Balay 1972d61bbb3SSatish Balay PetscFunctionBegin; 1982d61bbb3SSatish Balay if (m) *m = a->m; 1992d61bbb3SSatish Balay if (n) *n = a->n; 2002d61bbb3SSatish Balay PetscFunctionReturn(0); 2012d61bbb3SSatish Balay } 2022d61bbb3SSatish Balay 2032d61bbb3SSatish Balay #undef __FUNC__ 2042d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 2052d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2062d61bbb3SSatish Balay { 2072d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2082d61bbb3SSatish Balay 2092d61bbb3SSatish Balay PetscFunctionBegin; 2102d61bbb3SSatish Balay *m = 0; *n = a->m; 2112d61bbb3SSatish Balay PetscFunctionReturn(0); 2122d61bbb3SSatish Balay } 2132d61bbb3SSatish Balay 2142d61bbb3SSatish Balay #undef __FUNC__ 2152d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 2162d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2172d61bbb3SSatish Balay { 2182d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2192d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2203f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2213f1db9ecSBarry Smith Scalar *v_i; 2222d61bbb3SSatish Balay 2232d61bbb3SSatish Balay PetscFunctionBegin; 2242d61bbb3SSatish Balay bs = a->bs; 2252d61bbb3SSatish Balay ai = a->i; 2262d61bbb3SSatish Balay aj = a->j; 2272d61bbb3SSatish Balay aa = a->a; 2282d61bbb3SSatish Balay bs2 = a->bs2; 2292d61bbb3SSatish Balay 2302d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2312d61bbb3SSatish Balay 2322d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2332d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2342d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2352d61bbb3SSatish Balay *nz = bs*M; 2362d61bbb3SSatish Balay 2372d61bbb3SSatish Balay if (v) { 2382d61bbb3SSatish Balay *v = 0; 2392d61bbb3SSatish Balay if (*nz) { 2402d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v); 2412d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2422d61bbb3SSatish Balay v_i = *v + i*bs; 2432d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2442d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2452d61bbb3SSatish Balay } 2462d61bbb3SSatish Balay } 2472d61bbb3SSatish Balay } 2482d61bbb3SSatish Balay 2492d61bbb3SSatish Balay if (idx) { 2502d61bbb3SSatish Balay *idx = 0; 2512d61bbb3SSatish Balay if (*nz) { 2522d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx); 2532d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2542d61bbb3SSatish Balay idx_i = *idx + i*bs; 2552d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2562d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2572d61bbb3SSatish Balay } 2582d61bbb3SSatish Balay } 2592d61bbb3SSatish Balay } 2602d61bbb3SSatish Balay PetscFunctionReturn(0); 2612d61bbb3SSatish Balay } 2622d61bbb3SSatish Balay 2632d61bbb3SSatish Balay #undef __FUNC__ 2642d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2652d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2662d61bbb3SSatish Balay { 2672d61bbb3SSatish Balay PetscFunctionBegin; 2682d61bbb3SSatish Balay if (idx) {if (*idx) PetscFree(*idx);} 2692d61bbb3SSatish Balay if (v) {if (*v) PetscFree(*v);} 2702d61bbb3SSatish Balay PetscFunctionReturn(0); 2712d61bbb3SSatish Balay } 2722d61bbb3SSatish Balay 2732d61bbb3SSatish Balay #undef __FUNC__ 2742d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2752d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2762d61bbb3SSatish Balay { 2772d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2782d61bbb3SSatish Balay Mat C; 2792d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2802d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2813f1db9ecSBarry Smith MatScalar *array = a->a; 2822d61bbb3SSatish Balay 2832d61bbb3SSatish Balay PetscFunctionBegin; 2842d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2852d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 286549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2872d61bbb3SSatish Balay 2882d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2892d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 2902d61bbb3SSatish Balay PetscFree(col); 2912d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 2922d61bbb3SSatish Balay cols = rows + bs; 2932d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2942d61bbb3SSatish Balay cols[0] = i*bs; 2952d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 2962d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 2972d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 2982d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 2992d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 3002d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3012d61bbb3SSatish Balay array += bs2; 3022d61bbb3SSatish Balay } 3032d61bbb3SSatish Balay } 3042d61bbb3SSatish Balay PetscFree(rows); 3052d61bbb3SSatish Balay 3062d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3072d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3082d61bbb3SSatish Balay 3092d61bbb3SSatish Balay if (B != PETSC_NULL) { 3102d61bbb3SSatish Balay *B = C; 3112d61bbb3SSatish Balay } else { 312f830108cSBarry Smith PetscOps *Abops; 313cc2dc46cSBarry Smith MatOps Aops; 314f830108cSBarry Smith 3152d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 3162d61bbb3SSatish Balay PetscFree(a->a); 3172d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 3182d61bbb3SSatish Balay if (a->diag) PetscFree(a->diag); 3192d61bbb3SSatish Balay if (a->ilen) PetscFree(a->ilen); 3202d61bbb3SSatish Balay if (a->imax) PetscFree(a->imax); 3212d61bbb3SSatish Balay if (a->solve_work) PetscFree(a->solve_work); 3222d61bbb3SSatish Balay PetscFree(a); 323f830108cSBarry Smith 3247413bad6SBarry Smith 3257413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3267413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3277413bad6SBarry Smith 328f830108cSBarry Smith /* 329f830108cSBarry Smith This is horrible, horrible code. We need to keep the 330f830108cSBarry Smith A pointers for the bops and ops but copy everything 331f830108cSBarry Smith else from C. 332f830108cSBarry Smith */ 333f830108cSBarry Smith Abops = A->bops; 334f830108cSBarry Smith Aops = A->ops; 335549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 336f830108cSBarry Smith A->bops = Abops; 337f830108cSBarry Smith A->ops = Aops; 33827a8da17SBarry Smith A->qlist = 0; 339f830108cSBarry Smith 3402d61bbb3SSatish Balay PetscHeaderDestroy(C); 3412d61bbb3SSatish Balay } 3422d61bbb3SSatish Balay PetscFunctionReturn(0); 3432d61bbb3SSatish Balay } 3442d61bbb3SSatish Balay 3455615d1e5SSatish Balay #undef __FUNC__ 346d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 347b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3482593348eSBarry Smith { 349b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3509df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 351b6490206SBarry Smith Scalar *aa; 352ce6f0cecSBarry Smith FILE *file; 3532593348eSBarry Smith 3543a40ed3dSBarry Smith PetscFunctionBegin; 35590ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3562593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3572593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3583b2fbd54SBarry Smith 3592593348eSBarry Smith col_lens[1] = a->m; 3602593348eSBarry Smith col_lens[2] = a->n; 3617e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3622593348eSBarry Smith 3632593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 364b6490206SBarry Smith count = 0; 365b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 366b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 367b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 368b6490206SBarry Smith } 3692593348eSBarry Smith } 3700752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 3712593348eSBarry Smith PetscFree(col_lens); 3722593348eSBarry Smith 3732593348eSBarry Smith /* store column indices (zero start index) */ 37466963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj); 375b6490206SBarry Smith count = 0; 376b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 377b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 378b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 379b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 380b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3812593348eSBarry Smith } 3822593348eSBarry Smith } 383b6490206SBarry Smith } 384b6490206SBarry Smith } 3850752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 386b6490206SBarry Smith PetscFree(jj); 3872593348eSBarry Smith 3882593348eSBarry Smith /* store nonzero values */ 38966963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 390b6490206SBarry Smith count = 0; 391b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 392b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 393b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 394b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 3957e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 396b6490206SBarry Smith } 397b6490206SBarry Smith } 398b6490206SBarry Smith } 399b6490206SBarry Smith } 4000752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 401b6490206SBarry Smith PetscFree(aa); 402ce6f0cecSBarry Smith 403ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 404ce6f0cecSBarry Smith if (file) { 405ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 406ce6f0cecSBarry Smith } 4073a40ed3dSBarry Smith PetscFunctionReturn(0); 4082593348eSBarry Smith } 4092593348eSBarry Smith 4105615d1e5SSatish Balay #undef __FUNC__ 411d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 412b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4132593348eSBarry Smith { 414b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4159df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4162593348eSBarry Smith FILE *fd; 4172593348eSBarry Smith char *outputname; 4182593348eSBarry Smith 4193a40ed3dSBarry Smith PetscFunctionBegin; 42090ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 42177ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 42290ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 423639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 424d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4250ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 4267b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket format not supported"); 4270ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 42844cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 42944cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 43044cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 43144cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 43244cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 4333a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 4340ef38995SBarry Smith if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 43544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 436e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 4370ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 438766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 439e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 4400ef38995SBarry Smith } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 441e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 4420ef38995SBarry Smith } 44344cd7ae7SLois Curfman McInnes #else 4440ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 44544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 4460ef38995SBarry Smith } 44744cd7ae7SLois Curfman McInnes #endif 44844cd7ae7SLois Curfman McInnes } 44944cd7ae7SLois Curfman McInnes } 45044cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 45144cd7ae7SLois Curfman McInnes } 45244cd7ae7SLois Curfman McInnes } 4530ef38995SBarry Smith } else { 454b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 455b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 456b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 457b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 458b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 460e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 46188685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 462e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 4630ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 464766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 465e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 4660ef38995SBarry Smith } else { 467e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 46888685aaeSLois Curfman McInnes } 46988685aaeSLois Curfman McInnes #else 4707e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 47188685aaeSLois Curfman McInnes #endif 4722593348eSBarry Smith } 4732593348eSBarry Smith } 4742593348eSBarry Smith fprintf(fd,"\n"); 4752593348eSBarry Smith } 4762593348eSBarry Smith } 477b6490206SBarry Smith } 4782593348eSBarry Smith fflush(fd); 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 4802593348eSBarry Smith } 4812593348eSBarry Smith 4825615d1e5SSatish Balay #undef __FUNC__ 48377ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 48477ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4853270192aSSatish Balay { 48677ed5343SBarry Smith Mat A = (Mat) Aa; 4873270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4887dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 489fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 4903f1db9ecSBarry Smith MatScalar *aa; 49177ed5343SBarry Smith MPI_Comm comm; 49277ed5343SBarry Smith Viewer viewer; 4933270192aSSatish Balay 4943a40ed3dSBarry Smith PetscFunctionBegin; 49577ed5343SBarry Smith /* 49677ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 49777ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 49877ed5343SBarry Smith rest should return immediately. 49977ed5343SBarry Smith */ 50077ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 501d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 50277ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 5033270192aSSatish Balay 50477ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 50577ed5343SBarry Smith 50677ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 50777ed5343SBarry Smith 5083270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5093270192aSSatish Balay color = DRAW_BLUE; 5103270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5113270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5123270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5133270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5143270192aSSatish Balay aa = a->a + j*bs2; 5153270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5163270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5173270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 518b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5193270192aSSatish Balay } 5203270192aSSatish Balay } 5213270192aSSatish Balay } 5223270192aSSatish Balay } 5233270192aSSatish Balay color = DRAW_CYAN; 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 5383270192aSSatish Balay color = DRAW_RED; 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 } 55277ed5343SBarry Smith PetscFunctionReturn(0); 55377ed5343SBarry Smith } 5543270192aSSatish Balay 55577ed5343SBarry Smith #undef __FUNC__ 55677ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 55777ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 55877ed5343SBarry Smith { 55977ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5607dce120fSSatish Balay int ierr; 5617dce120fSSatish Balay double xl,yl,xr,yr,w,h; 56277ed5343SBarry Smith Draw draw; 56377ed5343SBarry Smith PetscTruth isnull; 5643270192aSSatish Balay 56577ed5343SBarry Smith PetscFunctionBegin; 56677ed5343SBarry Smith 56777ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 56877ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 56977ed5343SBarry Smith 57077ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 57177ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 57277ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5733270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 57477ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 57577ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5763a40ed3dSBarry Smith PetscFunctionReturn(0); 5773270192aSSatish Balay } 5783270192aSSatish Balay 5795615d1e5SSatish Balay #undef __FUNC__ 580d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 581e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5822593348eSBarry Smith { 58319bcc07fSBarry Smith ViewerType vtype; 58419bcc07fSBarry Smith int ierr; 5852593348eSBarry Smith 5863a40ed3dSBarry Smith PetscFunctionBegin; 58719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5887b2a1423SBarry Smith if (PetscTypeCompare(vtype,SOCKET_VIEWER)) { 5897b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 5903f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){ 5913a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 5923f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5933a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 5943f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 5953a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 5965cd90555SBarry Smith } else { 5975cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 5982593348eSBarry Smith } 5993a40ed3dSBarry Smith PetscFunctionReturn(0); 6002593348eSBarry Smith } 601b6490206SBarry Smith 602cd0e1443SSatish Balay 6035615d1e5SSatish Balay #undef __FUNC__ 6042d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6052d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 606cd0e1443SSatish Balay { 607cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6082d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6092d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6102d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6113f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 612cd0e1443SSatish Balay 6133a40ed3dSBarry Smith PetscFunctionBegin; 6142d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 615cd0e1443SSatish Balay row = im[k]; brow = row/bs; 616a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 617a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6182d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6192c3acbe9SBarry Smith nrow = ailen[brow]; 6202d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 621a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 622a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6232d61bbb3SSatish Balay col = in[l] ; 6242d61bbb3SSatish Balay bcol = col/bs; 6252d61bbb3SSatish Balay cidx = col%bs; 6262d61bbb3SSatish Balay ridx = row%bs; 6272d61bbb3SSatish Balay high = nrow; 6282d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6292d61bbb3SSatish Balay while (high-low > 5) { 630cd0e1443SSatish Balay t = (low+high)/2; 631cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 632cd0e1443SSatish Balay else low = t; 633cd0e1443SSatish Balay } 634cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 635cd0e1443SSatish Balay if (rp[i] > bcol) break; 636cd0e1443SSatish Balay if (rp[i] == bcol) { 6372d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6382d61bbb3SSatish Balay goto finished; 639cd0e1443SSatish Balay } 640cd0e1443SSatish Balay } 6412d61bbb3SSatish Balay *v++ = zero; 6422d61bbb3SSatish Balay finished:; 643cd0e1443SSatish Balay } 644cd0e1443SSatish Balay } 6453a40ed3dSBarry Smith PetscFunctionReturn(0); 646cd0e1443SSatish Balay } 647cd0e1443SSatish Balay 6482d61bbb3SSatish Balay 6495615d1e5SSatish Balay #undef __FUNC__ 65005a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 65192c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 65292c4ed94SBarry Smith { 65392c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6548a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 655dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 656549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6573f1db9ecSBarry Smith Scalar *value = v; 6583f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 65992c4ed94SBarry Smith 6603a40ed3dSBarry Smith PetscFunctionBegin; 6610e324ae4SSatish Balay if (roworiented) { 6620e324ae4SSatish Balay stepval = (n-1)*bs; 6630e324ae4SSatish Balay } else { 6640e324ae4SSatish Balay stepval = (m-1)*bs; 6650e324ae4SSatish Balay } 66692c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 66792c4ed94SBarry Smith row = im[k]; 6685ef9f2a5SBarry Smith if (row < 0) continue; 6693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 670a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 67192c4ed94SBarry Smith #endif 67292c4ed94SBarry Smith rp = aj + ai[row]; 67392c4ed94SBarry Smith ap = aa + bs2*ai[row]; 67492c4ed94SBarry Smith rmax = imax[row]; 67592c4ed94SBarry Smith nrow = ailen[row]; 67692c4ed94SBarry Smith low = 0; 67792c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6785ef9f2a5SBarry Smith if (in[l] < 0) continue; 6793a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 680a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 68192c4ed94SBarry Smith #endif 68292c4ed94SBarry Smith col = in[l]; 68392c4ed94SBarry Smith if (roworiented) { 6840e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6850e324ae4SSatish Balay } else { 6860e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 68792c4ed94SBarry Smith } 68892c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 68992c4ed94SBarry Smith while (high-low > 7) { 69092c4ed94SBarry Smith t = (low+high)/2; 69192c4ed94SBarry Smith if (rp[t] > col) high = t; 69292c4ed94SBarry Smith else low = t; 69392c4ed94SBarry Smith } 69492c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 69592c4ed94SBarry Smith if (rp[i] > col) break; 69692c4ed94SBarry Smith if (rp[i] == col) { 6978a84c255SSatish Balay bap = ap + bs2*i; 6980e324ae4SSatish Balay if (roworiented) { 6998a84c255SSatish Balay if (is == ADD_VALUES) { 700dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 701dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7028a84c255SSatish Balay bap[jj] += *value++; 703dd9472c6SBarry Smith } 704dd9472c6SBarry Smith } 7050e324ae4SSatish Balay } else { 706dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 707dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7080e324ae4SSatish Balay bap[jj] = *value++; 7098a84c255SSatish Balay } 710dd9472c6SBarry Smith } 711dd9472c6SBarry Smith } 7120e324ae4SSatish Balay } else { 7130e324ae4SSatish Balay if (is == ADD_VALUES) { 714dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 715dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7160e324ae4SSatish Balay *bap++ += *value++; 717dd9472c6SBarry Smith } 718dd9472c6SBarry Smith } 7190e324ae4SSatish Balay } else { 720dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 721dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7220e324ae4SSatish Balay *bap++ = *value++; 7230e324ae4SSatish Balay } 7248a84c255SSatish Balay } 725dd9472c6SBarry Smith } 726dd9472c6SBarry Smith } 727f1241b54SBarry Smith goto noinsert2; 72892c4ed94SBarry Smith } 72992c4ed94SBarry Smith } 73089280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 731a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 73292c4ed94SBarry Smith if (nrow >= rmax) { 73392c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 73492c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7353f1db9ecSBarry Smith MatScalar *new_a; 73692c4ed94SBarry Smith 737a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 73889280ab3SLois Curfman McInnes 73992c4ed94SBarry Smith /* malloc new storage space */ 7403f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7413f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 74292c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 74392c4ed94SBarry Smith new_i = new_j + new_nz; 74492c4ed94SBarry Smith 74592c4ed94SBarry Smith /* copy over old data into new slots */ 74692c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 74792c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 748549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 74992c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 750549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 751549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 752549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 753549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 75492c4ed94SBarry Smith /* free up old matrix storage */ 75592c4ed94SBarry Smith PetscFree(a->a); 75692c4ed94SBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 75792c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 75892c4ed94SBarry Smith a->singlemalloc = 1; 75992c4ed94SBarry Smith 76092c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 76192c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7623f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 76392c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 76492c4ed94SBarry Smith a->reallocs++; 76592c4ed94SBarry Smith a->nz++; 76692c4ed94SBarry Smith } 76792c4ed94SBarry Smith N = nrow++ - 1; 76892c4ed94SBarry Smith /* shift up all the later entries in this row */ 76992c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 77092c4ed94SBarry Smith rp[ii+1] = rp[ii]; 771549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 77292c4ed94SBarry Smith } 773549d3d68SSatish Balay if (N >= i) { 774549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 775549d3d68SSatish Balay } 77692c4ed94SBarry Smith rp[i] = col; 7778a84c255SSatish Balay bap = ap + bs2*i; 7780e324ae4SSatish Balay if (roworiented) { 779dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 780dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7810e324ae4SSatish Balay bap[jj] = *value++; 782dd9472c6SBarry Smith } 783dd9472c6SBarry Smith } 7840e324ae4SSatish Balay } else { 785dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 786dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7870e324ae4SSatish Balay *bap++ = *value++; 7880e324ae4SSatish Balay } 789dd9472c6SBarry Smith } 790dd9472c6SBarry Smith } 791f1241b54SBarry Smith noinsert2:; 79292c4ed94SBarry Smith low = i; 79392c4ed94SBarry Smith } 79492c4ed94SBarry Smith ailen[row] = nrow; 79592c4ed94SBarry Smith } 7963a40ed3dSBarry Smith PetscFunctionReturn(0); 79792c4ed94SBarry Smith } 79892c4ed94SBarry Smith 799f2501298SSatish Balay 8005615d1e5SSatish Balay #undef __FUNC__ 8015615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8028f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 803584200bdSSatish Balay { 804584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 805584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 806a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 807549d3d68SSatish Balay int mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr; 8083f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 809584200bdSSatish Balay 8103a40ed3dSBarry Smith PetscFunctionBegin; 8113a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 812584200bdSSatish Balay 81343ee02c3SBarry Smith if (m) rmax = ailen[0]; 814584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 815584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 816584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 817d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 818584200bdSSatish Balay if (fshift) { 819a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 820584200bdSSatish Balay N = ailen[i]; 821584200bdSSatish Balay for ( j=0; j<N; j++ ) { 822584200bdSSatish Balay ip[j-fshift] = ip[j]; 823549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 824584200bdSSatish Balay } 825584200bdSSatish Balay } 826584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 827584200bdSSatish Balay } 828584200bdSSatish Balay if (mbs) { 829584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 830584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 831584200bdSSatish Balay } 832584200bdSSatish Balay /* reset ilen and imax for each row */ 833584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 834584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 835584200bdSSatish Balay } 836a7c10996SSatish Balay a->nz = ai[mbs]; 837584200bdSSatish Balay 838584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 839584200bdSSatish Balay if (fshift && a->diag) { 840584200bdSSatish Balay PetscFree(a->diag); 841584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 842584200bdSSatish Balay a->diag = 0; 843584200bdSSatish Balay } 8443d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 845219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8463d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 847584200bdSSatish Balay a->reallocs); 848d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 849e2f3b5e9SSatish Balay a->reallocs = 0; 8504e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8514e220ebcSLois Curfman McInnes 8523a40ed3dSBarry Smith PetscFunctionReturn(0); 853584200bdSSatish Balay } 854584200bdSSatish Balay 8555a838052SSatish Balay 856bea157c4SSatish Balay 857bea157c4SSatish Balay /* 858bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 859bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 860bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 861bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 862bea157c4SSatish Balay */ 8635615d1e5SSatish Balay #undef __FUNC__ 864bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks" 865bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 866d9b7c43dSSatish Balay { 867bea157c4SSatish Balay int i,j,k,row; 868bea157c4SSatish Balay PetscTruth flg; 8693a40ed3dSBarry Smith 870bea157c4SSatish Balay /* PetscFunctionBegin;*/ 871bea157c4SSatish Balay for ( i=0,j=0; i<n; j++ ) { 872bea157c4SSatish Balay row = idx[i]; 873bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 874bea157c4SSatish Balay sizes[j] = 1; 875bea157c4SSatish Balay i++; 876e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 877bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 878bea157c4SSatish Balay i++; 879bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 880bea157c4SSatish Balay flg = PETSC_TRUE; 881bea157c4SSatish Balay for ( k=1; k<bs; k++ ) { 882bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 883bea157c4SSatish Balay flg = PETSC_FALSE; 884bea157c4SSatish Balay break; 885d9b7c43dSSatish Balay } 886bea157c4SSatish Balay } 887bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 888bea157c4SSatish Balay sizes[j] = bs; 889bea157c4SSatish Balay i+= bs; 890bea157c4SSatish Balay } else { 891bea157c4SSatish Balay sizes[j] = 1; 892bea157c4SSatish Balay i++; 893bea157c4SSatish Balay } 894bea157c4SSatish Balay } 895bea157c4SSatish Balay } 896bea157c4SSatish Balay *bs_max = j; 8973a40ed3dSBarry Smith PetscFunctionReturn(0); 898d9b7c43dSSatish Balay } 899d9b7c43dSSatish Balay 9005615d1e5SSatish Balay #undef __FUNC__ 9015615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 9028f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 903d9b7c43dSSatish Balay { 904d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 905b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 906bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9073f1db9ecSBarry Smith Scalar zero = 0.0; 9083f1db9ecSBarry Smith MatScalar *aa; 909d9b7c43dSSatish Balay 9103a40ed3dSBarry Smith PetscFunctionBegin; 911d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 912d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 913d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 914d9b7c43dSSatish Balay 915bea157c4SSatish Balay /* allocate memory for rows,sizes */ 916bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 917bea157c4SSatish Balay sizes = rows + is_n; 918bea157c4SSatish Balay 919bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 920bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 921bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 922bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 923b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 924bea157c4SSatish Balay 925bea157c4SSatish Balay for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) { 926bea157c4SSatish Balay row = rows[j]; 927b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 928bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 929bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 930bea157c4SSatish Balay if (sizes[i] == bs) { 931bea157c4SSatish Balay if (diag) { 932bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 933bea157c4SSatish Balay baij->ilen[row/bs] = 1; 934bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 935bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 936a07cd24cSSatish Balay } 937bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 938bea157c4SSatish Balay for ( k=0; k<bs; k++ ) { 939bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 940bea157c4SSatish Balay } 941bea157c4SSatish Balay } else { /* (!diag) */ 942bea157c4SSatish Balay baij->ilen[row/bs] = 0; 943bea157c4SSatish Balay } /* end (!diag) */ 944bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 945bea157c4SSatish Balay #if defined (USE_PETSC_DEBUG) 946bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 947bea157c4SSatish Balay #endif 948bea157c4SSatish Balay for ( k=0; k<count; k++ ) { 949d9b7c43dSSatish Balay aa[0] = zero; 950d9b7c43dSSatish Balay aa+=bs; 951d9b7c43dSSatish Balay } 952d9b7c43dSSatish Balay if (diag) { 953f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 954d9b7c43dSSatish Balay } 955d9b7c43dSSatish Balay } 956bea157c4SSatish Balay } 957bea157c4SSatish Balay 958bea157c4SSatish Balay PetscFree(rows); 9599a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9603a40ed3dSBarry Smith PetscFunctionReturn(0); 961d9b7c43dSSatish Balay } 9621c351548SSatish Balay 9635615d1e5SSatish Balay #undef __FUNC__ 9642d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9652d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9662d61bbb3SSatish Balay { 9672d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9682d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9692d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9702d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 971549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 9723f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9732d61bbb3SSatish Balay 9742d61bbb3SSatish Balay PetscFunctionBegin; 9752d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9762d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9775ef9f2a5SBarry Smith if (row < 0) continue; 9782d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 97932e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 9802d61bbb3SSatish Balay #endif 9812d61bbb3SSatish Balay rp = aj + ai[brow]; 9822d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9832d61bbb3SSatish Balay rmax = imax[brow]; 9842d61bbb3SSatish Balay nrow = ailen[brow]; 9852d61bbb3SSatish Balay low = 0; 9862d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9875ef9f2a5SBarry Smith if (in[l] < 0) continue; 9882d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g) 98932e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 9902d61bbb3SSatish Balay #endif 9912d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 9922d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 9932d61bbb3SSatish Balay if (roworiented) { 9945ef9f2a5SBarry Smith value = v[l + k*n]; 9952d61bbb3SSatish Balay } else { 9962d61bbb3SSatish Balay value = v[k + l*m]; 9972d61bbb3SSatish Balay } 9982d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 9992d61bbb3SSatish Balay while (high-low > 7) { 10002d61bbb3SSatish Balay t = (low+high)/2; 10012d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10022d61bbb3SSatish Balay else low = t; 10032d61bbb3SSatish Balay } 10042d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 10052d61bbb3SSatish Balay if (rp[i] > bcol) break; 10062d61bbb3SSatish Balay if (rp[i] == bcol) { 10072d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10082d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10092d61bbb3SSatish Balay else *bap = value; 10102d61bbb3SSatish Balay goto noinsert1; 10112d61bbb3SSatish Balay } 10122d61bbb3SSatish Balay } 10132d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10142d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10152d61bbb3SSatish Balay if (nrow >= rmax) { 10162d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10172d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10183f1db9ecSBarry Smith MatScalar *new_a; 10192d61bbb3SSatish Balay 10202d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10212d61bbb3SSatish Balay 10222d61bbb3SSatish Balay /* Malloc new storage space */ 10233f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10243f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 10252d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 10262d61bbb3SSatish Balay new_i = new_j + new_nz; 10272d61bbb3SSatish Balay 10282d61bbb3SSatish Balay /* copy over old data into new slots */ 10292d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10302d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1031549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10322d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1033549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1034549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1035549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1036549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10372d61bbb3SSatish Balay /* free up old matrix storage */ 10382d61bbb3SSatish Balay PetscFree(a->a); 10392d61bbb3SSatish Balay if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 10402d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 10412d61bbb3SSatish Balay a->singlemalloc = 1; 10422d61bbb3SSatish Balay 10432d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10442d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10453f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10462d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10472d61bbb3SSatish Balay a->reallocs++; 10482d61bbb3SSatish Balay a->nz++; 10492d61bbb3SSatish Balay } 10502d61bbb3SSatish Balay N = nrow++ - 1; 10512d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10522d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10532d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1054549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10552d61bbb3SSatish Balay } 1056549d3d68SSatish Balay if (N>=i) { 1057549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1058549d3d68SSatish Balay } 10592d61bbb3SSatish Balay rp[i] = bcol; 10602d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10612d61bbb3SSatish Balay noinsert1:; 10622d61bbb3SSatish Balay low = i; 10632d61bbb3SSatish Balay } 10642d61bbb3SSatish Balay ailen[brow] = nrow; 10652d61bbb3SSatish Balay } 10662d61bbb3SSatish Balay PetscFunctionReturn(0); 10672d61bbb3SSatish Balay } 10682d61bbb3SSatish Balay 10692d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10702d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10712d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 10727b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 10737b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 10742d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10752d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10762d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10772d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10782d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10792d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10802d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10812d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10822d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10832d61bbb3SSatish Balay 10842d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10852d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 10862d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 10872d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 10882d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 10892d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 109015091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 10912d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 10922d61bbb3SSatish Balay 10932d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 10942d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 10952d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 10962d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 10972d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 10982d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 109915091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11002d61bbb3SSatish Balay 11012d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11022d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11032d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11042d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11052d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 110615091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11072d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11082d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11092d61bbb3SSatish Balay 11102d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11112d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11122d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11132d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11142d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 111515091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11162d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11172d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11182d61bbb3SSatish Balay 11192d61bbb3SSatish Balay #undef __FUNC__ 11202d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 11215ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11222d61bbb3SSatish Balay { 11232d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 11242d61bbb3SSatish Balay Mat outA; 11252d61bbb3SSatish Balay int ierr; 1126667159a5SBarry Smith PetscTruth row_identity, col_identity; 11272d61bbb3SSatish Balay 11282d61bbb3SSatish Balay PetscFunctionBegin; 1129b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1130667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1131667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1132667159a5SBarry Smith if (!row_identity || !col_identity) { 1133b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1134667159a5SBarry Smith } 11352d61bbb3SSatish Balay 11362d61bbb3SSatish Balay outA = inA; 11372d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11382d61bbb3SSatish Balay a->row = row; 11392d61bbb3SSatish Balay a->col = col; 11402d61bbb3SSatish Balay 1141e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1142e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr); 11431e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1144e51c0b9cSSatish Balay 11452d61bbb3SSatish Balay if (!a->solve_work) { 11462d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 11472d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 11482d61bbb3SSatish Balay } 11492d61bbb3SSatish Balay 11502d61bbb3SSatish Balay if (!a->diag) { 11512d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA);CHKERRQ(ierr); 11522d61bbb3SSatish Balay } 11532d61bbb3SSatish Balay /* 115415091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 115515091d37SBarry Smith for ILU(0) factorization with natural ordering 11562d61bbb3SSatish Balay */ 115715091d37SBarry Smith switch (a->bs) { 115815091d37SBarry Smith case 2: 115915091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 116015091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 116115091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 116215091d37SBarry Smith break; 116315091d37SBarry Smith case 3: 116415091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 116515091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 116615091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 116715091d37SBarry Smith break; 116815091d37SBarry Smith case 4: 1169667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1170f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 117115091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 117215091d37SBarry Smith break; 117315091d37SBarry Smith case 5: 1174667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1175667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 117615091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 117715091d37SBarry Smith break; 117815091d37SBarry Smith case 6: 117915091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 118015091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 118115091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 118215091d37SBarry Smith break; 118315091d37SBarry Smith case 7: 118415091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 118515091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 118615091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 118715091d37SBarry Smith break; 11882d61bbb3SSatish Balay } 11892d61bbb3SSatish Balay 1190667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1191667159a5SBarry Smith 11922d61bbb3SSatish Balay PetscFunctionReturn(0); 11932d61bbb3SSatish Balay } 11942d61bbb3SSatish Balay #undef __FUNC__ 1195d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1196ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1197ba4ca20aSSatish Balay { 1198ba4ca20aSSatish Balay static int called = 0; 1199ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1200d132466eSBarry Smith int ierr; 1201ba4ca20aSSatish Balay 12023a40ed3dSBarry Smith PetscFunctionBegin; 12033a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 1204d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1205d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12063a40ed3dSBarry Smith PetscFunctionReturn(0); 1207ba4ca20aSSatish Balay } 1208d9b7c43dSSatish Balay 1209fb2e594dSBarry Smith EXTERN_C_BEGIN 121027a8da17SBarry Smith #undef __FUNC__ 121127a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 121227a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 121327a8da17SBarry Smith { 121427a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 121527a8da17SBarry Smith int i,nz,n; 121627a8da17SBarry Smith 121727a8da17SBarry Smith PetscFunctionBegin; 121827a8da17SBarry Smith nz = baij->maxnz; 121927a8da17SBarry Smith n = baij->n; 122027a8da17SBarry Smith for (i=0; i<nz; i++) { 122127a8da17SBarry Smith baij->j[i] = indices[i]; 122227a8da17SBarry Smith } 122327a8da17SBarry Smith baij->nz = nz; 122427a8da17SBarry Smith for ( i=0; i<n; i++ ) { 122527a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 122627a8da17SBarry Smith } 122727a8da17SBarry Smith 122827a8da17SBarry Smith PetscFunctionReturn(0); 122927a8da17SBarry Smith } 1230fb2e594dSBarry Smith EXTERN_C_END 123127a8da17SBarry Smith 123227a8da17SBarry Smith #undef __FUNC__ 123327a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 123427a8da17SBarry Smith /*@ 123527a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 123627a8da17SBarry Smith in the matrix. 123727a8da17SBarry Smith 123827a8da17SBarry Smith Input Parameters: 123927a8da17SBarry Smith + mat - the SeqBAIJ matrix 124027a8da17SBarry Smith - indices - the column indices 124127a8da17SBarry Smith 124215091d37SBarry Smith Level: advanced 124315091d37SBarry Smith 124427a8da17SBarry Smith Notes: 124527a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 124627a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 124727a8da17SBarry Smith of the MatSetValues() operation. 124827a8da17SBarry Smith 124927a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 125027a8da17SBarry Smith MatCreateSeqBAIJ(). 125127a8da17SBarry Smith 125227a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 125327a8da17SBarry Smith 125427a8da17SBarry Smith @*/ 125527a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 125627a8da17SBarry Smith { 125727a8da17SBarry Smith int ierr,(*f)(Mat,int *); 125827a8da17SBarry Smith 125927a8da17SBarry Smith PetscFunctionBegin; 126027a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1261549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 126227a8da17SBarry Smith if (f) { 126327a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 126427a8da17SBarry Smith } else { 126527a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 126627a8da17SBarry Smith } 126727a8da17SBarry Smith PetscFunctionReturn(0); 126827a8da17SBarry Smith } 126927a8da17SBarry Smith 12702593348eSBarry Smith /* -------------------------------------------------------------------*/ 1271cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1272cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1273cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1274cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1275cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1276cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1277cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1278cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1279cc2dc46cSBarry Smith 0, 1280cc2dc46cSBarry Smith 0, 1281cc2dc46cSBarry Smith 0, 1282cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1283cc2dc46cSBarry Smith 0, 1284b6490206SBarry Smith 0, 1285f2501298SSatish Balay MatTranspose_SeqBAIJ, 1286cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1287cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1288cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1289cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1290cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1291b6490206SBarry Smith 0, 1292cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1293cc2dc46cSBarry Smith 0, 1294cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1295cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1296cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1297cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1298cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1299cc2dc46cSBarry Smith 0, 1300cc2dc46cSBarry Smith 0, 1301cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1302cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1303cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1304cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1305cc2dc46cSBarry Smith 0, 1306cc2dc46cSBarry Smith 0, 1307cc2dc46cSBarry Smith 0, 13082e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1309cc2dc46cSBarry Smith 0, 1310cc2dc46cSBarry Smith 0, 1311cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1312cc2dc46cSBarry Smith 0, 1313cc2dc46cSBarry Smith 0, 1314cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1315cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1316cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1317cc2dc46cSBarry Smith 0, 1318cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1319cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1320cc2dc46cSBarry Smith 0, 1321cc2dc46cSBarry Smith 0, 1322cc2dc46cSBarry Smith 0, 1323cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 13243b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 132592c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1326cc2dc46cSBarry Smith 0, 1327cc2dc46cSBarry Smith 0, 1328cc2dc46cSBarry Smith 0, 1329cc2dc46cSBarry Smith 0, 1330cc2dc46cSBarry Smith 0, 1331cc2dc46cSBarry Smith 0, 1332d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1333cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1334cc2dc46cSBarry Smith 0, 1335cc2dc46cSBarry Smith 0, 1336cc2dc46cSBarry Smith MatGetMaps_Petsc}; 13372593348eSBarry Smith 13383e90b805SBarry Smith EXTERN_C_BEGIN 13393e90b805SBarry Smith #undef __FUNC__ 13403e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ" 13413e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 13423e90b805SBarry Smith { 13433e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 13443e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1345d132466eSBarry Smith int ierr; 13463e90b805SBarry Smith 13473e90b805SBarry Smith PetscFunctionBegin; 13483e90b805SBarry Smith if (aij->nonew != 1) { 13493e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13503e90b805SBarry Smith } 13513e90b805SBarry Smith 13523e90b805SBarry Smith /* allocate space for values if not already there */ 13533e90b805SBarry Smith if (!aij->saved_values) { 13543e90b805SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 13553e90b805SBarry Smith } 13563e90b805SBarry Smith 13573e90b805SBarry Smith /* copy values over */ 1358d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 13593e90b805SBarry Smith PetscFunctionReturn(0); 13603e90b805SBarry Smith } 13613e90b805SBarry Smith EXTERN_C_END 13623e90b805SBarry Smith 13633e90b805SBarry Smith EXTERN_C_BEGIN 13643e90b805SBarry Smith #undef __FUNC__ 136532e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ" 136632e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 13673e90b805SBarry Smith { 13683e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1369549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 13703e90b805SBarry Smith 13713e90b805SBarry Smith PetscFunctionBegin; 13723e90b805SBarry Smith if (aij->nonew != 1) { 13733e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13743e90b805SBarry Smith } 13753e90b805SBarry Smith if (!aij->saved_values) { 13763e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 13773e90b805SBarry Smith } 13783e90b805SBarry Smith 13793e90b805SBarry Smith /* copy values over */ 1380549d3d68SSatish Balay ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 13813e90b805SBarry Smith PetscFunctionReturn(0); 13823e90b805SBarry Smith } 13833e90b805SBarry Smith EXTERN_C_END 13843e90b805SBarry Smith 13855615d1e5SSatish Balay #undef __FUNC__ 13865615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 13872593348eSBarry Smith /*@C 138844cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 138944cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 139044cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 13917fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 13922bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 13932593348eSBarry Smith 1394db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1395db81eaa0SLois Curfman McInnes 13962593348eSBarry Smith Input Parameters: 1397db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1398b6490206SBarry Smith . bs - size of block 13992593348eSBarry Smith . m - number of rows 14002593348eSBarry Smith . n - number of columns 1401b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14027fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14032bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14042593348eSBarry Smith 14052593348eSBarry Smith Output Parameter: 14062593348eSBarry Smith . A - the matrix 14072593348eSBarry Smith 14080513a670SBarry Smith Options Database Keys: 1409db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1410db81eaa0SLois Curfman McInnes block calculations (much slower) 1411db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14120513a670SBarry Smith 141315091d37SBarry Smith Level: intermediate 141415091d37SBarry Smith 14152593348eSBarry Smith Notes: 141644cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 14172593348eSBarry Smith storage. That is, the stored row and column indices can begin at 141844cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 14192593348eSBarry Smith 14202593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 14212593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14222593348eSBarry Smith allocation. For additional details, see the users manual chapter on 14236da5968aSLois Curfman McInnes matrices. 14242593348eSBarry Smith 1425db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 14262593348eSBarry Smith @*/ 1427b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 14282593348eSBarry Smith { 14292593348eSBarry Smith Mat B; 1430b6490206SBarry Smith Mat_SeqBAIJ *b; 1431*302e33e3SBarry Smith int i,len,ierr,flg,mbs,nbs,bs2,size; 14323b2fbd54SBarry Smith 14333a40ed3dSBarry Smith PetscFunctionBegin; 1434d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1435a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1436b6490206SBarry Smith 14370513a670SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 1438*302e33e3SBarry Smith mbs = m/bs; 1439*302e33e3SBarry Smith nbs = n/bs; 1440*302e33e3SBarry Smith bs2 = bs*bs; 14410513a670SBarry Smith 14423a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1443a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 14443a40ed3dSBarry Smith } 14452593348eSBarry Smith 14462593348eSBarry Smith *A = 0; 14473f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 14482593348eSBarry Smith PLogObjectCreate(B); 1449b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1450549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1451549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 14521a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 14531a6a6d98SBarry Smith if (!flg) { 14547fc0212eSBarry Smith switch (bs) { 14557fc0212eSBarry Smith case 1: 1456f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1457f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1458f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1459f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 14607fc0212eSBarry Smith break; 14614eeb42bcSBarry Smith case 2: 1462f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1463f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1464f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1465f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 14666d84be18SBarry Smith break; 1467f327dd97SBarry Smith case 3: 1468f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1469f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1470f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1471f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 14724eeb42bcSBarry Smith break; 14736d84be18SBarry Smith case 4: 1474f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1475f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1476f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1477f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 14786d84be18SBarry Smith break; 14796d84be18SBarry Smith case 5: 1480f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1481f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1482f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1483f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 14846d84be18SBarry Smith break; 148515091d37SBarry Smith case 6: 148615091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 148715091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 148815091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 148915091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 149015091d37SBarry Smith break; 149148e9ddb2SSatish Balay case 7: 1492f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1493f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1494f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 149548e9ddb2SSatish Balay break; 14967fc0212eSBarry Smith } 14971a6a6d98SBarry Smith } 1498e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1499e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15002593348eSBarry Smith B->factor = 0; 15012593348eSBarry Smith B->lupivotthreshold = 1.0; 150290f02eecSBarry Smith B->mapping = 0; 15032593348eSBarry Smith b->row = 0; 15042593348eSBarry Smith b->col = 0; 1505e51c0b9cSSatish Balay b->icol = 0; 15062593348eSBarry Smith b->reallocs = 0; 15073e90b805SBarry Smith b->saved_values = 0; 15082593348eSBarry Smith 150944cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 151044cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1511a5ae1ecdSBarry Smith 15127413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 15137413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1514a5ae1ecdSBarry Smith 1515b6490206SBarry Smith b->mbs = mbs; 1516f2501298SSatish Balay b->nbs = nbs; 1517b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax); 15182593348eSBarry Smith if (nnz == PETSC_NULL) { 1519119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 15202593348eSBarry Smith else if (nz <= 0) nz = 1; 1521b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1522b6490206SBarry Smith nz = nz*mbs; 15233a40ed3dSBarry Smith } else { 15242593348eSBarry Smith nz = 0; 1525b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 15262593348eSBarry Smith } 15272593348eSBarry Smith 15282593348eSBarry Smith /* allocate the matrix space */ 15293f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 15303f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a); 1531549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15327e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 1533549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15342593348eSBarry Smith b->i = b->j + nz; 15352593348eSBarry Smith b->singlemalloc = 1; 15362593348eSBarry Smith 1537de6a44a3SBarry Smith b->i[0] = 0; 1538b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 15392593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 15402593348eSBarry Smith } 15412593348eSBarry Smith 1542b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1543b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1544f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1545b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 15462593348eSBarry Smith 1547b6490206SBarry Smith b->bs = bs; 15489df24120SSatish Balay b->bs2 = bs2; 1549b6490206SBarry Smith b->mbs = mbs; 15502593348eSBarry Smith b->nz = 0; 155118e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 15522593348eSBarry Smith b->sorted = 0; 15532593348eSBarry Smith b->roworiented = 1; 15542593348eSBarry Smith b->nonew = 0; 15552593348eSBarry Smith b->diag = 0; 15562593348eSBarry Smith b->solve_work = 0; 1557de6a44a3SBarry Smith b->mult_work = 0; 15582593348eSBarry Smith b->spptr = 0; 15594e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 15604e220ebcSLois Curfman McInnes 15612593348eSBarry Smith *A = B; 15622593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 15632593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 156427a8da17SBarry Smith 15653e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 15663e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 15673e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 15683e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 15693e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 15703e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 157127a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 157227a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 157327a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 15743a40ed3dSBarry Smith PetscFunctionReturn(0); 15752593348eSBarry Smith } 15762593348eSBarry Smith 15775615d1e5SSatish Balay #undef __FUNC__ 15782e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 15792e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 15802593348eSBarry Smith { 15812593348eSBarry Smith Mat C; 1582b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 158327a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1584de6a44a3SBarry Smith 15853a40ed3dSBarry Smith PetscFunctionBegin; 1586a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 15872593348eSBarry Smith 15882593348eSBarry Smith *B = 0; 15893f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 15902593348eSBarry Smith PLogObjectCreate(C); 1591b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1592549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1593e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1594e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 15952593348eSBarry Smith C->factor = A->factor; 15962593348eSBarry Smith c->row = 0; 15972593348eSBarry Smith c->col = 0; 1598cac97260SSatish Balay c->icol = 0; 159932e87ba3SBarry Smith c->saved_values = 0; 16002593348eSBarry Smith C->assembled = PETSC_TRUE; 16012593348eSBarry Smith 160244cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 160344cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 160444cd7ae7SLois Curfman McInnes C->M = a->m; 160544cd7ae7SLois Curfman McInnes C->N = a->n; 160644cd7ae7SLois Curfman McInnes 1607b6490206SBarry Smith c->bs = a->bs; 16089df24120SSatish Balay c->bs2 = a->bs2; 1609b6490206SBarry Smith c->mbs = a->mbs; 1610b6490206SBarry Smith c->nbs = a->nbs; 16112593348eSBarry Smith 1612b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1613b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1614b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16152593348eSBarry Smith c->imax[i] = a->imax[i]; 16162593348eSBarry Smith c->ilen[i] = a->ilen[i]; 16172593348eSBarry Smith } 16182593348eSBarry Smith 16192593348eSBarry Smith /* allocate the matrix space */ 16202593348eSBarry Smith c->singlemalloc = 1; 16213f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 16223f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a); 16237e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1624de6a44a3SBarry Smith c->i = c->j + nz; 1625549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1626b6490206SBarry Smith if (mbs > 0) { 1627549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16282e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1629549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16302e8a6d31SBarry Smith } else { 1631549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16322593348eSBarry Smith } 16332593348eSBarry Smith } 16342593348eSBarry Smith 1635f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16362593348eSBarry Smith c->sorted = a->sorted; 16372593348eSBarry Smith c->roworiented = a->roworiented; 16382593348eSBarry Smith c->nonew = a->nonew; 16392593348eSBarry Smith 16402593348eSBarry Smith if (a->diag) { 1641b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag); 1642b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1643b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16442593348eSBarry Smith c->diag[i] = a->diag[i]; 16452593348eSBarry Smith } 164698305bb5SBarry Smith } else c->diag = 0; 16472593348eSBarry Smith c->nz = a->nz; 16482593348eSBarry Smith c->maxnz = a->maxnz; 16492593348eSBarry Smith c->solve_work = 0; 16502593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 16517fc0212eSBarry Smith c->mult_work = 0; 16522593348eSBarry Smith *B = C; 16537b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 16543a40ed3dSBarry Smith PetscFunctionReturn(0); 16552593348eSBarry Smith } 16562593348eSBarry Smith 16575615d1e5SSatish Balay #undef __FUNC__ 16585615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 165919bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 16602593348eSBarry Smith { 1661b6490206SBarry Smith Mat_SeqBAIJ *a; 16622593348eSBarry Smith Mat B; 1663de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1664b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 166535aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1666a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1667b6490206SBarry Smith Scalar *aa; 166819bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 16692593348eSBarry Smith 16703a40ed3dSBarry Smith PetscFunctionBegin; 16711a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1672de6a44a3SBarry Smith bs2 = bs*bs; 1673b6490206SBarry Smith 1674d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1675cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 167690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16770752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1678a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 16792593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 16802593348eSBarry Smith 1681d64ed03dSBarry Smith if (header[3] < 0) { 1682a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1683d64ed03dSBarry Smith } 1684d64ed03dSBarry Smith 1685a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 168635aab85fSBarry Smith 168735aab85fSBarry Smith /* 168835aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 168935aab85fSBarry Smith divisible by the blocksize 169035aab85fSBarry Smith */ 1691b6490206SBarry Smith mbs = M/bs; 169235aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 169335aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 169435aab85fSBarry Smith else mbs++; 169535aab85fSBarry Smith if (extra_rows) { 1696537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 169735aab85fSBarry Smith } 1698b6490206SBarry Smith 16992593348eSBarry Smith /* read in row lengths */ 170035aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 17010752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 170235aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 17032593348eSBarry Smith 1704b6490206SBarry Smith /* read in column indices */ 170535aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj); 17060752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 170735aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1708b6490206SBarry Smith 1709b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1710b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1711549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 171235aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask); 1713549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 171435aab85fSBarry Smith masked = mask + mbs; 1715b6490206SBarry Smith rowcount = 0; nzcount = 0; 1716b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 171735aab85fSBarry Smith nmask = 0; 1718b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1719b6490206SBarry Smith kmax = rowlengths[rowcount]; 1720b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 172135aab85fSBarry Smith tmp = jj[nzcount++]/bs; 172235aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1723b6490206SBarry Smith } 1724b6490206SBarry Smith rowcount++; 1725b6490206SBarry Smith } 172635aab85fSBarry Smith browlengths[i] += nmask; 172735aab85fSBarry Smith /* zero out the mask elements we set */ 172835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1729b6490206SBarry Smith } 1730b6490206SBarry Smith 17312593348eSBarry Smith /* create our matrix */ 17323eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17332593348eSBarry Smith B = *A; 1734b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 17352593348eSBarry Smith 17362593348eSBarry Smith /* set matrix "i" values */ 1737de6a44a3SBarry Smith a->i[0] = 0; 1738b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1739b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1740b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17412593348eSBarry Smith } 1742b6490206SBarry Smith a->nz = 0; 1743b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 17442593348eSBarry Smith 1745b6490206SBarry Smith /* read in nonzero values */ 174635aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 17470752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 174835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1749b6490206SBarry Smith 1750b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1751b6490206SBarry Smith nzcount = 0; jcount = 0; 1752b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1753b6490206SBarry Smith nzcountb = nzcount; 175435aab85fSBarry Smith nmask = 0; 1755b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1756b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1757b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 175835aab85fSBarry Smith tmp = jj[nzcount++]/bs; 175935aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1760b6490206SBarry Smith } 1761b6490206SBarry Smith rowcount++; 1762b6490206SBarry Smith } 1763de6a44a3SBarry Smith /* sort the masked values */ 176477c4ece6SBarry Smith PetscSortInt(nmask,masked); 1765de6a44a3SBarry Smith 1766b6490206SBarry Smith /* set "j" values into matrix */ 1767b6490206SBarry Smith maskcount = 1; 176835aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 176935aab85fSBarry Smith a->j[jcount++] = masked[j]; 1770de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1771b6490206SBarry Smith } 1772b6490206SBarry Smith /* set "a" values into matrix */ 1773de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1774b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1775b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1776b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1777de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1778de6a44a3SBarry Smith block = mask[tmp] - 1; 1779de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1780de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1781b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1782b6490206SBarry Smith } 1783b6490206SBarry Smith } 178435aab85fSBarry Smith /* zero out the mask elements we set */ 178535aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1786b6490206SBarry Smith } 1787a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1788b6490206SBarry Smith 1789b6490206SBarry Smith PetscFree(rowlengths); 1790b6490206SBarry Smith PetscFree(browlengths); 1791b6490206SBarry Smith PetscFree(aa); 1792b6490206SBarry Smith PetscFree(jj); 1793b6490206SBarry Smith PetscFree(mask); 1794b6490206SBarry Smith 1795b6490206SBarry Smith B->assembled = PETSC_TRUE; 1796de6a44a3SBarry Smith 17979c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 17983a40ed3dSBarry Smith PetscFunctionReturn(0); 17992593348eSBarry Smith } 18002593348eSBarry Smith 18012593348eSBarry Smith 18022593348eSBarry Smith 1803