1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*faf3f1fcSBarry Smith static char vcid[] = "$Id: baij.c,v 1.181 1999/09/15 02:03:56 bsmith Exp bsmith $"; 32593348eSBarry Smith #endif 42593348eSBarry Smith 52593348eSBarry Smith /* 6b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 72593348eSBarry Smith matrix storage format. 82593348eSBarry Smith */ 93369ce9aSBarry Smith #include "sys.h" 1070f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 111a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 121a6a6d98SBarry Smith #include "src/inline/spops.h" 133b547af2SSatish Balay 142d61bbb3SSatish Balay #define CHUNKSIZE 10 15de6a44a3SBarry Smith 16be5855fcSBarry Smith /* 17be5855fcSBarry Smith Checks for missing diagonals 18be5855fcSBarry Smith */ 19be5855fcSBarry Smith #undef __FUNC__ 20be5855fcSBarry Smith #define __FUNC__ "MatMissingDiag_SeqBAIJ" 21be5855fcSBarry Smith int MatMissingDiag_SeqBAIJ(Mat A) 22be5855fcSBarry Smith { 23be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 24be5855fcSBarry Smith int *diag = a->diag, *jj = a->j,i; 25be5855fcSBarry Smith 26be5855fcSBarry Smith PetscFunctionBegin; 270e8e8aceSBarry Smith for ( i=0; i<a->mbs; i++ ) { 28be5855fcSBarry Smith if (jj[diag[i]] != i) { 29be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 30be5855fcSBarry Smith } 31be5855fcSBarry Smith } 32be5855fcSBarry Smith PetscFunctionReturn(0); 33be5855fcSBarry Smith } 34be5855fcSBarry Smith 355615d1e5SSatish Balay #undef __FUNC__ 365615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ" 37de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 38de6a44a3SBarry Smith { 39de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 407fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 41de6a44a3SBarry Smith 423a40ed3dSBarry Smith PetscFunctionBegin; 43de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 44de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 457fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 46ecc615b2SBarry Smith diag[i] = a->i[i+1]; 47de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 48de6a44a3SBarry Smith if (a->j[j] == i) { 49de6a44a3SBarry Smith diag[i] = j; 50de6a44a3SBarry Smith break; 51de6a44a3SBarry Smith } 52de6a44a3SBarry Smith } 53de6a44a3SBarry Smith } 54de6a44a3SBarry Smith a->diag = diag; 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56de6a44a3SBarry Smith } 572593348eSBarry Smith 582593348eSBarry Smith 593b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 603b2fbd54SBarry Smith 615615d1e5SSatish Balay #undef __FUNC__ 625615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 633b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 643b2fbd54SBarry Smith PetscTruth *done) 653b2fbd54SBarry Smith { 663b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 673b2fbd54SBarry Smith int ierr, n = a->mbs,i; 683b2fbd54SBarry Smith 693a40ed3dSBarry Smith PetscFunctionBegin; 703b2fbd54SBarry Smith *nn = n; 713a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 723b2fbd54SBarry Smith if (symmetric) { 733b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 743b2fbd54SBarry Smith } else if (oshift == 1) { 753b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 763b2fbd54SBarry Smith int nz = a->i[n] + 1; 773b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 783b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 793b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 803b2fbd54SBarry Smith } else { 813b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 823b2fbd54SBarry Smith } 833b2fbd54SBarry Smith 843a40ed3dSBarry Smith PetscFunctionReturn(0); 853b2fbd54SBarry Smith } 863b2fbd54SBarry Smith 875615d1e5SSatish Balay #undef __FUNC__ 88d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 893b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 903b2fbd54SBarry Smith PetscTruth *done) 913b2fbd54SBarry Smith { 923b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 93606d414cSSatish Balay int i,n = a->mbs,ierr; 943b2fbd54SBarry Smith 953a40ed3dSBarry Smith PetscFunctionBegin; 963a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 973b2fbd54SBarry Smith if (symmetric) { 98606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 99606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 100af5da2bfSBarry Smith } else if (oshift == 1) { 1013b2fbd54SBarry Smith int nz = a->i[n]; 1023b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 1033b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 1043b2fbd54SBarry Smith } 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 1063b2fbd54SBarry Smith } 1073b2fbd54SBarry Smith 1082d61bbb3SSatish Balay #undef __FUNC__ 1092d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 1102d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 1112d61bbb3SSatish Balay { 1122d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 1132d61bbb3SSatish Balay 1142d61bbb3SSatish Balay PetscFunctionBegin; 1152d61bbb3SSatish Balay *bs = baij->bs; 1162d61bbb3SSatish Balay PetscFunctionReturn(0); 1172d61bbb3SSatish Balay } 1182d61bbb3SSatish Balay 1192d61bbb3SSatish Balay 1202d61bbb3SSatish Balay #undef __FUNC__ 1212d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 122e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1232d61bbb3SSatish Balay { 1242d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 125e51c0b9cSSatish Balay int ierr; 1262d61bbb3SSatish Balay 12794d884c6SBarry Smith if (A->mapping) { 12894d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 12994d884c6SBarry Smith } 13094d884c6SBarry Smith if (A->bmapping) { 13194d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 13294d884c6SBarry Smith } 13361b13de0SBarry Smith if (A->rmap) { 13461b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 13561b13de0SBarry Smith } 13661b13de0SBarry Smith if (A->cmap) { 13761b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 13861b13de0SBarry Smith } 139aa482453SBarry Smith #if defined(PETSC_USE_LOG) 140e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1412d61bbb3SSatish Balay #endif 142606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 143606d414cSSatish Balay if (!a->singlemalloc) { 144606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 145606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 146606d414cSSatish Balay } 147606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 148606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 149606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 150606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 151606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 152e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 153606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 154606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1552d61bbb3SSatish Balay PLogObjectDestroy(A); 1562d61bbb3SSatish Balay PetscHeaderDestroy(A); 1572d61bbb3SSatish Balay PetscFunctionReturn(0); 1582d61bbb3SSatish Balay } 1592d61bbb3SSatish Balay 1602d61bbb3SSatish Balay #undef __FUNC__ 1612d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1622d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1632d61bbb3SSatish Balay { 1642d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1652d61bbb3SSatish Balay 1662d61bbb3SSatish Balay PetscFunctionBegin; 1672d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1682d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1692d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1702d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1712d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1724787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1734787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1742d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1752d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1762d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1772d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1782d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1792d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 180b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 181b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 182981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1832d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1842d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1852d61bbb3SSatish Balay } else { 1862d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1872d61bbb3SSatish Balay } 1882d61bbb3SSatish Balay PetscFunctionReturn(0); 1892d61bbb3SSatish Balay } 1902d61bbb3SSatish Balay 1912d61bbb3SSatish Balay 1922d61bbb3SSatish Balay #undef __FUNC__ 1932d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1942d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1952d61bbb3SSatish Balay { 1962d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1972d61bbb3SSatish Balay 1982d61bbb3SSatish Balay PetscFunctionBegin; 1992d61bbb3SSatish Balay if (m) *m = a->m; 2002d61bbb3SSatish Balay if (n) *n = a->n; 2012d61bbb3SSatish Balay PetscFunctionReturn(0); 2022d61bbb3SSatish Balay } 2032d61bbb3SSatish Balay 2042d61bbb3SSatish Balay #undef __FUNC__ 2052d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 2062d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2072d61bbb3SSatish Balay { 2082d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2092d61bbb3SSatish Balay 2102d61bbb3SSatish Balay PetscFunctionBegin; 2112d61bbb3SSatish Balay *m = 0; *n = a->m; 2122d61bbb3SSatish Balay PetscFunctionReturn(0); 2132d61bbb3SSatish Balay } 2142d61bbb3SSatish Balay 2152d61bbb3SSatish Balay #undef __FUNC__ 2162d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 2172d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2182d61bbb3SSatish Balay { 2192d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2202d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2213f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2223f1db9ecSBarry Smith Scalar *v_i; 2232d61bbb3SSatish Balay 2242d61bbb3SSatish Balay PetscFunctionBegin; 2252d61bbb3SSatish Balay bs = a->bs; 2262d61bbb3SSatish Balay ai = a->i; 2272d61bbb3SSatish Balay aj = a->j; 2282d61bbb3SSatish Balay aa = a->a; 2292d61bbb3SSatish Balay bs2 = a->bs2; 2302d61bbb3SSatish Balay 2312d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2322d61bbb3SSatish Balay 2332d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2342d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2352d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2362d61bbb3SSatish Balay *nz = bs*M; 2372d61bbb3SSatish Balay 2382d61bbb3SSatish Balay if (v) { 2392d61bbb3SSatish Balay *v = 0; 2402d61bbb3SSatish Balay if (*nz) { 2412d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v); 2422d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2432d61bbb3SSatish Balay v_i = *v + i*bs; 2442d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2452d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2462d61bbb3SSatish Balay } 2472d61bbb3SSatish Balay } 2482d61bbb3SSatish Balay } 2492d61bbb3SSatish Balay 2502d61bbb3SSatish Balay if (idx) { 2512d61bbb3SSatish Balay *idx = 0; 2522d61bbb3SSatish Balay if (*nz) { 2532d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx); 2542d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2552d61bbb3SSatish Balay idx_i = *idx + i*bs; 2562d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2572d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2582d61bbb3SSatish Balay } 2592d61bbb3SSatish Balay } 2602d61bbb3SSatish Balay } 2612d61bbb3SSatish Balay PetscFunctionReturn(0); 2622d61bbb3SSatish Balay } 2632d61bbb3SSatish Balay 2642d61bbb3SSatish Balay #undef __FUNC__ 2652d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2662d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2672d61bbb3SSatish Balay { 268606d414cSSatish Balay int ierr; 269606d414cSSatish Balay 2702d61bbb3SSatish Balay PetscFunctionBegin; 271606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 272606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2732d61bbb3SSatish Balay PetscFunctionReturn(0); 2742d61bbb3SSatish Balay } 2752d61bbb3SSatish Balay 2762d61bbb3SSatish Balay #undef __FUNC__ 2772d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2782d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2792d61bbb3SSatish Balay { 2802d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2812d61bbb3SSatish Balay Mat C; 2822d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2832d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2843f1db9ecSBarry Smith MatScalar *array = a->a; 2852d61bbb3SSatish Balay 2862d61bbb3SSatish Balay PetscFunctionBegin; 2872d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2882d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 289549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2902d61bbb3SSatish Balay 2912d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2922d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 293606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 2942d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 2952d61bbb3SSatish Balay cols = rows + bs; 2962d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2972d61bbb3SSatish Balay cols[0] = i*bs; 2982d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 2992d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3002d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 3012d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3022d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 3032d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3042d61bbb3SSatish Balay array += bs2; 3052d61bbb3SSatish Balay } 3062d61bbb3SSatish Balay } 307606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 3082d61bbb3SSatish Balay 3092d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3102d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3112d61bbb3SSatish Balay 3122d61bbb3SSatish Balay if (B != PETSC_NULL) { 3132d61bbb3SSatish Balay *B = C; 3142d61bbb3SSatish Balay } else { 315f830108cSBarry Smith PetscOps *Abops; 316cc2dc46cSBarry Smith MatOps Aops; 317f830108cSBarry Smith 3182d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 319606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 320606d414cSSatish Balay if (!a->singlemalloc) { 321606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 322606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 323606d414cSSatish Balay } 324606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 325606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 326606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 327606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 328606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 329f830108cSBarry Smith 3307413bad6SBarry Smith 3317413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3327413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3337413bad6SBarry Smith 334f830108cSBarry Smith /* 335f830108cSBarry Smith This is horrible, horrible code. We need to keep the 336f830108cSBarry Smith A pointers for the bops and ops but copy everything 337f830108cSBarry Smith else from C. 338f830108cSBarry Smith */ 339f830108cSBarry Smith Abops = A->bops; 340f830108cSBarry Smith Aops = A->ops; 341549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 342f830108cSBarry Smith A->bops = Abops; 343f830108cSBarry Smith A->ops = Aops; 34427a8da17SBarry Smith A->qlist = 0; 345f830108cSBarry Smith 3462d61bbb3SSatish Balay PetscHeaderDestroy(C); 3472d61bbb3SSatish Balay } 3482d61bbb3SSatish Balay PetscFunctionReturn(0); 3492d61bbb3SSatish Balay } 3502d61bbb3SSatish Balay 3515615d1e5SSatish Balay #undef __FUNC__ 352d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 353b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3542593348eSBarry Smith { 355b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3569df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 357b6490206SBarry Smith Scalar *aa; 358ce6f0cecSBarry Smith FILE *file; 3592593348eSBarry Smith 3603a40ed3dSBarry Smith PetscFunctionBegin; 36190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3622593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3632593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3643b2fbd54SBarry Smith 3652593348eSBarry Smith col_lens[1] = a->m; 3662593348eSBarry Smith col_lens[2] = a->n; 3677e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3682593348eSBarry Smith 3692593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 370b6490206SBarry Smith count = 0; 371b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 372b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 373b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 374b6490206SBarry Smith } 3752593348eSBarry Smith } 3760752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 377606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3782593348eSBarry Smith 3792593348eSBarry Smith /* store column indices (zero start index) */ 38066963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj); 381b6490206SBarry Smith count = 0; 382b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 383b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 384b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 385b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 386b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3872593348eSBarry Smith } 3882593348eSBarry Smith } 389b6490206SBarry Smith } 390b6490206SBarry Smith } 3910752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 392606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 3932593348eSBarry Smith 3942593348eSBarry Smith /* store nonzero values */ 39566963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 396b6490206SBarry Smith count = 0; 397b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 398b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 399b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 400b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4017e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 402b6490206SBarry Smith } 403b6490206SBarry Smith } 404b6490206SBarry Smith } 405b6490206SBarry Smith } 4060752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 407606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 408ce6f0cecSBarry Smith 409ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 410ce6f0cecSBarry Smith if (file) { 411ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 412ce6f0cecSBarry Smith } 4133a40ed3dSBarry Smith PetscFunctionReturn(0); 4142593348eSBarry Smith } 4152593348eSBarry Smith 4165615d1e5SSatish Balay #undef __FUNC__ 417d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 418b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4192593348eSBarry Smith { 420b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4219df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4222593348eSBarry Smith FILE *fd; 4232593348eSBarry Smith char *outputname; 4242593348eSBarry Smith 4253a40ed3dSBarry Smith PetscFunctionBegin; 42690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 42777ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 428888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 429639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 430d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4310ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 4327b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket format not supported"); 4330ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 43444cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 43544cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 43644cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 43744cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 43844cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 439aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4400ef38995SBarry Smith if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 44144cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 442e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 4430ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 444766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 445e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 4460ef38995SBarry Smith } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 447e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 4480ef38995SBarry Smith } 44944cd7ae7SLois Curfman McInnes #else 4500ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 45144cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 4520ef38995SBarry Smith } 45344cd7ae7SLois Curfman McInnes #endif 45444cd7ae7SLois Curfman McInnes } 45544cd7ae7SLois Curfman McInnes } 45644cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 45744cd7ae7SLois Curfman McInnes } 45844cd7ae7SLois Curfman McInnes } 4590ef38995SBarry Smith } else { 460b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 461b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 462b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 463b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 464b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 465aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 466e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 46788685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 468e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 4690ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 470766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 471e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 4720ef38995SBarry Smith } else { 473e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 47488685aaeSLois Curfman McInnes } 47588685aaeSLois Curfman McInnes #else 4767e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 47788685aaeSLois Curfman McInnes #endif 4782593348eSBarry Smith } 4792593348eSBarry Smith } 4802593348eSBarry Smith fprintf(fd,"\n"); 4812593348eSBarry Smith } 4822593348eSBarry Smith } 483b6490206SBarry Smith } 4842593348eSBarry Smith fflush(fd); 4853a40ed3dSBarry Smith PetscFunctionReturn(0); 4862593348eSBarry Smith } 4872593348eSBarry Smith 4885615d1e5SSatish Balay #undef __FUNC__ 48977ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 49077ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4913270192aSSatish Balay { 49277ed5343SBarry Smith Mat A = (Mat) Aa; 4933270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4947dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 495fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 4963f1db9ecSBarry Smith MatScalar *aa; 49777ed5343SBarry Smith MPI_Comm comm; 49877ed5343SBarry Smith Viewer viewer; 4993270192aSSatish Balay 5003a40ed3dSBarry Smith PetscFunctionBegin; 50177ed5343SBarry Smith /* 50277ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 50377ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 50477ed5343SBarry Smith rest should return immediately. 50577ed5343SBarry Smith */ 50677ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 507d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 50877ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 5093270192aSSatish Balay 51077ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 51177ed5343SBarry Smith 51277ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 51377ed5343SBarry Smith 5143270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5153270192aSSatish Balay color = DRAW_BLUE; 5163270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5173270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5183270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5193270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5203270192aSSatish Balay aa = a->a + j*bs2; 5213270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5223270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5233270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 524b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5253270192aSSatish Balay } 5263270192aSSatish Balay } 5273270192aSSatish Balay } 5283270192aSSatish Balay } 5293270192aSSatish Balay color = DRAW_CYAN; 5303270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5313270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5323270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5333270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5343270192aSSatish Balay aa = a->a + j*bs2; 5353270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5363270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5373270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 538b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5393270192aSSatish Balay } 5403270192aSSatish Balay } 5413270192aSSatish Balay } 5423270192aSSatish Balay } 5433270192aSSatish Balay 5443270192aSSatish Balay color = DRAW_RED; 5453270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5463270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5473270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5483270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5493270192aSSatish Balay aa = a->a + j*bs2; 5503270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5513270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5523270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 553b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5543270192aSSatish Balay } 5553270192aSSatish Balay } 5563270192aSSatish Balay } 5573270192aSSatish Balay } 55877ed5343SBarry Smith PetscFunctionReturn(0); 55977ed5343SBarry Smith } 5603270192aSSatish Balay 56177ed5343SBarry Smith #undef __FUNC__ 56277ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 56377ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 56477ed5343SBarry Smith { 56577ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5667dce120fSSatish Balay int ierr; 5677dce120fSSatish Balay double xl,yl,xr,yr,w,h; 56877ed5343SBarry Smith Draw draw; 56977ed5343SBarry Smith PetscTruth isnull; 5703270192aSSatish Balay 57177ed5343SBarry Smith PetscFunctionBegin; 57277ed5343SBarry Smith 57377ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 57477ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 57577ed5343SBarry Smith 57677ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 57777ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 57877ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5793270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 58077ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 58177ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5823a40ed3dSBarry Smith PetscFunctionReturn(0); 5833270192aSSatish Balay } 5843270192aSSatish Balay 5855615d1e5SSatish Balay #undef __FUNC__ 586d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 587e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5882593348eSBarry Smith { 58919bcc07fSBarry Smith ViewerType vtype; 59019bcc07fSBarry Smith int ierr; 5912593348eSBarry Smith 5923a40ed3dSBarry Smith PetscFunctionBegin; 59319bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5947b2a1423SBarry Smith if (PetscTypeCompare(vtype,SOCKET_VIEWER)) { 5957b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 5963f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){ 5973a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 5983f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 5993a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6003f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 6013a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6025cd90555SBarry Smith } else { 6035cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 6042593348eSBarry Smith } 6053a40ed3dSBarry Smith PetscFunctionReturn(0); 6062593348eSBarry Smith } 607b6490206SBarry Smith 608cd0e1443SSatish Balay 6095615d1e5SSatish Balay #undef __FUNC__ 6102d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6112d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 612cd0e1443SSatish Balay { 613cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6142d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6152d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6162d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6173f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 618cd0e1443SSatish Balay 6193a40ed3dSBarry Smith PetscFunctionBegin; 6202d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 621cd0e1443SSatish Balay row = im[k]; brow = row/bs; 622a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 623a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6242d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6252c3acbe9SBarry Smith nrow = ailen[brow]; 6262d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 627a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 628a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6292d61bbb3SSatish Balay col = in[l] ; 6302d61bbb3SSatish Balay bcol = col/bs; 6312d61bbb3SSatish Balay cidx = col%bs; 6322d61bbb3SSatish Balay ridx = row%bs; 6332d61bbb3SSatish Balay high = nrow; 6342d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6352d61bbb3SSatish Balay while (high-low > 5) { 636cd0e1443SSatish Balay t = (low+high)/2; 637cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 638cd0e1443SSatish Balay else low = t; 639cd0e1443SSatish Balay } 640cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 641cd0e1443SSatish Balay if (rp[i] > bcol) break; 642cd0e1443SSatish Balay if (rp[i] == bcol) { 6432d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6442d61bbb3SSatish Balay goto finished; 645cd0e1443SSatish Balay } 646cd0e1443SSatish Balay } 6472d61bbb3SSatish Balay *v++ = zero; 6482d61bbb3SSatish Balay finished:; 649cd0e1443SSatish Balay } 650cd0e1443SSatish Balay } 6513a40ed3dSBarry Smith PetscFunctionReturn(0); 652cd0e1443SSatish Balay } 653cd0e1443SSatish Balay 6542d61bbb3SSatish Balay 6555615d1e5SSatish Balay #undef __FUNC__ 65605a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 65792c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 65892c4ed94SBarry Smith { 65992c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6608a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 661dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 662549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6633f1db9ecSBarry Smith Scalar *value = v; 6643f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 66592c4ed94SBarry Smith 6663a40ed3dSBarry Smith PetscFunctionBegin; 6670e324ae4SSatish Balay if (roworiented) { 6680e324ae4SSatish Balay stepval = (n-1)*bs; 6690e324ae4SSatish Balay } else { 6700e324ae4SSatish Balay stepval = (m-1)*bs; 6710e324ae4SSatish Balay } 67292c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 67392c4ed94SBarry Smith row = im[k]; 6745ef9f2a5SBarry Smith if (row < 0) continue; 675aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 676a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 67792c4ed94SBarry Smith #endif 67892c4ed94SBarry Smith rp = aj + ai[row]; 67992c4ed94SBarry Smith ap = aa + bs2*ai[row]; 68092c4ed94SBarry Smith rmax = imax[row]; 68192c4ed94SBarry Smith nrow = ailen[row]; 68292c4ed94SBarry Smith low = 0; 68392c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6845ef9f2a5SBarry Smith if (in[l] < 0) continue; 685aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 686a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 68792c4ed94SBarry Smith #endif 68892c4ed94SBarry Smith col = in[l]; 68992c4ed94SBarry Smith if (roworiented) { 6900e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6910e324ae4SSatish Balay } else { 6920e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 69392c4ed94SBarry Smith } 69492c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 69592c4ed94SBarry Smith while (high-low > 7) { 69692c4ed94SBarry Smith t = (low+high)/2; 69792c4ed94SBarry Smith if (rp[t] > col) high = t; 69892c4ed94SBarry Smith else low = t; 69992c4ed94SBarry Smith } 70092c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 70192c4ed94SBarry Smith if (rp[i] > col) break; 70292c4ed94SBarry Smith if (rp[i] == col) { 7038a84c255SSatish Balay bap = ap + bs2*i; 7040e324ae4SSatish Balay if (roworiented) { 7058a84c255SSatish Balay if (is == ADD_VALUES) { 706dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 707dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7088a84c255SSatish Balay bap[jj] += *value++; 709dd9472c6SBarry Smith } 710dd9472c6SBarry Smith } 7110e324ae4SSatish Balay } else { 712dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 713dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7140e324ae4SSatish Balay bap[jj] = *value++; 7158a84c255SSatish Balay } 716dd9472c6SBarry Smith } 717dd9472c6SBarry Smith } 7180e324ae4SSatish Balay } else { 7190e324ae4SSatish Balay if (is == ADD_VALUES) { 720dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 721dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7220e324ae4SSatish Balay *bap++ += *value++; 723dd9472c6SBarry Smith } 724dd9472c6SBarry Smith } 7250e324ae4SSatish Balay } else { 726dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 727dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7280e324ae4SSatish Balay *bap++ = *value++; 7290e324ae4SSatish Balay } 7308a84c255SSatish Balay } 731dd9472c6SBarry Smith } 732dd9472c6SBarry Smith } 733f1241b54SBarry Smith goto noinsert2; 73492c4ed94SBarry Smith } 73592c4ed94SBarry Smith } 73689280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 737a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 73892c4ed94SBarry Smith if (nrow >= rmax) { 73992c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 74092c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7413f1db9ecSBarry Smith MatScalar *new_a; 74292c4ed94SBarry Smith 743a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74489280ab3SLois Curfman McInnes 74592c4ed94SBarry Smith /* malloc new storage space */ 7463f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7473f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 74892c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 74992c4ed94SBarry Smith new_i = new_j + new_nz; 75092c4ed94SBarry Smith 75192c4ed94SBarry Smith /* copy over old data into new slots */ 75292c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 75392c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 754549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 75592c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 756549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 757549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 758549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 759549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 76092c4ed94SBarry Smith /* free up old matrix storage */ 761606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 762606d414cSSatish Balay if (!a->singlemalloc) { 763606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 764606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 765606d414cSSatish Balay } 76692c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 76792c4ed94SBarry Smith a->singlemalloc = 1; 76892c4ed94SBarry Smith 76992c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 77092c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7713f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 77292c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 77392c4ed94SBarry Smith a->reallocs++; 77492c4ed94SBarry Smith a->nz++; 77592c4ed94SBarry Smith } 77692c4ed94SBarry Smith N = nrow++ - 1; 77792c4ed94SBarry Smith /* shift up all the later entries in this row */ 77892c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 77992c4ed94SBarry Smith rp[ii+1] = rp[ii]; 780549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 78192c4ed94SBarry Smith } 782549d3d68SSatish Balay if (N >= i) { 783549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 784549d3d68SSatish Balay } 78592c4ed94SBarry Smith rp[i] = col; 7868a84c255SSatish Balay bap = ap + bs2*i; 7870e324ae4SSatish Balay if (roworiented) { 788dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 789dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7900e324ae4SSatish Balay bap[jj] = *value++; 791dd9472c6SBarry Smith } 792dd9472c6SBarry Smith } 7930e324ae4SSatish Balay } else { 794dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 795dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7960e324ae4SSatish Balay *bap++ = *value++; 7970e324ae4SSatish Balay } 798dd9472c6SBarry Smith } 799dd9472c6SBarry Smith } 800f1241b54SBarry Smith noinsert2:; 80192c4ed94SBarry Smith low = i; 80292c4ed94SBarry Smith } 80392c4ed94SBarry Smith ailen[row] = nrow; 80492c4ed94SBarry Smith } 8053a40ed3dSBarry Smith PetscFunctionReturn(0); 80692c4ed94SBarry Smith } 80792c4ed94SBarry Smith 808f2501298SSatish Balay 8095615d1e5SSatish Balay #undef __FUNC__ 8105615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8118f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 812584200bdSSatish Balay { 813584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 814584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 815a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 816549d3d68SSatish Balay int mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr; 8173f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 818584200bdSSatish Balay 8193a40ed3dSBarry Smith PetscFunctionBegin; 8203a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 821584200bdSSatish Balay 82243ee02c3SBarry Smith if (m) rmax = ailen[0]; 823584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 824584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 825584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 826d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 827584200bdSSatish Balay if (fshift) { 828a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 829584200bdSSatish Balay N = ailen[i]; 830584200bdSSatish Balay for ( j=0; j<N; j++ ) { 831584200bdSSatish Balay ip[j-fshift] = ip[j]; 832549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 833584200bdSSatish Balay } 834584200bdSSatish Balay } 835584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 836584200bdSSatish Balay } 837584200bdSSatish Balay if (mbs) { 838584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 839584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 840584200bdSSatish Balay } 841584200bdSSatish Balay /* reset ilen and imax for each row */ 842584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 843584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 844584200bdSSatish Balay } 845a7c10996SSatish Balay a->nz = ai[mbs]; 846584200bdSSatish Balay 847584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 848584200bdSSatish Balay if (fshift && a->diag) { 849606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 850584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 851584200bdSSatish Balay a->diag = 0; 852584200bdSSatish Balay } 8533d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 854219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8553d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 856584200bdSSatish Balay a->reallocs); 857d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 858e2f3b5e9SSatish Balay a->reallocs = 0; 8594e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8604e220ebcSLois Curfman McInnes 8613a40ed3dSBarry Smith PetscFunctionReturn(0); 862584200bdSSatish Balay } 863584200bdSSatish Balay 8645a838052SSatish Balay 865bea157c4SSatish Balay 866bea157c4SSatish Balay /* 867bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 868bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 869bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 870bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 871bea157c4SSatish Balay */ 8725615d1e5SSatish Balay #undef __FUNC__ 873bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks" 874bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 875d9b7c43dSSatish Balay { 876bea157c4SSatish Balay int i,j,k,row; 877bea157c4SSatish Balay PetscTruth flg; 8783a40ed3dSBarry Smith 879bea157c4SSatish Balay /* PetscFunctionBegin;*/ 880bea157c4SSatish Balay for ( i=0,j=0; i<n; j++ ) { 881bea157c4SSatish Balay row = idx[i]; 882bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 883bea157c4SSatish Balay sizes[j] = 1; 884bea157c4SSatish Balay i++; 885e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 886bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 887bea157c4SSatish Balay i++; 888bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 889bea157c4SSatish Balay flg = PETSC_TRUE; 890bea157c4SSatish Balay for ( k=1; k<bs; k++ ) { 891bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 892bea157c4SSatish Balay flg = PETSC_FALSE; 893bea157c4SSatish Balay break; 894d9b7c43dSSatish Balay } 895bea157c4SSatish Balay } 896bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 897bea157c4SSatish Balay sizes[j] = bs; 898bea157c4SSatish Balay i+= bs; 899bea157c4SSatish Balay } else { 900bea157c4SSatish Balay sizes[j] = 1; 901bea157c4SSatish Balay i++; 902bea157c4SSatish Balay } 903bea157c4SSatish Balay } 904bea157c4SSatish Balay } 905bea157c4SSatish Balay *bs_max = j; 9063a40ed3dSBarry Smith PetscFunctionReturn(0); 907d9b7c43dSSatish Balay } 908d9b7c43dSSatish Balay 9095615d1e5SSatish Balay #undef __FUNC__ 9105615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 9118f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 912d9b7c43dSSatish Balay { 913d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 914b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 915bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9163f1db9ecSBarry Smith Scalar zero = 0.0; 9173f1db9ecSBarry Smith MatScalar *aa; 918d9b7c43dSSatish Balay 9193a40ed3dSBarry Smith PetscFunctionBegin; 920d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 921d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 922d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 923d9b7c43dSSatish Balay 924bea157c4SSatish Balay /* allocate memory for rows,sizes */ 925bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 926bea157c4SSatish Balay sizes = rows + is_n; 927bea157c4SSatish Balay 928bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 929bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 930bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 931bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 932b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 933bea157c4SSatish Balay 934bea157c4SSatish Balay for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) { 935bea157c4SSatish Balay row = rows[j]; 936b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 937bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 938bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 939bea157c4SSatish Balay if (sizes[i] == bs) { 940bea157c4SSatish Balay if (diag) { 941bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 942bea157c4SSatish Balay baij->ilen[row/bs] = 1; 943bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 944bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 945a07cd24cSSatish Balay } 946bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 947bea157c4SSatish Balay for ( k=0; k<bs; k++ ) { 948bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 949bea157c4SSatish Balay } 950bea157c4SSatish Balay } else { /* (!diag) */ 951bea157c4SSatish Balay baij->ilen[row/bs] = 0; 952bea157c4SSatish Balay } /* end (!diag) */ 953bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 954aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 955bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 956bea157c4SSatish Balay #endif 957bea157c4SSatish Balay for ( k=0; k<count; k++ ) { 958d9b7c43dSSatish Balay aa[0] = zero; 959d9b7c43dSSatish Balay aa+=bs; 960d9b7c43dSSatish Balay } 961d9b7c43dSSatish Balay if (diag) { 962f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 963d9b7c43dSSatish Balay } 964d9b7c43dSSatish Balay } 965bea157c4SSatish Balay } 966bea157c4SSatish Balay 967606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 9689a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9693a40ed3dSBarry Smith PetscFunctionReturn(0); 970d9b7c43dSSatish Balay } 9711c351548SSatish Balay 9725615d1e5SSatish Balay #undef __FUNC__ 9732d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9742d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9752d61bbb3SSatish Balay { 9762d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9772d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9782d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9792d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 980549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 9813f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9822d61bbb3SSatish Balay 9832d61bbb3SSatish Balay PetscFunctionBegin; 9842d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9852d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9865ef9f2a5SBarry Smith if (row < 0) continue; 987aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 98832e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 9892d61bbb3SSatish Balay #endif 9902d61bbb3SSatish Balay rp = aj + ai[brow]; 9912d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9922d61bbb3SSatish Balay rmax = imax[brow]; 9932d61bbb3SSatish Balay nrow = ailen[brow]; 9942d61bbb3SSatish Balay low = 0; 9952d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9965ef9f2a5SBarry Smith if (in[l] < 0) continue; 997aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 99832e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 9992d61bbb3SSatish Balay #endif 10002d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10012d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10022d61bbb3SSatish Balay if (roworiented) { 10035ef9f2a5SBarry Smith value = v[l + k*n]; 10042d61bbb3SSatish Balay } else { 10052d61bbb3SSatish Balay value = v[k + l*m]; 10062d61bbb3SSatish Balay } 10072d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10082d61bbb3SSatish Balay while (high-low > 7) { 10092d61bbb3SSatish Balay t = (low+high)/2; 10102d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10112d61bbb3SSatish Balay else low = t; 10122d61bbb3SSatish Balay } 10132d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 10142d61bbb3SSatish Balay if (rp[i] > bcol) break; 10152d61bbb3SSatish Balay if (rp[i] == bcol) { 10162d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10172d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10182d61bbb3SSatish Balay else *bap = value; 10192d61bbb3SSatish Balay goto noinsert1; 10202d61bbb3SSatish Balay } 10212d61bbb3SSatish Balay } 10222d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10232d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10242d61bbb3SSatish Balay if (nrow >= rmax) { 10252d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10262d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10273f1db9ecSBarry Smith MatScalar *new_a; 10282d61bbb3SSatish Balay 10292d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10302d61bbb3SSatish Balay 10312d61bbb3SSatish Balay /* Malloc new storage space */ 10323f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10333f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 10342d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 10352d61bbb3SSatish Balay new_i = new_j + new_nz; 10362d61bbb3SSatish Balay 10372d61bbb3SSatish Balay /* copy over old data into new slots */ 10382d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10392d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1040549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10412d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1042549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1043549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1044549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1045549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10462d61bbb3SSatish Balay /* free up old matrix storage */ 1047606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1048606d414cSSatish Balay if (!a->singlemalloc) { 1049606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1050606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1051606d414cSSatish Balay } 10522d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 10532d61bbb3SSatish Balay a->singlemalloc = 1; 10542d61bbb3SSatish Balay 10552d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10562d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10573f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10582d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10592d61bbb3SSatish Balay a->reallocs++; 10602d61bbb3SSatish Balay a->nz++; 10612d61bbb3SSatish Balay } 10622d61bbb3SSatish Balay N = nrow++ - 1; 10632d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10642d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10652d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1066549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10672d61bbb3SSatish Balay } 1068549d3d68SSatish Balay if (N>=i) { 1069549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1070549d3d68SSatish Balay } 10712d61bbb3SSatish Balay rp[i] = bcol; 10722d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10732d61bbb3SSatish Balay noinsert1:; 10742d61bbb3SSatish Balay low = i; 10752d61bbb3SSatish Balay } 10762d61bbb3SSatish Balay ailen[brow] = nrow; 10772d61bbb3SSatish Balay } 10782d61bbb3SSatish Balay PetscFunctionReturn(0); 10792d61bbb3SSatish Balay } 10802d61bbb3SSatish Balay 10812d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10822d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10832d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 10847b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 10857b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 10862d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10872d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10882d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10892d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10902d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10912d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10922d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10932d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10942d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10952d61bbb3SSatish Balay 10962d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10972d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 10982d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 10992d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 11002d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 11012d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 110215091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 11032d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 11042d61bbb3SSatish Balay 11052d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 11062d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 11072d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 11082d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 11092d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 11102d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 111115091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11122d61bbb3SSatish Balay 11132d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11142d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11152d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11162d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11172d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 111815091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11192d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11202d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11212d61bbb3SSatish Balay 11222d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11232d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11242d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11252d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11262d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 112715091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11282d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11292d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11302d61bbb3SSatish Balay 11312d61bbb3SSatish Balay #undef __FUNC__ 11322d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 11335ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11342d61bbb3SSatish Balay { 11352d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 11362d61bbb3SSatish Balay Mat outA; 11372d61bbb3SSatish Balay int ierr; 1138667159a5SBarry Smith PetscTruth row_identity, col_identity; 11392d61bbb3SSatish Balay 11402d61bbb3SSatish Balay PetscFunctionBegin; 1141b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1142667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1143667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1144667159a5SBarry Smith if (!row_identity || !col_identity) { 1145b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1146667159a5SBarry Smith } 11472d61bbb3SSatish Balay 11482d61bbb3SSatish Balay outA = inA; 11492d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11502d61bbb3SSatish Balay a->row = row; 11512d61bbb3SSatish Balay a->col = col; 11522d61bbb3SSatish Balay 1153e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1154e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr); 11551e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1156e51c0b9cSSatish Balay 11572d61bbb3SSatish Balay if (!a->solve_work) { 11582d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 11592d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 11602d61bbb3SSatish Balay } 11612d61bbb3SSatish Balay 11622d61bbb3SSatish Balay if (!a->diag) { 11632d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA);CHKERRQ(ierr); 11642d61bbb3SSatish Balay } 11652d61bbb3SSatish Balay /* 116615091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 116715091d37SBarry Smith for ILU(0) factorization with natural ordering 11682d61bbb3SSatish Balay */ 116915091d37SBarry Smith switch (a->bs) { 117015091d37SBarry Smith case 2: 117115091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 117215091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 117315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 117415091d37SBarry Smith break; 117515091d37SBarry Smith case 3: 117615091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 117715091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 117815091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 117915091d37SBarry Smith break; 118015091d37SBarry Smith case 4: 1181667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1182f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 118315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 118415091d37SBarry Smith break; 118515091d37SBarry Smith case 5: 1186667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1187667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 118815091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 118915091d37SBarry Smith break; 119015091d37SBarry Smith case 6: 119115091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 119215091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 119315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 119415091d37SBarry Smith break; 119515091d37SBarry Smith case 7: 119615091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 119715091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 119815091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 119915091d37SBarry Smith break; 12002d61bbb3SSatish Balay } 12012d61bbb3SSatish Balay 1202667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1203667159a5SBarry Smith 12042d61bbb3SSatish Balay PetscFunctionReturn(0); 12052d61bbb3SSatish Balay } 12062d61bbb3SSatish Balay #undef __FUNC__ 1207d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1208ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1209ba4ca20aSSatish Balay { 1210ba4ca20aSSatish Balay static int called = 0; 1211ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1212d132466eSBarry Smith int ierr; 1213ba4ca20aSSatish Balay 12143a40ed3dSBarry Smith PetscFunctionBegin; 12153a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 1216d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1217d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12183a40ed3dSBarry Smith PetscFunctionReturn(0); 1219ba4ca20aSSatish Balay } 1220d9b7c43dSSatish Balay 1221fb2e594dSBarry Smith EXTERN_C_BEGIN 122227a8da17SBarry Smith #undef __FUNC__ 122327a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 122427a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 122527a8da17SBarry Smith { 122627a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 122727a8da17SBarry Smith int i,nz,n; 122827a8da17SBarry Smith 122927a8da17SBarry Smith PetscFunctionBegin; 123027a8da17SBarry Smith nz = baij->maxnz; 123127a8da17SBarry Smith n = baij->n; 123227a8da17SBarry Smith for (i=0; i<nz; i++) { 123327a8da17SBarry Smith baij->j[i] = indices[i]; 123427a8da17SBarry Smith } 123527a8da17SBarry Smith baij->nz = nz; 123627a8da17SBarry Smith for ( i=0; i<n; i++ ) { 123727a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 123827a8da17SBarry Smith } 123927a8da17SBarry Smith 124027a8da17SBarry Smith PetscFunctionReturn(0); 124127a8da17SBarry Smith } 1242fb2e594dSBarry Smith EXTERN_C_END 124327a8da17SBarry Smith 124427a8da17SBarry Smith #undef __FUNC__ 124527a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 124627a8da17SBarry Smith /*@ 124727a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 124827a8da17SBarry Smith in the matrix. 124927a8da17SBarry Smith 125027a8da17SBarry Smith Input Parameters: 125127a8da17SBarry Smith + mat - the SeqBAIJ matrix 125227a8da17SBarry Smith - indices - the column indices 125327a8da17SBarry Smith 125415091d37SBarry Smith Level: advanced 125515091d37SBarry Smith 125627a8da17SBarry Smith Notes: 125727a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 125827a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 125927a8da17SBarry Smith of the MatSetValues() operation. 126027a8da17SBarry Smith 126127a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 126227a8da17SBarry Smith MatCreateSeqBAIJ(). 126327a8da17SBarry Smith 126427a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 126527a8da17SBarry Smith 126627a8da17SBarry Smith @*/ 126727a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 126827a8da17SBarry Smith { 126927a8da17SBarry Smith int ierr,(*f)(Mat,int *); 127027a8da17SBarry Smith 127127a8da17SBarry Smith PetscFunctionBegin; 127227a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1273549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 127427a8da17SBarry Smith if (f) { 127527a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 127627a8da17SBarry Smith } else { 127727a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 127827a8da17SBarry Smith } 127927a8da17SBarry Smith PetscFunctionReturn(0); 128027a8da17SBarry Smith } 128127a8da17SBarry Smith 12822593348eSBarry Smith /* -------------------------------------------------------------------*/ 1283cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1284cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1285cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1286cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1287cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1288cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1289cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1290cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1291cc2dc46cSBarry Smith 0, 1292cc2dc46cSBarry Smith 0, 1293cc2dc46cSBarry Smith 0, 1294cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1295cc2dc46cSBarry Smith 0, 1296b6490206SBarry Smith 0, 1297f2501298SSatish Balay MatTranspose_SeqBAIJ, 1298cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1299cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1300cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1301cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1302cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1303b6490206SBarry Smith 0, 1304cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1305cc2dc46cSBarry Smith 0, 1306cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1307cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1308cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1309cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1310cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1311cc2dc46cSBarry Smith 0, 1312cc2dc46cSBarry Smith 0, 1313cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1314cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1315cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1316cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1317cc2dc46cSBarry Smith 0, 1318cc2dc46cSBarry Smith 0, 1319cc2dc46cSBarry Smith 0, 13202e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1321cc2dc46cSBarry Smith 0, 1322cc2dc46cSBarry Smith 0, 1323cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1324cc2dc46cSBarry Smith 0, 1325cc2dc46cSBarry Smith 0, 1326cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1327cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1328cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1329cc2dc46cSBarry Smith 0, 1330cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1331cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1332cc2dc46cSBarry Smith 0, 1333cc2dc46cSBarry Smith 0, 1334cc2dc46cSBarry Smith 0, 1335cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 13363b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 133792c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1338cc2dc46cSBarry Smith 0, 1339cc2dc46cSBarry Smith 0, 1340cc2dc46cSBarry Smith 0, 1341cc2dc46cSBarry Smith 0, 1342cc2dc46cSBarry Smith 0, 1343cc2dc46cSBarry Smith 0, 1344d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1345cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1346cc2dc46cSBarry Smith 0, 1347cc2dc46cSBarry Smith 0, 1348cc2dc46cSBarry Smith MatGetMaps_Petsc}; 13492593348eSBarry Smith 13503e90b805SBarry Smith EXTERN_C_BEGIN 13513e90b805SBarry Smith #undef __FUNC__ 13523e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ" 13533e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 13543e90b805SBarry Smith { 13553e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 13563e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1357d132466eSBarry Smith int ierr; 13583e90b805SBarry Smith 13593e90b805SBarry Smith PetscFunctionBegin; 13603e90b805SBarry Smith if (aij->nonew != 1) { 13613e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13623e90b805SBarry Smith } 13633e90b805SBarry Smith 13643e90b805SBarry Smith /* allocate space for values if not already there */ 13653e90b805SBarry Smith if (!aij->saved_values) { 13663e90b805SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 13673e90b805SBarry Smith } 13683e90b805SBarry Smith 13693e90b805SBarry Smith /* copy values over */ 1370d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 13713e90b805SBarry Smith PetscFunctionReturn(0); 13723e90b805SBarry Smith } 13733e90b805SBarry Smith EXTERN_C_END 13743e90b805SBarry Smith 13753e90b805SBarry Smith EXTERN_C_BEGIN 13763e90b805SBarry Smith #undef __FUNC__ 137732e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ" 137832e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 13793e90b805SBarry Smith { 13803e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1381549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 13823e90b805SBarry Smith 13833e90b805SBarry Smith PetscFunctionBegin; 13843e90b805SBarry Smith if (aij->nonew != 1) { 13853e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13863e90b805SBarry Smith } 13873e90b805SBarry Smith if (!aij->saved_values) { 13883e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 13893e90b805SBarry Smith } 13903e90b805SBarry Smith 13913e90b805SBarry Smith /* copy values over */ 1392549d3d68SSatish Balay ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 13933e90b805SBarry Smith PetscFunctionReturn(0); 13943e90b805SBarry Smith } 13953e90b805SBarry Smith EXTERN_C_END 13963e90b805SBarry Smith 13975615d1e5SSatish Balay #undef __FUNC__ 13985615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 13992593348eSBarry Smith /*@C 140044cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 140144cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 140244cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 14037fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 14042bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 14052593348eSBarry Smith 1406db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1407db81eaa0SLois Curfman McInnes 14082593348eSBarry Smith Input Parameters: 1409db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1410b6490206SBarry Smith . bs - size of block 14112593348eSBarry Smith . m - number of rows 14122593348eSBarry Smith . n - number of columns 1413b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14147fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14152bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14162593348eSBarry Smith 14172593348eSBarry Smith Output Parameter: 14182593348eSBarry Smith . A - the matrix 14192593348eSBarry Smith 14200513a670SBarry Smith Options Database Keys: 1421db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1422db81eaa0SLois Curfman McInnes block calculations (much slower) 1423db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14240513a670SBarry Smith 142515091d37SBarry Smith Level: intermediate 142615091d37SBarry Smith 14272593348eSBarry Smith Notes: 142844cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 14292593348eSBarry Smith storage. That is, the stored row and column indices can begin at 143044cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 14312593348eSBarry Smith 14322593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 14332593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14342593348eSBarry Smith allocation. For additional details, see the users manual chapter on 14356da5968aSLois Curfman McInnes matrices. 14362593348eSBarry Smith 1437db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 14382593348eSBarry Smith @*/ 1439b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 14402593348eSBarry Smith { 14412593348eSBarry Smith Mat B; 1442b6490206SBarry Smith Mat_SeqBAIJ *b; 1443302e33e3SBarry Smith int i,len,ierr,flg,mbs,nbs,bs2,size; 14443b2fbd54SBarry Smith 14453a40ed3dSBarry Smith PetscFunctionBegin; 1446d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1447a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1448b6490206SBarry Smith 1449962fb4a0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1450302e33e3SBarry Smith mbs = m/bs; 1451302e33e3SBarry Smith nbs = n/bs; 1452302e33e3SBarry Smith bs2 = bs*bs; 14530513a670SBarry Smith 14543a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1455a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 14563a40ed3dSBarry Smith } 14572593348eSBarry Smith 1458b73539f3SBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz); 1459b73539f3SBarry Smith if (nnz) { 1460*faf3f1fcSBarry Smith for (i=0; i<mbs; i++) { 1461b73539f3SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1462b73539f3SBarry Smith } 1463b73539f3SBarry Smith } 1464b73539f3SBarry Smith 14652593348eSBarry Smith *A = 0; 14663f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 14672593348eSBarry Smith PLogObjectCreate(B); 1468b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1469549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1470549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 14711a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 14721a6a6d98SBarry Smith if (!flg) { 14737fc0212eSBarry Smith switch (bs) { 14747fc0212eSBarry Smith case 1: 1475f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1476f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1477f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1478f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 14797fc0212eSBarry Smith break; 14804eeb42bcSBarry Smith case 2: 1481f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1482f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1483f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1484f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 14856d84be18SBarry Smith break; 1486f327dd97SBarry Smith case 3: 1487f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1488f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1489f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1490f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 14914eeb42bcSBarry Smith break; 14926d84be18SBarry Smith case 4: 1493f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1494f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1495f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1496f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 14976d84be18SBarry Smith break; 14986d84be18SBarry Smith case 5: 1499f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1500f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1501f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1502f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 15036d84be18SBarry Smith break; 150415091d37SBarry Smith case 6: 150515091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 150615091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 150715091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 150815091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 150915091d37SBarry Smith break; 151048e9ddb2SSatish Balay case 7: 1511f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1512f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1513f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 151448e9ddb2SSatish Balay break; 15157fc0212eSBarry Smith } 15161a6a6d98SBarry Smith } 1517e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1518e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15192593348eSBarry Smith B->factor = 0; 15202593348eSBarry Smith B->lupivotthreshold = 1.0; 152190f02eecSBarry Smith B->mapping = 0; 15222593348eSBarry Smith b->row = 0; 15232593348eSBarry Smith b->col = 0; 1524e51c0b9cSSatish Balay b->icol = 0; 15252593348eSBarry Smith b->reallocs = 0; 15263e90b805SBarry Smith b->saved_values = 0; 15272593348eSBarry Smith 152844cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 152944cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1530a5ae1ecdSBarry Smith 15317413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 15327413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1533a5ae1ecdSBarry Smith 1534b6490206SBarry Smith b->mbs = mbs; 1535f2501298SSatish Balay b->nbs = nbs; 1536b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax); 15372593348eSBarry Smith if (nnz == PETSC_NULL) { 1538119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 15392593348eSBarry Smith else if (nz <= 0) nz = 1; 1540b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1541b6490206SBarry Smith nz = nz*mbs; 15423a40ed3dSBarry Smith } else { 15432593348eSBarry Smith nz = 0; 1544b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 15452593348eSBarry Smith } 15462593348eSBarry Smith 15472593348eSBarry Smith /* allocate the matrix space */ 15483f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 15493f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a); 1550549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15517e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 1552549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15532593348eSBarry Smith b->i = b->j + nz; 15542593348eSBarry Smith b->singlemalloc = 1; 15552593348eSBarry Smith 1556de6a44a3SBarry Smith b->i[0] = 0; 1557b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 15582593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 15592593348eSBarry Smith } 15602593348eSBarry Smith 1561b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1562b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1563f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1564b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 15652593348eSBarry Smith 1566b6490206SBarry Smith b->bs = bs; 15679df24120SSatish Balay b->bs2 = bs2; 1568b6490206SBarry Smith b->mbs = mbs; 15692593348eSBarry Smith b->nz = 0; 157018e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 15712593348eSBarry Smith b->sorted = 0; 15722593348eSBarry Smith b->roworiented = 1; 15732593348eSBarry Smith b->nonew = 0; 15742593348eSBarry Smith b->diag = 0; 15752593348eSBarry Smith b->solve_work = 0; 1576de6a44a3SBarry Smith b->mult_work = 0; 15772593348eSBarry Smith b->spptr = 0; 15784e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 15794e220ebcSLois Curfman McInnes 15802593348eSBarry Smith *A = B; 15812593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 15822593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 158327a8da17SBarry Smith 15843e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 15853e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 15863e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 15873e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 15883e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 15893e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 159027a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 159127a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 159227a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 15933a40ed3dSBarry Smith PetscFunctionReturn(0); 15942593348eSBarry Smith } 15952593348eSBarry Smith 15965615d1e5SSatish Balay #undef __FUNC__ 15972e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 15982e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 15992593348eSBarry Smith { 16002593348eSBarry Smith Mat C; 1601b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 160227a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1603de6a44a3SBarry Smith 16043a40ed3dSBarry Smith PetscFunctionBegin; 1605a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 16062593348eSBarry Smith 16072593348eSBarry Smith *B = 0; 16083f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 16092593348eSBarry Smith PLogObjectCreate(C); 1610b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1611549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1612e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1613e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 16142593348eSBarry Smith C->factor = A->factor; 16152593348eSBarry Smith c->row = 0; 16162593348eSBarry Smith c->col = 0; 1617cac97260SSatish Balay c->icol = 0; 161832e87ba3SBarry Smith c->saved_values = 0; 16192593348eSBarry Smith C->assembled = PETSC_TRUE; 16202593348eSBarry Smith 162144cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 162244cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 162344cd7ae7SLois Curfman McInnes C->M = a->m; 162444cd7ae7SLois Curfman McInnes C->N = a->n; 162544cd7ae7SLois Curfman McInnes 1626b6490206SBarry Smith c->bs = a->bs; 16279df24120SSatish Balay c->bs2 = a->bs2; 1628b6490206SBarry Smith c->mbs = a->mbs; 1629b6490206SBarry Smith c->nbs = a->nbs; 16302593348eSBarry Smith 1631b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1632b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1633b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16342593348eSBarry Smith c->imax[i] = a->imax[i]; 16352593348eSBarry Smith c->ilen[i] = a->ilen[i]; 16362593348eSBarry Smith } 16372593348eSBarry Smith 16382593348eSBarry Smith /* allocate the matrix space */ 16392593348eSBarry Smith c->singlemalloc = 1; 16403f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 16413f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a); 16427e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1643de6a44a3SBarry Smith c->i = c->j + nz; 1644549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1645b6490206SBarry Smith if (mbs > 0) { 1646549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16472e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1648549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16492e8a6d31SBarry Smith } else { 1650549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16512593348eSBarry Smith } 16522593348eSBarry Smith } 16532593348eSBarry Smith 1654f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16552593348eSBarry Smith c->sorted = a->sorted; 16562593348eSBarry Smith c->roworiented = a->roworiented; 16572593348eSBarry Smith c->nonew = a->nonew; 16582593348eSBarry Smith 16592593348eSBarry Smith if (a->diag) { 1660b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag); 1661b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1662b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16632593348eSBarry Smith c->diag[i] = a->diag[i]; 16642593348eSBarry Smith } 166598305bb5SBarry Smith } else c->diag = 0; 16662593348eSBarry Smith c->nz = a->nz; 16672593348eSBarry Smith c->maxnz = a->maxnz; 16682593348eSBarry Smith c->solve_work = 0; 16692593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 16707fc0212eSBarry Smith c->mult_work = 0; 16712593348eSBarry Smith *B = C; 16727b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 16733a40ed3dSBarry Smith PetscFunctionReturn(0); 16742593348eSBarry Smith } 16752593348eSBarry Smith 16765615d1e5SSatish Balay #undef __FUNC__ 16775615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 167819bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 16792593348eSBarry Smith { 1680b6490206SBarry Smith Mat_SeqBAIJ *a; 16812593348eSBarry Smith Mat B; 1682de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1683b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 168435aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1685a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1686b6490206SBarry Smith Scalar *aa; 168719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 16882593348eSBarry Smith 16893a40ed3dSBarry Smith PetscFunctionBegin; 16901a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1691de6a44a3SBarry Smith bs2 = bs*bs; 1692b6490206SBarry Smith 1693d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1694cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 169590ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16960752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1697a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 16982593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 16992593348eSBarry Smith 1700d64ed03dSBarry Smith if (header[3] < 0) { 1701a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1702d64ed03dSBarry Smith } 1703d64ed03dSBarry Smith 1704a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 170535aab85fSBarry Smith 170635aab85fSBarry Smith /* 170735aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 170835aab85fSBarry Smith divisible by the blocksize 170935aab85fSBarry Smith */ 1710b6490206SBarry Smith mbs = M/bs; 171135aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 171235aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 171335aab85fSBarry Smith else mbs++; 171435aab85fSBarry Smith if (extra_rows) { 1715537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 171635aab85fSBarry Smith } 1717b6490206SBarry Smith 17182593348eSBarry Smith /* read in row lengths */ 171935aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 17200752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 172135aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 17222593348eSBarry Smith 1723b6490206SBarry Smith /* read in column indices */ 172435aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj); 17250752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 172635aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1727b6490206SBarry Smith 1728b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1729b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1730549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 173135aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask); 1732549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 173335aab85fSBarry Smith masked = mask + mbs; 1734b6490206SBarry Smith rowcount = 0; nzcount = 0; 1735b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 173635aab85fSBarry Smith nmask = 0; 1737b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1738b6490206SBarry Smith kmax = rowlengths[rowcount]; 1739b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 174035aab85fSBarry Smith tmp = jj[nzcount++]/bs; 174135aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1742b6490206SBarry Smith } 1743b6490206SBarry Smith rowcount++; 1744b6490206SBarry Smith } 174535aab85fSBarry Smith browlengths[i] += nmask; 174635aab85fSBarry Smith /* zero out the mask elements we set */ 174735aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1748b6490206SBarry Smith } 1749b6490206SBarry Smith 17502593348eSBarry Smith /* create our matrix */ 17513eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17522593348eSBarry Smith B = *A; 1753b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 17542593348eSBarry Smith 17552593348eSBarry Smith /* set matrix "i" values */ 1756de6a44a3SBarry Smith a->i[0] = 0; 1757b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1758b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1759b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17602593348eSBarry Smith } 1761b6490206SBarry Smith a->nz = 0; 1762b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 17632593348eSBarry Smith 1764b6490206SBarry Smith /* read in nonzero values */ 176535aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 17660752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 176735aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1768b6490206SBarry Smith 1769b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1770b6490206SBarry Smith nzcount = 0; jcount = 0; 1771b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1772b6490206SBarry Smith nzcountb = nzcount; 177335aab85fSBarry Smith nmask = 0; 1774b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1775b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1776b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 177735aab85fSBarry Smith tmp = jj[nzcount++]/bs; 177835aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1779b6490206SBarry Smith } 1780b6490206SBarry Smith rowcount++; 1781b6490206SBarry Smith } 1782de6a44a3SBarry Smith /* sort the masked values */ 178377c4ece6SBarry Smith PetscSortInt(nmask,masked); 1784de6a44a3SBarry Smith 1785b6490206SBarry Smith /* set "j" values into matrix */ 1786b6490206SBarry Smith maskcount = 1; 178735aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 178835aab85fSBarry Smith a->j[jcount++] = masked[j]; 1789de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1790b6490206SBarry Smith } 1791b6490206SBarry Smith /* set "a" values into matrix */ 1792de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1793b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1794b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1795b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1796de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1797de6a44a3SBarry Smith block = mask[tmp] - 1; 1798de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1799de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1800b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1801b6490206SBarry Smith } 1802b6490206SBarry Smith } 180335aab85fSBarry Smith /* zero out the mask elements we set */ 180435aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1805b6490206SBarry Smith } 1806a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1807b6490206SBarry Smith 1808606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1809606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1810606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1811606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1812606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1813b6490206SBarry Smith 1814b6490206SBarry Smith B->assembled = PETSC_TRUE; 1815de6a44a3SBarry Smith 18169c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 18173a40ed3dSBarry Smith PetscFunctionReturn(0); 18182593348eSBarry Smith } 18192593348eSBarry Smith 18202593348eSBarry Smith 18212593348eSBarry Smith 1822