1*95d5f7c2SBarry Smith /*$Id: baij.c,v 1.202 2000/04/09 04:36:19 bsmith Exp bsmith $*/ 22593348eSBarry Smith 32593348eSBarry Smith /* 4b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 52593348eSBarry Smith matrix storage format. 62593348eSBarry Smith */ 73369ce9aSBarry Smith #include "sys.h" 870f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 91a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 101a6a6d98SBarry Smith #include "src/inline/spops.h" 113b547af2SSatish Balay 12*95d5f7c2SBarry Smith /* UGLY, ugly, ugly 13*95d5f7c2SBarry Smith When MatScalar == Scalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does 14*95d5f7c2SBarry Smith not exist. Otherwise ..._MatScalar() takes matrix elements in single precision and 15*95d5f7c2SBarry Smith inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ() 16*95d5f7c2SBarry Smith converts the entries into single precision and then calls ..._MatScalar() to put them 17*95d5f7c2SBarry Smith into the single precision data structures. 18*95d5f7c2SBarry Smith */ 19*95d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 20*95d5f7c2SBarry Smith extern int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode); 21*95d5f7c2SBarry Smith #else 22*95d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ 23*95d5f7c2SBarry Smith #endif 24*95d5f7c2SBarry Smith 252d61bbb3SSatish Balay #define CHUNKSIZE 10 26de6a44a3SBarry Smith 27be5855fcSBarry Smith /* 28be5855fcSBarry Smith Checks for missing diagonals 29be5855fcSBarry Smith */ 30be5855fcSBarry Smith #undef __FUNC__ 31b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMissingDiagonal_SeqBAIJ" 32c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A) 33be5855fcSBarry Smith { 34be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 35883fce79SBarry Smith int *diag,*jj = a->j,i,ierr; 36be5855fcSBarry Smith 37be5855fcSBarry Smith PetscFunctionBegin; 38c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 39883fce79SBarry Smith diag = a->diag; 400e8e8aceSBarry Smith for (i=0; i<a->mbs; i++) { 41be5855fcSBarry Smith if (jj[diag[i]] != i) { 42be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 43be5855fcSBarry Smith } 44be5855fcSBarry Smith } 45be5855fcSBarry Smith PetscFunctionReturn(0); 46be5855fcSBarry Smith } 47be5855fcSBarry Smith 485615d1e5SSatish Balay #undef __FUNC__ 49b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatMarkDiagonal_SeqBAIJ" 50c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A) 51de6a44a3SBarry Smith { 52de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 537fc0212eSBarry Smith int i,j,*diag,m = a->mbs; 54de6a44a3SBarry Smith 553a40ed3dSBarry Smith PetscFunctionBegin; 56883fce79SBarry Smith if (a->diag) PetscFunctionReturn(0); 57883fce79SBarry Smith 58de6a44a3SBarry Smith diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag); 59de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 607fc0212eSBarry Smith for (i=0; i<m; i++) { 61ecc615b2SBarry Smith diag[i] = a->i[i+1]; 62de6a44a3SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 63de6a44a3SBarry Smith if (a->j[j] == i) { 64de6a44a3SBarry Smith diag[i] = j; 65de6a44a3SBarry Smith break; 66de6a44a3SBarry Smith } 67de6a44a3SBarry Smith } 68de6a44a3SBarry Smith } 69de6a44a3SBarry Smith a->diag = diag; 703a40ed3dSBarry Smith PetscFunctionReturn(0); 71de6a44a3SBarry Smith } 722593348eSBarry Smith 732593348eSBarry Smith 743b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 753b2fbd54SBarry Smith 765615d1e5SSatish Balay #undef __FUNC__ 77b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRowIJ_SeqBAIJ" 780e6d2581SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 793b2fbd54SBarry Smith { 803b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 813b2fbd54SBarry Smith int ierr,n = a->mbs,i; 823b2fbd54SBarry Smith 833a40ed3dSBarry Smith PetscFunctionBegin; 843b2fbd54SBarry Smith *nn = n; 853a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 863b2fbd54SBarry Smith if (symmetric) { 873b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 883b2fbd54SBarry Smith } else if (oshift == 1) { 893b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 903b2fbd54SBarry Smith int nz = a->i[n] + 1; 913b2fbd54SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 923b2fbd54SBarry Smith for (i=0; i<n+1; i++) a->i[i]++; 933b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 943b2fbd54SBarry Smith } else { 953b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 963b2fbd54SBarry Smith } 973b2fbd54SBarry Smith 983a40ed3dSBarry Smith PetscFunctionReturn(0); 993b2fbd54SBarry Smith } 1003b2fbd54SBarry Smith 1015615d1e5SSatish Balay #undef __FUNC__ 102b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRowIJ_SeqBAIJ" 1033b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 1043b2fbd54SBarry Smith PetscTruth *done) 1053b2fbd54SBarry Smith { 1063b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 107606d414cSSatish Balay int i,n = a->mbs,ierr; 1083b2fbd54SBarry Smith 1093a40ed3dSBarry Smith PetscFunctionBegin; 1103a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 1113b2fbd54SBarry Smith if (symmetric) { 112606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 113606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 114af5da2bfSBarry Smith } else if (oshift == 1) { 1153b2fbd54SBarry Smith int nz = a->i[n]; 1163b2fbd54SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 1173b2fbd54SBarry Smith for (i=0; i<n+1; i++) a->i[i]--; 1183b2fbd54SBarry Smith } 1193a40ed3dSBarry Smith PetscFunctionReturn(0); 1203b2fbd54SBarry Smith } 1213b2fbd54SBarry Smith 1222d61bbb3SSatish Balay #undef __FUNC__ 123b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetBlockSize_SeqBAIJ" 1242d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs) 1252d61bbb3SSatish Balay { 1262d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data; 1272d61bbb3SSatish Balay 1282d61bbb3SSatish Balay PetscFunctionBegin; 1292d61bbb3SSatish Balay *bs = baij->bs; 1302d61bbb3SSatish Balay PetscFunctionReturn(0); 1312d61bbb3SSatish Balay } 1322d61bbb3SSatish Balay 1332d61bbb3SSatish Balay 1342d61bbb3SSatish Balay #undef __FUNC__ 135b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDestroy_SeqBAIJ" 136e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1372d61bbb3SSatish Balay { 1382d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 139e51c0b9cSSatish Balay int ierr; 1402d61bbb3SSatish Balay 141433994e6SBarry Smith PetscFunctionBegin; 14294d884c6SBarry Smith if (A->mapping) { 14394d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 14494d884c6SBarry Smith } 14594d884c6SBarry Smith if (A->bmapping) { 14694d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 14794d884c6SBarry Smith } 14861b13de0SBarry Smith if (A->rmap) { 14961b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 15061b13de0SBarry Smith } 15161b13de0SBarry Smith if (A->cmap) { 15261b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 15361b13de0SBarry Smith } 154aa482453SBarry Smith #if defined(PETSC_USE_LOG) 155e1311b90SBarry Smith PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1562d61bbb3SSatish Balay #endif 157606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 158606d414cSSatish Balay if (!a->singlemalloc) { 159606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 160606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 161606d414cSSatish Balay } 162c38d4ed2SBarry Smith if (a->row) { 163c38d4ed2SBarry Smith ierr = ISDestroy(a->row);CHKERRQ(ierr); 164c38d4ed2SBarry Smith } 165c38d4ed2SBarry Smith if (a->col) { 166c38d4ed2SBarry Smith ierr = ISDestroy(a->col);CHKERRQ(ierr); 167c38d4ed2SBarry Smith } 168606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 169606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 170606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 171606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 172606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 173e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 174606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 175606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1762d61bbb3SSatish Balay PLogObjectDestroy(A); 1772d61bbb3SSatish Balay PetscHeaderDestroy(A); 1782d61bbb3SSatish Balay PetscFunctionReturn(0); 1792d61bbb3SSatish Balay } 1802d61bbb3SSatish Balay 1812d61bbb3SSatish Balay #undef __FUNC__ 182b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_SeqBAIJ" 1832d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1842d61bbb3SSatish Balay { 1852d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1862d61bbb3SSatish Balay 1872d61bbb3SSatish Balay PetscFunctionBegin; 188c4992f7dSBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = PETSC_TRUE; 189c4992f7dSBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = PETSC_FALSE; 190c4992f7dSBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = PETSC_TRUE; 191c4992f7dSBarry Smith else if (op == MAT_COLUMNS_UNSORTED) a->sorted = PETSC_FALSE; 192c4992f7dSBarry Smith else if (op == MAT_KEEP_ZEROED_ROWS) a->keepzeroedrows = PETSC_TRUE; 1932d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1944787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1954787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1962d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1972d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1982d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1992d61bbb3SSatish Balay op == MAT_SYMMETRIC || 2002d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 2012d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 202b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 203b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 204981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 2052d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 2062d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 2072d61bbb3SSatish Balay } else { 2082d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 2092d61bbb3SSatish Balay } 2102d61bbb3SSatish Balay PetscFunctionReturn(0); 2112d61bbb3SSatish Balay } 2122d61bbb3SSatish Balay 2132d61bbb3SSatish Balay 2142d61bbb3SSatish Balay #undef __FUNC__ 215b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetSize_SeqBAIJ" 2162d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 2172d61bbb3SSatish Balay { 2182d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2192d61bbb3SSatish Balay 2202d61bbb3SSatish Balay PetscFunctionBegin; 2212d61bbb3SSatish Balay if (m) *m = a->m; 2222d61bbb3SSatish Balay if (n) *n = a->n; 2232d61bbb3SSatish Balay PetscFunctionReturn(0); 2242d61bbb3SSatish Balay } 2252d61bbb3SSatish Balay 2262d61bbb3SSatish Balay #undef __FUNC__ 227b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_SeqBAIJ" 2282d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2292d61bbb3SSatish Balay { 2302d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2312d61bbb3SSatish Balay 2322d61bbb3SSatish Balay PetscFunctionBegin; 2334c49b128SBarry Smith if (m) *m = 0; 2344c49b128SBarry Smith if (n) *n = a->m; 2352d61bbb3SSatish Balay PetscFunctionReturn(0); 2362d61bbb3SSatish Balay } 2372d61bbb3SSatish Balay 2382d61bbb3SSatish Balay #undef __FUNC__ 239b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_SeqBAIJ" 2402d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2412d61bbb3SSatish Balay { 2422d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2432d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2443f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2453f1db9ecSBarry Smith Scalar *v_i; 2462d61bbb3SSatish Balay 2472d61bbb3SSatish Balay PetscFunctionBegin; 2482d61bbb3SSatish Balay bs = a->bs; 2492d61bbb3SSatish Balay ai = a->i; 2502d61bbb3SSatish Balay aj = a->j; 2512d61bbb3SSatish Balay aa = a->a; 2522d61bbb3SSatish Balay bs2 = a->bs2; 2532d61bbb3SSatish Balay 2542d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2552d61bbb3SSatish Balay 2562d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2572d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2582d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2592d61bbb3SSatish Balay *nz = bs*M; 2602d61bbb3SSatish Balay 2612d61bbb3SSatish Balay if (v) { 2622d61bbb3SSatish Balay *v = 0; 2632d61bbb3SSatish Balay if (*nz) { 2642d61bbb3SSatish Balay *v = (Scalar*)PetscMalloc((*nz)*sizeof(Scalar));CHKPTRQ(*v); 2652d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 2662d61bbb3SSatish Balay v_i = *v + i*bs; 2672d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2682d61bbb3SSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 2692d61bbb3SSatish Balay } 2702d61bbb3SSatish Balay } 2712d61bbb3SSatish Balay } 2722d61bbb3SSatish Balay 2732d61bbb3SSatish Balay if (idx) { 2742d61bbb3SSatish Balay *idx = 0; 2752d61bbb3SSatish Balay if (*nz) { 2762d61bbb3SSatish Balay *idx = (int*)PetscMalloc((*nz)*sizeof(int));CHKPTRQ(*idx); 2772d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 2782d61bbb3SSatish Balay idx_i = *idx + i*bs; 2792d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2802d61bbb3SSatish Balay for (j=0; j<bs; j++) {idx_i[j] = itmp++;} 2812d61bbb3SSatish Balay } 2822d61bbb3SSatish Balay } 2832d61bbb3SSatish Balay } 2842d61bbb3SSatish Balay PetscFunctionReturn(0); 2852d61bbb3SSatish Balay } 2862d61bbb3SSatish Balay 2872d61bbb3SSatish Balay #undef __FUNC__ 288b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_SeqBAIJ" 2892d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2902d61bbb3SSatish Balay { 291606d414cSSatish Balay int ierr; 292606d414cSSatish Balay 2932d61bbb3SSatish Balay PetscFunctionBegin; 294606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 295606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2962d61bbb3SSatish Balay PetscFunctionReturn(0); 2972d61bbb3SSatish Balay } 2982d61bbb3SSatish Balay 2992d61bbb3SSatish Balay #undef __FUNC__ 300b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_SeqBAIJ" 3012d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 3022d61bbb3SSatish Balay { 3032d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 3042d61bbb3SSatish Balay Mat C; 3052d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 3060e6d2581SBarry Smith int *rows,*cols,bs2=a->bs2,refct; 307375fe846SBarry Smith Scalar *array; 3082d61bbb3SSatish Balay 3092d61bbb3SSatish Balay PetscFunctionBegin; 3100e6d2581SBarry Smith if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 3112d61bbb3SSatish Balay col = (int*)PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 312549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 3132d61bbb3SSatish Balay 314375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 315375fe846SBarry Smith array = (Scalar *)PetscMalloc(a->bs2*a->nz*sizeof(Scalar));CHKPTRQ(array); 316375fe846SBarry Smith for (i=0; i<a->bs2*a->nz; i++) array[i] = (Scalar)a->a[i]; 317375fe846SBarry Smith #else 3183eda8832SBarry Smith array = a->a; 319375fe846SBarry Smith #endif 320375fe846SBarry Smith 3212d61bbb3SSatish Balay for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 3222d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 323606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 3242d61bbb3SSatish Balay rows = (int*)PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 3252d61bbb3SSatish Balay cols = rows + bs; 3262d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3272d61bbb3SSatish Balay cols[0] = i*bs; 3282d61bbb3SSatish Balay for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 3292d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3302d61bbb3SSatish Balay for (j=0; j<len; j++) { 3312d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3322d61bbb3SSatish Balay for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 3332d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3342d61bbb3SSatish Balay array += bs2; 3352d61bbb3SSatish Balay } 3362d61bbb3SSatish Balay } 337606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 338375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 339375fe846SBarry Smith ierr = PetscFree(array); 340375fe846SBarry Smith #endif 3412d61bbb3SSatish Balay 3422d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3432d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3442d61bbb3SSatish Balay 345c4992f7dSBarry Smith if (B) { 3462d61bbb3SSatish Balay *B = C; 3472d61bbb3SSatish Balay } else { 348f830108cSBarry Smith PetscOps *Abops; 349cc2dc46cSBarry Smith MatOps Aops; 350f830108cSBarry Smith 3512d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 352606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 353606d414cSSatish Balay if (!a->singlemalloc) { 354606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 355606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 356606d414cSSatish Balay } 357606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 358606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 359606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 360606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 361606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 362f830108cSBarry Smith 3637413bad6SBarry Smith 3647413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3657413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3667413bad6SBarry Smith 367f830108cSBarry Smith /* 368f830108cSBarry Smith This is horrible, horrible code. We need to keep the 369f830108cSBarry Smith A pointers for the bops and ops but copy everything 370f830108cSBarry Smith else from C. 371f830108cSBarry Smith */ 372f830108cSBarry Smith Abops = A->bops; 373f830108cSBarry Smith Aops = A->ops; 3740e6d2581SBarry Smith refct = A->refct; 375549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 376f830108cSBarry Smith A->bops = Abops; 377f830108cSBarry Smith A->ops = Aops; 37827a8da17SBarry Smith A->qlist = 0; 3790e6d2581SBarry Smith A->refct = refct; 380f830108cSBarry Smith 3812d61bbb3SSatish Balay PetscHeaderDestroy(C); 3822d61bbb3SSatish Balay } 3832d61bbb3SSatish Balay PetscFunctionReturn(0); 3842d61bbb3SSatish Balay } 3852d61bbb3SSatish Balay 3865615d1e5SSatish Balay #undef __FUNC__ 387b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Binary" 388b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3892593348eSBarry Smith { 390b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 3919df24120SSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 392b6490206SBarry Smith Scalar *aa; 393ce6f0cecSBarry Smith FILE *file; 3942593348eSBarry Smith 3953a40ed3dSBarry Smith PetscFunctionBegin; 39690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3972593348eSBarry Smith col_lens = (int*)PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3982593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3993b2fbd54SBarry Smith 4002593348eSBarry Smith col_lens[1] = a->m; 4012593348eSBarry Smith col_lens[2] = a->n; 4027e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 4032593348eSBarry Smith 4042593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 405b6490206SBarry Smith count = 0; 406b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 407b6490206SBarry Smith for (j=0; j<bs; j++) { 408b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 409b6490206SBarry Smith } 4102593348eSBarry Smith } 4110752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 412606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 4132593348eSBarry Smith 4142593348eSBarry Smith /* store column indices (zero start index) */ 41566963ce1SSatish Balay jj = (int*)PetscMalloc((a->nz+1)*bs2*sizeof(int));CHKPTRQ(jj); 416b6490206SBarry Smith count = 0; 417b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 418b6490206SBarry Smith for (j=0; j<bs; j++) { 419b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 420b6490206SBarry Smith for (l=0; l<bs; l++) { 421b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 4222593348eSBarry Smith } 4232593348eSBarry Smith } 424b6490206SBarry Smith } 425b6490206SBarry Smith } 4260752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 427606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 4282593348eSBarry Smith 4292593348eSBarry Smith /* store nonzero values */ 43066963ce1SSatish Balay aa = (Scalar*)PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 431b6490206SBarry Smith count = 0; 432b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 433b6490206SBarry Smith for (j=0; j<bs; j++) { 434b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 435b6490206SBarry Smith for (l=0; l<bs; l++) { 4367e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 437b6490206SBarry Smith } 438b6490206SBarry Smith } 439b6490206SBarry Smith } 440b6490206SBarry Smith } 4410752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 442606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 443ce6f0cecSBarry Smith 444ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 445ce6f0cecSBarry Smith if (file) { 446ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 447ce6f0cecSBarry Smith } 4483a40ed3dSBarry Smith PetscFunctionReturn(0); 4492593348eSBarry Smith } 4502593348eSBarry Smith 4515615d1e5SSatish Balay #undef __FUNC__ 452b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_ASCII" 453b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4542593348eSBarry Smith { 455b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 4569df24120SSatish Balay int ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4572593348eSBarry Smith char *outputname; 4582593348eSBarry Smith 4593a40ed3dSBarry Smith PetscFunctionBegin; 46077ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 461888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 462639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 463d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4640ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 4656831982aSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 4660ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 4676831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 46844cd7ae7SLois Curfman McInnes for (i=0; i<a->mbs; i++) { 46944cd7ae7SLois Curfman McInnes for (j=0; j<bs; j++) { 4706831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 47144cd7ae7SLois Curfman McInnes for (k=a->i[i]; k<a->i[i+1]; k++) { 47244cd7ae7SLois Curfman McInnes for (l=0; l<bs; l++) { 473aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4740e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 4756831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l, 4760e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4770e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 4786831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l, 4790e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4800e6d2581SBarry Smith } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 4810e6d2581SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4820ef38995SBarry Smith } 48344cd7ae7SLois Curfman McInnes #else 4840ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 4856831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 4860ef38995SBarry Smith } 48744cd7ae7SLois Curfman McInnes #endif 48844cd7ae7SLois Curfman McInnes } 48944cd7ae7SLois Curfman McInnes } 4906831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 49144cd7ae7SLois Curfman McInnes } 49244cd7ae7SLois Curfman McInnes } 4936831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 4940ef38995SBarry Smith } else { 4956831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 496b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 497b6490206SBarry Smith for (j=0; j<bs; j++) { 4986831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 499b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 500b6490206SBarry Smith for (l=0; l<bs; l++) { 501aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5020e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 5036831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 5040e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5050e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 5066831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 5070e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5080ef38995SBarry Smith } else { 5090e6d2581SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 51088685aaeSLois Curfman McInnes } 51188685aaeSLois Curfman McInnes #else 5126831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 51388685aaeSLois Curfman McInnes #endif 5142593348eSBarry Smith } 5152593348eSBarry Smith } 5166831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 5172593348eSBarry Smith } 5182593348eSBarry Smith } 5196831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 520b6490206SBarry Smith } 5216831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 5223a40ed3dSBarry Smith PetscFunctionReturn(0); 5232593348eSBarry Smith } 5242593348eSBarry Smith 5255615d1e5SSatish Balay #undef __FUNC__ 526b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Draw_Zoom" 52777ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 5283270192aSSatish Balay { 52977ed5343SBarry Smith Mat A = (Mat) Aa; 5303270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 531b65db4caSBarry Smith int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2; 5320e6d2581SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 5333f1db9ecSBarry Smith MatScalar *aa; 53477ed5343SBarry Smith Viewer viewer; 5353270192aSSatish Balay 5363a40ed3dSBarry Smith PetscFunctionBegin; 5373270192aSSatish Balay 538b65db4caSBarry Smith /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/ 53977ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 54077ed5343SBarry Smith 54177ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 54277ed5343SBarry Smith 5433270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5443270192aSSatish Balay color = DRAW_BLUE; 5453270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5463270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 5473270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5483270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5493270192aSSatish Balay aa = a->a + j*bs2; 5503270192aSSatish Balay for (k=0; k<bs; k++) { 5513270192aSSatish Balay for (l=0; l<bs; l++) { 5520e6d2581SBarry Smith if (PetscRealPart(*aa++) >= 0.) continue; 553433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5543270192aSSatish Balay } 5553270192aSSatish Balay } 5563270192aSSatish Balay } 5573270192aSSatish Balay } 5583270192aSSatish Balay color = DRAW_CYAN; 5593270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5603270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 5613270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5623270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5633270192aSSatish Balay aa = a->a + j*bs2; 5643270192aSSatish Balay for (k=0; k<bs; k++) { 5653270192aSSatish Balay for (l=0; l<bs; l++) { 5660e6d2581SBarry Smith if (PetscRealPart(*aa++) != 0.) continue; 567433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5683270192aSSatish Balay } 5693270192aSSatish Balay } 5703270192aSSatish Balay } 5713270192aSSatish Balay } 5723270192aSSatish Balay 5733270192aSSatish Balay color = DRAW_RED; 5743270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5753270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 5763270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5773270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5783270192aSSatish Balay aa = a->a + j*bs2; 5793270192aSSatish Balay for (k=0; k<bs; k++) { 5803270192aSSatish Balay for (l=0; l<bs; l++) { 5810e6d2581SBarry Smith if (PetscRealPart(*aa++) <= 0.) continue; 582433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5833270192aSSatish Balay } 5843270192aSSatish Balay } 5853270192aSSatish Balay } 5863270192aSSatish Balay } 58777ed5343SBarry Smith PetscFunctionReturn(0); 58877ed5343SBarry Smith } 5893270192aSSatish Balay 59077ed5343SBarry Smith #undef __FUNC__ 591b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Draw" 59277ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 59377ed5343SBarry Smith { 59477ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 5957dce120fSSatish Balay int ierr; 5960e6d2581SBarry Smith PetscReal xl,yl,xr,yr,w,h; 59777ed5343SBarry Smith Draw draw; 59877ed5343SBarry Smith PetscTruth isnull; 5993270192aSSatish Balay 60077ed5343SBarry Smith PetscFunctionBegin; 60177ed5343SBarry Smith 60277ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 60377ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 60477ed5343SBarry Smith 60577ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 60677ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 60777ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 6083270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 60977ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 61077ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6123270192aSSatish Balay } 6133270192aSSatish Balay 6145615d1e5SSatish Balay #undef __FUNC__ 615b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ" 616e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 6172593348eSBarry Smith { 61819bcc07fSBarry Smith int ierr; 6196831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 6202593348eSBarry Smith 6213a40ed3dSBarry Smith PetscFunctionBegin; 6226831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 6236831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 6246831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 6256831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 6260f5bd95cSBarry Smith if (issocket) { 6277b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 6280f5bd95cSBarry Smith } else if (isascii){ 6293a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6300f5bd95cSBarry Smith } else if (isbinary) { 6313a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6320f5bd95cSBarry Smith } else if (isdraw) { 6333a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6345cd90555SBarry Smith } else { 6350f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 6362593348eSBarry Smith } 6373a40ed3dSBarry Smith PetscFunctionReturn(0); 6382593348eSBarry Smith } 639b6490206SBarry Smith 640cd0e1443SSatish Balay 6415615d1e5SSatish Balay #undef __FUNC__ 642b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetValues_SeqBAIJ" 6432d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 644cd0e1443SSatish Balay { 645cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6462d61bbb3SSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 6472d61bbb3SSatish Balay int *ai = a->i,*ailen = a->ilen; 6482d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6493f1db9ecSBarry Smith MatScalar *ap,*aa = a->a,zero = 0.0; 650cd0e1443SSatish Balay 6513a40ed3dSBarry Smith PetscFunctionBegin; 6522d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 653cd0e1443SSatish Balay row = im[k]; brow = row/bs; 654a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 655a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6562d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6572c3acbe9SBarry Smith nrow = ailen[brow]; 6582d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 659a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 660a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6612d61bbb3SSatish Balay col = in[l] ; 6622d61bbb3SSatish Balay bcol = col/bs; 6632d61bbb3SSatish Balay cidx = col%bs; 6642d61bbb3SSatish Balay ridx = row%bs; 6652d61bbb3SSatish Balay high = nrow; 6662d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6672d61bbb3SSatish Balay while (high-low > 5) { 668cd0e1443SSatish Balay t = (low+high)/2; 669cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 670cd0e1443SSatish Balay else low = t; 671cd0e1443SSatish Balay } 672cd0e1443SSatish Balay for (i=low; i<high; i++) { 673cd0e1443SSatish Balay if (rp[i] > bcol) break; 674cd0e1443SSatish Balay if (rp[i] == bcol) { 6752d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6762d61bbb3SSatish Balay goto finished; 677cd0e1443SSatish Balay } 678cd0e1443SSatish Balay } 6792d61bbb3SSatish Balay *v++ = zero; 6802d61bbb3SSatish Balay finished:; 681cd0e1443SSatish Balay } 682cd0e1443SSatish Balay } 6833a40ed3dSBarry Smith PetscFunctionReturn(0); 684cd0e1443SSatish Balay } 685cd0e1443SSatish Balay 686*95d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 687*95d5f7c2SBarry Smith #undef __FUNC__ 688*95d5f7c2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValuesBlocked_SeqBAIJ" 689*95d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 690*95d5f7c2SBarry Smith { 691*95d5f7c2SBarry Smith int ierr,i,N = m*n*((Mat_SeqBAIJ*)mat->data)->bs2; 692*95d5f7c2SBarry Smith MatScalar *vsingle = (MatScalar*)v; 693*95d5f7c2SBarry Smith 694*95d5f7c2SBarry Smith PetscFunctionBegin; 695*95d5f7c2SBarry Smith for (i=0; i<N; i++) { 696*95d5f7c2SBarry Smith vsingle[i] = v[i]; 697*95d5f7c2SBarry Smith } 698*95d5f7c2SBarry Smith ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr); 699*95d5f7c2SBarry Smith PetscFunctionReturn(0); 700*95d5f7c2SBarry Smith } 701*95d5f7c2SBarry Smith #endif 702*95d5f7c2SBarry Smith 7032d61bbb3SSatish Balay 7045615d1e5SSatish Balay #undef __FUNC__ 705b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValuesBlocked_SeqBAIJ" 706*95d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is) 70792c4ed94SBarry Smith { 70892c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 7098a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 710dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 711549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 712*95d5f7c2SBarry Smith MatScalar *value = v,*ap,*aa = a->a,*bap; 71392c4ed94SBarry Smith 7143a40ed3dSBarry Smith PetscFunctionBegin; 7150e324ae4SSatish Balay if (roworiented) { 7160e324ae4SSatish Balay stepval = (n-1)*bs; 7170e324ae4SSatish Balay } else { 7180e324ae4SSatish Balay stepval = (m-1)*bs; 7190e324ae4SSatish Balay } 72092c4ed94SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 72192c4ed94SBarry Smith row = im[k]; 7225ef9f2a5SBarry Smith if (row < 0) continue; 723aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 724a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 72592c4ed94SBarry Smith #endif 72692c4ed94SBarry Smith rp = aj + ai[row]; 72792c4ed94SBarry Smith ap = aa + bs2*ai[row]; 72892c4ed94SBarry Smith rmax = imax[row]; 72992c4ed94SBarry Smith nrow = ailen[row]; 73092c4ed94SBarry Smith low = 0; 73192c4ed94SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 7325ef9f2a5SBarry Smith if (in[l] < 0) continue; 733aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 734a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 73592c4ed94SBarry Smith #endif 73692c4ed94SBarry Smith col = in[l]; 73792c4ed94SBarry Smith if (roworiented) { 7380e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 7390e324ae4SSatish Balay } else { 7400e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 74192c4ed94SBarry Smith } 74292c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 74392c4ed94SBarry Smith while (high-low > 7) { 74492c4ed94SBarry Smith t = (low+high)/2; 74592c4ed94SBarry Smith if (rp[t] > col) high = t; 74692c4ed94SBarry Smith else low = t; 74792c4ed94SBarry Smith } 74892c4ed94SBarry Smith for (i=low; i<high; i++) { 74992c4ed94SBarry Smith if (rp[i] > col) break; 75092c4ed94SBarry Smith if (rp[i] == col) { 7518a84c255SSatish Balay bap = ap + bs2*i; 7520e324ae4SSatish Balay if (roworiented) { 7538a84c255SSatish Balay if (is == ADD_VALUES) { 754dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 755dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7568a84c255SSatish Balay bap[jj] += *value++; 757dd9472c6SBarry Smith } 758dd9472c6SBarry Smith } 7590e324ae4SSatish Balay } else { 760dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 761dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7620e324ae4SSatish Balay bap[jj] = *value++; 7638a84c255SSatish Balay } 764dd9472c6SBarry Smith } 765dd9472c6SBarry Smith } 7660e324ae4SSatish Balay } else { 7670e324ae4SSatish Balay if (is == ADD_VALUES) { 768dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 769dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7700e324ae4SSatish Balay *bap++ += *value++; 771dd9472c6SBarry Smith } 772dd9472c6SBarry Smith } 7730e324ae4SSatish Balay } else { 774dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 775dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7760e324ae4SSatish Balay *bap++ = *value++; 7770e324ae4SSatish Balay } 7788a84c255SSatish Balay } 779dd9472c6SBarry Smith } 780dd9472c6SBarry Smith } 781f1241b54SBarry Smith goto noinsert2; 78292c4ed94SBarry Smith } 78392c4ed94SBarry Smith } 78489280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 785a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 78692c4ed94SBarry Smith if (nrow >= rmax) { 78792c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 78892c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7893f1db9ecSBarry Smith MatScalar *new_a; 79092c4ed94SBarry Smith 791a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 79289280ab3SLois Curfman McInnes 79392c4ed94SBarry Smith /* malloc new storage space */ 7943f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7953f1db9ecSBarry Smith new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 79692c4ed94SBarry Smith new_j = (int*)(new_a + bs2*new_nz); 79792c4ed94SBarry Smith new_i = new_j + new_nz; 79892c4ed94SBarry Smith 79992c4ed94SBarry Smith /* copy over old data into new slots */ 80092c4ed94SBarry Smith for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 80192c4ed94SBarry Smith for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 802549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 80392c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 804549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 805549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 806549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 807549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 80892c4ed94SBarry Smith /* free up old matrix storage */ 809606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 810606d414cSSatish Balay if (!a->singlemalloc) { 811606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 812606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 813606d414cSSatish Balay } 81492c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 815c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 81692c4ed94SBarry Smith 81792c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 81892c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 8193f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 82092c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 82192c4ed94SBarry Smith a->reallocs++; 82292c4ed94SBarry Smith a->nz++; 82392c4ed94SBarry Smith } 82492c4ed94SBarry Smith N = nrow++ - 1; 82592c4ed94SBarry Smith /* shift up all the later entries in this row */ 82692c4ed94SBarry Smith for (ii=N; ii>=i; ii--) { 82792c4ed94SBarry Smith rp[ii+1] = rp[ii]; 828549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 82992c4ed94SBarry Smith } 830549d3d68SSatish Balay if (N >= i) { 831549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 832549d3d68SSatish Balay } 83392c4ed94SBarry Smith rp[i] = col; 8348a84c255SSatish Balay bap = ap + bs2*i; 8350e324ae4SSatish Balay if (roworiented) { 836dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 837dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 8380e324ae4SSatish Balay bap[jj] = *value++; 839dd9472c6SBarry Smith } 840dd9472c6SBarry Smith } 8410e324ae4SSatish Balay } else { 842dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 843dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 8440e324ae4SSatish Balay *bap++ = *value++; 8450e324ae4SSatish Balay } 846dd9472c6SBarry Smith } 847dd9472c6SBarry Smith } 848f1241b54SBarry Smith noinsert2:; 84992c4ed94SBarry Smith low = i; 85092c4ed94SBarry Smith } 85192c4ed94SBarry Smith ailen[row] = nrow; 85292c4ed94SBarry Smith } 8533a40ed3dSBarry Smith PetscFunctionReturn(0); 85492c4ed94SBarry Smith } 85592c4ed94SBarry Smith 856f2501298SSatish Balay 8575615d1e5SSatish Balay #undef __FUNC__ 858b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatAssemblyEnd_SeqBAIJ" 8598f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 860584200bdSSatish Balay { 861584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 862584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 863a7c10996SSatish Balay int m = a->m,*ip,N,*ailen = a->ilen; 864549d3d68SSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 8653f1db9ecSBarry Smith MatScalar *aa = a->a,*ap; 866584200bdSSatish Balay 8673a40ed3dSBarry Smith PetscFunctionBegin; 8683a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 869584200bdSSatish Balay 87043ee02c3SBarry Smith if (m) rmax = ailen[0]; 871584200bdSSatish Balay for (i=1; i<mbs; i++) { 872584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 873584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 874d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 875584200bdSSatish Balay if (fshift) { 876a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 877584200bdSSatish Balay N = ailen[i]; 878584200bdSSatish Balay for (j=0; j<N; j++) { 879584200bdSSatish Balay ip[j-fshift] = ip[j]; 880549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 881584200bdSSatish Balay } 882584200bdSSatish Balay } 883584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 884584200bdSSatish Balay } 885584200bdSSatish Balay if (mbs) { 886584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 887584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 888584200bdSSatish Balay } 889584200bdSSatish Balay /* reset ilen and imax for each row */ 890584200bdSSatish Balay for (i=0; i<mbs; i++) { 891584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 892584200bdSSatish Balay } 893a7c10996SSatish Balay a->nz = ai[mbs]; 894584200bdSSatish Balay 895584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 896584200bdSSatish Balay if (fshift && a->diag) { 897606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 898584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 899584200bdSSatish Balay a->diag = 0; 900584200bdSSatish Balay } 9013d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 902219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 9033d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 904584200bdSSatish Balay a->reallocs); 905d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 906e2f3b5e9SSatish Balay a->reallocs = 0; 9070e6d2581SBarry Smith A->info.nz_unneeded = (PetscReal)fshift*bs2; 9084e220ebcSLois Curfman McInnes 9093a40ed3dSBarry Smith PetscFunctionReturn(0); 910584200bdSSatish Balay } 911584200bdSSatish Balay 9125a838052SSatish Balay 913bea157c4SSatish Balay 914bea157c4SSatish Balay /* 915bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 916bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 917bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 918bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 919bea157c4SSatish Balay */ 9205615d1e5SSatish Balay #undef __FUNC__ 921b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqBAIJ_Check_Blocks" 922bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 923d9b7c43dSSatish Balay { 924bea157c4SSatish Balay int i,j,k,row; 925bea157c4SSatish Balay PetscTruth flg; 9263a40ed3dSBarry Smith 927433994e6SBarry Smith PetscFunctionBegin; 928bea157c4SSatish Balay for (i=0,j=0; i<n; j++) { 929bea157c4SSatish Balay row = idx[i]; 930bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 931bea157c4SSatish Balay sizes[j] = 1; 932bea157c4SSatish Balay i++; 933e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 934bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 935bea157c4SSatish Balay i++; 936bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 937bea157c4SSatish Balay flg = PETSC_TRUE; 938bea157c4SSatish Balay for (k=1; k<bs; k++) { 939bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 940bea157c4SSatish Balay flg = PETSC_FALSE; 941bea157c4SSatish Balay break; 942d9b7c43dSSatish Balay } 943bea157c4SSatish Balay } 944bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 945bea157c4SSatish Balay sizes[j] = bs; 946bea157c4SSatish Balay i+= bs; 947bea157c4SSatish Balay } else { 948bea157c4SSatish Balay sizes[j] = 1; 949bea157c4SSatish Balay i++; 950bea157c4SSatish Balay } 951bea157c4SSatish Balay } 952bea157c4SSatish Balay } 953bea157c4SSatish Balay *bs_max = j; 9543a40ed3dSBarry Smith PetscFunctionReturn(0); 955d9b7c43dSSatish Balay } 956d9b7c43dSSatish Balay 9575615d1e5SSatish Balay #undef __FUNC__ 958b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqBAIJ" 9598f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is,Scalar *diag) 960d9b7c43dSSatish Balay { 961d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 962b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 963bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9643f1db9ecSBarry Smith Scalar zero = 0.0; 9653f1db9ecSBarry Smith MatScalar *aa; 966d9b7c43dSSatish Balay 9673a40ed3dSBarry Smith PetscFunctionBegin; 968d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 969d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 970d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 971d9b7c43dSSatish Balay 972bea157c4SSatish Balay /* allocate memory for rows,sizes */ 973bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 974bea157c4SSatish Balay sizes = rows + is_n; 975bea157c4SSatish Balay 976bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 977bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 978bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 979dffd3267SBarry Smith if (baij->keepzeroedrows) { 980dffd3267SBarry Smith for (i=0; i<is_n; i++) { sizes[i] = 1; } 981dffd3267SBarry Smith bs_max = is_n; 982dffd3267SBarry Smith } else { 983bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 984dffd3267SBarry Smith } 985b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 986bea157c4SSatish Balay 987bea157c4SSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 988bea157c4SSatish Balay row = rows[j]; 989b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 990bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 991bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 992dffd3267SBarry Smith if (sizes[i] == bs && !baij->keepzeroedrows) { 993bea157c4SSatish Balay if (diag) { 994bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 995bea157c4SSatish Balay baij->ilen[row/bs] = 1; 996bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 997bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 998a07cd24cSSatish Balay } 999bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 1000bea157c4SSatish Balay for (k=0; k<bs; k++) { 1001bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 1002bea157c4SSatish Balay } 1003bea157c4SSatish Balay } else { /* (!diag) */ 1004bea157c4SSatish Balay baij->ilen[row/bs] = 0; 1005bea157c4SSatish Balay } /* end (!diag) */ 1006bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 1007aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 1008bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 1009bea157c4SSatish Balay #endif 1010bea157c4SSatish Balay for (k=0; k<count; k++) { 1011d9b7c43dSSatish Balay aa[0] = zero; 1012d9b7c43dSSatish Balay aa+=bs; 1013d9b7c43dSSatish Balay } 1014d9b7c43dSSatish Balay if (diag) { 1015f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 1016d9b7c43dSSatish Balay } 1017d9b7c43dSSatish Balay } 1018bea157c4SSatish Balay } 1019bea157c4SSatish Balay 1020606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 10219a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10223a40ed3dSBarry Smith PetscFunctionReturn(0); 1023d9b7c43dSSatish Balay } 10241c351548SSatish Balay 10255615d1e5SSatish Balay #undef __FUNC__ 1026b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValues_SeqBAIJ" 10272d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 10282d61bbb3SSatish Balay { 10292d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 10302d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 10312d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 10322d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 1033549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 10343f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 10352d61bbb3SSatish Balay 10362d61bbb3SSatish Balay PetscFunctionBegin; 10372d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 10382d61bbb3SSatish Balay row = im[k]; brow = row/bs; 10395ef9f2a5SBarry Smith if (row < 0) continue; 1040aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 104132e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 10422d61bbb3SSatish Balay #endif 10432d61bbb3SSatish Balay rp = aj + ai[brow]; 10442d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 10452d61bbb3SSatish Balay rmax = imax[brow]; 10462d61bbb3SSatish Balay nrow = ailen[brow]; 10472d61bbb3SSatish Balay low = 0; 10482d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 10495ef9f2a5SBarry Smith if (in[l] < 0) continue; 1050aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 105132e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 10522d61bbb3SSatish Balay #endif 10532d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10542d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10552d61bbb3SSatish Balay if (roworiented) { 10565ef9f2a5SBarry Smith value = v[l + k*n]; 10572d61bbb3SSatish Balay } else { 10582d61bbb3SSatish Balay value = v[k + l*m]; 10592d61bbb3SSatish Balay } 10602d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10612d61bbb3SSatish Balay while (high-low > 7) { 10622d61bbb3SSatish Balay t = (low+high)/2; 10632d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10642d61bbb3SSatish Balay else low = t; 10652d61bbb3SSatish Balay } 10662d61bbb3SSatish Balay for (i=low; i<high; i++) { 10672d61bbb3SSatish Balay if (rp[i] > bcol) break; 10682d61bbb3SSatish Balay if (rp[i] == bcol) { 10692d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10702d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10712d61bbb3SSatish Balay else *bap = value; 10722d61bbb3SSatish Balay goto noinsert1; 10732d61bbb3SSatish Balay } 10742d61bbb3SSatish Balay } 10752d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10762d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10772d61bbb3SSatish Balay if (nrow >= rmax) { 10782d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10792d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10803f1db9ecSBarry Smith MatScalar *new_a; 10812d61bbb3SSatish Balay 10822d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10832d61bbb3SSatish Balay 10842d61bbb3SSatish Balay /* Malloc new storage space */ 10853f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10863f1db9ecSBarry Smith new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 10872d61bbb3SSatish Balay new_j = (int*)(new_a + bs2*new_nz); 10882d61bbb3SSatish Balay new_i = new_j + new_nz; 10892d61bbb3SSatish Balay 10902d61bbb3SSatish Balay /* copy over old data into new slots */ 10912d61bbb3SSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 10922d61bbb3SSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1093549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10942d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1095549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1096549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1097549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1098549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10992d61bbb3SSatish Balay /* free up old matrix storage */ 1100606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1101606d414cSSatish Balay if (!a->singlemalloc) { 1102606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1103606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1104606d414cSSatish Balay } 11052d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1106c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 11072d61bbb3SSatish Balay 11082d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 11092d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 11103f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 11112d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 11122d61bbb3SSatish Balay a->reallocs++; 11132d61bbb3SSatish Balay a->nz++; 11142d61bbb3SSatish Balay } 11152d61bbb3SSatish Balay N = nrow++ - 1; 11162d61bbb3SSatish Balay /* shift up all the later entries in this row */ 11172d61bbb3SSatish Balay for (ii=N; ii>=i; ii--) { 11182d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1119549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 11202d61bbb3SSatish Balay } 1121549d3d68SSatish Balay if (N>=i) { 1122549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1123549d3d68SSatish Balay } 11242d61bbb3SSatish Balay rp[i] = bcol; 11252d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 11262d61bbb3SSatish Balay noinsert1:; 11272d61bbb3SSatish Balay low = i; 11282d61bbb3SSatish Balay } 11292d61bbb3SSatish Balay ailen[brow] = nrow; 11302d61bbb3SSatish Balay } 11312d61bbb3SSatish Balay PetscFunctionReturn(0); 11322d61bbb3SSatish Balay } 11332d61bbb3SSatish Balay 11340e6d2581SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,PetscReal,Mat*); 11350e6d2581SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,PetscReal); 11362d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 11377b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 11387b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 11397c922b88SBarry Smith extern int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec); 11407c922b88SBarry Smith extern int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 11412d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 11420e6d2581SBarry Smith extern int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *); 11432d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*); 11442d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 11452d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 11462d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 11472d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 11482d61bbb3SSatish Balay 11492d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 11502d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 11512d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 11522d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 11532d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 11542d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 115515091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 11562d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 11577c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec); 11587c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec); 11597c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec); 11607c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec); 11617c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec); 11627c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec); 11637c922b88SBarry Smith extern int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec); 11642d61bbb3SSatish Balay 11652d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 11662d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 11672d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 11682d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 11692d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 11702d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 117115091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11722d61bbb3SSatish Balay 11732d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11742d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11752d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11762d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11772d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 117815091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11792d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11802d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11812d61bbb3SSatish Balay 11822d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11832d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11842d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11852d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11862d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 118715091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11882d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11892d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11902d61bbb3SSatish Balay 11912d61bbb3SSatish Balay #undef __FUNC__ 1192b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatILUFactor_SeqBAIJ" 11935ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11942d61bbb3SSatish Balay { 11952d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 11962d61bbb3SSatish Balay Mat outA; 11972d61bbb3SSatish Balay int ierr; 1198667159a5SBarry Smith PetscTruth row_identity,col_identity; 11992d61bbb3SSatish Balay 12002d61bbb3SSatish Balay PetscFunctionBegin; 1201b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1202667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1203667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1204667159a5SBarry Smith if (!row_identity || !col_identity) { 1205b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1206667159a5SBarry Smith } 12072d61bbb3SSatish Balay 12082d61bbb3SSatish Balay outA = inA; 12092d61bbb3SSatish Balay inA->factor = FACTOR_LU; 12102d61bbb3SSatish Balay 12112d61bbb3SSatish Balay if (!a->diag) { 1212c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 12132d61bbb3SSatish Balay } 12142d61bbb3SSatish Balay /* 121515091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 121615091d37SBarry Smith for ILU(0) factorization with natural ordering 12172d61bbb3SSatish Balay */ 121815091d37SBarry Smith switch (a->bs) { 1219f1af5d2fSBarry Smith case 1: 12207c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 1221f1af5d2fSBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 122215091d37SBarry Smith case 2: 122315091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 122415091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 12257c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 122615091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 122715091d37SBarry Smith break; 122815091d37SBarry Smith case 3: 122915091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 123015091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 12317c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 123215091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 123315091d37SBarry Smith break; 123415091d37SBarry Smith case 4: 1235667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1236f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 12377c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering; 123815091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 123915091d37SBarry Smith break; 124015091d37SBarry Smith case 5: 1241667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1242667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 12437c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering; 124415091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 124515091d37SBarry Smith break; 124615091d37SBarry Smith case 6: 124715091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 124815091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 12497c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering; 125015091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 125115091d37SBarry Smith break; 125215091d37SBarry Smith case 7: 125315091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 12547c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering; 125515091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 125615091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 125715091d37SBarry Smith break; 1258c38d4ed2SBarry Smith default: 1259c38d4ed2SBarry Smith a->row = row; 1260c38d4ed2SBarry Smith a->col = col; 1261c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1262c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1263c38d4ed2SBarry Smith 1264c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 12654c49b128SBarry Smith ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); 1266c38d4ed2SBarry Smith PLogObjectParent(inA,a->icol); 1267c38d4ed2SBarry Smith 1268c38d4ed2SBarry Smith if (!a->solve_work) { 1269c38d4ed2SBarry Smith a->solve_work = (Scalar*)PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1270c38d4ed2SBarry Smith PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 1271c38d4ed2SBarry Smith } 12722d61bbb3SSatish Balay } 12732d61bbb3SSatish Balay 1274667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1275667159a5SBarry Smith 12762d61bbb3SSatish Balay PetscFunctionReturn(0); 12772d61bbb3SSatish Balay } 12782d61bbb3SSatish Balay #undef __FUNC__ 1279b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatPrintHelp_SeqBAIJ" 1280ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1281ba4ca20aSSatish Balay { 1282c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1283ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1284d132466eSBarry Smith int ierr; 1285ba4ca20aSSatish Balay 12863a40ed3dSBarry Smith PetscFunctionBegin; 1287c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1288d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1289d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12903a40ed3dSBarry Smith PetscFunctionReturn(0); 1291ba4ca20aSSatish Balay } 1292d9b7c43dSSatish Balay 1293fb2e594dSBarry Smith EXTERN_C_BEGIN 129427a8da17SBarry Smith #undef __FUNC__ 1295b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetColumnIndices_SeqBAIJ" 129627a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 129727a8da17SBarry Smith { 129827a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 129927a8da17SBarry Smith int i,nz,n; 130027a8da17SBarry Smith 130127a8da17SBarry Smith PetscFunctionBegin; 130227a8da17SBarry Smith nz = baij->maxnz; 130327a8da17SBarry Smith n = baij->n; 130427a8da17SBarry Smith for (i=0; i<nz; i++) { 130527a8da17SBarry Smith baij->j[i] = indices[i]; 130627a8da17SBarry Smith } 130727a8da17SBarry Smith baij->nz = nz; 130827a8da17SBarry Smith for (i=0; i<n; i++) { 130927a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 131027a8da17SBarry Smith } 131127a8da17SBarry Smith 131227a8da17SBarry Smith PetscFunctionReturn(0); 131327a8da17SBarry Smith } 1314fb2e594dSBarry Smith EXTERN_C_END 131527a8da17SBarry Smith 131627a8da17SBarry Smith #undef __FUNC__ 1317b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetColumnIndices" 131827a8da17SBarry Smith /*@ 131927a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 132027a8da17SBarry Smith in the matrix. 132127a8da17SBarry Smith 132227a8da17SBarry Smith Input Parameters: 132327a8da17SBarry Smith + mat - the SeqBAIJ matrix 132427a8da17SBarry Smith - indices - the column indices 132527a8da17SBarry Smith 132615091d37SBarry Smith Level: advanced 132715091d37SBarry Smith 132827a8da17SBarry Smith Notes: 132927a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 133027a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 133127a8da17SBarry Smith of the MatSetValues() operation. 133227a8da17SBarry Smith 133327a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 133427a8da17SBarry Smith MatCreateSeqBAIJ(). 133527a8da17SBarry Smith 133627a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 133727a8da17SBarry Smith 133827a8da17SBarry Smith @*/ 133927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 134027a8da17SBarry Smith { 134127a8da17SBarry Smith int ierr,(*f)(Mat,int *); 134227a8da17SBarry Smith 134327a8da17SBarry Smith PetscFunctionBegin; 134427a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1345549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 134627a8da17SBarry Smith if (f) { 134727a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 134827a8da17SBarry Smith } else { 134927a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 135027a8da17SBarry Smith } 135127a8da17SBarry Smith PetscFunctionReturn(0); 135227a8da17SBarry Smith } 135327a8da17SBarry Smith 13542593348eSBarry Smith /* -------------------------------------------------------------------*/ 1355cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1356cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1357cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1358cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1359cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 13607c922b88SBarry Smith MatMultTranspose_SeqBAIJ, 13617c922b88SBarry Smith MatMultTransposeAdd_SeqBAIJ, 1362cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1363cc2dc46cSBarry Smith 0, 1364cc2dc46cSBarry Smith 0, 1365cc2dc46cSBarry Smith 0, 1366cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1367cc2dc46cSBarry Smith 0, 1368b6490206SBarry Smith 0, 1369f2501298SSatish Balay MatTranspose_SeqBAIJ, 1370cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1371cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1372cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1373cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1374cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1375b6490206SBarry Smith 0, 1376cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1377cc2dc46cSBarry Smith 0, 1378cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1379cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1380cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1381cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1382cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1383cc2dc46cSBarry Smith 0, 1384cc2dc46cSBarry Smith 0, 1385cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1386cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1387cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1388cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1389cc2dc46cSBarry Smith 0, 1390cc2dc46cSBarry Smith 0, 1391cc2dc46cSBarry Smith 0, 13922e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1393cc2dc46cSBarry Smith 0, 1394cc2dc46cSBarry Smith 0, 1395cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1396cc2dc46cSBarry Smith 0, 1397cc2dc46cSBarry Smith 0, 1398cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1399cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1400cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1401cc2dc46cSBarry Smith 0, 1402cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1403cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1404cc2dc46cSBarry Smith 0, 1405cc2dc46cSBarry Smith 0, 1406cc2dc46cSBarry Smith 0, 1407cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 14083b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 140992c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1410cc2dc46cSBarry Smith 0, 1411cc2dc46cSBarry Smith 0, 1412cc2dc46cSBarry Smith 0, 1413cc2dc46cSBarry Smith 0, 1414cc2dc46cSBarry Smith 0, 1415cc2dc46cSBarry Smith 0, 1416d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1417cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1418cc2dc46cSBarry Smith 0, 1419cc2dc46cSBarry Smith 0, 1420cc2dc46cSBarry Smith MatGetMaps_Petsc}; 14212593348eSBarry Smith 14223e90b805SBarry Smith EXTERN_C_BEGIN 14233e90b805SBarry Smith #undef __FUNC__ 1424b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatStoreValues_SeqBAIJ" 14253e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 14263e90b805SBarry Smith { 14273e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 14283e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1429d132466eSBarry Smith int ierr; 14303e90b805SBarry Smith 14313e90b805SBarry Smith PetscFunctionBegin; 14323e90b805SBarry Smith if (aij->nonew != 1) { 14333e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14343e90b805SBarry Smith } 14353e90b805SBarry Smith 14363e90b805SBarry Smith /* allocate space for values if not already there */ 14373e90b805SBarry Smith if (!aij->saved_values) { 14383e90b805SBarry Smith aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 14393e90b805SBarry Smith } 14403e90b805SBarry Smith 14413e90b805SBarry Smith /* copy values over */ 1442d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 14433e90b805SBarry Smith PetscFunctionReturn(0); 14443e90b805SBarry Smith } 14453e90b805SBarry Smith EXTERN_C_END 14463e90b805SBarry Smith 14473e90b805SBarry Smith EXTERN_C_BEGIN 14483e90b805SBarry Smith #undef __FUNC__ 1449b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRetrieveValues_SeqBAIJ" 145032e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 14513e90b805SBarry Smith { 14523e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1453549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 14543e90b805SBarry Smith 14553e90b805SBarry Smith PetscFunctionBegin; 14563e90b805SBarry Smith if (aij->nonew != 1) { 14573e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14583e90b805SBarry Smith } 14593e90b805SBarry Smith if (!aij->saved_values) { 14603e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 14613e90b805SBarry Smith } 14623e90b805SBarry Smith 14633e90b805SBarry Smith /* copy values over */ 1464549d3d68SSatish Balay ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 14653e90b805SBarry Smith PetscFunctionReturn(0); 14663e90b805SBarry Smith } 14673e90b805SBarry Smith EXTERN_C_END 14683e90b805SBarry Smith 14695615d1e5SSatish Balay #undef __FUNC__ 1470b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateSeqBAIJ" 14712593348eSBarry Smith /*@C 147244cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 147344cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 147444cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 14757fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 14762bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 14772593348eSBarry Smith 1478db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1479db81eaa0SLois Curfman McInnes 14802593348eSBarry Smith Input Parameters: 1481db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1482b6490206SBarry Smith . bs - size of block 14832593348eSBarry Smith . m - number of rows 14842593348eSBarry Smith . n - number of columns 1485b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14867fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14872bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14882593348eSBarry Smith 14892593348eSBarry Smith Output Parameter: 14902593348eSBarry Smith . A - the matrix 14912593348eSBarry Smith 14920513a670SBarry Smith Options Database Keys: 1493db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1494db81eaa0SLois Curfman McInnes block calculations (much slower) 1495db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14960513a670SBarry Smith 149715091d37SBarry Smith Level: intermediate 149815091d37SBarry Smith 14992593348eSBarry Smith Notes: 150044cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 15012593348eSBarry Smith storage. That is, the stored row and column indices can begin at 150244cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 15032593348eSBarry Smith 15042593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 15052593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 15062593348eSBarry Smith allocation. For additional details, see the users manual chapter on 15076da5968aSLois Curfman McInnes matrices. 15082593348eSBarry Smith 1509db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 15102593348eSBarry Smith @*/ 1511b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 15122593348eSBarry Smith { 15132593348eSBarry Smith Mat B; 1514b6490206SBarry Smith Mat_SeqBAIJ *b; 1515f1af5d2fSBarry Smith int i,len,ierr,mbs,nbs,bs2,size; 1516f1af5d2fSBarry Smith PetscTruth flg; 15173b2fbd54SBarry Smith 15183a40ed3dSBarry Smith PetscFunctionBegin; 1519d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1520a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1521b6490206SBarry Smith 1522962fb4a0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1523302e33e3SBarry Smith mbs = m/bs; 1524302e33e3SBarry Smith nbs = n/bs; 1525302e33e3SBarry Smith bs2 = bs*bs; 15260513a670SBarry Smith 15273a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1528a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 15293a40ed3dSBarry Smith } 15302593348eSBarry Smith 1531b73539f3SBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz); 1532b73539f3SBarry Smith if (nnz) { 1533faf3f1fcSBarry Smith for (i=0; i<mbs; i++) { 1534b73539f3SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1535b73539f3SBarry Smith } 1536b73539f3SBarry Smith } 1537b73539f3SBarry Smith 15382593348eSBarry Smith *A = 0; 15393f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 15402593348eSBarry Smith PLogObjectCreate(B); 1541b6490206SBarry Smith B->data = (void*)(b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1542549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1543549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 15441a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 15451a6a6d98SBarry Smith if (!flg) { 15467fc0212eSBarry Smith switch (bs) { 15477fc0212eSBarry Smith case 1: 1548f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1549f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 15507c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 1551f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1552f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 15537fc0212eSBarry Smith break; 15544eeb42bcSBarry Smith case 2: 1555f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1556f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 15577c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 1558f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1559f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 15606d84be18SBarry Smith break; 1561f327dd97SBarry Smith case 3: 1562f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1563f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 15647c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 1565f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1566f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 15674eeb42bcSBarry Smith break; 15686d84be18SBarry Smith case 4: 1569f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1570f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 15717c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4; 1572f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1573f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 15746d84be18SBarry Smith break; 15756d84be18SBarry Smith case 5: 1576f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1577f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 15787c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5; 1579f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1580f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 15816d84be18SBarry Smith break; 158215091d37SBarry Smith case 6: 158315091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 158415091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 15857c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6; 158615091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 158715091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 158815091d37SBarry Smith break; 158948e9ddb2SSatish Balay case 7: 1590f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1591f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 15927c922b88SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7; 1593f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 159448e9ddb2SSatish Balay break; 15957fc0212eSBarry Smith } 15961a6a6d98SBarry Smith } 1597e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1598e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15992593348eSBarry Smith B->factor = 0; 16002593348eSBarry Smith B->lupivotthreshold = 1.0; 160190f02eecSBarry Smith B->mapping = 0; 16022593348eSBarry Smith b->row = 0; 16032593348eSBarry Smith b->col = 0; 1604e51c0b9cSSatish Balay b->icol = 0; 16052593348eSBarry Smith b->reallocs = 0; 16063e90b805SBarry Smith b->saved_values = 0; 16072593348eSBarry Smith 160844cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 160944cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1610a5ae1ecdSBarry Smith 16117413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 16127413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1613a5ae1ecdSBarry Smith 1614b6490206SBarry Smith b->mbs = mbs; 1615f2501298SSatish Balay b->nbs = nbs; 1616b6490206SBarry Smith b->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax); 1617c4992f7dSBarry Smith if (!nnz) { 1618119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 16192593348eSBarry Smith else if (nz <= 0) nz = 1; 1620b6490206SBarry Smith for (i=0; i<mbs; i++) b->imax[i] = nz; 1621b6490206SBarry Smith nz = nz*mbs; 16223a40ed3dSBarry Smith } else { 16232593348eSBarry Smith nz = 0; 1624b6490206SBarry Smith for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 16252593348eSBarry Smith } 16262593348eSBarry Smith 16272593348eSBarry Smith /* allocate the matrix space */ 16283f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 16293f1db9ecSBarry Smith b->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a); 1630549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 16317e67e3f9SSatish Balay b->j = (int*)(b->a + nz*bs2); 1632549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 16332593348eSBarry Smith b->i = b->j + nz; 1634c4992f7dSBarry Smith b->singlemalloc = PETSC_TRUE; 16352593348eSBarry Smith 1636de6a44a3SBarry Smith b->i[0] = 0; 1637b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 16382593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 16392593348eSBarry Smith } 16402593348eSBarry Smith 1641b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1642b6490206SBarry Smith b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int)); 1643f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1644b6490206SBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 16452593348eSBarry Smith 1646b6490206SBarry Smith b->bs = bs; 16479df24120SSatish Balay b->bs2 = bs2; 1648b6490206SBarry Smith b->mbs = mbs; 16492593348eSBarry Smith b->nz = 0; 165018e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 1651c4992f7dSBarry Smith b->sorted = PETSC_FALSE; 1652c4992f7dSBarry Smith b->roworiented = PETSC_TRUE; 16532593348eSBarry Smith b->nonew = 0; 16542593348eSBarry Smith b->diag = 0; 16552593348eSBarry Smith b->solve_work = 0; 1656de6a44a3SBarry Smith b->mult_work = 0; 16572593348eSBarry Smith b->spptr = 0; 16580e6d2581SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 1659883fce79SBarry Smith b->keepzeroedrows = PETSC_FALSE; 16604e220ebcSLois Curfman McInnes 16612593348eSBarry Smith *A = B; 16622593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 16632593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 166427a8da17SBarry Smith 1665f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 16663e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 16673e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1668f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 16693e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 16703e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1671f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 167227a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 167327a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 16743a40ed3dSBarry Smith PetscFunctionReturn(0); 16752593348eSBarry Smith } 16762593348eSBarry Smith 16775615d1e5SSatish Balay #undef __FUNC__ 1678b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_SeqBAIJ" 16792e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 16802593348eSBarry Smith { 16812593348eSBarry Smith Mat C; 1682b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data; 168327a8da17SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1684de6a44a3SBarry Smith 16853a40ed3dSBarry Smith PetscFunctionBegin; 1686a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 16872593348eSBarry Smith 16882593348eSBarry Smith *B = 0; 16893f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 16902593348eSBarry Smith PLogObjectCreate(C); 1691b6490206SBarry Smith C->data = (void*)(c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1692549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1693e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1694e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 16952593348eSBarry Smith C->factor = A->factor; 16962593348eSBarry Smith c->row = 0; 16972593348eSBarry Smith c->col = 0; 1698cac97260SSatish Balay c->icol = 0; 169932e87ba3SBarry Smith c->saved_values = 0; 1700f1e2ffcdSBarry Smith c->keepzeroedrows = a->keepzeroedrows; 17012593348eSBarry Smith C->assembled = PETSC_TRUE; 17022593348eSBarry Smith 170344cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 170444cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 170544cd7ae7SLois Curfman McInnes C->M = a->m; 170644cd7ae7SLois Curfman McInnes C->N = a->n; 170744cd7ae7SLois Curfman McInnes 1708b6490206SBarry Smith c->bs = a->bs; 17099df24120SSatish Balay c->bs2 = a->bs2; 1710b6490206SBarry Smith c->mbs = a->mbs; 1711b6490206SBarry Smith c->nbs = a->nbs; 17122593348eSBarry Smith 1713b6490206SBarry Smith c->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1714b6490206SBarry Smith c->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1715b6490206SBarry Smith for (i=0; i<mbs; i++) { 17162593348eSBarry Smith c->imax[i] = a->imax[i]; 17172593348eSBarry Smith c->ilen[i] = a->ilen[i]; 17182593348eSBarry Smith } 17192593348eSBarry Smith 17202593348eSBarry Smith /* allocate the matrix space */ 1721c4992f7dSBarry Smith c->singlemalloc = PETSC_TRUE; 17223f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 17233f1db9ecSBarry Smith c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a); 17247e67e3f9SSatish Balay c->j = (int*)(c->a + nz*bs2); 1725de6a44a3SBarry Smith c->i = c->j + nz; 1726549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1727b6490206SBarry Smith if (mbs > 0) { 1728549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 17292e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1730549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 17312e8a6d31SBarry Smith } else { 1732549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 17332593348eSBarry Smith } 17342593348eSBarry Smith } 17352593348eSBarry Smith 1736f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 17372593348eSBarry Smith c->sorted = a->sorted; 17382593348eSBarry Smith c->roworiented = a->roworiented; 17392593348eSBarry Smith c->nonew = a->nonew; 17402593348eSBarry Smith 17412593348eSBarry Smith if (a->diag) { 1742b6490206SBarry Smith c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag); 1743b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1744b6490206SBarry Smith for (i=0; i<mbs; i++) { 17452593348eSBarry Smith c->diag[i] = a->diag[i]; 17462593348eSBarry Smith } 174798305bb5SBarry Smith } else c->diag = 0; 17482593348eSBarry Smith c->nz = a->nz; 17492593348eSBarry Smith c->maxnz = a->maxnz; 17502593348eSBarry Smith c->solve_work = 0; 17512593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 17527fc0212eSBarry Smith c->mult_work = 0; 17532593348eSBarry Smith *B = C; 17547b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 17553a40ed3dSBarry Smith PetscFunctionReturn(0); 17562593348eSBarry Smith } 17572593348eSBarry Smith 17585615d1e5SSatish Balay #undef __FUNC__ 1759b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_SeqBAIJ" 176019bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 17612593348eSBarry Smith { 1762b6490206SBarry Smith Mat_SeqBAIJ *a; 17632593348eSBarry Smith Mat B; 1764f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1765b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 176635aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1767a7c10996SSatish Balay int *masked,nmask,tmp,bs2,ishift; 1768b6490206SBarry Smith Scalar *aa; 176919bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 17702593348eSBarry Smith 17713a40ed3dSBarry Smith PetscFunctionBegin; 1772f1af5d2fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1773de6a44a3SBarry Smith bs2 = bs*bs; 1774b6490206SBarry Smith 1775d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1776cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 177790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 17780752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1779a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 17802593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 17812593348eSBarry Smith 1782d64ed03dSBarry Smith if (header[3] < 0) { 1783a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1784d64ed03dSBarry Smith } 1785d64ed03dSBarry Smith 1786a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 178735aab85fSBarry Smith 178835aab85fSBarry Smith /* 178935aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 179035aab85fSBarry Smith divisible by the blocksize 179135aab85fSBarry Smith */ 1792b6490206SBarry Smith mbs = M/bs; 179335aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 179435aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 179535aab85fSBarry Smith else mbs++; 179635aab85fSBarry Smith if (extra_rows) { 1797537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 179835aab85fSBarry Smith } 1799b6490206SBarry Smith 18002593348eSBarry Smith /* read in row lengths */ 180135aab85fSBarry Smith rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 18020752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 180335aab85fSBarry Smith for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 18042593348eSBarry Smith 1805b6490206SBarry Smith /* read in column indices */ 180635aab85fSBarry Smith jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj); 18070752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 180835aab85fSBarry Smith for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 1809b6490206SBarry Smith 1810b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1811b6490206SBarry Smith browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1812549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 181335aab85fSBarry Smith mask = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask); 1814549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 181535aab85fSBarry Smith masked = mask + mbs; 1816b6490206SBarry Smith rowcount = 0; nzcount = 0; 1817b6490206SBarry Smith for (i=0; i<mbs; i++) { 181835aab85fSBarry Smith nmask = 0; 1819b6490206SBarry Smith for (j=0; j<bs; j++) { 1820b6490206SBarry Smith kmax = rowlengths[rowcount]; 1821b6490206SBarry Smith for (k=0; k<kmax; k++) { 182235aab85fSBarry Smith tmp = jj[nzcount++]/bs; 182335aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1824b6490206SBarry Smith } 1825b6490206SBarry Smith rowcount++; 1826b6490206SBarry Smith } 182735aab85fSBarry Smith browlengths[i] += nmask; 182835aab85fSBarry Smith /* zero out the mask elements we set */ 182935aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1830b6490206SBarry Smith } 1831b6490206SBarry Smith 18322593348eSBarry Smith /* create our matrix */ 18333eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 18342593348eSBarry Smith B = *A; 1835b6490206SBarry Smith a = (Mat_SeqBAIJ*)B->data; 18362593348eSBarry Smith 18372593348eSBarry Smith /* set matrix "i" values */ 1838de6a44a3SBarry Smith a->i[0] = 0; 1839b6490206SBarry Smith for (i=1; i<= mbs; i++) { 1840b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1841b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 18422593348eSBarry Smith } 1843b6490206SBarry Smith a->nz = 0; 1844b6490206SBarry Smith for (i=0; i<mbs; i++) a->nz += browlengths[i]; 18452593348eSBarry Smith 1846b6490206SBarry Smith /* read in nonzero values */ 184735aab85fSBarry Smith aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 18480752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 184935aab85fSBarry Smith for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 1850b6490206SBarry Smith 1851b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1852b6490206SBarry Smith nzcount = 0; jcount = 0; 1853b6490206SBarry Smith for (i=0; i<mbs; i++) { 1854b6490206SBarry Smith nzcountb = nzcount; 185535aab85fSBarry Smith nmask = 0; 1856b6490206SBarry Smith for (j=0; j<bs; j++) { 1857b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1858b6490206SBarry Smith for (k=0; k<kmax; k++) { 185935aab85fSBarry Smith tmp = jj[nzcount++]/bs; 186035aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1861b6490206SBarry Smith } 1862b6490206SBarry Smith rowcount++; 1863b6490206SBarry Smith } 1864de6a44a3SBarry Smith /* sort the masked values */ 1865433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1866de6a44a3SBarry Smith 1867b6490206SBarry Smith /* set "j" values into matrix */ 1868b6490206SBarry Smith maskcount = 1; 186935aab85fSBarry Smith for (j=0; j<nmask; j++) { 187035aab85fSBarry Smith a->j[jcount++] = masked[j]; 1871de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1872b6490206SBarry Smith } 1873b6490206SBarry Smith /* set "a" values into matrix */ 1874de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1875b6490206SBarry Smith for (j=0; j<bs; j++) { 1876b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1877b6490206SBarry Smith for (k=0; k<kmax; k++) { 1878de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1879de6a44a3SBarry Smith block = mask[tmp] - 1; 1880de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1881de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1882375fe846SBarry Smith a->a[idx] = (MatScalar)aa[nzcountb++]; 1883b6490206SBarry Smith } 1884b6490206SBarry Smith } 188535aab85fSBarry Smith /* zero out the mask elements we set */ 188635aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1887b6490206SBarry Smith } 1888a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1889b6490206SBarry Smith 1890606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1891606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1892606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1893606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1894606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1895b6490206SBarry Smith 1896b6490206SBarry Smith B->assembled = PETSC_TRUE; 1897de6a44a3SBarry Smith 18989c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 18993a40ed3dSBarry Smith PetscFunctionReturn(0); 19002593348eSBarry Smith } 19012593348eSBarry Smith 19022593348eSBarry Smith 19032593348eSBarry Smith 1904