1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*606d414cSSatish Balay static char vcid[] = "$Id: baij.c,v 1.178 1999/06/08 22:56:09 balay Exp balay $"; 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; 93*606d414cSSatish Balay int i,n = a->mbs,ierr; 943b2fbd54SBarry Smith 953a40ed3dSBarry Smith PetscFunctionBegin; 963a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 973b2fbd54SBarry Smith if (symmetric) { 98*606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 99*606d414cSSatish 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->refct > 0) PetscFunctionReturn(0); 12894d884c6SBarry Smith 12994d884c6SBarry Smith if (A->mapping) { 13094d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 13194d884c6SBarry Smith } 13294d884c6SBarry Smith if (A->bmapping) { 13394d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 13494d884c6SBarry Smith } 13561b13de0SBarry Smith if (A->rmap) { 13661b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 13761b13de0SBarry Smith } 13861b13de0SBarry Smith if (A->cmap) { 13961b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 14061b13de0SBarry Smith } 141aa482453SBarry Smith #if defined(PETSC_USE_LOG) 142e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1432d61bbb3SSatish Balay #endif 144*606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 145*606d414cSSatish Balay if (!a->singlemalloc) { 146*606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 147*606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 148*606d414cSSatish Balay } 149*606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 150*606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 151*606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 152*606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 153*606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 154e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 155*606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 156*606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1572d61bbb3SSatish Balay PLogObjectDestroy(A); 1582d61bbb3SSatish Balay PetscHeaderDestroy(A); 1592d61bbb3SSatish Balay PetscFunctionReturn(0); 1602d61bbb3SSatish Balay } 1612d61bbb3SSatish Balay 1622d61bbb3SSatish Balay #undef __FUNC__ 1632d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1642d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1652d61bbb3SSatish Balay { 1662d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1672d61bbb3SSatish Balay 1682d61bbb3SSatish Balay PetscFunctionBegin; 1692d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1702d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1712d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1722d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1732d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1744787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1754787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1762d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1772d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1782d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1792d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1802d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1812d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 182b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 183b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 184981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1852d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1862d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1872d61bbb3SSatish Balay } else { 1882d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1892d61bbb3SSatish Balay } 1902d61bbb3SSatish Balay PetscFunctionReturn(0); 1912d61bbb3SSatish Balay } 1922d61bbb3SSatish Balay 1932d61bbb3SSatish Balay 1942d61bbb3SSatish Balay #undef __FUNC__ 1952d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1962d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1972d61bbb3SSatish Balay { 1982d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1992d61bbb3SSatish Balay 2002d61bbb3SSatish Balay PetscFunctionBegin; 2012d61bbb3SSatish Balay if (m) *m = a->m; 2022d61bbb3SSatish Balay if (n) *n = a->n; 2032d61bbb3SSatish Balay PetscFunctionReturn(0); 2042d61bbb3SSatish Balay } 2052d61bbb3SSatish Balay 2062d61bbb3SSatish Balay #undef __FUNC__ 2072d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 2082d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2092d61bbb3SSatish Balay { 2102d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2112d61bbb3SSatish Balay 2122d61bbb3SSatish Balay PetscFunctionBegin; 2132d61bbb3SSatish Balay *m = 0; *n = a->m; 2142d61bbb3SSatish Balay PetscFunctionReturn(0); 2152d61bbb3SSatish Balay } 2162d61bbb3SSatish Balay 2172d61bbb3SSatish Balay #undef __FUNC__ 2182d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 2192d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2202d61bbb3SSatish Balay { 2212d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2222d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2233f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2243f1db9ecSBarry Smith Scalar *v_i; 2252d61bbb3SSatish Balay 2262d61bbb3SSatish Balay PetscFunctionBegin; 2272d61bbb3SSatish Balay bs = a->bs; 2282d61bbb3SSatish Balay ai = a->i; 2292d61bbb3SSatish Balay aj = a->j; 2302d61bbb3SSatish Balay aa = a->a; 2312d61bbb3SSatish Balay bs2 = a->bs2; 2322d61bbb3SSatish Balay 2332d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2342d61bbb3SSatish Balay 2352d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2362d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2372d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2382d61bbb3SSatish Balay *nz = bs*M; 2392d61bbb3SSatish Balay 2402d61bbb3SSatish Balay if (v) { 2412d61bbb3SSatish Balay *v = 0; 2422d61bbb3SSatish Balay if (*nz) { 2432d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v); 2442d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2452d61bbb3SSatish Balay v_i = *v + i*bs; 2462d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2472d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2482d61bbb3SSatish Balay } 2492d61bbb3SSatish Balay } 2502d61bbb3SSatish Balay } 2512d61bbb3SSatish Balay 2522d61bbb3SSatish Balay if (idx) { 2532d61bbb3SSatish Balay *idx = 0; 2542d61bbb3SSatish Balay if (*nz) { 2552d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx); 2562d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2572d61bbb3SSatish Balay idx_i = *idx + i*bs; 2582d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2592d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2602d61bbb3SSatish Balay } 2612d61bbb3SSatish Balay } 2622d61bbb3SSatish Balay } 2632d61bbb3SSatish Balay PetscFunctionReturn(0); 2642d61bbb3SSatish Balay } 2652d61bbb3SSatish Balay 2662d61bbb3SSatish Balay #undef __FUNC__ 2672d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2682d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2692d61bbb3SSatish Balay { 270*606d414cSSatish Balay int ierr; 271*606d414cSSatish Balay 2722d61bbb3SSatish Balay PetscFunctionBegin; 273*606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 274*606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2752d61bbb3SSatish Balay PetscFunctionReturn(0); 2762d61bbb3SSatish Balay } 2772d61bbb3SSatish Balay 2782d61bbb3SSatish Balay #undef __FUNC__ 2792d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2802d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2812d61bbb3SSatish Balay { 2822d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2832d61bbb3SSatish Balay Mat C; 2842d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2852d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2863f1db9ecSBarry Smith MatScalar *array = a->a; 2872d61bbb3SSatish Balay 2882d61bbb3SSatish Balay PetscFunctionBegin; 2892d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2902d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 291549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2922d61bbb3SSatish Balay 2932d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2942d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 295*606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 2962d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 2972d61bbb3SSatish Balay cols = rows + bs; 2982d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2992d61bbb3SSatish Balay cols[0] = i*bs; 3002d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 3012d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3022d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 3032d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3042d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 3052d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3062d61bbb3SSatish Balay array += bs2; 3072d61bbb3SSatish Balay } 3082d61bbb3SSatish Balay } 309*606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 3102d61bbb3SSatish Balay 3112d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3122d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3132d61bbb3SSatish Balay 3142d61bbb3SSatish Balay if (B != PETSC_NULL) { 3152d61bbb3SSatish Balay *B = C; 3162d61bbb3SSatish Balay } else { 317f830108cSBarry Smith PetscOps *Abops; 318cc2dc46cSBarry Smith MatOps Aops; 319f830108cSBarry Smith 3202d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 321*606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 322*606d414cSSatish Balay if (!a->singlemalloc) { 323*606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 324*606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 325*606d414cSSatish Balay } 326*606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 327*606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 328*606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 329*606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 330*606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 331f830108cSBarry Smith 3327413bad6SBarry Smith 3337413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3347413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3357413bad6SBarry Smith 336f830108cSBarry Smith /* 337f830108cSBarry Smith This is horrible, horrible code. We need to keep the 338f830108cSBarry Smith A pointers for the bops and ops but copy everything 339f830108cSBarry Smith else from C. 340f830108cSBarry Smith */ 341f830108cSBarry Smith Abops = A->bops; 342f830108cSBarry Smith Aops = A->ops; 343549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 344f830108cSBarry Smith A->bops = Abops; 345f830108cSBarry Smith A->ops = Aops; 34627a8da17SBarry Smith A->qlist = 0; 347f830108cSBarry Smith 3482d61bbb3SSatish Balay PetscHeaderDestroy(C); 3492d61bbb3SSatish Balay } 3502d61bbb3SSatish Balay PetscFunctionReturn(0); 3512d61bbb3SSatish Balay } 3522d61bbb3SSatish Balay 3535615d1e5SSatish Balay #undef __FUNC__ 354d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 355b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3562593348eSBarry Smith { 357b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3589df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 359b6490206SBarry Smith Scalar *aa; 360ce6f0cecSBarry Smith FILE *file; 3612593348eSBarry Smith 3623a40ed3dSBarry Smith PetscFunctionBegin; 36390ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3642593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3652593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3663b2fbd54SBarry Smith 3672593348eSBarry Smith col_lens[1] = a->m; 3682593348eSBarry Smith col_lens[2] = a->n; 3697e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3702593348eSBarry Smith 3712593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 372b6490206SBarry Smith count = 0; 373b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 374b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 375b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 376b6490206SBarry Smith } 3772593348eSBarry Smith } 3780752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 379*606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3802593348eSBarry Smith 3812593348eSBarry Smith /* store column indices (zero start index) */ 38266963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj); 383b6490206SBarry Smith count = 0; 384b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 385b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 386b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 387b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 388b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3892593348eSBarry Smith } 3902593348eSBarry Smith } 391b6490206SBarry Smith } 392b6490206SBarry Smith } 3930752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 394*606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 3952593348eSBarry Smith 3962593348eSBarry Smith /* store nonzero values */ 39766963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 398b6490206SBarry Smith count = 0; 399b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 400b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 401b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 402b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4037e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 404b6490206SBarry Smith } 405b6490206SBarry Smith } 406b6490206SBarry Smith } 407b6490206SBarry Smith } 4080752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 409*606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 410ce6f0cecSBarry Smith 411ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 412ce6f0cecSBarry Smith if (file) { 413ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 414ce6f0cecSBarry Smith } 4153a40ed3dSBarry Smith PetscFunctionReturn(0); 4162593348eSBarry Smith } 4172593348eSBarry Smith 4185615d1e5SSatish Balay #undef __FUNC__ 419d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 420b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4212593348eSBarry Smith { 422b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4239df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4242593348eSBarry Smith FILE *fd; 4252593348eSBarry Smith char *outputname; 4262593348eSBarry Smith 4273a40ed3dSBarry Smith PetscFunctionBegin; 42890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 42977ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 430888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 431639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 432d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4330ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 4347b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket format not supported"); 4350ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 43644cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 43744cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 43844cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i*bs+j); 43944cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 44044cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 441aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4420ef38995SBarry Smith if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 44344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 444e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 4450ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 446766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 447e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 4480ef38995SBarry Smith } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 449e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 4500ef38995SBarry Smith } 45144cd7ae7SLois Curfman McInnes #else 4520ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 45344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 4540ef38995SBarry Smith } 45544cd7ae7SLois Curfman McInnes #endif 45644cd7ae7SLois Curfman McInnes } 45744cd7ae7SLois Curfman McInnes } 45844cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 45944cd7ae7SLois Curfman McInnes } 46044cd7ae7SLois Curfman McInnes } 4610ef38995SBarry Smith } else { 462b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 463b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 464b6490206SBarry Smith fprintf(fd,"row %d:",i*bs+j); 465b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 466b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 467aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 468e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 46988685aaeSLois Curfman McInnes fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 470e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j])); 4710ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 472766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",bs*a->j[k]+l, 473e20fef11SSatish Balay PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j])); 4740ef38995SBarry Smith } else { 475e20fef11SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j])); 47688685aaeSLois Curfman McInnes } 47788685aaeSLois Curfman McInnes #else 4787e67e3f9SSatish Balay fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 47988685aaeSLois Curfman McInnes #endif 4802593348eSBarry Smith } 4812593348eSBarry Smith } 4822593348eSBarry Smith fprintf(fd,"\n"); 4832593348eSBarry Smith } 4842593348eSBarry Smith } 485b6490206SBarry Smith } 4862593348eSBarry Smith fflush(fd); 4873a40ed3dSBarry Smith PetscFunctionReturn(0); 4882593348eSBarry Smith } 4892593348eSBarry Smith 4905615d1e5SSatish Balay #undef __FUNC__ 49177ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 49277ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4933270192aSSatish Balay { 49477ed5343SBarry Smith Mat A = (Mat) Aa; 4953270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4967dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 497fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 4983f1db9ecSBarry Smith MatScalar *aa; 49977ed5343SBarry Smith MPI_Comm comm; 50077ed5343SBarry Smith Viewer viewer; 5013270192aSSatish Balay 5023a40ed3dSBarry Smith PetscFunctionBegin; 50377ed5343SBarry Smith /* 50477ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 50577ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 50677ed5343SBarry Smith rest should return immediately. 50777ed5343SBarry Smith */ 50877ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 509d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 51077ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 5113270192aSSatish Balay 51277ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 51377ed5343SBarry Smith 51477ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 51577ed5343SBarry Smith 5163270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5173270192aSSatish Balay color = DRAW_BLUE; 5183270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5193270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5203270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5213270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5223270192aSSatish Balay aa = a->a + j*bs2; 5233270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5243270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5253270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 526b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5273270192aSSatish Balay } 5283270192aSSatish Balay } 5293270192aSSatish Balay } 5303270192aSSatish Balay } 5313270192aSSatish Balay color = DRAW_CYAN; 5323270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5333270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5343270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5353270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5363270192aSSatish Balay aa = a->a + j*bs2; 5373270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5383270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5393270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 540b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5413270192aSSatish Balay } 5423270192aSSatish Balay } 5433270192aSSatish Balay } 5443270192aSSatish Balay } 5453270192aSSatish Balay 5463270192aSSatish Balay color = DRAW_RED; 5473270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5483270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5493270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5503270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5513270192aSSatish Balay aa = a->a + j*bs2; 5523270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5533270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5543270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 555b34f160eSSatish Balay DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color); 5563270192aSSatish Balay } 5573270192aSSatish Balay } 5583270192aSSatish Balay } 5593270192aSSatish Balay } 56077ed5343SBarry Smith PetscFunctionReturn(0); 56177ed5343SBarry Smith } 5623270192aSSatish Balay 56377ed5343SBarry Smith #undef __FUNC__ 56477ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 56577ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 56677ed5343SBarry Smith { 56777ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5687dce120fSSatish Balay int ierr; 5697dce120fSSatish Balay double xl,yl,xr,yr,w,h; 57077ed5343SBarry Smith Draw draw; 57177ed5343SBarry Smith PetscTruth isnull; 5723270192aSSatish Balay 57377ed5343SBarry Smith PetscFunctionBegin; 57477ed5343SBarry Smith 57577ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 57677ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 57777ed5343SBarry Smith 57877ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 57977ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 58077ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5813270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 58277ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 58377ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5843a40ed3dSBarry Smith PetscFunctionReturn(0); 5853270192aSSatish Balay } 5863270192aSSatish Balay 5875615d1e5SSatish Balay #undef __FUNC__ 588d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 589e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5902593348eSBarry Smith { 59119bcc07fSBarry Smith ViewerType vtype; 59219bcc07fSBarry Smith int ierr; 5932593348eSBarry Smith 5943a40ed3dSBarry Smith PetscFunctionBegin; 59519bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 5967b2a1423SBarry Smith if (PetscTypeCompare(vtype,SOCKET_VIEWER)) { 5977b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 5983f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){ 5993a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6003f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 6013a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6023f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 6033a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6045cd90555SBarry Smith } else { 6055cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 6062593348eSBarry Smith } 6073a40ed3dSBarry Smith PetscFunctionReturn(0); 6082593348eSBarry Smith } 609b6490206SBarry Smith 610cd0e1443SSatish Balay 6115615d1e5SSatish Balay #undef __FUNC__ 6122d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6132d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 614cd0e1443SSatish Balay { 615cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6162d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6172d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6182d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6193f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 620cd0e1443SSatish Balay 6213a40ed3dSBarry Smith PetscFunctionBegin; 6222d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 623cd0e1443SSatish Balay row = im[k]; brow = row/bs; 624a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 625a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6262d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6272c3acbe9SBarry Smith nrow = ailen[brow]; 6282d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 629a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 630a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6312d61bbb3SSatish Balay col = in[l] ; 6322d61bbb3SSatish Balay bcol = col/bs; 6332d61bbb3SSatish Balay cidx = col%bs; 6342d61bbb3SSatish Balay ridx = row%bs; 6352d61bbb3SSatish Balay high = nrow; 6362d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6372d61bbb3SSatish Balay while (high-low > 5) { 638cd0e1443SSatish Balay t = (low+high)/2; 639cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 640cd0e1443SSatish Balay else low = t; 641cd0e1443SSatish Balay } 642cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 643cd0e1443SSatish Balay if (rp[i] > bcol) break; 644cd0e1443SSatish Balay if (rp[i] == bcol) { 6452d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6462d61bbb3SSatish Balay goto finished; 647cd0e1443SSatish Balay } 648cd0e1443SSatish Balay } 6492d61bbb3SSatish Balay *v++ = zero; 6502d61bbb3SSatish Balay finished:; 651cd0e1443SSatish Balay } 652cd0e1443SSatish Balay } 6533a40ed3dSBarry Smith PetscFunctionReturn(0); 654cd0e1443SSatish Balay } 655cd0e1443SSatish Balay 6562d61bbb3SSatish Balay 6575615d1e5SSatish Balay #undef __FUNC__ 65805a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 65992c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 66092c4ed94SBarry Smith { 66192c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6628a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 663dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 664549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6653f1db9ecSBarry Smith Scalar *value = v; 6663f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 66792c4ed94SBarry Smith 6683a40ed3dSBarry Smith PetscFunctionBegin; 6690e324ae4SSatish Balay if (roworiented) { 6700e324ae4SSatish Balay stepval = (n-1)*bs; 6710e324ae4SSatish Balay } else { 6720e324ae4SSatish Balay stepval = (m-1)*bs; 6730e324ae4SSatish Balay } 67492c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 67592c4ed94SBarry Smith row = im[k]; 6765ef9f2a5SBarry Smith if (row < 0) continue; 677aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 678a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 67992c4ed94SBarry Smith #endif 68092c4ed94SBarry Smith rp = aj + ai[row]; 68192c4ed94SBarry Smith ap = aa + bs2*ai[row]; 68292c4ed94SBarry Smith rmax = imax[row]; 68392c4ed94SBarry Smith nrow = ailen[row]; 68492c4ed94SBarry Smith low = 0; 68592c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6865ef9f2a5SBarry Smith if (in[l] < 0) continue; 687aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 688a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 68992c4ed94SBarry Smith #endif 69092c4ed94SBarry Smith col = in[l]; 69192c4ed94SBarry Smith if (roworiented) { 6920e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6930e324ae4SSatish Balay } else { 6940e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 69592c4ed94SBarry Smith } 69692c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 69792c4ed94SBarry Smith while (high-low > 7) { 69892c4ed94SBarry Smith t = (low+high)/2; 69992c4ed94SBarry Smith if (rp[t] > col) high = t; 70092c4ed94SBarry Smith else low = t; 70192c4ed94SBarry Smith } 70292c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 70392c4ed94SBarry Smith if (rp[i] > col) break; 70492c4ed94SBarry Smith if (rp[i] == col) { 7058a84c255SSatish Balay bap = ap + bs2*i; 7060e324ae4SSatish Balay if (roworiented) { 7078a84c255SSatish Balay if (is == ADD_VALUES) { 708dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 709dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7108a84c255SSatish Balay bap[jj] += *value++; 711dd9472c6SBarry Smith } 712dd9472c6SBarry Smith } 7130e324ae4SSatish Balay } else { 714dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 715dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7160e324ae4SSatish Balay bap[jj] = *value++; 7178a84c255SSatish Balay } 718dd9472c6SBarry Smith } 719dd9472c6SBarry Smith } 7200e324ae4SSatish Balay } else { 7210e324ae4SSatish Balay if (is == ADD_VALUES) { 722dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 723dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7240e324ae4SSatish Balay *bap++ += *value++; 725dd9472c6SBarry Smith } 726dd9472c6SBarry Smith } 7270e324ae4SSatish Balay } else { 728dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 729dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7300e324ae4SSatish Balay *bap++ = *value++; 7310e324ae4SSatish Balay } 7328a84c255SSatish Balay } 733dd9472c6SBarry Smith } 734dd9472c6SBarry Smith } 735f1241b54SBarry Smith goto noinsert2; 73692c4ed94SBarry Smith } 73792c4ed94SBarry Smith } 73889280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 739a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74092c4ed94SBarry Smith if (nrow >= rmax) { 74192c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 74292c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7433f1db9ecSBarry Smith MatScalar *new_a; 74492c4ed94SBarry Smith 745a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74689280ab3SLois Curfman McInnes 74792c4ed94SBarry Smith /* malloc new storage space */ 7483f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7493f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 75092c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 75192c4ed94SBarry Smith new_i = new_j + new_nz; 75292c4ed94SBarry Smith 75392c4ed94SBarry Smith /* copy over old data into new slots */ 75492c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 75592c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 756549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 75792c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 758549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 759549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 760549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 761549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 76292c4ed94SBarry Smith /* free up old matrix storage */ 763*606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 764*606d414cSSatish Balay if (!a->singlemalloc) { 765*606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 766*606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 767*606d414cSSatish Balay } 76892c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 76992c4ed94SBarry Smith a->singlemalloc = 1; 77092c4ed94SBarry Smith 77192c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 77292c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7733f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 77492c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 77592c4ed94SBarry Smith a->reallocs++; 77692c4ed94SBarry Smith a->nz++; 77792c4ed94SBarry Smith } 77892c4ed94SBarry Smith N = nrow++ - 1; 77992c4ed94SBarry Smith /* shift up all the later entries in this row */ 78092c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 78192c4ed94SBarry Smith rp[ii+1] = rp[ii]; 782549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 78392c4ed94SBarry Smith } 784549d3d68SSatish Balay if (N >= i) { 785549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 786549d3d68SSatish Balay } 78792c4ed94SBarry Smith rp[i] = col; 7888a84c255SSatish Balay bap = ap + bs2*i; 7890e324ae4SSatish Balay if (roworiented) { 790dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 791dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7920e324ae4SSatish Balay bap[jj] = *value++; 793dd9472c6SBarry Smith } 794dd9472c6SBarry Smith } 7950e324ae4SSatish Balay } else { 796dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 797dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7980e324ae4SSatish Balay *bap++ = *value++; 7990e324ae4SSatish Balay } 800dd9472c6SBarry Smith } 801dd9472c6SBarry Smith } 802f1241b54SBarry Smith noinsert2:; 80392c4ed94SBarry Smith low = i; 80492c4ed94SBarry Smith } 80592c4ed94SBarry Smith ailen[row] = nrow; 80692c4ed94SBarry Smith } 8073a40ed3dSBarry Smith PetscFunctionReturn(0); 80892c4ed94SBarry Smith } 80992c4ed94SBarry Smith 810f2501298SSatish Balay 8115615d1e5SSatish Balay #undef __FUNC__ 8125615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8138f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 814584200bdSSatish Balay { 815584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 816584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 817a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 818549d3d68SSatish Balay int mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr; 8193f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 820584200bdSSatish Balay 8213a40ed3dSBarry Smith PetscFunctionBegin; 8223a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 823584200bdSSatish Balay 82443ee02c3SBarry Smith if (m) rmax = ailen[0]; 825584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 826584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 827584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 828d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 829584200bdSSatish Balay if (fshift) { 830a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 831584200bdSSatish Balay N = ailen[i]; 832584200bdSSatish Balay for ( j=0; j<N; j++ ) { 833584200bdSSatish Balay ip[j-fshift] = ip[j]; 834549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 835584200bdSSatish Balay } 836584200bdSSatish Balay } 837584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 838584200bdSSatish Balay } 839584200bdSSatish Balay if (mbs) { 840584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 841584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 842584200bdSSatish Balay } 843584200bdSSatish Balay /* reset ilen and imax for each row */ 844584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 845584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 846584200bdSSatish Balay } 847a7c10996SSatish Balay a->nz = ai[mbs]; 848584200bdSSatish Balay 849584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 850584200bdSSatish Balay if (fshift && a->diag) { 851*606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 852584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 853584200bdSSatish Balay a->diag = 0; 854584200bdSSatish Balay } 8553d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 856219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8573d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 858584200bdSSatish Balay a->reallocs); 859d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 860e2f3b5e9SSatish Balay a->reallocs = 0; 8614e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8624e220ebcSLois Curfman McInnes 8633a40ed3dSBarry Smith PetscFunctionReturn(0); 864584200bdSSatish Balay } 865584200bdSSatish Balay 8665a838052SSatish Balay 867bea157c4SSatish Balay 868bea157c4SSatish Balay /* 869bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 870bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 871bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 872bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 873bea157c4SSatish Balay */ 8745615d1e5SSatish Balay #undef __FUNC__ 875bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks" 876bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 877d9b7c43dSSatish Balay { 878bea157c4SSatish Balay int i,j,k,row; 879bea157c4SSatish Balay PetscTruth flg; 8803a40ed3dSBarry Smith 881bea157c4SSatish Balay /* PetscFunctionBegin;*/ 882bea157c4SSatish Balay for ( i=0,j=0; i<n; j++ ) { 883bea157c4SSatish Balay row = idx[i]; 884bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 885bea157c4SSatish Balay sizes[j] = 1; 886bea157c4SSatish Balay i++; 887e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 888bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 889bea157c4SSatish Balay i++; 890bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 891bea157c4SSatish Balay flg = PETSC_TRUE; 892bea157c4SSatish Balay for ( k=1; k<bs; k++ ) { 893bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 894bea157c4SSatish Balay flg = PETSC_FALSE; 895bea157c4SSatish Balay break; 896d9b7c43dSSatish Balay } 897bea157c4SSatish Balay } 898bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 899bea157c4SSatish Balay sizes[j] = bs; 900bea157c4SSatish Balay i+= bs; 901bea157c4SSatish Balay } else { 902bea157c4SSatish Balay sizes[j] = 1; 903bea157c4SSatish Balay i++; 904bea157c4SSatish Balay } 905bea157c4SSatish Balay } 906bea157c4SSatish Balay } 907bea157c4SSatish Balay *bs_max = j; 9083a40ed3dSBarry Smith PetscFunctionReturn(0); 909d9b7c43dSSatish Balay } 910d9b7c43dSSatish Balay 9115615d1e5SSatish Balay #undef __FUNC__ 9125615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 9138f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 914d9b7c43dSSatish Balay { 915d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 916b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 917bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9183f1db9ecSBarry Smith Scalar zero = 0.0; 9193f1db9ecSBarry Smith MatScalar *aa; 920d9b7c43dSSatish Balay 9213a40ed3dSBarry Smith PetscFunctionBegin; 922d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 923d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 924d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 925d9b7c43dSSatish Balay 926bea157c4SSatish Balay /* allocate memory for rows,sizes */ 927bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 928bea157c4SSatish Balay sizes = rows + is_n; 929bea157c4SSatish Balay 930bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 931bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 932bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 933bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 934b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 935bea157c4SSatish Balay 936bea157c4SSatish Balay for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) { 937bea157c4SSatish Balay row = rows[j]; 938b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 939bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 940bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 941bea157c4SSatish Balay if (sizes[i] == bs) { 942bea157c4SSatish Balay if (diag) { 943bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 944bea157c4SSatish Balay baij->ilen[row/bs] = 1; 945bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 946bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 947a07cd24cSSatish Balay } 948bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 949bea157c4SSatish Balay for ( k=0; k<bs; k++ ) { 950bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 951bea157c4SSatish Balay } 952bea157c4SSatish Balay } else { /* (!diag) */ 953bea157c4SSatish Balay baij->ilen[row/bs] = 0; 954bea157c4SSatish Balay } /* end (!diag) */ 955bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 956aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 957bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 958bea157c4SSatish Balay #endif 959bea157c4SSatish Balay for ( k=0; k<count; k++ ) { 960d9b7c43dSSatish Balay aa[0] = zero; 961d9b7c43dSSatish Balay aa+=bs; 962d9b7c43dSSatish Balay } 963d9b7c43dSSatish Balay if (diag) { 964f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 965d9b7c43dSSatish Balay } 966d9b7c43dSSatish Balay } 967bea157c4SSatish Balay } 968bea157c4SSatish Balay 969*606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 9709a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9713a40ed3dSBarry Smith PetscFunctionReturn(0); 972d9b7c43dSSatish Balay } 9731c351548SSatish Balay 9745615d1e5SSatish Balay #undef __FUNC__ 9752d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9762d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9772d61bbb3SSatish Balay { 9782d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9792d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9802d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9812d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 982549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 9833f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9842d61bbb3SSatish Balay 9852d61bbb3SSatish Balay PetscFunctionBegin; 9862d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9872d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9885ef9f2a5SBarry Smith if (row < 0) continue; 989aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 99032e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 9912d61bbb3SSatish Balay #endif 9922d61bbb3SSatish Balay rp = aj + ai[brow]; 9932d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9942d61bbb3SSatish Balay rmax = imax[brow]; 9952d61bbb3SSatish Balay nrow = ailen[brow]; 9962d61bbb3SSatish Balay low = 0; 9972d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 9985ef9f2a5SBarry Smith if (in[l] < 0) continue; 999aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 100032e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 10012d61bbb3SSatish Balay #endif 10022d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10032d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10042d61bbb3SSatish Balay if (roworiented) { 10055ef9f2a5SBarry Smith value = v[l + k*n]; 10062d61bbb3SSatish Balay } else { 10072d61bbb3SSatish Balay value = v[k + l*m]; 10082d61bbb3SSatish Balay } 10092d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10102d61bbb3SSatish Balay while (high-low > 7) { 10112d61bbb3SSatish Balay t = (low+high)/2; 10122d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10132d61bbb3SSatish Balay else low = t; 10142d61bbb3SSatish Balay } 10152d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 10162d61bbb3SSatish Balay if (rp[i] > bcol) break; 10172d61bbb3SSatish Balay if (rp[i] == bcol) { 10182d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10192d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10202d61bbb3SSatish Balay else *bap = value; 10212d61bbb3SSatish Balay goto noinsert1; 10222d61bbb3SSatish Balay } 10232d61bbb3SSatish Balay } 10242d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10252d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10262d61bbb3SSatish Balay if (nrow >= rmax) { 10272d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10282d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10293f1db9ecSBarry Smith MatScalar *new_a; 10302d61bbb3SSatish Balay 10312d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10322d61bbb3SSatish Balay 10332d61bbb3SSatish Balay /* Malloc new storage space */ 10343f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10353f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 10362d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 10372d61bbb3SSatish Balay new_i = new_j + new_nz; 10382d61bbb3SSatish Balay 10392d61bbb3SSatish Balay /* copy over old data into new slots */ 10402d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10412d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1042549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10432d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1044549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1045549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1046549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1047549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10482d61bbb3SSatish Balay /* free up old matrix storage */ 1049*606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1050*606d414cSSatish Balay if (!a->singlemalloc) { 1051*606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1052*606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1053*606d414cSSatish Balay } 10542d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 10552d61bbb3SSatish Balay a->singlemalloc = 1; 10562d61bbb3SSatish Balay 10572d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10582d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10593f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10602d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10612d61bbb3SSatish Balay a->reallocs++; 10622d61bbb3SSatish Balay a->nz++; 10632d61bbb3SSatish Balay } 10642d61bbb3SSatish Balay N = nrow++ - 1; 10652d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10662d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10672d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1068549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10692d61bbb3SSatish Balay } 1070549d3d68SSatish Balay if (N>=i) { 1071549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1072549d3d68SSatish Balay } 10732d61bbb3SSatish Balay rp[i] = bcol; 10742d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10752d61bbb3SSatish Balay noinsert1:; 10762d61bbb3SSatish Balay low = i; 10772d61bbb3SSatish Balay } 10782d61bbb3SSatish Balay ailen[brow] = nrow; 10792d61bbb3SSatish Balay } 10802d61bbb3SSatish Balay PetscFunctionReturn(0); 10812d61bbb3SSatish Balay } 10822d61bbb3SSatish Balay 10832d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10842d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10852d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 10867b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 10877b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 10882d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10892d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10902d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10912d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10922d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10932d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10942d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10952d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 10962d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 10972d61bbb3SSatish Balay 10982d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 10992d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 11002d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 11012d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 11022d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 11032d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 110415091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 11052d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 11062d61bbb3SSatish Balay 11072d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 11082d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 11092d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 11102d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 11112d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 11122d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 111315091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11142d61bbb3SSatish Balay 11152d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11162d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11172d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11182d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11192d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 112015091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11212d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11222d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11232d61bbb3SSatish Balay 11242d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11252d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11262d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11272d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11282d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 112915091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11302d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11312d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11322d61bbb3SSatish Balay 11332d61bbb3SSatish Balay #undef __FUNC__ 11342d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 11355ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11362d61bbb3SSatish Balay { 11372d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 11382d61bbb3SSatish Balay Mat outA; 11392d61bbb3SSatish Balay int ierr; 1140667159a5SBarry Smith PetscTruth row_identity, col_identity; 11412d61bbb3SSatish Balay 11422d61bbb3SSatish Balay PetscFunctionBegin; 1143b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1144667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1145667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1146667159a5SBarry Smith if (!row_identity || !col_identity) { 1147b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1148667159a5SBarry Smith } 11492d61bbb3SSatish Balay 11502d61bbb3SSatish Balay outA = inA; 11512d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11522d61bbb3SSatish Balay a->row = row; 11532d61bbb3SSatish Balay a->col = col; 11542d61bbb3SSatish Balay 1155e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1156e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr); 11571e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1158e51c0b9cSSatish Balay 11592d61bbb3SSatish Balay if (!a->solve_work) { 11602d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 11612d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 11622d61bbb3SSatish Balay } 11632d61bbb3SSatish Balay 11642d61bbb3SSatish Balay if (!a->diag) { 11652d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA);CHKERRQ(ierr); 11662d61bbb3SSatish Balay } 11672d61bbb3SSatish Balay /* 116815091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 116915091d37SBarry Smith for ILU(0) factorization with natural ordering 11702d61bbb3SSatish Balay */ 117115091d37SBarry Smith switch (a->bs) { 117215091d37SBarry Smith case 2: 117315091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 117415091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 117515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 117615091d37SBarry Smith break; 117715091d37SBarry Smith case 3: 117815091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 117915091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 118015091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 118115091d37SBarry Smith break; 118215091d37SBarry Smith case 4: 1183667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1184f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 118515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 118615091d37SBarry Smith break; 118715091d37SBarry Smith case 5: 1188667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1189667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 119015091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 119115091d37SBarry Smith break; 119215091d37SBarry Smith case 6: 119315091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 119415091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 119515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 119615091d37SBarry Smith break; 119715091d37SBarry Smith case 7: 119815091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 119915091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 120015091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 120115091d37SBarry Smith break; 12022d61bbb3SSatish Balay } 12032d61bbb3SSatish Balay 1204667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1205667159a5SBarry Smith 12062d61bbb3SSatish Balay PetscFunctionReturn(0); 12072d61bbb3SSatish Balay } 12082d61bbb3SSatish Balay #undef __FUNC__ 1209d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1210ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1211ba4ca20aSSatish Balay { 1212ba4ca20aSSatish Balay static int called = 0; 1213ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1214d132466eSBarry Smith int ierr; 1215ba4ca20aSSatish Balay 12163a40ed3dSBarry Smith PetscFunctionBegin; 12173a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 1218d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1219d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12203a40ed3dSBarry Smith PetscFunctionReturn(0); 1221ba4ca20aSSatish Balay } 1222d9b7c43dSSatish Balay 1223fb2e594dSBarry Smith EXTERN_C_BEGIN 122427a8da17SBarry Smith #undef __FUNC__ 122527a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 122627a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 122727a8da17SBarry Smith { 122827a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 122927a8da17SBarry Smith int i,nz,n; 123027a8da17SBarry Smith 123127a8da17SBarry Smith PetscFunctionBegin; 123227a8da17SBarry Smith nz = baij->maxnz; 123327a8da17SBarry Smith n = baij->n; 123427a8da17SBarry Smith for (i=0; i<nz; i++) { 123527a8da17SBarry Smith baij->j[i] = indices[i]; 123627a8da17SBarry Smith } 123727a8da17SBarry Smith baij->nz = nz; 123827a8da17SBarry Smith for ( i=0; i<n; i++ ) { 123927a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 124027a8da17SBarry Smith } 124127a8da17SBarry Smith 124227a8da17SBarry Smith PetscFunctionReturn(0); 124327a8da17SBarry Smith } 1244fb2e594dSBarry Smith EXTERN_C_END 124527a8da17SBarry Smith 124627a8da17SBarry Smith #undef __FUNC__ 124727a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 124827a8da17SBarry Smith /*@ 124927a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 125027a8da17SBarry Smith in the matrix. 125127a8da17SBarry Smith 125227a8da17SBarry Smith Input Parameters: 125327a8da17SBarry Smith + mat - the SeqBAIJ matrix 125427a8da17SBarry Smith - indices - the column indices 125527a8da17SBarry Smith 125615091d37SBarry Smith Level: advanced 125715091d37SBarry Smith 125827a8da17SBarry Smith Notes: 125927a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 126027a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 126127a8da17SBarry Smith of the MatSetValues() operation. 126227a8da17SBarry Smith 126327a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 126427a8da17SBarry Smith MatCreateSeqBAIJ(). 126527a8da17SBarry Smith 126627a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 126727a8da17SBarry Smith 126827a8da17SBarry Smith @*/ 126927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 127027a8da17SBarry Smith { 127127a8da17SBarry Smith int ierr,(*f)(Mat,int *); 127227a8da17SBarry Smith 127327a8da17SBarry Smith PetscFunctionBegin; 127427a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1275549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 127627a8da17SBarry Smith if (f) { 127727a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 127827a8da17SBarry Smith } else { 127927a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 128027a8da17SBarry Smith } 128127a8da17SBarry Smith PetscFunctionReturn(0); 128227a8da17SBarry Smith } 128327a8da17SBarry Smith 12842593348eSBarry Smith /* -------------------------------------------------------------------*/ 1285cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1286cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1287cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1288cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1289cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1290cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1291cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1292cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1293cc2dc46cSBarry Smith 0, 1294cc2dc46cSBarry Smith 0, 1295cc2dc46cSBarry Smith 0, 1296cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1297cc2dc46cSBarry Smith 0, 1298b6490206SBarry Smith 0, 1299f2501298SSatish Balay MatTranspose_SeqBAIJ, 1300cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1301cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1302cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1303cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1304cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1305b6490206SBarry Smith 0, 1306cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1307cc2dc46cSBarry Smith 0, 1308cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1309cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1310cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1311cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1312cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1313cc2dc46cSBarry Smith 0, 1314cc2dc46cSBarry Smith 0, 1315cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1316cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1317cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1318cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1319cc2dc46cSBarry Smith 0, 1320cc2dc46cSBarry Smith 0, 1321cc2dc46cSBarry Smith 0, 13222e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1323cc2dc46cSBarry Smith 0, 1324cc2dc46cSBarry Smith 0, 1325cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1326cc2dc46cSBarry Smith 0, 1327cc2dc46cSBarry Smith 0, 1328cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1329cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1330cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1331cc2dc46cSBarry Smith 0, 1332cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1333cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1334cc2dc46cSBarry Smith 0, 1335cc2dc46cSBarry Smith 0, 1336cc2dc46cSBarry Smith 0, 1337cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 13383b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 133992c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1340cc2dc46cSBarry Smith 0, 1341cc2dc46cSBarry Smith 0, 1342cc2dc46cSBarry Smith 0, 1343cc2dc46cSBarry Smith 0, 1344cc2dc46cSBarry Smith 0, 1345cc2dc46cSBarry Smith 0, 1346d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1347cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1348cc2dc46cSBarry Smith 0, 1349cc2dc46cSBarry Smith 0, 1350cc2dc46cSBarry Smith MatGetMaps_Petsc}; 13512593348eSBarry Smith 13523e90b805SBarry Smith EXTERN_C_BEGIN 13533e90b805SBarry Smith #undef __FUNC__ 13543e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ" 13553e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 13563e90b805SBarry Smith { 13573e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 13583e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1359d132466eSBarry Smith int ierr; 13603e90b805SBarry Smith 13613e90b805SBarry Smith PetscFunctionBegin; 13623e90b805SBarry Smith if (aij->nonew != 1) { 13633e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13643e90b805SBarry Smith } 13653e90b805SBarry Smith 13663e90b805SBarry Smith /* allocate space for values if not already there */ 13673e90b805SBarry Smith if (!aij->saved_values) { 13683e90b805SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 13693e90b805SBarry Smith } 13703e90b805SBarry Smith 13713e90b805SBarry Smith /* copy values over */ 1372d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 13733e90b805SBarry Smith PetscFunctionReturn(0); 13743e90b805SBarry Smith } 13753e90b805SBarry Smith EXTERN_C_END 13763e90b805SBarry Smith 13773e90b805SBarry Smith EXTERN_C_BEGIN 13783e90b805SBarry Smith #undef __FUNC__ 137932e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ" 138032e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 13813e90b805SBarry Smith { 13823e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1383549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 13843e90b805SBarry Smith 13853e90b805SBarry Smith PetscFunctionBegin; 13863e90b805SBarry Smith if (aij->nonew != 1) { 13873e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13883e90b805SBarry Smith } 13893e90b805SBarry Smith if (!aij->saved_values) { 13903e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 13913e90b805SBarry Smith } 13923e90b805SBarry Smith 13933e90b805SBarry Smith /* copy values over */ 1394549d3d68SSatish Balay ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 13953e90b805SBarry Smith PetscFunctionReturn(0); 13963e90b805SBarry Smith } 13973e90b805SBarry Smith EXTERN_C_END 13983e90b805SBarry Smith 13995615d1e5SSatish Balay #undef __FUNC__ 14005615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 14012593348eSBarry Smith /*@C 140244cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 140344cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 140444cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 14057fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 14062bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 14072593348eSBarry Smith 1408db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1409db81eaa0SLois Curfman McInnes 14102593348eSBarry Smith Input Parameters: 1411db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1412b6490206SBarry Smith . bs - size of block 14132593348eSBarry Smith . m - number of rows 14142593348eSBarry Smith . n - number of columns 1415b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14167fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14172bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14182593348eSBarry Smith 14192593348eSBarry Smith Output Parameter: 14202593348eSBarry Smith . A - the matrix 14212593348eSBarry Smith 14220513a670SBarry Smith Options Database Keys: 1423db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1424db81eaa0SLois Curfman McInnes block calculations (much slower) 1425db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14260513a670SBarry Smith 142715091d37SBarry Smith Level: intermediate 142815091d37SBarry Smith 14292593348eSBarry Smith Notes: 143044cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 14312593348eSBarry Smith storage. That is, the stored row and column indices can begin at 143244cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 14332593348eSBarry Smith 14342593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 14352593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14362593348eSBarry Smith allocation. For additional details, see the users manual chapter on 14376da5968aSLois Curfman McInnes matrices. 14382593348eSBarry Smith 1439db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 14402593348eSBarry Smith @*/ 1441b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 14422593348eSBarry Smith { 14432593348eSBarry Smith Mat B; 1444b6490206SBarry Smith Mat_SeqBAIJ *b; 1445302e33e3SBarry Smith int i,len,ierr,flg,mbs,nbs,bs2,size; 14463b2fbd54SBarry Smith 14473a40ed3dSBarry Smith PetscFunctionBegin; 1448d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1449a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1450b6490206SBarry Smith 1451962fb4a0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1452302e33e3SBarry Smith mbs = m/bs; 1453302e33e3SBarry Smith nbs = n/bs; 1454302e33e3SBarry Smith bs2 = bs*bs; 14550513a670SBarry Smith 14563a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1457a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 14583a40ed3dSBarry Smith } 14592593348eSBarry Smith 14602593348eSBarry Smith *A = 0; 14613f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 14622593348eSBarry Smith PLogObjectCreate(B); 1463b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1464549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1465549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 14661a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 14671a6a6d98SBarry Smith if (!flg) { 14687fc0212eSBarry Smith switch (bs) { 14697fc0212eSBarry Smith case 1: 1470f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1471f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1472f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1473f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 14747fc0212eSBarry Smith break; 14754eeb42bcSBarry Smith case 2: 1476f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1477f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1478f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1479f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 14806d84be18SBarry Smith break; 1481f327dd97SBarry Smith case 3: 1482f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1483f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1484f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1485f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 14864eeb42bcSBarry Smith break; 14876d84be18SBarry Smith case 4: 1488f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1489f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1490f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1491f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 14926d84be18SBarry Smith break; 14936d84be18SBarry Smith case 5: 1494f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1495f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1496f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1497f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 14986d84be18SBarry Smith break; 149915091d37SBarry Smith case 6: 150015091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 150115091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 150215091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 150315091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 150415091d37SBarry Smith break; 150548e9ddb2SSatish Balay case 7: 1506f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1507f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1508f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 150948e9ddb2SSatish Balay break; 15107fc0212eSBarry Smith } 15111a6a6d98SBarry Smith } 1512e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1513e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15142593348eSBarry Smith B->factor = 0; 15152593348eSBarry Smith B->lupivotthreshold = 1.0; 151690f02eecSBarry Smith B->mapping = 0; 15172593348eSBarry Smith b->row = 0; 15182593348eSBarry Smith b->col = 0; 1519e51c0b9cSSatish Balay b->icol = 0; 15202593348eSBarry Smith b->reallocs = 0; 15213e90b805SBarry Smith b->saved_values = 0; 15222593348eSBarry Smith 152344cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 152444cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1525a5ae1ecdSBarry Smith 15267413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 15277413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1528a5ae1ecdSBarry Smith 1529b6490206SBarry Smith b->mbs = mbs; 1530f2501298SSatish Balay b->nbs = nbs; 1531b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax); 15322593348eSBarry Smith if (nnz == PETSC_NULL) { 1533119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 15342593348eSBarry Smith else if (nz <= 0) nz = 1; 1535b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1536b6490206SBarry Smith nz = nz*mbs; 15373a40ed3dSBarry Smith } else { 15382593348eSBarry Smith nz = 0; 1539b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 15402593348eSBarry Smith } 15412593348eSBarry Smith 15422593348eSBarry Smith /* allocate the matrix space */ 15433f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 15443f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a); 1545549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15467e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 1547549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15482593348eSBarry Smith b->i = b->j + nz; 15492593348eSBarry Smith b->singlemalloc = 1; 15502593348eSBarry Smith 1551de6a44a3SBarry Smith b->i[0] = 0; 1552b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 15532593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 15542593348eSBarry Smith } 15552593348eSBarry Smith 1556b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1557b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1558f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1559b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 15602593348eSBarry Smith 1561b6490206SBarry Smith b->bs = bs; 15629df24120SSatish Balay b->bs2 = bs2; 1563b6490206SBarry Smith b->mbs = mbs; 15642593348eSBarry Smith b->nz = 0; 156518e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 15662593348eSBarry Smith b->sorted = 0; 15672593348eSBarry Smith b->roworiented = 1; 15682593348eSBarry Smith b->nonew = 0; 15692593348eSBarry Smith b->diag = 0; 15702593348eSBarry Smith b->solve_work = 0; 1571de6a44a3SBarry Smith b->mult_work = 0; 15722593348eSBarry Smith b->spptr = 0; 15734e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 15744e220ebcSLois Curfman McInnes 15752593348eSBarry Smith *A = B; 15762593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 15772593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 157827a8da17SBarry Smith 15793e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 15803e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 15813e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 15823e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 15833e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 15843e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 158527a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 158627a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 158727a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 15883a40ed3dSBarry Smith PetscFunctionReturn(0); 15892593348eSBarry Smith } 15902593348eSBarry Smith 15915615d1e5SSatish Balay #undef __FUNC__ 15922e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 15932e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 15942593348eSBarry Smith { 15952593348eSBarry Smith Mat C; 1596b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 159727a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1598de6a44a3SBarry Smith 15993a40ed3dSBarry Smith PetscFunctionBegin; 1600a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 16012593348eSBarry Smith 16022593348eSBarry Smith *B = 0; 16033f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 16042593348eSBarry Smith PLogObjectCreate(C); 1605b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1606549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1607e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1608e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 16092593348eSBarry Smith C->factor = A->factor; 16102593348eSBarry Smith c->row = 0; 16112593348eSBarry Smith c->col = 0; 1612cac97260SSatish Balay c->icol = 0; 161332e87ba3SBarry Smith c->saved_values = 0; 16142593348eSBarry Smith C->assembled = PETSC_TRUE; 16152593348eSBarry Smith 161644cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 161744cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 161844cd7ae7SLois Curfman McInnes C->M = a->m; 161944cd7ae7SLois Curfman McInnes C->N = a->n; 162044cd7ae7SLois Curfman McInnes 1621b6490206SBarry Smith c->bs = a->bs; 16229df24120SSatish Balay c->bs2 = a->bs2; 1623b6490206SBarry Smith c->mbs = a->mbs; 1624b6490206SBarry Smith c->nbs = a->nbs; 16252593348eSBarry Smith 1626b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1627b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1628b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16292593348eSBarry Smith c->imax[i] = a->imax[i]; 16302593348eSBarry Smith c->ilen[i] = a->ilen[i]; 16312593348eSBarry Smith } 16322593348eSBarry Smith 16332593348eSBarry Smith /* allocate the matrix space */ 16342593348eSBarry Smith c->singlemalloc = 1; 16353f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 16363f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a); 16377e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1638de6a44a3SBarry Smith c->i = c->j + nz; 1639549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1640b6490206SBarry Smith if (mbs > 0) { 1641549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16422e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1643549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16442e8a6d31SBarry Smith } else { 1645549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16462593348eSBarry Smith } 16472593348eSBarry Smith } 16482593348eSBarry Smith 1649f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16502593348eSBarry Smith c->sorted = a->sorted; 16512593348eSBarry Smith c->roworiented = a->roworiented; 16522593348eSBarry Smith c->nonew = a->nonew; 16532593348eSBarry Smith 16542593348eSBarry Smith if (a->diag) { 1655b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag); 1656b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1657b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16582593348eSBarry Smith c->diag[i] = a->diag[i]; 16592593348eSBarry Smith } 166098305bb5SBarry Smith } else c->diag = 0; 16612593348eSBarry Smith c->nz = a->nz; 16622593348eSBarry Smith c->maxnz = a->maxnz; 16632593348eSBarry Smith c->solve_work = 0; 16642593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 16657fc0212eSBarry Smith c->mult_work = 0; 16662593348eSBarry Smith *B = C; 16677b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 16683a40ed3dSBarry Smith PetscFunctionReturn(0); 16692593348eSBarry Smith } 16702593348eSBarry Smith 16715615d1e5SSatish Balay #undef __FUNC__ 16725615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 167319bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 16742593348eSBarry Smith { 1675b6490206SBarry Smith Mat_SeqBAIJ *a; 16762593348eSBarry Smith Mat B; 1677de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1678b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 167935aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1680a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1681b6490206SBarry Smith Scalar *aa; 168219bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 16832593348eSBarry Smith 16843a40ed3dSBarry Smith PetscFunctionBegin; 16851a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1686de6a44a3SBarry Smith bs2 = bs*bs; 1687b6490206SBarry Smith 1688d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1689cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 169090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16910752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1692a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 16932593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 16942593348eSBarry Smith 1695d64ed03dSBarry Smith if (header[3] < 0) { 1696a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1697d64ed03dSBarry Smith } 1698d64ed03dSBarry Smith 1699a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 170035aab85fSBarry Smith 170135aab85fSBarry Smith /* 170235aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 170335aab85fSBarry Smith divisible by the blocksize 170435aab85fSBarry Smith */ 1705b6490206SBarry Smith mbs = M/bs; 170635aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 170735aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 170835aab85fSBarry Smith else mbs++; 170935aab85fSBarry Smith if (extra_rows) { 1710537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 171135aab85fSBarry Smith } 1712b6490206SBarry Smith 17132593348eSBarry Smith /* read in row lengths */ 171435aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 17150752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 171635aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 17172593348eSBarry Smith 1718b6490206SBarry Smith /* read in column indices */ 171935aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj); 17200752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 172135aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1722b6490206SBarry Smith 1723b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1724b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1725549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 172635aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask); 1727549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 172835aab85fSBarry Smith masked = mask + mbs; 1729b6490206SBarry Smith rowcount = 0; nzcount = 0; 1730b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 173135aab85fSBarry Smith nmask = 0; 1732b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1733b6490206SBarry Smith kmax = rowlengths[rowcount]; 1734b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 173535aab85fSBarry Smith tmp = jj[nzcount++]/bs; 173635aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1737b6490206SBarry Smith } 1738b6490206SBarry Smith rowcount++; 1739b6490206SBarry Smith } 174035aab85fSBarry Smith browlengths[i] += nmask; 174135aab85fSBarry Smith /* zero out the mask elements we set */ 174235aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1743b6490206SBarry Smith } 1744b6490206SBarry Smith 17452593348eSBarry Smith /* create our matrix */ 17463eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17472593348eSBarry Smith B = *A; 1748b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 17492593348eSBarry Smith 17502593348eSBarry Smith /* set matrix "i" values */ 1751de6a44a3SBarry Smith a->i[0] = 0; 1752b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1753b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1754b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17552593348eSBarry Smith } 1756b6490206SBarry Smith a->nz = 0; 1757b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 17582593348eSBarry Smith 1759b6490206SBarry Smith /* read in nonzero values */ 176035aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 17610752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 176235aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1763b6490206SBarry Smith 1764b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1765b6490206SBarry Smith nzcount = 0; jcount = 0; 1766b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1767b6490206SBarry Smith nzcountb = nzcount; 176835aab85fSBarry Smith nmask = 0; 1769b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1770b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1771b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 177235aab85fSBarry Smith tmp = jj[nzcount++]/bs; 177335aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1774b6490206SBarry Smith } 1775b6490206SBarry Smith rowcount++; 1776b6490206SBarry Smith } 1777de6a44a3SBarry Smith /* sort the masked values */ 177877c4ece6SBarry Smith PetscSortInt(nmask,masked); 1779de6a44a3SBarry Smith 1780b6490206SBarry Smith /* set "j" values into matrix */ 1781b6490206SBarry Smith maskcount = 1; 178235aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 178335aab85fSBarry Smith a->j[jcount++] = masked[j]; 1784de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1785b6490206SBarry Smith } 1786b6490206SBarry Smith /* set "a" values into matrix */ 1787de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1788b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1789b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1790b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1791de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1792de6a44a3SBarry Smith block = mask[tmp] - 1; 1793de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1794de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1795b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1796b6490206SBarry Smith } 1797b6490206SBarry Smith } 179835aab85fSBarry Smith /* zero out the mask elements we set */ 179935aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1800b6490206SBarry Smith } 1801a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1802b6490206SBarry Smith 1803*606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1804*606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1805*606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1806*606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1807*606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1808b6490206SBarry Smith 1809b6490206SBarry Smith B->assembled = PETSC_TRUE; 1810de6a44a3SBarry Smith 18119c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 18123a40ed3dSBarry Smith PetscFunctionReturn(0); 18132593348eSBarry Smith } 18142593348eSBarry Smith 18152593348eSBarry Smith 18162593348eSBarry Smith 1817