1*c38d4ed2SBarry Smith /*$Id: baij.c,v 1.188 1999/11/05 14:45:32 bsmith Exp bsmith $*/ 22593348eSBarry Smith 32593348eSBarry Smith /* 4b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 52593348eSBarry Smith matrix storage format. 62593348eSBarry Smith */ 73369ce9aSBarry Smith #include "sys.h" 870f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 91a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 101a6a6d98SBarry Smith #include "src/inline/spops.h" 113b547af2SSatish Balay 122d61bbb3SSatish Balay #define CHUNKSIZE 10 13de6a44a3SBarry Smith 14be5855fcSBarry Smith /* 15be5855fcSBarry Smith Checks for missing diagonals 16be5855fcSBarry Smith */ 17be5855fcSBarry Smith #undef __FUNC__ 18be5855fcSBarry Smith #define __FUNC__ "MatMissingDiag_SeqBAIJ" 19be5855fcSBarry Smith int MatMissingDiag_SeqBAIJ(Mat A) 20be5855fcSBarry Smith { 21be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 22be5855fcSBarry Smith int *diag = a->diag, *jj = a->j,i; 23be5855fcSBarry Smith 24be5855fcSBarry Smith PetscFunctionBegin; 250e8e8aceSBarry Smith for ( i=0; i<a->mbs; i++ ) { 26be5855fcSBarry Smith if (jj[diag[i]] != i) { 27be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 28be5855fcSBarry Smith } 29be5855fcSBarry Smith } 30be5855fcSBarry Smith PetscFunctionReturn(0); 31be5855fcSBarry Smith } 32be5855fcSBarry Smith 335615d1e5SSatish Balay #undef __FUNC__ 345615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ" 35de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 36de6a44a3SBarry Smith { 37de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 387fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 39de6a44a3SBarry Smith 403a40ed3dSBarry Smith PetscFunctionBegin; 41de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 42de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 437fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 44ecc615b2SBarry Smith diag[i] = a->i[i+1]; 45de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 46de6a44a3SBarry Smith if (a->j[j] == i) { 47de6a44a3SBarry Smith diag[i] = j; 48de6a44a3SBarry Smith break; 49de6a44a3SBarry Smith } 50de6a44a3SBarry Smith } 51de6a44a3SBarry Smith } 52de6a44a3SBarry Smith a->diag = diag; 533a40ed3dSBarry Smith PetscFunctionReturn(0); 54de6a44a3SBarry Smith } 552593348eSBarry Smith 562593348eSBarry Smith 573b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 583b2fbd54SBarry Smith 595615d1e5SSatish Balay #undef __FUNC__ 605615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 613b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 623b2fbd54SBarry Smith PetscTruth *done) 633b2fbd54SBarry Smith { 643b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 653b2fbd54SBarry Smith int ierr, n = a->mbs,i; 663b2fbd54SBarry Smith 673a40ed3dSBarry Smith PetscFunctionBegin; 683b2fbd54SBarry Smith *nn = n; 693a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 703b2fbd54SBarry Smith if (symmetric) { 713b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 723b2fbd54SBarry Smith } else if (oshift == 1) { 733b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 743b2fbd54SBarry Smith int nz = a->i[n] + 1; 753b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 763b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 773b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 783b2fbd54SBarry Smith } else { 793b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 803b2fbd54SBarry Smith } 813b2fbd54SBarry Smith 823a40ed3dSBarry Smith PetscFunctionReturn(0); 833b2fbd54SBarry Smith } 843b2fbd54SBarry Smith 855615d1e5SSatish Balay #undef __FUNC__ 86d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 873b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 883b2fbd54SBarry Smith PetscTruth *done) 893b2fbd54SBarry Smith { 903b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 91606d414cSSatish Balay int i,n = a->mbs,ierr; 923b2fbd54SBarry Smith 933a40ed3dSBarry Smith PetscFunctionBegin; 943a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 953b2fbd54SBarry Smith if (symmetric) { 96606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 97606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 98af5da2bfSBarry Smith } else if (oshift == 1) { 993b2fbd54SBarry Smith int nz = a->i[n]; 1003b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 1013b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 1023b2fbd54SBarry Smith } 1033a40ed3dSBarry Smith PetscFunctionReturn(0); 1043b2fbd54SBarry Smith } 1053b2fbd54SBarry Smith 1062d61bbb3SSatish Balay #undef __FUNC__ 1072d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 1082d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 1092d61bbb3SSatish Balay { 1102d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 1112d61bbb3SSatish Balay 1122d61bbb3SSatish Balay PetscFunctionBegin; 1132d61bbb3SSatish Balay *bs = baij->bs; 1142d61bbb3SSatish Balay PetscFunctionReturn(0); 1152d61bbb3SSatish Balay } 1162d61bbb3SSatish Balay 1172d61bbb3SSatish Balay 1182d61bbb3SSatish Balay #undef __FUNC__ 1192d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 120e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1212d61bbb3SSatish Balay { 1222d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 123e51c0b9cSSatish Balay int ierr; 1242d61bbb3SSatish Balay 125433994e6SBarry Smith PetscFunctionBegin; 12694d884c6SBarry Smith if (A->mapping) { 12794d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 12894d884c6SBarry Smith } 12994d884c6SBarry Smith if (A->bmapping) { 13094d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 13194d884c6SBarry Smith } 13261b13de0SBarry Smith if (A->rmap) { 13361b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 13461b13de0SBarry Smith } 13561b13de0SBarry Smith if (A->cmap) { 13661b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 13761b13de0SBarry Smith } 138aa482453SBarry Smith #if defined(PETSC_USE_LOG) 139e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1402d61bbb3SSatish Balay #endif 141606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 142606d414cSSatish Balay if (!a->singlemalloc) { 143606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 144606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 145606d414cSSatish Balay } 146*c38d4ed2SBarry Smith if (a->row) { 147*c38d4ed2SBarry Smith ierr = ISDestroy(a->row);CHKERRQ(ierr); 148*c38d4ed2SBarry Smith } 149*c38d4ed2SBarry Smith if (a->col) { 150*c38d4ed2SBarry Smith ierr = ISDestroy(a->col);CHKERRQ(ierr); 151*c38d4ed2SBarry Smith } 152606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 153606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 154606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 155606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 156606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 157e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 158606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 159606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1602d61bbb3SSatish Balay PLogObjectDestroy(A); 1612d61bbb3SSatish Balay PetscHeaderDestroy(A); 1622d61bbb3SSatish Balay PetscFunctionReturn(0); 1632d61bbb3SSatish Balay } 1642d61bbb3SSatish Balay 1652d61bbb3SSatish Balay #undef __FUNC__ 1662d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1672d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1682d61bbb3SSatish Balay { 1692d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1702d61bbb3SSatish Balay 1712d61bbb3SSatish Balay PetscFunctionBegin; 1722d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1732d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1742d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1752d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1762d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1774787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1784787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1792d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1802d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1812d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1822d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1832d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1842d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 185b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 186b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 187981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1882d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1892d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1902d61bbb3SSatish Balay } else { 1912d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1922d61bbb3SSatish Balay } 1932d61bbb3SSatish Balay PetscFunctionReturn(0); 1942d61bbb3SSatish Balay } 1952d61bbb3SSatish Balay 1962d61bbb3SSatish Balay 1972d61bbb3SSatish Balay #undef __FUNC__ 1982d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1992d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 2002d61bbb3SSatish Balay { 2012d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2022d61bbb3SSatish Balay 2032d61bbb3SSatish Balay PetscFunctionBegin; 2042d61bbb3SSatish Balay if (m) *m = a->m; 2052d61bbb3SSatish Balay if (n) *n = a->n; 2062d61bbb3SSatish Balay PetscFunctionReturn(0); 2072d61bbb3SSatish Balay } 2082d61bbb3SSatish Balay 2092d61bbb3SSatish Balay #undef __FUNC__ 2102d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 2112d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2122d61bbb3SSatish Balay { 2132d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2142d61bbb3SSatish Balay 2152d61bbb3SSatish Balay PetscFunctionBegin; 2162d61bbb3SSatish Balay *m = 0; *n = a->m; 2172d61bbb3SSatish Balay PetscFunctionReturn(0); 2182d61bbb3SSatish Balay } 2192d61bbb3SSatish Balay 2202d61bbb3SSatish Balay #undef __FUNC__ 2212d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 2222d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2232d61bbb3SSatish Balay { 2242d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2252d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2263f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2273f1db9ecSBarry Smith Scalar *v_i; 2282d61bbb3SSatish Balay 2292d61bbb3SSatish Balay PetscFunctionBegin; 2302d61bbb3SSatish Balay bs = a->bs; 2312d61bbb3SSatish Balay ai = a->i; 2322d61bbb3SSatish Balay aj = a->j; 2332d61bbb3SSatish Balay aa = a->a; 2342d61bbb3SSatish Balay bs2 = a->bs2; 2352d61bbb3SSatish Balay 2362d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2372d61bbb3SSatish Balay 2382d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2392d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2402d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2412d61bbb3SSatish Balay *nz = bs*M; 2422d61bbb3SSatish Balay 2432d61bbb3SSatish Balay if (v) { 2442d61bbb3SSatish Balay *v = 0; 2452d61bbb3SSatish Balay if (*nz) { 2462d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v); 2472d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2482d61bbb3SSatish Balay v_i = *v + i*bs; 2492d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2502d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2512d61bbb3SSatish Balay } 2522d61bbb3SSatish Balay } 2532d61bbb3SSatish Balay } 2542d61bbb3SSatish Balay 2552d61bbb3SSatish Balay if (idx) { 2562d61bbb3SSatish Balay *idx = 0; 2572d61bbb3SSatish Balay if (*nz) { 2582d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx); 2592d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2602d61bbb3SSatish Balay idx_i = *idx + i*bs; 2612d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2622d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2632d61bbb3SSatish Balay } 2642d61bbb3SSatish Balay } 2652d61bbb3SSatish Balay } 2662d61bbb3SSatish Balay PetscFunctionReturn(0); 2672d61bbb3SSatish Balay } 2682d61bbb3SSatish Balay 2692d61bbb3SSatish Balay #undef __FUNC__ 2702d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2712d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2722d61bbb3SSatish Balay { 273606d414cSSatish Balay int ierr; 274606d414cSSatish Balay 2752d61bbb3SSatish Balay PetscFunctionBegin; 276606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 277606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2782d61bbb3SSatish Balay PetscFunctionReturn(0); 2792d61bbb3SSatish Balay } 2802d61bbb3SSatish Balay 2812d61bbb3SSatish Balay #undef __FUNC__ 2822d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2832d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2842d61bbb3SSatish Balay { 2852d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2862d61bbb3SSatish Balay Mat C; 2872d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2882d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2893f1db9ecSBarry Smith MatScalar *array = a->a; 2902d61bbb3SSatish Balay 2912d61bbb3SSatish Balay PetscFunctionBegin; 2922d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2932d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 294549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2952d61bbb3SSatish Balay 2962d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2972d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 298606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 2992d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 3002d61bbb3SSatish Balay cols = rows + bs; 3012d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 3022d61bbb3SSatish Balay cols[0] = i*bs; 3032d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 3042d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3052d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 3062d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3072d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 3082d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3092d61bbb3SSatish Balay array += bs2; 3102d61bbb3SSatish Balay } 3112d61bbb3SSatish Balay } 312606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 3132d61bbb3SSatish Balay 3142d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3152d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3162d61bbb3SSatish Balay 3172d61bbb3SSatish Balay if (B != PETSC_NULL) { 3182d61bbb3SSatish Balay *B = C; 3192d61bbb3SSatish Balay } else { 320f830108cSBarry Smith PetscOps *Abops; 321cc2dc46cSBarry Smith MatOps Aops; 322f830108cSBarry Smith 3232d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 324606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 325606d414cSSatish Balay if (!a->singlemalloc) { 326606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 327606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 328606d414cSSatish Balay } 329606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 330606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 331606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 332606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 333606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 334f830108cSBarry Smith 3357413bad6SBarry Smith 3367413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3377413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3387413bad6SBarry Smith 339f830108cSBarry Smith /* 340f830108cSBarry Smith This is horrible, horrible code. We need to keep the 341f830108cSBarry Smith A pointers for the bops and ops but copy everything 342f830108cSBarry Smith else from C. 343f830108cSBarry Smith */ 344f830108cSBarry Smith Abops = A->bops; 345f830108cSBarry Smith Aops = A->ops; 346549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 347f830108cSBarry Smith A->bops = Abops; 348f830108cSBarry Smith A->ops = Aops; 34927a8da17SBarry Smith A->qlist = 0; 350f830108cSBarry Smith 3512d61bbb3SSatish Balay PetscHeaderDestroy(C); 3522d61bbb3SSatish Balay } 3532d61bbb3SSatish Balay PetscFunctionReturn(0); 3542d61bbb3SSatish Balay } 3552d61bbb3SSatish Balay 3565615d1e5SSatish Balay #undef __FUNC__ 357d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 358b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3592593348eSBarry Smith { 360b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3619df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 362b6490206SBarry Smith Scalar *aa; 363ce6f0cecSBarry Smith FILE *file; 3642593348eSBarry Smith 3653a40ed3dSBarry Smith PetscFunctionBegin; 36690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3672593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3682593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3693b2fbd54SBarry Smith 3702593348eSBarry Smith col_lens[1] = a->m; 3712593348eSBarry Smith col_lens[2] = a->n; 3727e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3732593348eSBarry Smith 3742593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 375b6490206SBarry Smith count = 0; 376b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 377b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 378b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 379b6490206SBarry Smith } 3802593348eSBarry Smith } 3810752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 382606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3832593348eSBarry Smith 3842593348eSBarry Smith /* store column indices (zero start index) */ 38566963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj); 386b6490206SBarry Smith count = 0; 387b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 388b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 389b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 390b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 391b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3922593348eSBarry Smith } 3932593348eSBarry Smith } 394b6490206SBarry Smith } 395b6490206SBarry Smith } 3960752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 397606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 3982593348eSBarry Smith 3992593348eSBarry Smith /* store nonzero values */ 40066963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 401b6490206SBarry Smith count = 0; 402b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 403b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 404b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 405b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4067e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 407b6490206SBarry Smith } 408b6490206SBarry Smith } 409b6490206SBarry Smith } 410b6490206SBarry Smith } 4110752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 412606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 413ce6f0cecSBarry Smith 414ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 415ce6f0cecSBarry Smith if (file) { 416ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 417ce6f0cecSBarry Smith } 4183a40ed3dSBarry Smith PetscFunctionReturn(0); 4192593348eSBarry Smith } 4202593348eSBarry Smith 4215615d1e5SSatish Balay #undef __FUNC__ 422d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 423b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4242593348eSBarry Smith { 425b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4269df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4272593348eSBarry Smith char *outputname; 4282593348eSBarry Smith 4293a40ed3dSBarry Smith PetscFunctionBegin; 43077ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 431888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 432639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 433d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4340ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 4356831982aSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 4360ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 4376831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 43844cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 43944cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 4406831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 44144cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 44244cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 443aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4440ef38995SBarry Smith if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 4456831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 4466831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4470ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 4486831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 4496831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4500ef38995SBarry Smith } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 4516831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4520ef38995SBarry Smith } 45344cd7ae7SLois Curfman McInnes #else 4540ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 4556831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 4560ef38995SBarry Smith } 45744cd7ae7SLois Curfman McInnes #endif 45844cd7ae7SLois Curfman McInnes } 45944cd7ae7SLois Curfman McInnes } 4606831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 46144cd7ae7SLois Curfman McInnes } 46244cd7ae7SLois Curfman McInnes } 4636831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 4640ef38995SBarry Smith } else { 4656831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 466b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 467b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 4686831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 469b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 470b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 471aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 472e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 4736831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 4746831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4750ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 4766831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 4776831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4780ef38995SBarry Smith } else { 4796831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 48088685aaeSLois Curfman McInnes } 48188685aaeSLois Curfman McInnes #else 4826831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 48388685aaeSLois Curfman McInnes #endif 4842593348eSBarry Smith } 4852593348eSBarry Smith } 4866831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 4872593348eSBarry Smith } 4882593348eSBarry Smith } 4896831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 490b6490206SBarry Smith } 4916831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 4923a40ed3dSBarry Smith PetscFunctionReturn(0); 4932593348eSBarry Smith } 4942593348eSBarry Smith 4955615d1e5SSatish Balay #undef __FUNC__ 49677ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 49777ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4983270192aSSatish Balay { 49977ed5343SBarry Smith Mat A = (Mat) Aa; 5003270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5017dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 502fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 5033f1db9ecSBarry Smith MatScalar *aa; 50477ed5343SBarry Smith MPI_Comm comm; 50577ed5343SBarry Smith Viewer viewer; 5063270192aSSatish Balay 5073a40ed3dSBarry Smith PetscFunctionBegin; 50877ed5343SBarry Smith /* 50977ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 51077ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 51177ed5343SBarry Smith rest should return immediately. 51277ed5343SBarry Smith */ 51377ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 514d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 51577ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 5163270192aSSatish Balay 51777ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 51877ed5343SBarry Smith 51977ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 52077ed5343SBarry Smith 5213270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5223270192aSSatish Balay color = DRAW_BLUE; 5233270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5243270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5253270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5263270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5273270192aSSatish Balay aa = a->a + j*bs2; 5283270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5293270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5303270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 531433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5323270192aSSatish Balay } 5333270192aSSatish Balay } 5343270192aSSatish Balay } 5353270192aSSatish Balay } 5363270192aSSatish Balay color = DRAW_CYAN; 5373270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5383270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5393270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5403270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5413270192aSSatish Balay aa = a->a + j*bs2; 5423270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5433270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5443270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 545433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5463270192aSSatish Balay } 5473270192aSSatish Balay } 5483270192aSSatish Balay } 5493270192aSSatish Balay } 5503270192aSSatish Balay 5513270192aSSatish Balay color = DRAW_RED; 5523270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5533270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5543270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5553270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5563270192aSSatish Balay aa = a->a + j*bs2; 5573270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5583270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5593270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 560433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5613270192aSSatish Balay } 5623270192aSSatish Balay } 5633270192aSSatish Balay } 5643270192aSSatish Balay } 56577ed5343SBarry Smith PetscFunctionReturn(0); 56677ed5343SBarry Smith } 5673270192aSSatish Balay 56877ed5343SBarry Smith #undef __FUNC__ 56977ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 57077ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 57177ed5343SBarry Smith { 57277ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5737dce120fSSatish Balay int ierr; 5747dce120fSSatish Balay double xl,yl,xr,yr,w,h; 57577ed5343SBarry Smith Draw draw; 57677ed5343SBarry Smith PetscTruth isnull; 5773270192aSSatish Balay 57877ed5343SBarry Smith PetscFunctionBegin; 57977ed5343SBarry Smith 58077ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 58177ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 58277ed5343SBarry Smith 58377ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 58477ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 58577ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5863270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 58777ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 58877ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5893a40ed3dSBarry Smith PetscFunctionReturn(0); 5903270192aSSatish Balay } 5913270192aSSatish Balay 5925615d1e5SSatish Balay #undef __FUNC__ 593d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 594e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5952593348eSBarry Smith { 59619bcc07fSBarry Smith int ierr; 5976831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 5982593348eSBarry Smith 5993a40ed3dSBarry Smith PetscFunctionBegin; 6006831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 6016831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 6026831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 6036831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 6040f5bd95cSBarry Smith if (issocket) { 6057b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 6060f5bd95cSBarry Smith } else if (isascii){ 6073a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6080f5bd95cSBarry Smith } else if (isbinary) { 6093a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6100f5bd95cSBarry Smith } else if (isdraw) { 6113a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6125cd90555SBarry Smith } else { 6130f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 6142593348eSBarry Smith } 6153a40ed3dSBarry Smith PetscFunctionReturn(0); 6162593348eSBarry Smith } 617b6490206SBarry Smith 618cd0e1443SSatish Balay 6195615d1e5SSatish Balay #undef __FUNC__ 6202d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6212d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 622cd0e1443SSatish Balay { 623cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6242d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6252d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6262d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6273f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 628cd0e1443SSatish Balay 6293a40ed3dSBarry Smith PetscFunctionBegin; 6302d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 631cd0e1443SSatish Balay row = im[k]; brow = row/bs; 632a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 633a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6342d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6352c3acbe9SBarry Smith nrow = ailen[brow]; 6362d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 637a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 638a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6392d61bbb3SSatish Balay col = in[l] ; 6402d61bbb3SSatish Balay bcol = col/bs; 6412d61bbb3SSatish Balay cidx = col%bs; 6422d61bbb3SSatish Balay ridx = row%bs; 6432d61bbb3SSatish Balay high = nrow; 6442d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6452d61bbb3SSatish Balay while (high-low > 5) { 646cd0e1443SSatish Balay t = (low+high)/2; 647cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 648cd0e1443SSatish Balay else low = t; 649cd0e1443SSatish Balay } 650cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 651cd0e1443SSatish Balay if (rp[i] > bcol) break; 652cd0e1443SSatish Balay if (rp[i] == bcol) { 6532d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6542d61bbb3SSatish Balay goto finished; 655cd0e1443SSatish Balay } 656cd0e1443SSatish Balay } 6572d61bbb3SSatish Balay *v++ = zero; 6582d61bbb3SSatish Balay finished:; 659cd0e1443SSatish Balay } 660cd0e1443SSatish Balay } 6613a40ed3dSBarry Smith PetscFunctionReturn(0); 662cd0e1443SSatish Balay } 663cd0e1443SSatish Balay 6642d61bbb3SSatish Balay 6655615d1e5SSatish Balay #undef __FUNC__ 66605a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 66792c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 66892c4ed94SBarry Smith { 66992c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6708a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 671dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 672549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6733f1db9ecSBarry Smith Scalar *value = v; 6743f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 67592c4ed94SBarry Smith 6763a40ed3dSBarry Smith PetscFunctionBegin; 6770e324ae4SSatish Balay if (roworiented) { 6780e324ae4SSatish Balay stepval = (n-1)*bs; 6790e324ae4SSatish Balay } else { 6800e324ae4SSatish Balay stepval = (m-1)*bs; 6810e324ae4SSatish Balay } 68292c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 68392c4ed94SBarry Smith row = im[k]; 6845ef9f2a5SBarry Smith if (row < 0) continue; 685aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 686a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 68792c4ed94SBarry Smith #endif 68892c4ed94SBarry Smith rp = aj + ai[row]; 68992c4ed94SBarry Smith ap = aa + bs2*ai[row]; 69092c4ed94SBarry Smith rmax = imax[row]; 69192c4ed94SBarry Smith nrow = ailen[row]; 69292c4ed94SBarry Smith low = 0; 69392c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6945ef9f2a5SBarry Smith if (in[l] < 0) continue; 695aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 696a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 69792c4ed94SBarry Smith #endif 69892c4ed94SBarry Smith col = in[l]; 69992c4ed94SBarry Smith if (roworiented) { 7000e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 7010e324ae4SSatish Balay } else { 7020e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 70392c4ed94SBarry Smith } 70492c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 70592c4ed94SBarry Smith while (high-low > 7) { 70692c4ed94SBarry Smith t = (low+high)/2; 70792c4ed94SBarry Smith if (rp[t] > col) high = t; 70892c4ed94SBarry Smith else low = t; 70992c4ed94SBarry Smith } 71092c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 71192c4ed94SBarry Smith if (rp[i] > col) break; 71292c4ed94SBarry Smith if (rp[i] == col) { 7138a84c255SSatish Balay bap = ap + bs2*i; 7140e324ae4SSatish Balay if (roworiented) { 7158a84c255SSatish Balay if (is == ADD_VALUES) { 716dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 717dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7188a84c255SSatish Balay bap[jj] += *value++; 719dd9472c6SBarry Smith } 720dd9472c6SBarry Smith } 7210e324ae4SSatish Balay } else { 722dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 723dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7240e324ae4SSatish Balay bap[jj] = *value++; 7258a84c255SSatish Balay } 726dd9472c6SBarry Smith } 727dd9472c6SBarry Smith } 7280e324ae4SSatish Balay } else { 7290e324ae4SSatish Balay if (is == ADD_VALUES) { 730dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 731dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7320e324ae4SSatish Balay *bap++ += *value++; 733dd9472c6SBarry Smith } 734dd9472c6SBarry Smith } 7350e324ae4SSatish Balay } else { 736dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 737dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7380e324ae4SSatish Balay *bap++ = *value++; 7390e324ae4SSatish Balay } 7408a84c255SSatish Balay } 741dd9472c6SBarry Smith } 742dd9472c6SBarry Smith } 743f1241b54SBarry Smith goto noinsert2; 74492c4ed94SBarry Smith } 74592c4ed94SBarry Smith } 74689280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 747a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74892c4ed94SBarry Smith if (nrow >= rmax) { 74992c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 75092c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7513f1db9ecSBarry Smith MatScalar *new_a; 75292c4ed94SBarry Smith 753a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 75489280ab3SLois Curfman McInnes 75592c4ed94SBarry Smith /* malloc new storage space */ 7563f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7573f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 75892c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 75992c4ed94SBarry Smith new_i = new_j + new_nz; 76092c4ed94SBarry Smith 76192c4ed94SBarry Smith /* copy over old data into new slots */ 76292c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 76392c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 764549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 76592c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 766549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 767549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 768549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 769549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 77092c4ed94SBarry Smith /* free up old matrix storage */ 771606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 772606d414cSSatish Balay if (!a->singlemalloc) { 773606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 774606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 775606d414cSSatish Balay } 77692c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 77792c4ed94SBarry Smith a->singlemalloc = 1; 77892c4ed94SBarry Smith 77992c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 78092c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7813f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 78292c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 78392c4ed94SBarry Smith a->reallocs++; 78492c4ed94SBarry Smith a->nz++; 78592c4ed94SBarry Smith } 78692c4ed94SBarry Smith N = nrow++ - 1; 78792c4ed94SBarry Smith /* shift up all the later entries in this row */ 78892c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 78992c4ed94SBarry Smith rp[ii+1] = rp[ii]; 790549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 79192c4ed94SBarry Smith } 792549d3d68SSatish Balay if (N >= i) { 793549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 794549d3d68SSatish Balay } 79592c4ed94SBarry Smith rp[i] = col; 7968a84c255SSatish Balay bap = ap + bs2*i; 7970e324ae4SSatish Balay if (roworiented) { 798dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 799dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 8000e324ae4SSatish Balay bap[jj] = *value++; 801dd9472c6SBarry Smith } 802dd9472c6SBarry Smith } 8030e324ae4SSatish Balay } else { 804dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 805dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 8060e324ae4SSatish Balay *bap++ = *value++; 8070e324ae4SSatish Balay } 808dd9472c6SBarry Smith } 809dd9472c6SBarry Smith } 810f1241b54SBarry Smith noinsert2:; 81192c4ed94SBarry Smith low = i; 81292c4ed94SBarry Smith } 81392c4ed94SBarry Smith ailen[row] = nrow; 81492c4ed94SBarry Smith } 8153a40ed3dSBarry Smith PetscFunctionReturn(0); 81692c4ed94SBarry Smith } 81792c4ed94SBarry Smith 818f2501298SSatish Balay 8195615d1e5SSatish Balay #undef __FUNC__ 8205615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8218f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 822584200bdSSatish Balay { 823584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 824584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 825a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 826549d3d68SSatish Balay int mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr; 8273f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 828584200bdSSatish Balay 8293a40ed3dSBarry Smith PetscFunctionBegin; 8303a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 831584200bdSSatish Balay 83243ee02c3SBarry Smith if (m) rmax = ailen[0]; 833584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 834584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 835584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 836d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 837584200bdSSatish Balay if (fshift) { 838a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 839584200bdSSatish Balay N = ailen[i]; 840584200bdSSatish Balay for ( j=0; j<N; j++ ) { 841584200bdSSatish Balay ip[j-fshift] = ip[j]; 842549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 843584200bdSSatish Balay } 844584200bdSSatish Balay } 845584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 846584200bdSSatish Balay } 847584200bdSSatish Balay if (mbs) { 848584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 849584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 850584200bdSSatish Balay } 851584200bdSSatish Balay /* reset ilen and imax for each row */ 852584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 853584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 854584200bdSSatish Balay } 855a7c10996SSatish Balay a->nz = ai[mbs]; 856584200bdSSatish Balay 857584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 858584200bdSSatish Balay if (fshift && a->diag) { 859606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 860584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 861584200bdSSatish Balay a->diag = 0; 862584200bdSSatish Balay } 8633d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 864219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8653d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 866584200bdSSatish Balay a->reallocs); 867d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 868e2f3b5e9SSatish Balay a->reallocs = 0; 8694e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8704e220ebcSLois Curfman McInnes 8713a40ed3dSBarry Smith PetscFunctionReturn(0); 872584200bdSSatish Balay } 873584200bdSSatish Balay 8745a838052SSatish Balay 875bea157c4SSatish Balay 876bea157c4SSatish Balay /* 877bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 878bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 879bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 880bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 881bea157c4SSatish Balay */ 8825615d1e5SSatish Balay #undef __FUNC__ 883bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks" 884bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 885d9b7c43dSSatish Balay { 886bea157c4SSatish Balay int i,j,k,row; 887bea157c4SSatish Balay PetscTruth flg; 8883a40ed3dSBarry Smith 889433994e6SBarry Smith PetscFunctionBegin; 890bea157c4SSatish Balay for ( i=0,j=0; i<n; j++ ) { 891bea157c4SSatish Balay row = idx[i]; 892bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 893bea157c4SSatish Balay sizes[j] = 1; 894bea157c4SSatish Balay i++; 895e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 896bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 897bea157c4SSatish Balay i++; 898bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 899bea157c4SSatish Balay flg = PETSC_TRUE; 900bea157c4SSatish Balay for ( k=1; k<bs; k++ ) { 901bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 902bea157c4SSatish Balay flg = PETSC_FALSE; 903bea157c4SSatish Balay break; 904d9b7c43dSSatish Balay } 905bea157c4SSatish Balay } 906bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 907bea157c4SSatish Balay sizes[j] = bs; 908bea157c4SSatish Balay i+= bs; 909bea157c4SSatish Balay } else { 910bea157c4SSatish Balay sizes[j] = 1; 911bea157c4SSatish Balay i++; 912bea157c4SSatish Balay } 913bea157c4SSatish Balay } 914bea157c4SSatish Balay } 915bea157c4SSatish Balay *bs_max = j; 9163a40ed3dSBarry Smith PetscFunctionReturn(0); 917d9b7c43dSSatish Balay } 918d9b7c43dSSatish Balay 9195615d1e5SSatish Balay #undef __FUNC__ 9205615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 9218f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 922d9b7c43dSSatish Balay { 923d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 924b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 925bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9263f1db9ecSBarry Smith Scalar zero = 0.0; 9273f1db9ecSBarry Smith MatScalar *aa; 928d9b7c43dSSatish Balay 9293a40ed3dSBarry Smith PetscFunctionBegin; 930d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 931d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 932d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 933d9b7c43dSSatish Balay 934bea157c4SSatish Balay /* allocate memory for rows,sizes */ 935bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 936bea157c4SSatish Balay sizes = rows + is_n; 937bea157c4SSatish Balay 938bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 939bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 940bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 941bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 942b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 943bea157c4SSatish Balay 944bea157c4SSatish Balay for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) { 945bea157c4SSatish Balay row = rows[j]; 946b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 947bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 948bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 949bea157c4SSatish Balay if (sizes[i] == bs) { 950bea157c4SSatish Balay if (diag) { 951bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 952bea157c4SSatish Balay baij->ilen[row/bs] = 1; 953bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 954bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 955a07cd24cSSatish Balay } 956bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 957bea157c4SSatish Balay for ( k=0; k<bs; k++ ) { 958bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 959bea157c4SSatish Balay } 960bea157c4SSatish Balay } else { /* (!diag) */ 961bea157c4SSatish Balay baij->ilen[row/bs] = 0; 962bea157c4SSatish Balay } /* end (!diag) */ 963bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 964aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 965bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 966bea157c4SSatish Balay #endif 967bea157c4SSatish Balay for ( k=0; k<count; k++ ) { 968d9b7c43dSSatish Balay aa[0] = zero; 969d9b7c43dSSatish Balay aa+=bs; 970d9b7c43dSSatish Balay } 971d9b7c43dSSatish Balay if (diag) { 972f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 973d9b7c43dSSatish Balay } 974d9b7c43dSSatish Balay } 975bea157c4SSatish Balay } 976bea157c4SSatish Balay 977606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 9789a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9793a40ed3dSBarry Smith PetscFunctionReturn(0); 980d9b7c43dSSatish Balay } 9811c351548SSatish Balay 9825615d1e5SSatish Balay #undef __FUNC__ 9832d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9842d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9852d61bbb3SSatish Balay { 9862d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9872d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9882d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9892d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 990549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 9913f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9922d61bbb3SSatish Balay 9932d61bbb3SSatish Balay PetscFunctionBegin; 9942d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9952d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9965ef9f2a5SBarry Smith if (row < 0) continue; 997aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 99832e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 9992d61bbb3SSatish Balay #endif 10002d61bbb3SSatish Balay rp = aj + ai[brow]; 10012d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 10022d61bbb3SSatish Balay rmax = imax[brow]; 10032d61bbb3SSatish Balay nrow = ailen[brow]; 10042d61bbb3SSatish Balay low = 0; 10052d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 10065ef9f2a5SBarry Smith if (in[l] < 0) continue; 1007aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 100832e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 10092d61bbb3SSatish Balay #endif 10102d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10112d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10122d61bbb3SSatish Balay if (roworiented) { 10135ef9f2a5SBarry Smith value = v[l + k*n]; 10142d61bbb3SSatish Balay } else { 10152d61bbb3SSatish Balay value = v[k + l*m]; 10162d61bbb3SSatish Balay } 10172d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10182d61bbb3SSatish Balay while (high-low > 7) { 10192d61bbb3SSatish Balay t = (low+high)/2; 10202d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10212d61bbb3SSatish Balay else low = t; 10222d61bbb3SSatish Balay } 10232d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 10242d61bbb3SSatish Balay if (rp[i] > bcol) break; 10252d61bbb3SSatish Balay if (rp[i] == bcol) { 10262d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10272d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10282d61bbb3SSatish Balay else *bap = value; 10292d61bbb3SSatish Balay goto noinsert1; 10302d61bbb3SSatish Balay } 10312d61bbb3SSatish Balay } 10322d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10332d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10342d61bbb3SSatish Balay if (nrow >= rmax) { 10352d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10362d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10373f1db9ecSBarry Smith MatScalar *new_a; 10382d61bbb3SSatish Balay 10392d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10402d61bbb3SSatish Balay 10412d61bbb3SSatish Balay /* Malloc new storage space */ 10423f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10433f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 10442d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 10452d61bbb3SSatish Balay new_i = new_j + new_nz; 10462d61bbb3SSatish Balay 10472d61bbb3SSatish Balay /* copy over old data into new slots */ 10482d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10492d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1050549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10512d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1052549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1053549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1054549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1055549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10562d61bbb3SSatish Balay /* free up old matrix storage */ 1057606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1058606d414cSSatish Balay if (!a->singlemalloc) { 1059606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1060606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1061606d414cSSatish Balay } 10622d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 10632d61bbb3SSatish Balay a->singlemalloc = 1; 10642d61bbb3SSatish Balay 10652d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10662d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10673f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10682d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10692d61bbb3SSatish Balay a->reallocs++; 10702d61bbb3SSatish Balay a->nz++; 10712d61bbb3SSatish Balay } 10722d61bbb3SSatish Balay N = nrow++ - 1; 10732d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10742d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10752d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1076549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10772d61bbb3SSatish Balay } 1078549d3d68SSatish Balay if (N>=i) { 1079549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1080549d3d68SSatish Balay } 10812d61bbb3SSatish Balay rp[i] = bcol; 10822d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10832d61bbb3SSatish Balay noinsert1:; 10842d61bbb3SSatish Balay low = i; 10852d61bbb3SSatish Balay } 10862d61bbb3SSatish Balay ailen[brow] = nrow; 10872d61bbb3SSatish Balay } 10882d61bbb3SSatish Balay PetscFunctionReturn(0); 10892d61bbb3SSatish Balay } 10902d61bbb3SSatish Balay 10912d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10922d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10932d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 10947b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 10957b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 10962d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10972d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10982d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10992d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 11002d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 11012d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 11022d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 11032d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 11042d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 11052d61bbb3SSatish Balay 11062d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 11072d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 11082d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 11092d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 11102d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 11112d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 111215091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 11132d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 1114f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_7(Mat,Vec,Vec); 1115f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_6(Mat,Vec,Vec); 1116f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_5(Mat,Vec,Vec); 1117f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_4(Mat,Vec,Vec); 1118f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_3(Mat,Vec,Vec); 1119f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_2(Mat,Vec,Vec); 1120f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_1(Mat,Vec,Vec); 11212d61bbb3SSatish Balay 11222d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 11232d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 11242d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 11252d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 11262d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 11272d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 112815091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11292d61bbb3SSatish Balay 11302d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11312d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11322d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11332d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11342d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 113515091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11362d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11372d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11382d61bbb3SSatish Balay 11392d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11402d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11412d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11422d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11432d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 114415091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11452d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11462d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11472d61bbb3SSatish Balay 11482d61bbb3SSatish Balay #undef __FUNC__ 11492d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 11505ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11512d61bbb3SSatish Balay { 11522d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 11532d61bbb3SSatish Balay Mat outA; 11542d61bbb3SSatish Balay int ierr; 1155667159a5SBarry Smith PetscTruth row_identity, col_identity; 11562d61bbb3SSatish Balay 11572d61bbb3SSatish Balay PetscFunctionBegin; 1158b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1159667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1160667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1161667159a5SBarry Smith if (!row_identity || !col_identity) { 1162b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1163667159a5SBarry Smith } 11642d61bbb3SSatish Balay 11652d61bbb3SSatish Balay outA = inA; 11662d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11672d61bbb3SSatish Balay 11682d61bbb3SSatish Balay if (!a->diag) { 11692d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA);CHKERRQ(ierr); 11702d61bbb3SSatish Balay } 11712d61bbb3SSatish Balay /* 117215091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 117315091d37SBarry Smith for ILU(0) factorization with natural ordering 11742d61bbb3SSatish Balay */ 117515091d37SBarry Smith switch (a->bs) { 1176f1af5d2fSBarry Smith case 1: 1177f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_2_NaturalOrdering; 1178f1af5d2fSBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 117915091d37SBarry Smith case 2: 118015091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 118115091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 1182f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_2_NaturalOrdering; 118315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 118415091d37SBarry Smith break; 118515091d37SBarry Smith case 3: 118615091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 118715091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 1188f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_3_NaturalOrdering; 118915091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 119015091d37SBarry Smith break; 119115091d37SBarry Smith case 4: 1192667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1193f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 1194f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_4_NaturalOrdering; 119515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 119615091d37SBarry Smith break; 119715091d37SBarry Smith case 5: 1198667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1199667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 1200f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_5_NaturalOrdering; 120115091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 120215091d37SBarry Smith break; 120315091d37SBarry Smith case 6: 120415091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 120515091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 1206f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_6_NaturalOrdering; 120715091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 120815091d37SBarry Smith break; 120915091d37SBarry Smith case 7: 121015091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 1211f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_7_NaturalOrdering; 121215091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 121315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 121415091d37SBarry Smith break; 1215*c38d4ed2SBarry Smith default: 1216*c38d4ed2SBarry Smith a->row = row; 1217*c38d4ed2SBarry Smith a->col = col; 1218*c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1219*c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1220*c38d4ed2SBarry Smith 1221*c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1222*c38d4ed2SBarry Smith ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr); 1223*c38d4ed2SBarry Smith PLogObjectParent(inA,a->icol); 1224*c38d4ed2SBarry Smith 1225*c38d4ed2SBarry Smith if (!a->solve_work) { 1226*c38d4ed2SBarry Smith a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1227*c38d4ed2SBarry Smith PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 1228*c38d4ed2SBarry Smith } 12292d61bbb3SSatish Balay } 12302d61bbb3SSatish Balay 1231667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1232667159a5SBarry Smith 12332d61bbb3SSatish Balay PetscFunctionReturn(0); 12342d61bbb3SSatish Balay } 12352d61bbb3SSatish Balay #undef __FUNC__ 1236d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1237ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1238ba4ca20aSSatish Balay { 1239*c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1240ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1241d132466eSBarry Smith int ierr; 1242ba4ca20aSSatish Balay 12433a40ed3dSBarry Smith PetscFunctionBegin; 1244*c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1245d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1246d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12473a40ed3dSBarry Smith PetscFunctionReturn(0); 1248ba4ca20aSSatish Balay } 1249d9b7c43dSSatish Balay 1250fb2e594dSBarry Smith EXTERN_C_BEGIN 125127a8da17SBarry Smith #undef __FUNC__ 125227a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 125327a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 125427a8da17SBarry Smith { 125527a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 125627a8da17SBarry Smith int i,nz,n; 125727a8da17SBarry Smith 125827a8da17SBarry Smith PetscFunctionBegin; 125927a8da17SBarry Smith nz = baij->maxnz; 126027a8da17SBarry Smith n = baij->n; 126127a8da17SBarry Smith for (i=0; i<nz; i++) { 126227a8da17SBarry Smith baij->j[i] = indices[i]; 126327a8da17SBarry Smith } 126427a8da17SBarry Smith baij->nz = nz; 126527a8da17SBarry Smith for ( i=0; i<n; i++ ) { 126627a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 126727a8da17SBarry Smith } 126827a8da17SBarry Smith 126927a8da17SBarry Smith PetscFunctionReturn(0); 127027a8da17SBarry Smith } 1271fb2e594dSBarry Smith EXTERN_C_END 127227a8da17SBarry Smith 127327a8da17SBarry Smith #undef __FUNC__ 127427a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 127527a8da17SBarry Smith /*@ 127627a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 127727a8da17SBarry Smith in the matrix. 127827a8da17SBarry Smith 127927a8da17SBarry Smith Input Parameters: 128027a8da17SBarry Smith + mat - the SeqBAIJ matrix 128127a8da17SBarry Smith - indices - the column indices 128227a8da17SBarry Smith 128315091d37SBarry Smith Level: advanced 128415091d37SBarry Smith 128527a8da17SBarry Smith Notes: 128627a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 128727a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 128827a8da17SBarry Smith of the MatSetValues() operation. 128927a8da17SBarry Smith 129027a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 129127a8da17SBarry Smith MatCreateSeqBAIJ(). 129227a8da17SBarry Smith 129327a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 129427a8da17SBarry Smith 129527a8da17SBarry Smith @*/ 129627a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 129727a8da17SBarry Smith { 129827a8da17SBarry Smith int ierr,(*f)(Mat,int *); 129927a8da17SBarry Smith 130027a8da17SBarry Smith PetscFunctionBegin; 130127a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1302549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 130327a8da17SBarry Smith if (f) { 130427a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 130527a8da17SBarry Smith } else { 130627a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 130727a8da17SBarry Smith } 130827a8da17SBarry Smith PetscFunctionReturn(0); 130927a8da17SBarry Smith } 131027a8da17SBarry Smith 13112593348eSBarry Smith /* -------------------------------------------------------------------*/ 1312cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1313cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1314cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1315cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1316cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1317cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1318cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1319cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1320cc2dc46cSBarry Smith 0, 1321cc2dc46cSBarry Smith 0, 1322cc2dc46cSBarry Smith 0, 1323cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1324cc2dc46cSBarry Smith 0, 1325b6490206SBarry Smith 0, 1326f2501298SSatish Balay MatTranspose_SeqBAIJ, 1327cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1328cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1329cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1330cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1331cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1332b6490206SBarry Smith 0, 1333cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1334cc2dc46cSBarry Smith 0, 1335cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1336cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1337cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1338cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1339cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1340cc2dc46cSBarry Smith 0, 1341cc2dc46cSBarry Smith 0, 1342cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1343cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1344cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1345cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1346cc2dc46cSBarry Smith 0, 1347cc2dc46cSBarry Smith 0, 1348cc2dc46cSBarry Smith 0, 13492e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1350cc2dc46cSBarry Smith 0, 1351cc2dc46cSBarry Smith 0, 1352cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1353cc2dc46cSBarry Smith 0, 1354cc2dc46cSBarry Smith 0, 1355cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1356cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1357cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1358cc2dc46cSBarry Smith 0, 1359cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1360cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1361cc2dc46cSBarry Smith 0, 1362cc2dc46cSBarry Smith 0, 1363cc2dc46cSBarry Smith 0, 1364cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 13653b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 136692c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1367cc2dc46cSBarry Smith 0, 1368cc2dc46cSBarry Smith 0, 1369cc2dc46cSBarry Smith 0, 1370cc2dc46cSBarry Smith 0, 1371cc2dc46cSBarry Smith 0, 1372cc2dc46cSBarry Smith 0, 1373d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1374cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1375cc2dc46cSBarry Smith 0, 1376cc2dc46cSBarry Smith 0, 1377cc2dc46cSBarry Smith MatGetMaps_Petsc}; 13782593348eSBarry Smith 13793e90b805SBarry Smith EXTERN_C_BEGIN 13803e90b805SBarry Smith #undef __FUNC__ 13813e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ" 13823e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 13833e90b805SBarry Smith { 13843e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 13853e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1386d132466eSBarry Smith int ierr; 13873e90b805SBarry Smith 13883e90b805SBarry Smith PetscFunctionBegin; 13893e90b805SBarry Smith if (aij->nonew != 1) { 13903e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13913e90b805SBarry Smith } 13923e90b805SBarry Smith 13933e90b805SBarry Smith /* allocate space for values if not already there */ 13943e90b805SBarry Smith if (!aij->saved_values) { 13953e90b805SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 13963e90b805SBarry Smith } 13973e90b805SBarry Smith 13983e90b805SBarry Smith /* copy values over */ 1399d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 14003e90b805SBarry Smith PetscFunctionReturn(0); 14013e90b805SBarry Smith } 14023e90b805SBarry Smith EXTERN_C_END 14033e90b805SBarry Smith 14043e90b805SBarry Smith EXTERN_C_BEGIN 14053e90b805SBarry Smith #undef __FUNC__ 140632e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ" 140732e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 14083e90b805SBarry Smith { 14093e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1410549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 14113e90b805SBarry Smith 14123e90b805SBarry Smith PetscFunctionBegin; 14133e90b805SBarry Smith if (aij->nonew != 1) { 14143e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14153e90b805SBarry Smith } 14163e90b805SBarry Smith if (!aij->saved_values) { 14173e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 14183e90b805SBarry Smith } 14193e90b805SBarry Smith 14203e90b805SBarry Smith /* copy values over */ 1421549d3d68SSatish Balay ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 14223e90b805SBarry Smith PetscFunctionReturn(0); 14233e90b805SBarry Smith } 14243e90b805SBarry Smith EXTERN_C_END 14253e90b805SBarry Smith 14265615d1e5SSatish Balay #undef __FUNC__ 14275615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 14282593348eSBarry Smith /*@C 142944cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 143044cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 143144cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 14327fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 14332bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 14342593348eSBarry Smith 1435db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1436db81eaa0SLois Curfman McInnes 14372593348eSBarry Smith Input Parameters: 1438db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1439b6490206SBarry Smith . bs - size of block 14402593348eSBarry Smith . m - number of rows 14412593348eSBarry Smith . n - number of columns 1442b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14437fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14442bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14452593348eSBarry Smith 14462593348eSBarry Smith Output Parameter: 14472593348eSBarry Smith . A - the matrix 14482593348eSBarry Smith 14490513a670SBarry Smith Options Database Keys: 1450db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1451db81eaa0SLois Curfman McInnes block calculations (much slower) 1452db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14530513a670SBarry Smith 145415091d37SBarry Smith Level: intermediate 145515091d37SBarry Smith 14562593348eSBarry Smith Notes: 145744cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 14582593348eSBarry Smith storage. That is, the stored row and column indices can begin at 145944cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 14602593348eSBarry Smith 14612593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 14622593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14632593348eSBarry Smith allocation. For additional details, see the users manual chapter on 14646da5968aSLois Curfman McInnes matrices. 14652593348eSBarry Smith 1466db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 14672593348eSBarry Smith @*/ 1468b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 14692593348eSBarry Smith { 14702593348eSBarry Smith Mat B; 1471b6490206SBarry Smith Mat_SeqBAIJ *b; 1472f1af5d2fSBarry Smith int i,len,ierr,mbs,nbs,bs2,size; 1473f1af5d2fSBarry Smith PetscTruth flg; 14743b2fbd54SBarry Smith 14753a40ed3dSBarry Smith PetscFunctionBegin; 1476d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1477a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1478b6490206SBarry Smith 1479962fb4a0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1480302e33e3SBarry Smith mbs = m/bs; 1481302e33e3SBarry Smith nbs = n/bs; 1482302e33e3SBarry Smith bs2 = bs*bs; 14830513a670SBarry Smith 14843a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1485a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 14863a40ed3dSBarry Smith } 14872593348eSBarry Smith 1488b73539f3SBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz); 1489b73539f3SBarry Smith if (nnz) { 1490faf3f1fcSBarry Smith for (i=0; i<mbs; i++) { 1491b73539f3SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1492b73539f3SBarry Smith } 1493b73539f3SBarry Smith } 1494b73539f3SBarry Smith 14952593348eSBarry Smith *A = 0; 14963f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 14972593348eSBarry Smith PLogObjectCreate(B); 1498b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1499549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1500549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 15011a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 15021a6a6d98SBarry Smith if (!flg) { 15037fc0212eSBarry Smith switch (bs) { 15047fc0212eSBarry Smith case 1: 1505f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1506f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1507f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_1; 1508f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1509f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 15107fc0212eSBarry Smith break; 15114eeb42bcSBarry Smith case 2: 1512f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1513f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1514f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_2; 1515f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1516f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 15176d84be18SBarry Smith break; 1518f327dd97SBarry Smith case 3: 1519f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1520f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1521f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_3; 1522f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1523f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 15244eeb42bcSBarry Smith break; 15256d84be18SBarry Smith case 4: 1526f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1527f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1528f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_4; 1529f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1530f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 15316d84be18SBarry Smith break; 15326d84be18SBarry Smith case 5: 1533f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1534f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1535f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_5; 1536f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1537f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 15386d84be18SBarry Smith break; 153915091d37SBarry Smith case 6: 154015091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 154115091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 1542f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_6; 154315091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 154415091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 154515091d37SBarry Smith break; 154648e9ddb2SSatish Balay case 7: 1547f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1548f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1549f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_7; 1550f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 155148e9ddb2SSatish Balay break; 15527fc0212eSBarry Smith } 15531a6a6d98SBarry Smith } 1554e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1555e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15562593348eSBarry Smith B->factor = 0; 15572593348eSBarry Smith B->lupivotthreshold = 1.0; 155890f02eecSBarry Smith B->mapping = 0; 15592593348eSBarry Smith b->row = 0; 15602593348eSBarry Smith b->col = 0; 1561e51c0b9cSSatish Balay b->icol = 0; 15622593348eSBarry Smith b->reallocs = 0; 15633e90b805SBarry Smith b->saved_values = 0; 15642593348eSBarry Smith 156544cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 156644cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1567a5ae1ecdSBarry Smith 15687413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 15697413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1570a5ae1ecdSBarry Smith 1571b6490206SBarry Smith b->mbs = mbs; 1572f2501298SSatish Balay b->nbs = nbs; 1573b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax); 15742593348eSBarry Smith if (nnz == PETSC_NULL) { 1575119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 15762593348eSBarry Smith else if (nz <= 0) nz = 1; 1577b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1578b6490206SBarry Smith nz = nz*mbs; 15793a40ed3dSBarry Smith } else { 15802593348eSBarry Smith nz = 0; 1581b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 15822593348eSBarry Smith } 15832593348eSBarry Smith 15842593348eSBarry Smith /* allocate the matrix space */ 15853f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 15863f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a); 1587549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15887e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 1589549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15902593348eSBarry Smith b->i = b->j + nz; 15912593348eSBarry Smith b->singlemalloc = 1; 15922593348eSBarry Smith 1593de6a44a3SBarry Smith b->i[0] = 0; 1594b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 15952593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 15962593348eSBarry Smith } 15972593348eSBarry Smith 1598b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1599b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1600f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1601b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 16022593348eSBarry Smith 1603b6490206SBarry Smith b->bs = bs; 16049df24120SSatish Balay b->bs2 = bs2; 1605b6490206SBarry Smith b->mbs = mbs; 16062593348eSBarry Smith b->nz = 0; 160718e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 16082593348eSBarry Smith b->sorted = 0; 16092593348eSBarry Smith b->roworiented = 1; 16102593348eSBarry Smith b->nonew = 0; 16112593348eSBarry Smith b->diag = 0; 16122593348eSBarry Smith b->solve_work = 0; 1613de6a44a3SBarry Smith b->mult_work = 0; 16142593348eSBarry Smith b->spptr = 0; 16154e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 16164e220ebcSLois Curfman McInnes 16172593348eSBarry Smith *A = B; 16182593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 16192593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 162027a8da17SBarry Smith 1621f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 16223e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 16233e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1624f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 16253e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 16263e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1627f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 162827a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 162927a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 16303a40ed3dSBarry Smith PetscFunctionReturn(0); 16312593348eSBarry Smith } 16322593348eSBarry Smith 16335615d1e5SSatish Balay #undef __FUNC__ 16342e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 16352e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 16362593348eSBarry Smith { 16372593348eSBarry Smith Mat C; 1638b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 163927a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1640de6a44a3SBarry Smith 16413a40ed3dSBarry Smith PetscFunctionBegin; 1642a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 16432593348eSBarry Smith 16442593348eSBarry Smith *B = 0; 16453f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 16462593348eSBarry Smith PLogObjectCreate(C); 1647b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1648549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1649e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1650e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 16512593348eSBarry Smith C->factor = A->factor; 16522593348eSBarry Smith c->row = 0; 16532593348eSBarry Smith c->col = 0; 1654cac97260SSatish Balay c->icol = 0; 165532e87ba3SBarry Smith c->saved_values = 0; 16562593348eSBarry Smith C->assembled = PETSC_TRUE; 16572593348eSBarry Smith 165844cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 165944cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 166044cd7ae7SLois Curfman McInnes C->M = a->m; 166144cd7ae7SLois Curfman McInnes C->N = a->n; 166244cd7ae7SLois Curfman McInnes 1663b6490206SBarry Smith c->bs = a->bs; 16649df24120SSatish Balay c->bs2 = a->bs2; 1665b6490206SBarry Smith c->mbs = a->mbs; 1666b6490206SBarry Smith c->nbs = a->nbs; 16672593348eSBarry Smith 1668b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1669b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1670b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16712593348eSBarry Smith c->imax[i] = a->imax[i]; 16722593348eSBarry Smith c->ilen[i] = a->ilen[i]; 16732593348eSBarry Smith } 16742593348eSBarry Smith 16752593348eSBarry Smith /* allocate the matrix space */ 16762593348eSBarry Smith c->singlemalloc = 1; 16773f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 16783f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a); 16797e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1680de6a44a3SBarry Smith c->i = c->j + nz; 1681549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1682b6490206SBarry Smith if (mbs > 0) { 1683549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16842e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1685549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16862e8a6d31SBarry Smith } else { 1687549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16882593348eSBarry Smith } 16892593348eSBarry Smith } 16902593348eSBarry Smith 1691f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16922593348eSBarry Smith c->sorted = a->sorted; 16932593348eSBarry Smith c->roworiented = a->roworiented; 16942593348eSBarry Smith c->nonew = a->nonew; 16952593348eSBarry Smith 16962593348eSBarry Smith if (a->diag) { 1697b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag); 1698b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1699b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 17002593348eSBarry Smith c->diag[i] = a->diag[i]; 17012593348eSBarry Smith } 170298305bb5SBarry Smith } else c->diag = 0; 17032593348eSBarry Smith c->nz = a->nz; 17042593348eSBarry Smith c->maxnz = a->maxnz; 17052593348eSBarry Smith c->solve_work = 0; 17062593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 17077fc0212eSBarry Smith c->mult_work = 0; 17082593348eSBarry Smith *B = C; 17097b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 17103a40ed3dSBarry Smith PetscFunctionReturn(0); 17112593348eSBarry Smith } 17122593348eSBarry Smith 17135615d1e5SSatish Balay #undef __FUNC__ 17145615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 171519bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 17162593348eSBarry Smith { 1717b6490206SBarry Smith Mat_SeqBAIJ *a; 17182593348eSBarry Smith Mat B; 1719f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1720b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 172135aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1722a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1723b6490206SBarry Smith Scalar *aa; 172419bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 17252593348eSBarry Smith 17263a40ed3dSBarry Smith PetscFunctionBegin; 1727f1af5d2fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1728de6a44a3SBarry Smith bs2 = bs*bs; 1729b6490206SBarry Smith 1730d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1731cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 173290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 17330752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1734a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 17352593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 17362593348eSBarry Smith 1737d64ed03dSBarry Smith if (header[3] < 0) { 1738a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1739d64ed03dSBarry Smith } 1740d64ed03dSBarry Smith 1741a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 174235aab85fSBarry Smith 174335aab85fSBarry Smith /* 174435aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 174535aab85fSBarry Smith divisible by the blocksize 174635aab85fSBarry Smith */ 1747b6490206SBarry Smith mbs = M/bs; 174835aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 174935aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 175035aab85fSBarry Smith else mbs++; 175135aab85fSBarry Smith if (extra_rows) { 1752537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 175335aab85fSBarry Smith } 1754b6490206SBarry Smith 17552593348eSBarry Smith /* read in row lengths */ 175635aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 17570752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 175835aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 17592593348eSBarry Smith 1760b6490206SBarry Smith /* read in column indices */ 176135aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj); 17620752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 176335aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1764b6490206SBarry Smith 1765b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1766b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1767549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 176835aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask); 1769549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 177035aab85fSBarry Smith masked = mask + mbs; 1771b6490206SBarry Smith rowcount = 0; nzcount = 0; 1772b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 177335aab85fSBarry Smith nmask = 0; 1774b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1775b6490206SBarry Smith kmax = rowlengths[rowcount]; 1776b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 177735aab85fSBarry Smith tmp = jj[nzcount++]/bs; 177835aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1779b6490206SBarry Smith } 1780b6490206SBarry Smith rowcount++; 1781b6490206SBarry Smith } 178235aab85fSBarry Smith browlengths[i] += nmask; 178335aab85fSBarry Smith /* zero out the mask elements we set */ 178435aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1785b6490206SBarry Smith } 1786b6490206SBarry Smith 17872593348eSBarry Smith /* create our matrix */ 17883eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17892593348eSBarry Smith B = *A; 1790b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 17912593348eSBarry Smith 17922593348eSBarry Smith /* set matrix "i" values */ 1793de6a44a3SBarry Smith a->i[0] = 0; 1794b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1795b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1796b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17972593348eSBarry Smith } 1798b6490206SBarry Smith a->nz = 0; 1799b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 18002593348eSBarry Smith 1801b6490206SBarry Smith /* read in nonzero values */ 180235aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 18030752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 180435aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1805b6490206SBarry Smith 1806b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1807b6490206SBarry Smith nzcount = 0; jcount = 0; 1808b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1809b6490206SBarry Smith nzcountb = nzcount; 181035aab85fSBarry Smith nmask = 0; 1811b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1812b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1813b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 181435aab85fSBarry Smith tmp = jj[nzcount++]/bs; 181535aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1816b6490206SBarry Smith } 1817b6490206SBarry Smith rowcount++; 1818b6490206SBarry Smith } 1819de6a44a3SBarry Smith /* sort the masked values */ 1820433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1821de6a44a3SBarry Smith 1822b6490206SBarry Smith /* set "j" values into matrix */ 1823b6490206SBarry Smith maskcount = 1; 182435aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 182535aab85fSBarry Smith a->j[jcount++] = masked[j]; 1826de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1827b6490206SBarry Smith } 1828b6490206SBarry Smith /* set "a" values into matrix */ 1829de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1830b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1831b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1832b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1833de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1834de6a44a3SBarry Smith block = mask[tmp] - 1; 1835de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1836de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1837b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1838b6490206SBarry Smith } 1839b6490206SBarry Smith } 184035aab85fSBarry Smith /* zero out the mask elements we set */ 184135aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1842b6490206SBarry Smith } 1843a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1844b6490206SBarry Smith 1845606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1846606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1847606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1848606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1849606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1850b6490206SBarry Smith 1851b6490206SBarry Smith B->assembled = PETSC_TRUE; 1852de6a44a3SBarry Smith 18539c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 18543a40ed3dSBarry Smith PetscFunctionReturn(0); 18552593348eSBarry Smith } 18562593348eSBarry Smith 18572593348eSBarry Smith 18582593348eSBarry Smith 1859