1*35d8aa7fSBarry Smith /*$Id: baij.c,v 1.213 2000/10/24 20:25:52 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 */ 7e090d566SSatish Balay #include "petscsys.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 1295d5f7c2SBarry Smith /* UGLY, ugly, ugly 1395d5f7c2SBarry Smith When MatScalar == Scalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does 1495d5f7c2SBarry Smith not exist. Otherwise ..._MatScalar() takes matrix elements in single precision and 1595d5f7c2SBarry Smith inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ() 1695d5f7c2SBarry Smith converts the entries into single precision and then calls ..._MatScalar() to put them 1795d5f7c2SBarry Smith into the single precision data structures. 1895d5f7c2SBarry Smith */ 1995d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 20ca44d042SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode); 2195d5f7c2SBarry Smith #else 2295d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ 2395d5f7c2SBarry Smith #endif 2495d5f7c2SBarry 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) { 4229bbc08cSBarry Smith SETERRQ1(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 74ca44d042SBarry 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; 142aa482453SBarry Smith #if defined(PETSC_USE_LOG) 143273d9f13SBarry Smith PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz); 1442d61bbb3SSatish Balay #endif 145606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 146606d414cSSatish Balay if (!a->singlemalloc) { 147606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 148606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 149606d414cSSatish Balay } 150c38d4ed2SBarry Smith if (a->row) { 151c38d4ed2SBarry Smith ierr = ISDestroy(a->row);CHKERRQ(ierr); 152c38d4ed2SBarry Smith } 153c38d4ed2SBarry Smith if (a->col) { 154c38d4ed2SBarry Smith ierr = ISDestroy(a->col);CHKERRQ(ierr); 155c38d4ed2SBarry Smith } 156606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 157606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 158606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 159606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 160606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 161e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 162606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 163563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 164563b5814SBarry Smith if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);} 165563b5814SBarry Smith #endif 166606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1672d61bbb3SSatish Balay PetscFunctionReturn(0); 1682d61bbb3SSatish Balay } 1692d61bbb3SSatish Balay 1702d61bbb3SSatish Balay #undef __FUNC__ 171b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetOption_SeqBAIJ" 1722d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1732d61bbb3SSatish Balay { 1742d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1752d61bbb3SSatish Balay 1762d61bbb3SSatish Balay PetscFunctionBegin; 177c4992f7dSBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = PETSC_TRUE; 178c4992f7dSBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = PETSC_FALSE; 179c4992f7dSBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = PETSC_TRUE; 180c4992f7dSBarry Smith else if (op == MAT_COLUMNS_UNSORTED) a->sorted = PETSC_FALSE; 181c4992f7dSBarry Smith else if (op == MAT_KEEP_ZEROED_ROWS) a->keepzeroedrows = PETSC_TRUE; 1822d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1834787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1844787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1852d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1862d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1872d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1882d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1892d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1902d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 191b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 192b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 193981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1942d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 19529bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 1962d61bbb3SSatish Balay } else { 19729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 1982d61bbb3SSatish Balay } 1992d61bbb3SSatish Balay PetscFunctionReturn(0); 2002d61bbb3SSatish Balay } 2012d61bbb3SSatish Balay 2022d61bbb3SSatish Balay #undef __FUNC__ 203b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetOwnershipRange_SeqBAIJ" 2042d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2052d61bbb3SSatish Balay { 2062d61bbb3SSatish Balay PetscFunctionBegin; 2074c49b128SBarry Smith if (m) *m = 0; 208273d9f13SBarry Smith if (n) *n = A->m; 2092d61bbb3SSatish Balay PetscFunctionReturn(0); 2102d61bbb3SSatish Balay } 2112d61bbb3SSatish Balay 2122d61bbb3SSatish Balay #undef __FUNC__ 213b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetRow_SeqBAIJ" 2142d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2152d61bbb3SSatish Balay { 2162d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2172d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2183f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2193f1db9ecSBarry Smith Scalar *v_i; 2202d61bbb3SSatish Balay 2212d61bbb3SSatish Balay PetscFunctionBegin; 2222d61bbb3SSatish Balay bs = a->bs; 2232d61bbb3SSatish Balay ai = a->i; 2242d61bbb3SSatish Balay aj = a->j; 2252d61bbb3SSatish Balay aa = a->a; 2262d61bbb3SSatish Balay bs2 = a->bs2; 2272d61bbb3SSatish Balay 228273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 2292d61bbb3SSatish Balay 2302d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2312d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2322d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2332d61bbb3SSatish Balay *nz = bs*M; 2342d61bbb3SSatish Balay 2352d61bbb3SSatish Balay if (v) { 2362d61bbb3SSatish Balay *v = 0; 2372d61bbb3SSatish Balay if (*nz) { 2382d61bbb3SSatish Balay *v = (Scalar*)PetscMalloc((*nz)*sizeof(Scalar));CHKPTRQ(*v); 2392d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 2402d61bbb3SSatish Balay v_i = *v + i*bs; 2412d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2422d61bbb3SSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 2432d61bbb3SSatish Balay } 2442d61bbb3SSatish Balay } 2452d61bbb3SSatish Balay } 2462d61bbb3SSatish Balay 2472d61bbb3SSatish Balay if (idx) { 2482d61bbb3SSatish Balay *idx = 0; 2492d61bbb3SSatish Balay if (*nz) { 2502d61bbb3SSatish Balay *idx = (int*)PetscMalloc((*nz)*sizeof(int));CHKPTRQ(*idx); 2512d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 2522d61bbb3SSatish Balay idx_i = *idx + i*bs; 2532d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2542d61bbb3SSatish Balay for (j=0; j<bs; j++) {idx_i[j] = itmp++;} 2552d61bbb3SSatish Balay } 2562d61bbb3SSatish Balay } 2572d61bbb3SSatish Balay } 2582d61bbb3SSatish Balay PetscFunctionReturn(0); 2592d61bbb3SSatish Balay } 2602d61bbb3SSatish Balay 2612d61bbb3SSatish Balay #undef __FUNC__ 262b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRestoreRow_SeqBAIJ" 2632d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2642d61bbb3SSatish Balay { 265606d414cSSatish Balay int ierr; 266606d414cSSatish Balay 2672d61bbb3SSatish Balay PetscFunctionBegin; 268606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 269606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2702d61bbb3SSatish Balay PetscFunctionReturn(0); 2712d61bbb3SSatish Balay } 2722d61bbb3SSatish Balay 2732d61bbb3SSatish Balay #undef __FUNC__ 274b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatTranspose_SeqBAIJ" 2752d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2762d61bbb3SSatish Balay { 2772d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2782d61bbb3SSatish Balay Mat C; 2792d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 280273d9f13SBarry Smith int *rows,*cols,bs2=a->bs2; 281375fe846SBarry Smith Scalar *array; 2822d61bbb3SSatish Balay 2832d61bbb3SSatish Balay PetscFunctionBegin; 28429bbc08cSBarry Smith if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place"); 2852d61bbb3SSatish Balay col = (int*)PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 286549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2872d61bbb3SSatish Balay 288375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 289375fe846SBarry Smith array = (Scalar *)PetscMalloc(a->bs2*a->nz*sizeof(Scalar));CHKPTRQ(array); 290375fe846SBarry Smith for (i=0; i<a->bs2*a->nz; i++) array[i] = (Scalar)a->a[i]; 291375fe846SBarry Smith #else 2923eda8832SBarry Smith array = a->a; 293375fe846SBarry Smith #endif 294375fe846SBarry Smith 2952d61bbb3SSatish Balay for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 296273d9f13SBarry Smith ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 297606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 2982d61bbb3SSatish Balay rows = (int*)PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 2992d61bbb3SSatish Balay cols = rows + bs; 3002d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3012d61bbb3SSatish Balay cols[0] = i*bs; 3022d61bbb3SSatish Balay for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 3032d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3042d61bbb3SSatish Balay for (j=0; j<len; j++) { 3052d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3062d61bbb3SSatish Balay for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 3072d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3082d61bbb3SSatish Balay array += bs2; 3092d61bbb3SSatish Balay } 3102d61bbb3SSatish Balay } 311606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 312375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 313375fe846SBarry Smith ierr = PetscFree(array); 314375fe846SBarry Smith #endif 3152d61bbb3SSatish Balay 3162d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3172d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3182d61bbb3SSatish Balay 319c4992f7dSBarry Smith if (B) { 3202d61bbb3SSatish Balay *B = C; 3212d61bbb3SSatish Balay } else { 322273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 3232d61bbb3SSatish Balay } 3242d61bbb3SSatish Balay PetscFunctionReturn(0); 3252d61bbb3SSatish Balay } 3262d61bbb3SSatish Balay 3275615d1e5SSatish Balay #undef __FUNC__ 328b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Binary" 329b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3302593348eSBarry Smith { 331b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 3329df24120SSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 333b6490206SBarry Smith Scalar *aa; 334ce6f0cecSBarry Smith FILE *file; 3352593348eSBarry Smith 3363a40ed3dSBarry Smith PetscFunctionBegin; 33790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 338273d9f13SBarry Smith col_lens = (int*)PetscMalloc((4+A->m)*sizeof(int));CHKPTRQ(col_lens); 3392593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3403b2fbd54SBarry Smith 341273d9f13SBarry Smith col_lens[1] = A->m; 342273d9f13SBarry Smith col_lens[2] = A->n; 3437e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3442593348eSBarry Smith 3452593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 346b6490206SBarry Smith count = 0; 347b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 348b6490206SBarry Smith for (j=0; j<bs; j++) { 349b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 350b6490206SBarry Smith } 3512593348eSBarry Smith } 352273d9f13SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 353606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3542593348eSBarry Smith 3552593348eSBarry Smith /* store column indices (zero start index) */ 35666963ce1SSatish Balay jj = (int*)PetscMalloc((a->nz+1)*bs2*sizeof(int));CHKPTRQ(jj); 357b6490206SBarry Smith count = 0; 358b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 359b6490206SBarry Smith for (j=0; j<bs; j++) { 360b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 361b6490206SBarry Smith for (l=0; l<bs; l++) { 362b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3632593348eSBarry Smith } 3642593348eSBarry Smith } 365b6490206SBarry Smith } 366b6490206SBarry Smith } 3670752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 368606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 3692593348eSBarry Smith 3702593348eSBarry Smith /* store nonzero values */ 37166963ce1SSatish Balay aa = (Scalar*)PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 372b6490206SBarry Smith count = 0; 373b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 374b6490206SBarry Smith for (j=0; j<bs; j++) { 375b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 376b6490206SBarry Smith for (l=0; l<bs; l++) { 3777e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 378b6490206SBarry Smith } 379b6490206SBarry Smith } 380b6490206SBarry Smith } 381b6490206SBarry Smith } 3820752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 383606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 384ce6f0cecSBarry Smith 385ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 386ce6f0cecSBarry Smith if (file) { 387ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 388ce6f0cecSBarry Smith } 3893a40ed3dSBarry Smith PetscFunctionReturn(0); 3902593348eSBarry Smith } 3912593348eSBarry Smith 3925615d1e5SSatish Balay #undef __FUNC__ 393b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_ASCII" 394b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 3952593348eSBarry Smith { 396b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 3979df24120SSatish Balay int ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2; 3982593348eSBarry Smith char *outputname; 3992593348eSBarry Smith 4003a40ed3dSBarry Smith PetscFunctionBegin; 40177ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 402888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 403639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 404d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4050ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 40629bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 4070ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 4086831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 40944cd7ae7SLois Curfman McInnes for (i=0; i<a->mbs; i++) { 41044cd7ae7SLois Curfman McInnes for (j=0; j<bs; j++) { 4116831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 41244cd7ae7SLois Curfman McInnes for (k=a->i[i]; k<a->i[i+1]; k++) { 41344cd7ae7SLois Curfman McInnes for (l=0; l<bs; l++) { 414aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4150e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 4166831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l, 4170e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4180e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 4196831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l, 4200e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4210e6d2581SBarry Smith } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 4220e6d2581SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4230ef38995SBarry Smith } 42444cd7ae7SLois Curfman McInnes #else 4250ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 4266831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 4270ef38995SBarry Smith } 42844cd7ae7SLois Curfman McInnes #endif 42944cd7ae7SLois Curfman McInnes } 43044cd7ae7SLois Curfman McInnes } 4316831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 43244cd7ae7SLois Curfman McInnes } 43344cd7ae7SLois Curfman McInnes } 4346831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 4350ef38995SBarry Smith } else { 4366831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 437b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 438b6490206SBarry Smith for (j=0; j<bs; j++) { 4396831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 440b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 441b6490206SBarry Smith for (l=0; l<bs; l++) { 442aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4430e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 4446831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 4450e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4460e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 4476831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 4480e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4490ef38995SBarry Smith } else { 4500e6d2581SBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 45188685aaeSLois Curfman McInnes } 45288685aaeSLois Curfman McInnes #else 4536831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 45488685aaeSLois Curfman McInnes #endif 4552593348eSBarry Smith } 4562593348eSBarry Smith } 4576831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 4582593348eSBarry Smith } 4592593348eSBarry Smith } 4606831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 461b6490206SBarry Smith } 4626831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 4633a40ed3dSBarry Smith PetscFunctionReturn(0); 4642593348eSBarry Smith } 4652593348eSBarry Smith 4665615d1e5SSatish Balay #undef __FUNC__ 467b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Draw_Zoom" 46877ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4693270192aSSatish Balay { 47077ed5343SBarry Smith Mat A = (Mat) Aa; 4713270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 472b65db4caSBarry Smith int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2; 4730e6d2581SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 4743f1db9ecSBarry Smith MatScalar *aa; 47577ed5343SBarry Smith Viewer viewer; 4763270192aSSatish Balay 4773a40ed3dSBarry Smith PetscFunctionBegin; 4783270192aSSatish Balay 479b65db4caSBarry Smith /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/ 48077ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 48177ed5343SBarry Smith 48277ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 48377ed5343SBarry Smith 4843270192aSSatish Balay /* loop over matrix elements drawing boxes */ 4853270192aSSatish Balay color = DRAW_BLUE; 4863270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 4873270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 488273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 4893270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 4903270192aSSatish Balay aa = a->a + j*bs2; 4913270192aSSatish Balay for (k=0; k<bs; k++) { 4923270192aSSatish Balay for (l=0; l<bs; l++) { 4930e6d2581SBarry Smith if (PetscRealPart(*aa++) >= 0.) continue; 494433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 4953270192aSSatish Balay } 4963270192aSSatish Balay } 4973270192aSSatish Balay } 4983270192aSSatish Balay } 4993270192aSSatish Balay color = DRAW_CYAN; 5003270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5013270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 502273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 5033270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5043270192aSSatish Balay aa = a->a + j*bs2; 5053270192aSSatish Balay for (k=0; k<bs; k++) { 5063270192aSSatish Balay for (l=0; l<bs; l++) { 5070e6d2581SBarry Smith if (PetscRealPart(*aa++) != 0.) continue; 508433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5093270192aSSatish Balay } 5103270192aSSatish Balay } 5113270192aSSatish Balay } 5123270192aSSatish Balay } 5133270192aSSatish Balay 5143270192aSSatish Balay color = DRAW_RED; 5153270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5163270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 517273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 5183270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5193270192aSSatish Balay aa = a->a + j*bs2; 5203270192aSSatish Balay for (k=0; k<bs; k++) { 5213270192aSSatish Balay for (l=0; l<bs; l++) { 5220e6d2581SBarry Smith if (PetscRealPart(*aa++) <= 0.) continue; 523433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5243270192aSSatish Balay } 5253270192aSSatish Balay } 5263270192aSSatish Balay } 5273270192aSSatish Balay } 52877ed5343SBarry Smith PetscFunctionReturn(0); 52977ed5343SBarry Smith } 5303270192aSSatish Balay 53177ed5343SBarry Smith #undef __FUNC__ 532b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ_Draw" 53377ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 53477ed5343SBarry Smith { 5357dce120fSSatish Balay int ierr; 5360e6d2581SBarry Smith PetscReal xl,yl,xr,yr,w,h; 53777ed5343SBarry Smith Draw draw; 53877ed5343SBarry Smith PetscTruth isnull; 5393270192aSSatish Balay 54077ed5343SBarry Smith PetscFunctionBegin; 54177ed5343SBarry Smith 54277ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 54377ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 54477ed5343SBarry Smith 54577ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 546273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 54777ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5483270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 54977ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 55077ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5513a40ed3dSBarry Smith PetscFunctionReturn(0); 5523270192aSSatish Balay } 5533270192aSSatish Balay 5545615d1e5SSatish Balay #undef __FUNC__ 555b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatView_SeqBAIJ" 556e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5572593348eSBarry Smith { 55819bcc07fSBarry Smith int ierr; 5596831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 5602593348eSBarry Smith 5613a40ed3dSBarry Smith PetscFunctionBegin; 5626831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 5636831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 5646831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 5656831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 5660f5bd95cSBarry Smith if (issocket) { 56729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 5680f5bd95cSBarry Smith } else if (isascii){ 5693a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 5700f5bd95cSBarry Smith } else if (isbinary) { 5713a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 5720f5bd95cSBarry Smith } else if (isdraw) { 5733a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 5745cd90555SBarry Smith } else { 57529bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 5762593348eSBarry Smith } 5773a40ed3dSBarry Smith PetscFunctionReturn(0); 5782593348eSBarry Smith } 579b6490206SBarry Smith 580cd0e1443SSatish Balay 5815615d1e5SSatish Balay #undef __FUNC__ 582b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatGetValues_SeqBAIJ" 5832d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 584cd0e1443SSatish Balay { 585cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 5862d61bbb3SSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 5872d61bbb3SSatish Balay int *ai = a->i,*ailen = a->ilen; 5882d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 5893f1db9ecSBarry Smith MatScalar *ap,*aa = a->a,zero = 0.0; 590cd0e1443SSatish Balay 5913a40ed3dSBarry Smith PetscFunctionBegin; 5922d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 593cd0e1443SSatish Balay row = im[k]; brow = row/bs; 59429bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 595273d9f13SBarry Smith if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 5962d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 5972c3acbe9SBarry Smith nrow = ailen[brow]; 5982d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 59929bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 600273d9f13SBarry Smith if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 6012d61bbb3SSatish Balay col = in[l] ; 6022d61bbb3SSatish Balay bcol = col/bs; 6032d61bbb3SSatish Balay cidx = col%bs; 6042d61bbb3SSatish Balay ridx = row%bs; 6052d61bbb3SSatish Balay high = nrow; 6062d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6072d61bbb3SSatish Balay while (high-low > 5) { 608cd0e1443SSatish Balay t = (low+high)/2; 609cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 610cd0e1443SSatish Balay else low = t; 611cd0e1443SSatish Balay } 612cd0e1443SSatish Balay for (i=low; i<high; i++) { 613cd0e1443SSatish Balay if (rp[i] > bcol) break; 614cd0e1443SSatish Balay if (rp[i] == bcol) { 6152d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6162d61bbb3SSatish Balay goto finished; 617cd0e1443SSatish Balay } 618cd0e1443SSatish Balay } 6192d61bbb3SSatish Balay *v++ = zero; 6202d61bbb3SSatish Balay finished:; 621cd0e1443SSatish Balay } 622cd0e1443SSatish Balay } 6233a40ed3dSBarry Smith PetscFunctionReturn(0); 624cd0e1443SSatish Balay } 625cd0e1443SSatish Balay 62695d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 62795d5f7c2SBarry Smith #undef __FUNC__ 62895d5f7c2SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValuesBlocked_SeqBAIJ" 62995d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 63095d5f7c2SBarry Smith { 631563b5814SBarry Smith Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data; 632563b5814SBarry Smith int ierr,i,N = m*n*b->bs2; 633563b5814SBarry Smith MatScalar *vsingle; 63495d5f7c2SBarry Smith 63595d5f7c2SBarry Smith PetscFunctionBegin; 636563b5814SBarry Smith if (N > b->setvalueslen) { 637563b5814SBarry Smith if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);} 638563b5814SBarry Smith b->setvaluescopy = (MatScalar*)PetscMalloc(N*sizeof(MatScalar));CHKPTRQ(b->setvaluescopy); 639563b5814SBarry Smith b->setvalueslen = N; 640563b5814SBarry Smith } 641563b5814SBarry Smith vsingle = b->setvaluescopy; 64295d5f7c2SBarry Smith for (i=0; i<N; i++) { 64395d5f7c2SBarry Smith vsingle[i] = v[i]; 64495d5f7c2SBarry Smith } 64595d5f7c2SBarry Smith ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr); 64695d5f7c2SBarry Smith PetscFunctionReturn(0); 64795d5f7c2SBarry Smith } 64895d5f7c2SBarry Smith #endif 64995d5f7c2SBarry Smith 6502d61bbb3SSatish Balay 6515615d1e5SSatish Balay #undef __FUNC__ 652b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValuesBlocked_SeqBAIJ" 65395d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is) 65492c4ed94SBarry Smith { 65592c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6568a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 657273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 658549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 659273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 66095d5f7c2SBarry Smith MatScalar *value = v,*ap,*aa = a->a,*bap; 66192c4ed94SBarry Smith 6623a40ed3dSBarry Smith PetscFunctionBegin; 6630e324ae4SSatish Balay if (roworiented) { 6640e324ae4SSatish Balay stepval = (n-1)*bs; 6650e324ae4SSatish Balay } else { 6660e324ae4SSatish Balay stepval = (m-1)*bs; 6670e324ae4SSatish Balay } 66892c4ed94SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 66992c4ed94SBarry Smith row = im[k]; 6705ef9f2a5SBarry Smith if (row < 0) continue; 671aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 67229bbc08cSBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 67392c4ed94SBarry Smith #endif 67492c4ed94SBarry Smith rp = aj + ai[row]; 67592c4ed94SBarry Smith ap = aa + bs2*ai[row]; 67692c4ed94SBarry Smith rmax = imax[row]; 67792c4ed94SBarry Smith nrow = ailen[row]; 67892c4ed94SBarry Smith low = 0; 67992c4ed94SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 6805ef9f2a5SBarry Smith if (in[l] < 0) continue; 681aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 68229bbc08cSBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 68392c4ed94SBarry Smith #endif 68492c4ed94SBarry Smith col = in[l]; 68592c4ed94SBarry Smith if (roworiented) { 6860e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6870e324ae4SSatish Balay } else { 6880e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 68992c4ed94SBarry Smith } 69092c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 69192c4ed94SBarry Smith while (high-low > 7) { 69292c4ed94SBarry Smith t = (low+high)/2; 69392c4ed94SBarry Smith if (rp[t] > col) high = t; 69492c4ed94SBarry Smith else low = t; 69592c4ed94SBarry Smith } 69692c4ed94SBarry Smith for (i=low; i<high; i++) { 69792c4ed94SBarry Smith if (rp[i] > col) break; 69892c4ed94SBarry Smith if (rp[i] == col) { 6998a84c255SSatish Balay bap = ap + bs2*i; 7000e324ae4SSatish Balay if (roworiented) { 7018a84c255SSatish Balay if (is == ADD_VALUES) { 702dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 703dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7048a84c255SSatish Balay bap[jj] += *value++; 705dd9472c6SBarry Smith } 706dd9472c6SBarry Smith } 7070e324ae4SSatish Balay } else { 708dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 709dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7100e324ae4SSatish Balay bap[jj] = *value++; 7118a84c255SSatish Balay } 712dd9472c6SBarry Smith } 713dd9472c6SBarry Smith } 7140e324ae4SSatish Balay } else { 7150e324ae4SSatish Balay if (is == ADD_VALUES) { 716dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 717dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7180e324ae4SSatish Balay *bap++ += *value++; 719dd9472c6SBarry Smith } 720dd9472c6SBarry Smith } 7210e324ae4SSatish Balay } else { 722dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 723dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7240e324ae4SSatish Balay *bap++ = *value++; 7250e324ae4SSatish Balay } 7268a84c255SSatish Balay } 727dd9472c6SBarry Smith } 728dd9472c6SBarry Smith } 729f1241b54SBarry Smith goto noinsert2; 73092c4ed94SBarry Smith } 73192c4ed94SBarry Smith } 73289280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 73329bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 73492c4ed94SBarry Smith if (nrow >= rmax) { 73592c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 73692c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7373f1db9ecSBarry Smith MatScalar *new_a; 73892c4ed94SBarry Smith 73929bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 74089280ab3SLois Curfman McInnes 74192c4ed94SBarry Smith /* malloc new storage space */ 7423f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7433f1db9ecSBarry Smith new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 74492c4ed94SBarry Smith new_j = (int*)(new_a + bs2*new_nz); 74592c4ed94SBarry Smith new_i = new_j + new_nz; 74692c4ed94SBarry Smith 74792c4ed94SBarry Smith /* copy over old data into new slots */ 74892c4ed94SBarry Smith for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 74992c4ed94SBarry Smith for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 750549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 75192c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 752549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 753549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 754549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 755549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 75692c4ed94SBarry Smith /* free up old matrix storage */ 757606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 758606d414cSSatish Balay if (!a->singlemalloc) { 759606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 760606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 761606d414cSSatish Balay } 76292c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 763c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 76492c4ed94SBarry Smith 76592c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 76692c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7673f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 76892c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 76992c4ed94SBarry Smith a->reallocs++; 77092c4ed94SBarry Smith a->nz++; 77192c4ed94SBarry Smith } 77292c4ed94SBarry Smith N = nrow++ - 1; 77392c4ed94SBarry Smith /* shift up all the later entries in this row */ 77492c4ed94SBarry Smith for (ii=N; ii>=i; ii--) { 77592c4ed94SBarry Smith rp[ii+1] = rp[ii]; 776549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 77792c4ed94SBarry Smith } 778549d3d68SSatish Balay if (N >= i) { 779549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 780549d3d68SSatish Balay } 78192c4ed94SBarry Smith rp[i] = col; 7828a84c255SSatish Balay bap = ap + bs2*i; 7830e324ae4SSatish Balay if (roworiented) { 784dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 785dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7860e324ae4SSatish Balay bap[jj] = *value++; 787dd9472c6SBarry Smith } 788dd9472c6SBarry Smith } 7890e324ae4SSatish Balay } else { 790dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 791dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7920e324ae4SSatish Balay *bap++ = *value++; 7930e324ae4SSatish Balay } 794dd9472c6SBarry Smith } 795dd9472c6SBarry Smith } 796f1241b54SBarry Smith noinsert2:; 79792c4ed94SBarry Smith low = i; 79892c4ed94SBarry Smith } 79992c4ed94SBarry Smith ailen[row] = nrow; 80092c4ed94SBarry Smith } 8013a40ed3dSBarry Smith PetscFunctionReturn(0); 80292c4ed94SBarry Smith } 80392c4ed94SBarry Smith 804f2501298SSatish Balay 8055615d1e5SSatish Balay #undef __FUNC__ 806b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatAssemblyEnd_SeqBAIJ" 8078f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 808584200bdSSatish Balay { 809584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 810584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 811273d9f13SBarry Smith int m = A->m,*ip,N,*ailen = a->ilen; 812549d3d68SSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 8133f1db9ecSBarry Smith MatScalar *aa = a->a,*ap; 814584200bdSSatish Balay 8153a40ed3dSBarry Smith PetscFunctionBegin; 8163a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 817584200bdSSatish Balay 81843ee02c3SBarry Smith if (m) rmax = ailen[0]; 819584200bdSSatish Balay for (i=1; i<mbs; i++) { 820584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 821584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 822d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 823584200bdSSatish Balay if (fshift) { 824a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 825584200bdSSatish Balay N = ailen[i]; 826584200bdSSatish Balay for (j=0; j<N; j++) { 827584200bdSSatish Balay ip[j-fshift] = ip[j]; 828549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 829584200bdSSatish Balay } 830584200bdSSatish Balay } 831584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 832584200bdSSatish Balay } 833584200bdSSatish Balay if (mbs) { 834584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 835584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 836584200bdSSatish Balay } 837584200bdSSatish Balay /* reset ilen and imax for each row */ 838584200bdSSatish Balay for (i=0; i<mbs; i++) { 839584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 840584200bdSSatish Balay } 841a7c10996SSatish Balay a->nz = ai[mbs]; 842584200bdSSatish Balay 843584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 844584200bdSSatish Balay if (fshift && a->diag) { 845606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 846584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 847584200bdSSatish Balay a->diag = 0; 848584200bdSSatish Balay } 8493d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 850273d9f13SBarry Smith m,A->n,a->bs,fshift*bs2,a->nz*bs2); 8513d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 852584200bdSSatish Balay a->reallocs); 853d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 854e2f3b5e9SSatish Balay a->reallocs = 0; 8550e6d2581SBarry Smith A->info.nz_unneeded = (PetscReal)fshift*bs2; 8564e220ebcSLois Curfman McInnes 8573a40ed3dSBarry Smith PetscFunctionReturn(0); 858584200bdSSatish Balay } 859584200bdSSatish Balay 8605a838052SSatish Balay 861bea157c4SSatish Balay 862bea157c4SSatish Balay /* 863bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 864bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 865bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 866bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 867bea157c4SSatish Balay */ 8685615d1e5SSatish Balay #undef __FUNC__ 869b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqBAIJ_Check_Blocks" 870bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 871d9b7c43dSSatish Balay { 872bea157c4SSatish Balay int i,j,k,row; 873bea157c4SSatish Balay PetscTruth flg; 8743a40ed3dSBarry Smith 875433994e6SBarry Smith PetscFunctionBegin; 876bea157c4SSatish Balay for (i=0,j=0; i<n; j++) { 877bea157c4SSatish Balay row = idx[i]; 878bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 879bea157c4SSatish Balay sizes[j] = 1; 880bea157c4SSatish Balay i++; 881e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 882bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 883bea157c4SSatish Balay i++; 884bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 885bea157c4SSatish Balay flg = PETSC_TRUE; 886bea157c4SSatish Balay for (k=1; k<bs; k++) { 887bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 888bea157c4SSatish Balay flg = PETSC_FALSE; 889bea157c4SSatish Balay break; 890d9b7c43dSSatish Balay } 891bea157c4SSatish Balay } 892bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 893bea157c4SSatish Balay sizes[j] = bs; 894bea157c4SSatish Balay i+= bs; 895bea157c4SSatish Balay } else { 896bea157c4SSatish Balay sizes[j] = 1; 897bea157c4SSatish Balay i++; 898bea157c4SSatish Balay } 899bea157c4SSatish Balay } 900bea157c4SSatish Balay } 901bea157c4SSatish Balay *bs_max = j; 9023a40ed3dSBarry Smith PetscFunctionReturn(0); 903d9b7c43dSSatish Balay } 904d9b7c43dSSatish Balay 9055615d1e5SSatish Balay #undef __FUNC__ 906b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatZeroRows_SeqBAIJ" 9078f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is,Scalar *diag) 908d9b7c43dSSatish Balay { 909d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 910b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 911bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9123f1db9ecSBarry Smith Scalar zero = 0.0; 9133f1db9ecSBarry Smith MatScalar *aa; 914d9b7c43dSSatish Balay 9153a40ed3dSBarry Smith PetscFunctionBegin; 916d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 917b9b97703SBarry Smith ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr); 918d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 919d9b7c43dSSatish Balay 920bea157c4SSatish Balay /* allocate memory for rows,sizes */ 921bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 922bea157c4SSatish Balay sizes = rows + is_n; 923bea157c4SSatish Balay 924563b5814SBarry Smith /* copy IS values to rows, and sort them */ 925bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 926bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 927dffd3267SBarry Smith if (baij->keepzeroedrows) { 928dffd3267SBarry Smith for (i=0; i<is_n; i++) { sizes[i] = 1; } 929dffd3267SBarry Smith bs_max = is_n; 930dffd3267SBarry Smith } else { 931bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 932dffd3267SBarry Smith } 933b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 934bea157c4SSatish Balay 935bea157c4SSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 936bea157c4SSatish Balay row = rows[j]; 937273d9f13SBarry Smith if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 938bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 939bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 940dffd3267SBarry Smith if (sizes[i] == bs && !baij->keepzeroedrows) { 941bea157c4SSatish Balay if (diag) { 942bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 943bea157c4SSatish Balay baij->ilen[row/bs] = 1; 944bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 945bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 946a07cd24cSSatish Balay } 947563b5814SBarry Smith /* Now insert all the diagonal values for this bs */ 948bea157c4SSatish Balay for (k=0; k<bs; k++) { 949bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 950bea157c4SSatish Balay } 951bea157c4SSatish Balay } else { /* (!diag) */ 952bea157c4SSatish Balay baij->ilen[row/bs] = 0; 953bea157c4SSatish Balay } /* end (!diag) */ 954bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 955aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 95629bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 957bea157c4SSatish Balay #endif 958bea157c4SSatish Balay for (k=0; k<count; k++) { 959d9b7c43dSSatish Balay aa[0] = zero; 960d9b7c43dSSatish Balay aa += bs; 961d9b7c43dSSatish Balay } 962d9b7c43dSSatish Balay if (diag) { 963f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 964d9b7c43dSSatish Balay } 965d9b7c43dSSatish Balay } 966bea157c4SSatish Balay } 967bea157c4SSatish Balay 968606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 9699a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9703a40ed3dSBarry Smith PetscFunctionReturn(0); 971d9b7c43dSSatish Balay } 9721c351548SSatish Balay 9735615d1e5SSatish Balay #undef __FUNC__ 974b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSetValues_SeqBAIJ" 9752d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9762d61bbb3SSatish Balay { 9772d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 9782d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 979273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 9802d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 981549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 982273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 9833f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9842d61bbb3SSatish Balay 9852d61bbb3SSatish Balay PetscFunctionBegin; 9862d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 9872d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9885ef9f2a5SBarry Smith if (row < 0) continue; 989aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 990273d9f13SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 9912d61bbb3SSatish Balay #endif 9922d61bbb3SSatish Balay rp = aj + ai[brow]; 9932d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9942d61bbb3SSatish Balay rmax = imax[brow]; 9952d61bbb3SSatish Balay nrow = ailen[brow]; 9962d61bbb3SSatish Balay low = 0; 9972d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 9985ef9f2a5SBarry Smith if (in[l] < 0) continue; 999aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1000273d9f13SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n); 10012d61bbb3SSatish Balay #endif 10022d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10032d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10042d61bbb3SSatish Balay if (roworiented) { 10055ef9f2a5SBarry Smith value = v[l + k*n]; 10062d61bbb3SSatish Balay } else { 10072d61bbb3SSatish Balay value = v[k + l*m]; 10082d61bbb3SSatish Balay } 10092d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10102d61bbb3SSatish Balay while (high-low > 7) { 10112d61bbb3SSatish Balay t = (low+high)/2; 10122d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10132d61bbb3SSatish Balay else low = t; 10142d61bbb3SSatish Balay } 10152d61bbb3SSatish Balay for (i=low; i<high; i++) { 10162d61bbb3SSatish Balay if (rp[i] > bcol) break; 10172d61bbb3SSatish Balay if (rp[i] == bcol) { 10182d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10192d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10202d61bbb3SSatish Balay else *bap = value; 10212d61bbb3SSatish Balay goto noinsert1; 10222d61bbb3SSatish Balay } 10232d61bbb3SSatish Balay } 10242d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 102529bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 10262d61bbb3SSatish Balay if (nrow >= rmax) { 10272d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10282d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10293f1db9ecSBarry Smith MatScalar *new_a; 10302d61bbb3SSatish Balay 103129bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 10322d61bbb3SSatish Balay 10332d61bbb3SSatish Balay /* Malloc new storage space */ 10343f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10353f1db9ecSBarry Smith new_a = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a); 10362d61bbb3SSatish Balay new_j = (int*)(new_a + bs2*new_nz); 10372d61bbb3SSatish Balay new_i = new_j + new_nz; 10382d61bbb3SSatish Balay 10392d61bbb3SSatish Balay /* copy over old data into new slots */ 10402d61bbb3SSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 10412d61bbb3SSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1042549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10432d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1044549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1045549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1046549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1047549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10482d61bbb3SSatish Balay /* free up old matrix storage */ 1049606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1050606d414cSSatish Balay if (!a->singlemalloc) { 1051606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1052606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1053606d414cSSatish Balay } 10542d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1055c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 10562d61bbb3SSatish Balay 10572d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10582d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10593f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10602d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10612d61bbb3SSatish Balay a->reallocs++; 10622d61bbb3SSatish Balay a->nz++; 10632d61bbb3SSatish Balay } 10642d61bbb3SSatish Balay N = nrow++ - 1; 10652d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10662d61bbb3SSatish Balay for (ii=N; ii>=i; ii--) { 10672d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1068549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10692d61bbb3SSatish Balay } 1070549d3d68SSatish Balay if (N>=i) { 1071549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1072549d3d68SSatish Balay } 10732d61bbb3SSatish Balay rp[i] = bcol; 10742d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10752d61bbb3SSatish Balay noinsert1:; 10762d61bbb3SSatish Balay low = i; 10772d61bbb3SSatish Balay } 10782d61bbb3SSatish Balay ailen[brow] = nrow; 10792d61bbb3SSatish Balay } 10802d61bbb3SSatish Balay PetscFunctionReturn(0); 10812d61bbb3SSatish Balay } 10822d61bbb3SSatish Balay 1083b9b97703SBarry Smith EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatLUInfo*,Mat*); 1084b9b97703SBarry Smith EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatLUInfo*); 1085ca44d042SBarry Smith EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1086ca44d042SBarry Smith EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 1087ca44d042SBarry Smith EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 1088ca44d042SBarry Smith EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec); 1089ca44d042SBarry Smith EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 1090ca44d042SBarry Smith EXTERN int MatScale_SeqBAIJ(Scalar*,Mat); 1091ca44d042SBarry Smith EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *); 1092ca44d042SBarry Smith EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*); 1093ca44d042SBarry Smith EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec); 1094ca44d042SBarry Smith EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 1095ca44d042SBarry Smith EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 1096ca44d042SBarry Smith EXTERN int MatZeroEntries_SeqBAIJ(Mat); 10972d61bbb3SSatish Balay 1098ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 1099ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 1100ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 1101ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 1102ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 1103ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 1104ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 1105ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 1106ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec); 1107ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec); 1108ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec); 1109ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec); 1110ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec); 1111ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec); 1112ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec); 11132d61bbb3SSatish Balay 1114ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 1115ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 1116ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 1117ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 1118ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 1119ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 1120ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11212d61bbb3SSatish Balay 1122ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 1123ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 1124ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 1125ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 1126ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 1127ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 1128ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 1129ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11302d61bbb3SSatish Balay 1131ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 1132ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 1133ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 1134ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 1135ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 1136ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 1137ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 1138ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11392d61bbb3SSatish Balay 11402d61bbb3SSatish Balay #undef __FUNC__ 1141b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatILUFactor_SeqBAIJ" 11425ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11432d61bbb3SSatish Balay { 11442d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 11452d61bbb3SSatish Balay Mat outA; 11462d61bbb3SSatish Balay int ierr; 1147667159a5SBarry Smith PetscTruth row_identity,col_identity; 11482d61bbb3SSatish Balay 11492d61bbb3SSatish Balay PetscFunctionBegin; 115029bbc08cSBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 1151667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1152667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1153667159a5SBarry Smith if (!row_identity || !col_identity) { 115429bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ILU"); 1155667159a5SBarry Smith } 11562d61bbb3SSatish Balay 11572d61bbb3SSatish Balay outA = inA; 11582d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11592d61bbb3SSatish Balay 11602d61bbb3SSatish Balay if (!a->diag) { 1161c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 11622d61bbb3SSatish Balay } 11632d61bbb3SSatish Balay /* 116415091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 116515091d37SBarry Smith for ILU(0) factorization with natural ordering 11662d61bbb3SSatish Balay */ 116715091d37SBarry Smith switch (a->bs) { 1168f1af5d2fSBarry Smith case 1: 1169e37c484bSBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1170e37c484bSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering; 1171e37c484bSBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering; 1172e37c484bSBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=1\n"); 117315091d37SBarry Smith case 2: 117415091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 117515091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 11767c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 117715091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 117815091d37SBarry Smith break; 117915091d37SBarry Smith case 3: 118015091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 118115091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 11827c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 118315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 118415091d37SBarry Smith break; 118515091d37SBarry Smith case 4: 1186667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1187f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 11887c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering; 118915091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 119015091d37SBarry Smith break; 119115091d37SBarry Smith case 5: 1192667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1193667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 11947c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering; 119515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 119615091d37SBarry Smith break; 119715091d37SBarry Smith case 6: 119815091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 119915091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 12007c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering; 120115091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 120215091d37SBarry Smith break; 120315091d37SBarry Smith case 7: 120415091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 12057c922b88SBarry Smith inA->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering; 120615091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 120715091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 120815091d37SBarry Smith break; 1209c38d4ed2SBarry Smith default: 1210c38d4ed2SBarry Smith a->row = row; 1211c38d4ed2SBarry Smith a->col = col; 1212c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1213c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1214c38d4ed2SBarry Smith 1215c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 12164c49b128SBarry Smith ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); 1217c38d4ed2SBarry Smith PLogObjectParent(inA,a->icol); 1218c38d4ed2SBarry Smith 1219c38d4ed2SBarry Smith if (!a->solve_work) { 1220273d9f13SBarry Smith a->solve_work = (Scalar*)PetscMalloc((inA->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1221273d9f13SBarry Smith PLogObjectMemory(inA,(inA->m+a->bs)*sizeof(Scalar)); 1222c38d4ed2SBarry Smith } 12232d61bbb3SSatish Balay } 12242d61bbb3SSatish Balay 1225667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1226667159a5SBarry Smith 12272d61bbb3SSatish Balay PetscFunctionReturn(0); 12282d61bbb3SSatish Balay } 12292d61bbb3SSatish Balay #undef __FUNC__ 1230b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatPrintHelp_SeqBAIJ" 1231ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1232ba4ca20aSSatish Balay { 1233c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1234ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1235d132466eSBarry Smith int ierr; 1236ba4ca20aSSatish Balay 12373a40ed3dSBarry Smith PetscFunctionBegin; 1238c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1239d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1240d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12413a40ed3dSBarry Smith PetscFunctionReturn(0); 1242ba4ca20aSSatish Balay } 1243d9b7c43dSSatish Balay 1244fb2e594dSBarry Smith EXTERN_C_BEGIN 124527a8da17SBarry Smith #undef __FUNC__ 1246b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetColumnIndices_SeqBAIJ" 124727a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 124827a8da17SBarry Smith { 124927a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 125027a8da17SBarry Smith int i,nz,n; 125127a8da17SBarry Smith 125227a8da17SBarry Smith PetscFunctionBegin; 125327a8da17SBarry Smith nz = baij->maxnz; 1254273d9f13SBarry Smith n = mat->n; 125527a8da17SBarry Smith for (i=0; i<nz; i++) { 125627a8da17SBarry Smith baij->j[i] = indices[i]; 125727a8da17SBarry Smith } 125827a8da17SBarry Smith baij->nz = nz; 125927a8da17SBarry Smith for (i=0; i<n; i++) { 126027a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 126127a8da17SBarry Smith } 126227a8da17SBarry Smith 126327a8da17SBarry Smith PetscFunctionReturn(0); 126427a8da17SBarry Smith } 1265fb2e594dSBarry Smith EXTERN_C_END 126627a8da17SBarry Smith 126727a8da17SBarry Smith #undef __FUNC__ 1268b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetColumnIndices" 126927a8da17SBarry Smith /*@ 127027a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 127127a8da17SBarry Smith in the matrix. 127227a8da17SBarry Smith 127327a8da17SBarry Smith Input Parameters: 127427a8da17SBarry Smith + mat - the SeqBAIJ matrix 127527a8da17SBarry Smith - indices - the column indices 127627a8da17SBarry Smith 127715091d37SBarry Smith Level: advanced 127815091d37SBarry Smith 127927a8da17SBarry Smith Notes: 128027a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 128127a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 128227a8da17SBarry Smith of the MatSetValues() operation. 128327a8da17SBarry Smith 128427a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 128527a8da17SBarry Smith MatCreateSeqBAIJ(). 128627a8da17SBarry Smith 128727a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 128827a8da17SBarry Smith 128927a8da17SBarry Smith @*/ 129027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 129127a8da17SBarry Smith { 129227a8da17SBarry Smith int ierr,(*f)(Mat,int *); 129327a8da17SBarry Smith 129427a8da17SBarry Smith PetscFunctionBegin; 129527a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1296549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 129727a8da17SBarry Smith if (f) { 129827a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 129927a8da17SBarry Smith } else { 130029bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 130127a8da17SBarry Smith } 130227a8da17SBarry Smith PetscFunctionReturn(0); 130327a8da17SBarry Smith } 130427a8da17SBarry Smith 1305273d9f13SBarry Smith #undef __FUNC__ 1306273d9f13SBarry Smith #define __FUNC__ "MatGetRowMax_SeqBAIJ" 1307273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v) 1308273d9f13SBarry Smith { 1309273d9f13SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1310273d9f13SBarry Smith int ierr,i,j,n,row,bs,*ai,*aj,mbs; 1311273d9f13SBarry Smith PetscReal atmp; 1312273d9f13SBarry Smith Scalar *x,zero = 0.0; 1313273d9f13SBarry Smith MatScalar *aa; 1314273d9f13SBarry Smith int ncols,brow,krow,kcol; 1315273d9f13SBarry Smith 1316273d9f13SBarry Smith PetscFunctionBegin; 1317273d9f13SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1318273d9f13SBarry Smith bs = a->bs; 1319273d9f13SBarry Smith aa = a->a; 1320273d9f13SBarry Smith ai = a->i; 1321273d9f13SBarry Smith aj = a->j; 1322273d9f13SBarry Smith mbs = a->mbs; 1323273d9f13SBarry Smith 1324273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 1325273d9f13SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1326273d9f13SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1327273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1328273d9f13SBarry Smith for (i=0; i<mbs; i++) { 1329273d9f13SBarry Smith ncols = ai[1] - ai[0]; ai++; 1330273d9f13SBarry Smith brow = bs*i; 1331273d9f13SBarry Smith for (j=0; j<ncols; j++){ 1332273d9f13SBarry Smith /* bcol = bs*(*aj); */ 1333273d9f13SBarry Smith for (kcol=0; kcol<bs; kcol++){ 1334273d9f13SBarry Smith for (krow=0; krow<bs; krow++){ 1335273d9f13SBarry Smith atmp = PetscAbsScalar(*aa); aa++; 1336273d9f13SBarry Smith row = brow + krow; /* row index */ 1337273d9f13SBarry Smith /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */ 1338273d9f13SBarry Smith if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp; 1339273d9f13SBarry Smith } 1340273d9f13SBarry Smith } 1341273d9f13SBarry Smith aj++; 1342273d9f13SBarry Smith } 1343273d9f13SBarry Smith } 1344273d9f13SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1345273d9f13SBarry Smith PetscFunctionReturn(0); 1346273d9f13SBarry Smith } 1347273d9f13SBarry Smith 1348273d9f13SBarry Smith #undef __FUNC__ 1349273d9f13SBarry Smith #define __FUNC__ /*<a name="MatSetUpPreallocation_SeqBAIJ"></a>*/"MatSetUpPreallocation_SeqBAIJ" 1350273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A) 1351273d9f13SBarry Smith { 1352273d9f13SBarry Smith int ierr; 1353273d9f13SBarry Smith 1354273d9f13SBarry Smith PetscFunctionBegin; 1355273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1356273d9f13SBarry Smith PetscFunctionReturn(0); 1357273d9f13SBarry Smith } 1358273d9f13SBarry Smith 13592593348eSBarry Smith /* -------------------------------------------------------------------*/ 1360cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1361cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1362cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1363cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1364cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 13657c922b88SBarry Smith MatMultTranspose_SeqBAIJ, 13667c922b88SBarry Smith MatMultTransposeAdd_SeqBAIJ, 1367cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1368cc2dc46cSBarry Smith 0, 1369cc2dc46cSBarry Smith 0, 1370cc2dc46cSBarry Smith 0, 1371cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1372cc2dc46cSBarry Smith 0, 1373b6490206SBarry Smith 0, 1374f2501298SSatish Balay MatTranspose_SeqBAIJ, 1375cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1376cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1377cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1378cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1379cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1380b6490206SBarry Smith 0, 1381cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1382cc2dc46cSBarry Smith 0, 1383cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1384cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1385cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1386cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1387cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1388cc2dc46cSBarry Smith 0, 1389cc2dc46cSBarry Smith 0, 1390273d9f13SBarry Smith MatSetUpPreallocation_SeqBAIJ, 1391273d9f13SBarry Smith 0, 1392cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1393cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1394cc2dc46cSBarry Smith 0, 1395cc2dc46cSBarry Smith 0, 1396cc2dc46cSBarry Smith 0, 13972e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1398cc2dc46cSBarry Smith 0, 1399cc2dc46cSBarry Smith 0, 1400cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1401cc2dc46cSBarry Smith 0, 1402cc2dc46cSBarry Smith 0, 1403cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1404cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1405cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1406cc2dc46cSBarry Smith 0, 1407cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1408cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1409cc2dc46cSBarry Smith 0, 1410cc2dc46cSBarry Smith 0, 1411cc2dc46cSBarry Smith 0, 1412cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 14133b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 141492c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1415cc2dc46cSBarry Smith 0, 1416cc2dc46cSBarry Smith 0, 1417cc2dc46cSBarry Smith 0, 1418cc2dc46cSBarry Smith 0, 1419cc2dc46cSBarry Smith 0, 1420cc2dc46cSBarry Smith 0, 1421d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1422cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1423b9b97703SBarry Smith MatDestroy_SeqBAIJ, 1424b9b97703SBarry Smith MatView_SeqBAIJ, 1425273d9f13SBarry Smith MatGetMaps_Petsc, 1426273d9f13SBarry Smith 0, 1427273d9f13SBarry Smith 0, 1428273d9f13SBarry Smith 0, 1429273d9f13SBarry Smith 0, 1430273d9f13SBarry Smith 0, 1431273d9f13SBarry Smith 0, 1432273d9f13SBarry Smith MatGetRowMax_SeqBAIJ, 1433273d9f13SBarry Smith MatConvert_Basic}; 14342593348eSBarry Smith 14353e90b805SBarry Smith EXTERN_C_BEGIN 14363e90b805SBarry Smith #undef __FUNC__ 1437b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatStoreValues_SeqBAIJ" 14383e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 14393e90b805SBarry Smith { 14403e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1441273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2; 1442d132466eSBarry Smith int ierr; 14433e90b805SBarry Smith 14443e90b805SBarry Smith PetscFunctionBegin; 14453e90b805SBarry Smith if (aij->nonew != 1) { 144629bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14473e90b805SBarry Smith } 14483e90b805SBarry Smith 14493e90b805SBarry Smith /* allocate space for values if not already there */ 14503e90b805SBarry Smith if (!aij->saved_values) { 14513e90b805SBarry Smith aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 14523e90b805SBarry Smith } 14533e90b805SBarry Smith 14543e90b805SBarry Smith /* copy values over */ 1455d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 14563e90b805SBarry Smith PetscFunctionReturn(0); 14573e90b805SBarry Smith } 14583e90b805SBarry Smith EXTERN_C_END 14593e90b805SBarry Smith 14603e90b805SBarry Smith EXTERN_C_BEGIN 14613e90b805SBarry Smith #undef __FUNC__ 1462b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatRetrieveValues_SeqBAIJ" 146332e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 14643e90b805SBarry Smith { 14653e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1466273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 14673e90b805SBarry Smith 14683e90b805SBarry Smith PetscFunctionBegin; 14693e90b805SBarry Smith if (aij->nonew != 1) { 147029bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14713e90b805SBarry Smith } 14723e90b805SBarry Smith if (!aij->saved_values) { 147329bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 14743e90b805SBarry Smith } 14753e90b805SBarry Smith 14763e90b805SBarry Smith /* copy values over */ 1477549d3d68SSatish Balay ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 14783e90b805SBarry Smith PetscFunctionReturn(0); 14793e90b805SBarry Smith } 14803e90b805SBarry Smith EXTERN_C_END 14813e90b805SBarry Smith 1482273d9f13SBarry Smith EXTERN_C_BEGIN 1483273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*); 1484273d9f13SBarry Smith EXTERN_C_END 1485273d9f13SBarry Smith 1486273d9f13SBarry Smith EXTERN_C_BEGIN 14875615d1e5SSatish Balay #undef __FUNC__ 1488273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreate_SeqBAIJ" 1489273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B) 14902593348eSBarry Smith { 1491273d9f13SBarry Smith int ierr,size; 1492b6490206SBarry Smith Mat_SeqBAIJ *b; 14933b2fbd54SBarry Smith 14943a40ed3dSBarry Smith PetscFunctionBegin; 1495273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 149629bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1497b6490206SBarry Smith 1498273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1499273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1500b6490206SBarry Smith B->data = (void*)(b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1501549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1502549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 15032593348eSBarry Smith B->factor = 0; 15042593348eSBarry Smith B->lupivotthreshold = 1.0; 150590f02eecSBarry Smith B->mapping = 0; 15062593348eSBarry Smith b->row = 0; 15072593348eSBarry Smith b->col = 0; 1508e51c0b9cSSatish Balay b->icol = 0; 15092593348eSBarry Smith b->reallocs = 0; 15103e90b805SBarry Smith b->saved_values = 0; 1511563b5814SBarry Smith #if defined(PEYSC_USE_MAT_SINGLE) 1512563b5814SBarry Smith b->setvalueslen = 0; 1513563b5814SBarry Smith b->setvaluescopy = PETSC_NULL; 1514563b5814SBarry Smith #endif 15152593348eSBarry Smith 1516273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 1517273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1518a5ae1ecdSBarry Smith 1519c4992f7dSBarry Smith b->sorted = PETSC_FALSE; 1520c4992f7dSBarry Smith b->roworiented = PETSC_TRUE; 15212593348eSBarry Smith b->nonew = 0; 15222593348eSBarry Smith b->diag = 0; 15232593348eSBarry Smith b->solve_work = 0; 1524de6a44a3SBarry Smith b->mult_work = 0; 15252593348eSBarry Smith b->spptr = 0; 15260e6d2581SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 1527883fce79SBarry Smith b->keepzeroedrows = PETSC_FALSE; 15284e220ebcSLois Curfman McInnes 1529f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 15303e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 1531bc4b532fSSatish Balay MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1532f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 15333e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 1534bc4b532fSSatish Balay MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1535f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 153627a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1537bc4b532fSSatish Balay MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 1538273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_SeqBAIJ_SeqAIJ_C", 1539273d9f13SBarry Smith "MatConvert_SeqBAIJ_SeqAIJ", 1540273d9f13SBarry Smith MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr); 15413a40ed3dSBarry Smith PetscFunctionReturn(0); 15422593348eSBarry Smith } 1543273d9f13SBarry Smith EXTERN_C_END 15442593348eSBarry Smith 15455615d1e5SSatish Balay #undef __FUNC__ 1546b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatDuplicate_SeqBAIJ" 15472e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 15482593348eSBarry Smith { 15492593348eSBarry Smith Mat C; 1550b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data; 155127a8da17SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1552de6a44a3SBarry Smith 15533a40ed3dSBarry Smith PetscFunctionBegin; 155429bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 15552593348eSBarry Smith 15562593348eSBarry Smith *B = 0; 1557273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1558273d9f13SBarry Smith ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr); 1559273d9f13SBarry Smith c = (Mat_SeqBAIJ*)C->data; 156044cd7ae7SLois Curfman McInnes 1561b6490206SBarry Smith c->bs = a->bs; 15629df24120SSatish Balay c->bs2 = a->bs2; 1563b6490206SBarry Smith c->mbs = a->mbs; 1564b6490206SBarry Smith c->nbs = a->nbs; 1565*35d8aa7fSBarry Smith ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 15662593348eSBarry Smith 1567b6490206SBarry Smith c->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1568b6490206SBarry Smith c->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1569b6490206SBarry Smith for (i=0; i<mbs; i++) { 15702593348eSBarry Smith c->imax[i] = a->imax[i]; 15712593348eSBarry Smith c->ilen[i] = a->ilen[i]; 15722593348eSBarry Smith } 15732593348eSBarry Smith 15742593348eSBarry Smith /* allocate the matrix space */ 1575c4992f7dSBarry Smith c->singlemalloc = PETSC_TRUE; 15763f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 15773f1db9ecSBarry Smith c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a); 15787e67e3f9SSatish Balay c->j = (int*)(c->a + nz*bs2); 1579de6a44a3SBarry Smith c->i = c->j + nz; 1580549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1581b6490206SBarry Smith if (mbs > 0) { 1582549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 15832e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1584549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 15852e8a6d31SBarry Smith } else { 1586549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 15872593348eSBarry Smith } 15882593348eSBarry Smith } 15892593348eSBarry Smith 1590f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 15912593348eSBarry Smith c->sorted = a->sorted; 15922593348eSBarry Smith c->roworiented = a->roworiented; 15932593348eSBarry Smith c->nonew = a->nonew; 15942593348eSBarry Smith 15952593348eSBarry Smith if (a->diag) { 1596b6490206SBarry Smith c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag); 1597b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1598b6490206SBarry Smith for (i=0; i<mbs; i++) { 15992593348eSBarry Smith c->diag[i] = a->diag[i]; 16002593348eSBarry Smith } 160198305bb5SBarry Smith } else c->diag = 0; 16022593348eSBarry Smith c->nz = a->nz; 16032593348eSBarry Smith c->maxnz = a->maxnz; 16042593348eSBarry Smith c->solve_work = 0; 16052593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 16067fc0212eSBarry Smith c->mult_work = 0; 1607273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 1608273d9f13SBarry Smith C->assembled = PETSC_TRUE; 16092593348eSBarry Smith *B = C; 16107b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 16113a40ed3dSBarry Smith PetscFunctionReturn(0); 16122593348eSBarry Smith } 16132593348eSBarry Smith 1614273d9f13SBarry Smith EXTERN_C_BEGIN 16155615d1e5SSatish Balay #undef __FUNC__ 1616b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatLoad_SeqBAIJ" 161719bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 16182593348eSBarry Smith { 1619b6490206SBarry Smith Mat_SeqBAIJ *a; 16202593348eSBarry Smith Mat B; 1621f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1622b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 162335aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1624a7c10996SSatish Balay int *masked,nmask,tmp,bs2,ishift; 1625b6490206SBarry Smith Scalar *aa; 162619bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 16272593348eSBarry Smith 16283a40ed3dSBarry Smith PetscFunctionBegin; 1629f1af5d2fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1630de6a44a3SBarry Smith bs2 = bs*bs; 1631b6490206SBarry Smith 1632d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 163329bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 163490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16350752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 163629bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 16372593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 16382593348eSBarry Smith 1639d64ed03dSBarry Smith if (header[3] < 0) { 164029bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ"); 1641d64ed03dSBarry Smith } 1642d64ed03dSBarry Smith 164329bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 164435aab85fSBarry Smith 164535aab85fSBarry Smith /* 164635aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 164735aab85fSBarry Smith divisible by the blocksize 164835aab85fSBarry Smith */ 1649b6490206SBarry Smith mbs = M/bs; 165035aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 165135aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 165235aab85fSBarry Smith else mbs++; 165335aab85fSBarry Smith if (extra_rows) { 1654537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 165535aab85fSBarry Smith } 1656b6490206SBarry Smith 16572593348eSBarry Smith /* read in row lengths */ 165835aab85fSBarry Smith rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 16590752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 166035aab85fSBarry Smith for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 16612593348eSBarry Smith 1662b6490206SBarry Smith /* read in column indices */ 166335aab85fSBarry Smith jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj); 16640752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 166535aab85fSBarry Smith for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 1666b6490206SBarry Smith 1667b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1668b6490206SBarry Smith browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1669549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 167035aab85fSBarry Smith mask = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask); 1671549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 167235aab85fSBarry Smith masked = mask + mbs; 1673b6490206SBarry Smith rowcount = 0; nzcount = 0; 1674b6490206SBarry Smith for (i=0; i<mbs; i++) { 167535aab85fSBarry Smith nmask = 0; 1676b6490206SBarry Smith for (j=0; j<bs; j++) { 1677b6490206SBarry Smith kmax = rowlengths[rowcount]; 1678b6490206SBarry Smith for (k=0; k<kmax; k++) { 167935aab85fSBarry Smith tmp = jj[nzcount++]/bs; 168035aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1681b6490206SBarry Smith } 1682b6490206SBarry Smith rowcount++; 1683b6490206SBarry Smith } 168435aab85fSBarry Smith browlengths[i] += nmask; 168535aab85fSBarry Smith /* zero out the mask elements we set */ 168635aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1687b6490206SBarry Smith } 1688b6490206SBarry Smith 16892593348eSBarry Smith /* create our matrix */ 16903eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 16912593348eSBarry Smith B = *A; 1692b6490206SBarry Smith a = (Mat_SeqBAIJ*)B->data; 16932593348eSBarry Smith 16942593348eSBarry Smith /* set matrix "i" values */ 1695de6a44a3SBarry Smith a->i[0] = 0; 1696b6490206SBarry Smith for (i=1; i<= mbs; i++) { 1697b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1698b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 16992593348eSBarry Smith } 1700b6490206SBarry Smith a->nz = 0; 1701b6490206SBarry Smith for (i=0; i<mbs; i++) a->nz += browlengths[i]; 17022593348eSBarry Smith 1703b6490206SBarry Smith /* read in nonzero values */ 170435aab85fSBarry Smith aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 17050752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 170635aab85fSBarry Smith for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 1707b6490206SBarry Smith 1708b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1709b6490206SBarry Smith nzcount = 0; jcount = 0; 1710b6490206SBarry Smith for (i=0; i<mbs; i++) { 1711b6490206SBarry Smith nzcountb = nzcount; 171235aab85fSBarry Smith nmask = 0; 1713b6490206SBarry Smith for (j=0; j<bs; j++) { 1714b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1715b6490206SBarry Smith for (k=0; k<kmax; k++) { 171635aab85fSBarry Smith tmp = jj[nzcount++]/bs; 171735aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1718b6490206SBarry Smith } 1719b6490206SBarry Smith rowcount++; 1720b6490206SBarry Smith } 1721de6a44a3SBarry Smith /* sort the masked values */ 1722433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1723de6a44a3SBarry Smith 1724b6490206SBarry Smith /* set "j" values into matrix */ 1725b6490206SBarry Smith maskcount = 1; 172635aab85fSBarry Smith for (j=0; j<nmask; j++) { 172735aab85fSBarry Smith a->j[jcount++] = masked[j]; 1728de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1729b6490206SBarry Smith } 1730b6490206SBarry Smith /* set "a" values into matrix */ 1731de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1732b6490206SBarry Smith for (j=0; j<bs; j++) { 1733b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1734b6490206SBarry Smith for (k=0; k<kmax; k++) { 1735de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1736de6a44a3SBarry Smith block = mask[tmp] - 1; 1737de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1738de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1739375fe846SBarry Smith a->a[idx] = (MatScalar)aa[nzcountb++]; 1740b6490206SBarry Smith } 1741b6490206SBarry Smith } 174235aab85fSBarry Smith /* zero out the mask elements we set */ 174335aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1744b6490206SBarry Smith } 174529bbc08cSBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 1746b6490206SBarry Smith 1747606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1748606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1749606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1750606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1751606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1752b6490206SBarry Smith 1753b6490206SBarry Smith B->assembled = PETSC_TRUE; 1754de6a44a3SBarry Smith 17559c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 17563a40ed3dSBarry Smith PetscFunctionReturn(0); 17572593348eSBarry Smith } 1758273d9f13SBarry Smith EXTERN_C_END 17592593348eSBarry Smith 1760273d9f13SBarry Smith #undef __FUNC__ 1761273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatCreateSeqBAIJ" 1762273d9f13SBarry Smith /*@C 1763273d9f13SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 1764273d9f13SBarry Smith compressed row) format. For good matrix assembly performance the 1765273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 1766273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 1767273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 17682593348eSBarry Smith 1769273d9f13SBarry Smith Collective on MPI_Comm 1770273d9f13SBarry Smith 1771273d9f13SBarry Smith Input Parameters: 1772273d9f13SBarry Smith + comm - MPI communicator, set to PETSC_COMM_SELF 1773273d9f13SBarry Smith . bs - size of block 1774273d9f13SBarry Smith . m - number of rows 1775273d9f13SBarry Smith . n - number of columns 1776*35d8aa7fSBarry Smith . nz - number of nonzero blocks per block row (same for all rows) 1777*35d8aa7fSBarry Smith - nnz - array containing the number of nonzero blocks in the various block rows 1778273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 1779273d9f13SBarry Smith 1780273d9f13SBarry Smith Output Parameter: 1781273d9f13SBarry Smith . A - the matrix 1782273d9f13SBarry Smith 1783273d9f13SBarry Smith Options Database Keys: 1784273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 1785273d9f13SBarry Smith block calculations (much slower) 1786273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 1787273d9f13SBarry Smith 1788273d9f13SBarry Smith Level: intermediate 1789273d9f13SBarry Smith 1790273d9f13SBarry Smith Notes: 1791*35d8aa7fSBarry Smith A nonzero block is any block that as 1 or more nonzeros in it 1792*35d8aa7fSBarry Smith 1793273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 1794273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 1795273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 1796273d9f13SBarry Smith 1797273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1798273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1799273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 1800273d9f13SBarry Smith matrices. 1801273d9f13SBarry Smith 1802273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 1803273d9f13SBarry Smith @*/ 1804273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1805273d9f13SBarry Smith { 1806273d9f13SBarry Smith int ierr; 1807273d9f13SBarry Smith 1808273d9f13SBarry Smith PetscFunctionBegin; 1809273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1810273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr); 1811273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1812273d9f13SBarry Smith PetscFunctionReturn(0); 1813273d9f13SBarry Smith } 1814273d9f13SBarry Smith 1815273d9f13SBarry Smith #undef __FUNC__ 1816273d9f13SBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatSeqBAIJSetPreallocation" 1817273d9f13SBarry Smith /*@C 1818273d9f13SBarry Smith MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros 1819273d9f13SBarry Smith per row in the matrix. For good matrix assembly performance the 1820273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 1821273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 1822273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 1823273d9f13SBarry Smith 1824273d9f13SBarry Smith Collective on MPI_Comm 1825273d9f13SBarry Smith 1826273d9f13SBarry Smith Input Parameters: 1827273d9f13SBarry Smith + A - the matrix 1828273d9f13SBarry Smith . bs - size of block 1829273d9f13SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 1830273d9f13SBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 1831273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 1832273d9f13SBarry Smith 1833273d9f13SBarry Smith Options Database Keys: 1834273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 1835273d9f13SBarry Smith block calculations (much slower) 1836273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 1837273d9f13SBarry Smith 1838273d9f13SBarry Smith Level: intermediate 1839273d9f13SBarry Smith 1840273d9f13SBarry Smith Notes: 1841273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 1842273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 1843273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 1844273d9f13SBarry Smith 1845273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1846273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1847273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 1848273d9f13SBarry Smith matrices. 1849273d9f13SBarry Smith 1850273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 1851273d9f13SBarry Smith @*/ 1852273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 1853273d9f13SBarry Smith { 1854273d9f13SBarry Smith Mat_SeqBAIJ *b; 1855273d9f13SBarry Smith int i,len,ierr,mbs,nbs,bs2,newbs = bs; 1856273d9f13SBarry Smith PetscTruth flg; 1857273d9f13SBarry Smith 1858273d9f13SBarry Smith PetscFunctionBegin; 1859273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr); 1860273d9f13SBarry Smith if (!flg) PetscFunctionReturn(0); 1861273d9f13SBarry Smith 1862273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1863273d9f13SBarry Smith ierr = OptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr); 1864273d9f13SBarry Smith if (nnz && newbs != bs) { 1865273d9f13SBarry Smith SETERRQ(1,"Cannot change blocksize from command line if setting nnz"); 1866273d9f13SBarry Smith } 1867273d9f13SBarry Smith bs = newbs; 1868273d9f13SBarry Smith 1869273d9f13SBarry Smith mbs = B->m/bs; 1870273d9f13SBarry Smith nbs = B->n/bs; 1871273d9f13SBarry Smith bs2 = bs*bs; 1872273d9f13SBarry Smith 1873273d9f13SBarry Smith if (mbs*bs!=B->m || nbs*bs!=B->n) { 1874*35d8aa7fSBarry Smith SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs); 1875273d9f13SBarry Smith } 1876273d9f13SBarry Smith 1877273d9f13SBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than -2: value %d",nz); 1878273d9f13SBarry Smith if (nnz) { 1879273d9f13SBarry Smith for (i=0; i<mbs; i++) { 1880273d9f13SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1881273d9f13SBarry Smith } 1882273d9f13SBarry Smith } 1883273d9f13SBarry Smith 1884273d9f13SBarry Smith b = (Mat_SeqBAIJ*)B->data; 1885273d9f13SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 1886273d9f13SBarry Smith if (!flg) { 1887273d9f13SBarry Smith switch (bs) { 1888273d9f13SBarry Smith case 1: 1889273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1890273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1891273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 1892273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1893273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 1894273d9f13SBarry Smith break; 1895273d9f13SBarry Smith case 2: 1896273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1897273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1898273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 1899273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1900273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 1901273d9f13SBarry Smith break; 1902273d9f13SBarry Smith case 3: 1903273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1904273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1905273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 1906273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1907273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 1908273d9f13SBarry Smith break; 1909273d9f13SBarry Smith case 4: 1910273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1911273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1912273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4; 1913273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1914273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 1915273d9f13SBarry Smith break; 1916273d9f13SBarry Smith case 5: 1917273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1918273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1919273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5; 1920273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1921273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 1922273d9f13SBarry Smith break; 1923273d9f13SBarry Smith case 6: 1924273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 1925273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 1926273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6; 1927273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 1928273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 1929273d9f13SBarry Smith break; 1930273d9f13SBarry Smith case 7: 1931273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1932273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1933273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7; 1934273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 1935273d9f13SBarry Smith break; 1936273d9f13SBarry Smith default: 1937273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_N; 1938273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_N; 1939273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_N; 1940273d9f13SBarry Smith break; 1941273d9f13SBarry Smith } 1942273d9f13SBarry Smith } 1943273d9f13SBarry Smith b->bs = bs; 1944273d9f13SBarry Smith b->mbs = mbs; 1945273d9f13SBarry Smith b->nbs = nbs; 1946273d9f13SBarry Smith b->imax = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax); 1947273d9f13SBarry Smith if (!nnz) { 1948273d9f13SBarry Smith if (nz == PETSC_DEFAULT) nz = 5; 1949273d9f13SBarry Smith else if (nz <= 0) nz = 1; 1950273d9f13SBarry Smith for (i=0; i<mbs; i++) b->imax[i] = nz; 1951273d9f13SBarry Smith nz = nz*mbs; 1952273d9f13SBarry Smith } else { 1953273d9f13SBarry Smith nz = 0; 1954273d9f13SBarry Smith for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 1955273d9f13SBarry Smith } 1956273d9f13SBarry Smith 1957273d9f13SBarry Smith /* allocate the matrix space */ 1958273d9f13SBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 1959273d9f13SBarry Smith b->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a); 1960273d9f13SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1961273d9f13SBarry Smith b->j = (int*)(b->a + nz*bs2); 1962273d9f13SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 1963273d9f13SBarry Smith b->i = b->j + nz; 1964273d9f13SBarry Smith b->singlemalloc = PETSC_TRUE; 1965273d9f13SBarry Smith 1966273d9f13SBarry Smith b->i[0] = 0; 1967273d9f13SBarry Smith for (i=1; i<mbs+1; i++) { 1968273d9f13SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 1969273d9f13SBarry Smith } 1970273d9f13SBarry Smith 1971273d9f13SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1972273d9f13SBarry Smith b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int)); 1973273d9f13SBarry Smith PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1974273d9f13SBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 1975273d9f13SBarry Smith 1976273d9f13SBarry Smith b->bs = bs; 1977273d9f13SBarry Smith b->bs2 = bs2; 1978273d9f13SBarry Smith b->mbs = mbs; 1979273d9f13SBarry Smith b->nz = 0; 1980273d9f13SBarry Smith b->maxnz = nz*bs2; 1981273d9f13SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 1982273d9f13SBarry Smith PetscFunctionReturn(0); 1983273d9f13SBarry Smith } 19842593348eSBarry Smith 1985