1bba1ac68SSatish Balay /*$Id: baij.c,v 1.245 2001/08/07 03:02:55 balay 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 */ 770f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 81a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 91a6a6d98SBarry Smith #include "src/inline/spops.h" 10325e03aeSBarry Smith #include "petscsys.h" /*I "petscmat.h" I*/ 113b547af2SSatish Balay 12af674e45SBarry Smith /* 13af674e45SBarry Smith Special version for Fun3d sequential benchmark 14af674e45SBarry Smith */ 15af674e45SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS) 16af674e45SBarry Smith #define matsetvaluesblocked4_ MATSETVALUESBLOCKED4 17af674e45SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 18af674e45SBarry Smith #define matsetvaluesblocked4_ matsetvaluesblocked4 19af674e45SBarry Smith #endif 20af674e45SBarry Smith 219c8c1041SBarry Smith EXTERN_C_BEGIN 22af674e45SBarry Smith #undef __FUNCT__ 23af674e45SBarry Smith #define __FUNCT__ "matsetvaluesblocked4_" 24f15d580aSBarry Smith void matsetvaluesblocked4_(Mat *AA,int *mm,const int im[],int *nn,const int in[],const PetscScalar v[]) 25af674e45SBarry Smith { 26af674e45SBarry Smith Mat A = *AA; 27af674e45SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 28a037b02bSBarry Smith int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,N,m = *mm,n = *nn; 2938fde467SHong Zhang int *ai=a->i,*ailen=a->ilen; 30a037b02bSBarry Smith int *aj=a->j,stepval; 31f15d580aSBarry Smith const PetscScalar *value = v; 324bb09213Spetsc MatScalar *ap,*aa = a->a,*bap; 33af674e45SBarry Smith 34af674e45SBarry Smith PetscFunctionBegin; 35af674e45SBarry Smith stepval = (n-1)*4; 36af674e45SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 37af674e45SBarry Smith row = im[k]; 38af674e45SBarry Smith rp = aj + ai[row]; 39af674e45SBarry Smith ap = aa + 16*ai[row]; 40af674e45SBarry Smith nrow = ailen[row]; 41af674e45SBarry Smith low = 0; 42af674e45SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 43af674e45SBarry Smith col = in[l]; 44af674e45SBarry Smith value = v + k*(stepval+4)*4 + l*4; 45af674e45SBarry Smith low = 0; high = nrow; 46af674e45SBarry Smith while (high-low > 7) { 47af674e45SBarry Smith t = (low+high)/2; 48af674e45SBarry Smith if (rp[t] > col) high = t; 49af674e45SBarry Smith else low = t; 50af674e45SBarry Smith } 51af674e45SBarry Smith for (i=low; i<high; i++) { 52af674e45SBarry Smith if (rp[i] > col) break; 53af674e45SBarry Smith if (rp[i] == col) { 54af674e45SBarry Smith bap = ap + 16*i; 55af674e45SBarry Smith for (ii=0; ii<4; ii++,value+=stepval) { 56af674e45SBarry Smith for (jj=ii; jj<16; jj+=4) { 57af674e45SBarry Smith bap[jj] += *value++; 58af674e45SBarry Smith } 59af674e45SBarry Smith } 60af674e45SBarry Smith goto noinsert2; 61af674e45SBarry Smith } 62af674e45SBarry Smith } 63af674e45SBarry Smith N = nrow++ - 1; 64af674e45SBarry Smith /* shift up all the later entries in this row */ 65af674e45SBarry Smith for (ii=N; ii>=i; ii--) { 66af674e45SBarry Smith rp[ii+1] = rp[ii]; 67a037b02bSBarry Smith PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar)); 68af674e45SBarry Smith } 69af674e45SBarry Smith if (N >= i) { 70a037b02bSBarry Smith PetscMemzero(ap+16*i,16*sizeof(MatScalar)); 71af674e45SBarry Smith } 72af674e45SBarry Smith rp[i] = col; 73af674e45SBarry Smith bap = ap + 16*i; 74af674e45SBarry Smith for (ii=0; ii<4; ii++,value+=stepval) { 75af674e45SBarry Smith for (jj=ii; jj<16; jj+=4) { 76af674e45SBarry Smith bap[jj] = *value++; 77af674e45SBarry Smith } 78af674e45SBarry Smith } 79af674e45SBarry Smith noinsert2:; 80af674e45SBarry Smith low = i; 81af674e45SBarry Smith } 82af674e45SBarry Smith ailen[row] = nrow; 83af674e45SBarry Smith } 84af674e45SBarry Smith } 859c8c1041SBarry Smith EXTERN_C_END 86af674e45SBarry Smith 87af674e45SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS) 88af674e45SBarry Smith #define matsetvalues4_ MATSETVALUES4 89af674e45SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 90af674e45SBarry Smith #define matsetvalues4_ matsetvalues4 91af674e45SBarry Smith #endif 92af674e45SBarry Smith 939c8c1041SBarry Smith EXTERN_C_BEGIN 94af674e45SBarry Smith #undef __FUNCT__ 95af674e45SBarry Smith #define __FUNCT__ "MatSetValues4_" 964bb09213Spetsc void matsetvalues4_(Mat *AA,int *mm,int *im,int *nn,int *in,PetscScalar *v) 97af674e45SBarry Smith { 98af674e45SBarry Smith Mat A = *AA; 99af674e45SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 100a037b02bSBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,N,n = *nn,m = *mm; 10138fde467SHong Zhang int *ai=a->i,*ailen=a->ilen; 102af674e45SBarry Smith int *aj=a->j,brow,bcol; 10338fde467SHong Zhang int ridx,cidx; 104af674e45SBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 105af674e45SBarry Smith 106af674e45SBarry Smith PetscFunctionBegin; 107af674e45SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 108af674e45SBarry Smith row = im[k]; brow = row/4; 109af674e45SBarry Smith rp = aj + ai[brow]; 110af674e45SBarry Smith ap = aa + 16*ai[brow]; 111af674e45SBarry Smith nrow = ailen[brow]; 112af674e45SBarry Smith low = 0; 113af674e45SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 114af674e45SBarry Smith col = in[l]; bcol = col/4; 115af674e45SBarry Smith ridx = row % 4; cidx = col % 4; 116af674e45SBarry Smith value = v[l + k*n]; 117af674e45SBarry Smith low = 0; high = nrow; 118af674e45SBarry Smith while (high-low > 7) { 119af674e45SBarry Smith t = (low+high)/2; 120af674e45SBarry Smith if (rp[t] > bcol) high = t; 121af674e45SBarry Smith else low = t; 122af674e45SBarry Smith } 123af674e45SBarry Smith for (i=low; i<high; i++) { 124af674e45SBarry Smith if (rp[i] > bcol) break; 125af674e45SBarry Smith if (rp[i] == bcol) { 126af674e45SBarry Smith bap = ap + 16*i + 4*cidx + ridx; 127af674e45SBarry Smith *bap += value; 128af674e45SBarry Smith goto noinsert1; 129af674e45SBarry Smith } 130af674e45SBarry Smith } 131af674e45SBarry Smith N = nrow++ - 1; 132af674e45SBarry Smith /* shift up all the later entries in this row */ 133af674e45SBarry Smith for (ii=N; ii>=i; ii--) { 134af674e45SBarry Smith rp[ii+1] = rp[ii]; 135a037b02bSBarry Smith PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar)); 136af674e45SBarry Smith } 137af674e45SBarry Smith if (N>=i) { 138a037b02bSBarry Smith PetscMemzero(ap+16*i,16*sizeof(MatScalar)); 139af674e45SBarry Smith } 140af674e45SBarry Smith rp[i] = bcol; 141af674e45SBarry Smith ap[16*i + 4*cidx + ridx] = value; 142af674e45SBarry Smith noinsert1:; 143af674e45SBarry Smith low = i; 144af674e45SBarry Smith } 145af674e45SBarry Smith ailen[brow] = nrow; 146af674e45SBarry Smith } 147af674e45SBarry Smith } 1489c8c1041SBarry Smith EXTERN_C_END 149af674e45SBarry Smith 15095d5f7c2SBarry Smith /* UGLY, ugly, ugly 15187828ca2SBarry Smith When MatScalar == PetscScalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does 1523477af42SKris Buschelman not exist. Otherwise ..._MatScalar() takes matrix dlements in single precision and 15395d5f7c2SBarry Smith inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ() 15495d5f7c2SBarry Smith converts the entries into single precision and then calls ..._MatScalar() to put them 15595d5f7c2SBarry Smith into the single precision data structures. 15695d5f7c2SBarry Smith */ 15795d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 158f15d580aSBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,const int[],int,const int[],const MatScalar[],InsertMode); 15995d5f7c2SBarry Smith #else 16095d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ 16195d5f7c2SBarry Smith #endif 16295d5f7c2SBarry Smith 1632d61bbb3SSatish Balay #define CHUNKSIZE 10 164de6a44a3SBarry Smith 165be5855fcSBarry Smith /* 166be5855fcSBarry Smith Checks for missing diagonals 167be5855fcSBarry Smith */ 1684a2ae208SSatish Balay #undef __FUNCT__ 1694a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ" 170c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A) 171be5855fcSBarry Smith { 172be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 173883fce79SBarry Smith int *diag,*jj = a->j,i,ierr; 174be5855fcSBarry Smith 175be5855fcSBarry Smith PetscFunctionBegin; 176c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 177883fce79SBarry Smith diag = a->diag; 1780e8e8aceSBarry Smith for (i=0; i<a->mbs; i++) { 179be5855fcSBarry Smith if (jj[diag[i]] != i) { 18029bbc08cSBarry Smith SETERRQ1(1,"Matrix is missing diagonal number %d",i); 181be5855fcSBarry Smith } 182be5855fcSBarry Smith } 183be5855fcSBarry Smith PetscFunctionReturn(0); 184be5855fcSBarry Smith } 185be5855fcSBarry Smith 1864a2ae208SSatish Balay #undef __FUNCT__ 1874a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ" 188c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A) 189de6a44a3SBarry Smith { 190de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 19182502324SSatish Balay int i,j,*diag,m = a->mbs,ierr; 192de6a44a3SBarry Smith 1933a40ed3dSBarry Smith PetscFunctionBegin; 194883fce79SBarry Smith if (a->diag) PetscFunctionReturn(0); 195883fce79SBarry Smith 196b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 197b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 1987fc0212eSBarry Smith for (i=0; i<m; i++) { 199ecc615b2SBarry Smith diag[i] = a->i[i+1]; 200de6a44a3SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 201de6a44a3SBarry Smith if (a->j[j] == i) { 202de6a44a3SBarry Smith diag[i] = j; 203de6a44a3SBarry Smith break; 204de6a44a3SBarry Smith } 205de6a44a3SBarry Smith } 206de6a44a3SBarry Smith } 207de6a44a3SBarry Smith a->diag = diag; 2083a40ed3dSBarry Smith PetscFunctionReturn(0); 209de6a44a3SBarry Smith } 2102593348eSBarry Smith 2112593348eSBarry Smith 212ca44d042SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 2133b2fbd54SBarry Smith 2144a2ae208SSatish Balay #undef __FUNCT__ 2154a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqBAIJ" 2164e7234bfSBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done) 2173b2fbd54SBarry Smith { 2183b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2193b2fbd54SBarry Smith int ierr,n = a->mbs,i; 2203b2fbd54SBarry Smith 2213a40ed3dSBarry Smith PetscFunctionBegin; 2223b2fbd54SBarry Smith *nn = n; 2233a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 2243b2fbd54SBarry Smith if (symmetric) { 2253b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 2263b2fbd54SBarry Smith } else if (oshift == 1) { 2273b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 228435da068SBarry Smith int nz = a->i[n]; 2293b2fbd54SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 2303b2fbd54SBarry Smith for (i=0; i<n+1; i++) a->i[i]++; 2313b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 2323b2fbd54SBarry Smith } else { 2333b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 2343b2fbd54SBarry Smith } 2353b2fbd54SBarry Smith 2363a40ed3dSBarry Smith PetscFunctionReturn(0); 2373b2fbd54SBarry Smith } 2383b2fbd54SBarry Smith 2394a2ae208SSatish Balay #undef __FUNCT__ 2404a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ" 2414e7234bfSBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done) 2423b2fbd54SBarry Smith { 2433b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 244606d414cSSatish Balay int i,n = a->mbs,ierr; 2453b2fbd54SBarry Smith 2463a40ed3dSBarry Smith PetscFunctionBegin; 2473a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 2483b2fbd54SBarry Smith if (symmetric) { 249606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 250606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 251af5da2bfSBarry Smith } else if (oshift == 1) { 252435da068SBarry Smith int nz = a->i[n]-1; 2533b2fbd54SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 2543b2fbd54SBarry Smith for (i=0; i<n+1; i++) a->i[i]--; 2553b2fbd54SBarry Smith } 2563a40ed3dSBarry Smith PetscFunctionReturn(0); 2573b2fbd54SBarry Smith } 2583b2fbd54SBarry Smith 2594a2ae208SSatish Balay #undef __FUNCT__ 2604a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqBAIJ" 2612d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs) 2622d61bbb3SSatish Balay { 2632d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data; 2642d61bbb3SSatish Balay 2652d61bbb3SSatish Balay PetscFunctionBegin; 2662d61bbb3SSatish Balay *bs = baij->bs; 2672d61bbb3SSatish Balay PetscFunctionReturn(0); 2682d61bbb3SSatish Balay } 2692d61bbb3SSatish Balay 2702d61bbb3SSatish Balay 2714a2ae208SSatish Balay #undef __FUNCT__ 2724a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqBAIJ" 273e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 2742d61bbb3SSatish Balay { 2752d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 276e51c0b9cSSatish Balay int ierr; 2772d61bbb3SSatish Balay 278433994e6SBarry Smith PetscFunctionBegin; 279aa482453SBarry Smith #if defined(PETSC_USE_LOG) 280b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz); 2812d61bbb3SSatish Balay #endif 282606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 283606d414cSSatish Balay if (!a->singlemalloc) { 284606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 285606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 286606d414cSSatish Balay } 287c38d4ed2SBarry Smith if (a->row) { 288c38d4ed2SBarry Smith ierr = ISDestroy(a->row);CHKERRQ(ierr); 289c38d4ed2SBarry Smith } 290c38d4ed2SBarry Smith if (a->col) { 291c38d4ed2SBarry Smith ierr = ISDestroy(a->col);CHKERRQ(ierr); 292c38d4ed2SBarry Smith } 293606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 294606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 295606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 296606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 297606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 298e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 299606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 300563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 301563b5814SBarry Smith if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);} 302563b5814SBarry Smith #endif 303c4319e64SHong Zhang if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);} 304c4319e64SHong Zhang 305606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 3062d61bbb3SSatish Balay PetscFunctionReturn(0); 3072d61bbb3SSatish Balay } 3082d61bbb3SSatish Balay 3094a2ae208SSatish Balay #undef __FUNCT__ 3104a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqBAIJ" 3112d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 3122d61bbb3SSatish Balay { 3132d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 3142d61bbb3SSatish Balay 3152d61bbb3SSatish Balay PetscFunctionBegin; 316aa275fccSKris Buschelman switch (op) { 317aa275fccSKris Buschelman case MAT_ROW_ORIENTED: 318aa275fccSKris Buschelman a->roworiented = PETSC_TRUE; 319aa275fccSKris Buschelman break; 320aa275fccSKris Buschelman case MAT_COLUMN_ORIENTED: 321aa275fccSKris Buschelman a->roworiented = PETSC_FALSE; 322aa275fccSKris Buschelman break; 323aa275fccSKris Buschelman case MAT_COLUMNS_SORTED: 324aa275fccSKris Buschelman a->sorted = PETSC_TRUE; 325aa275fccSKris Buschelman break; 326aa275fccSKris Buschelman case MAT_COLUMNS_UNSORTED: 327aa275fccSKris Buschelman a->sorted = PETSC_FALSE; 328aa275fccSKris Buschelman break; 329aa275fccSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 330aa275fccSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 331aa275fccSKris Buschelman break; 332aa275fccSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 333aa275fccSKris Buschelman a->nonew = 1; 334aa275fccSKris Buschelman break; 335aa275fccSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 336aa275fccSKris Buschelman a->nonew = -1; 337aa275fccSKris Buschelman break; 338aa275fccSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 339aa275fccSKris Buschelman a->nonew = -2; 340aa275fccSKris Buschelman break; 341aa275fccSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 342aa275fccSKris Buschelman a->nonew = 0; 343aa275fccSKris Buschelman break; 344aa275fccSKris Buschelman case MAT_ROWS_SORTED: 345aa275fccSKris Buschelman case MAT_ROWS_UNSORTED: 346aa275fccSKris Buschelman case MAT_YES_NEW_DIAGONALS: 347aa275fccSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 348aa275fccSKris Buschelman case MAT_USE_HASH_TABLE: 349b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 350aa275fccSKris Buschelman break; 351aa275fccSKris Buschelman case MAT_NO_NEW_DIAGONALS: 35229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 353*77e54ba9SKris Buschelman case MAT_SYMMETRIC: 354*77e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 355*77e54ba9SKris Buschelman break; 356aa275fccSKris Buschelman default: 35729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 3582d61bbb3SSatish Balay } 3592d61bbb3SSatish Balay PetscFunctionReturn(0); 3602d61bbb3SSatish Balay } 3612d61bbb3SSatish Balay 3624a2ae208SSatish Balay #undef __FUNCT__ 3634a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ" 36487828ca2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 3652d61bbb3SSatish Balay { 3662d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 36782502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr; 3683f1db9ecSBarry Smith MatScalar *aa,*aa_i; 36987828ca2SBarry Smith PetscScalar *v_i; 3702d61bbb3SSatish Balay 3712d61bbb3SSatish Balay PetscFunctionBegin; 3722d61bbb3SSatish Balay bs = a->bs; 3732d61bbb3SSatish Balay ai = a->i; 3742d61bbb3SSatish Balay aj = a->j; 3752d61bbb3SSatish Balay aa = a->a; 3762d61bbb3SSatish Balay bs2 = a->bs2; 3772d61bbb3SSatish Balay 378a45adfd6SMatthew Knepley if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d out of range", row); 3792d61bbb3SSatish Balay 3802d61bbb3SSatish Balay bn = row/bs; /* Block number */ 3812d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 3822d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 3832d61bbb3SSatish Balay *nz = bs*M; 3842d61bbb3SSatish Balay 3852d61bbb3SSatish Balay if (v) { 3862d61bbb3SSatish Balay *v = 0; 3872d61bbb3SSatish Balay if (*nz) { 38887828ca2SBarry Smith ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr); 3892d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 3902d61bbb3SSatish Balay v_i = *v + i*bs; 3912d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 3922d61bbb3SSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 3932d61bbb3SSatish Balay } 3942d61bbb3SSatish Balay } 3952d61bbb3SSatish Balay } 3962d61bbb3SSatish Balay 3972d61bbb3SSatish Balay if (idx) { 3982d61bbb3SSatish Balay *idx = 0; 3992d61bbb3SSatish Balay if (*nz) { 400b0a32e0cSBarry Smith ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr); 4012d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 4022d61bbb3SSatish Balay idx_i = *idx + i*bs; 4032d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 4042d61bbb3SSatish Balay for (j=0; j<bs; j++) {idx_i[j] = itmp++;} 4052d61bbb3SSatish Balay } 4062d61bbb3SSatish Balay } 4072d61bbb3SSatish Balay } 4082d61bbb3SSatish Balay PetscFunctionReturn(0); 4092d61bbb3SSatish Balay } 4102d61bbb3SSatish Balay 4114a2ae208SSatish Balay #undef __FUNCT__ 4124a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ" 41387828ca2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 4142d61bbb3SSatish Balay { 415606d414cSSatish Balay int ierr; 416606d414cSSatish Balay 4172d61bbb3SSatish Balay PetscFunctionBegin; 418606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 419606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 4202d61bbb3SSatish Balay PetscFunctionReturn(0); 4212d61bbb3SSatish Balay } 4222d61bbb3SSatish Balay 4234a2ae208SSatish Balay #undef __FUNCT__ 4244a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ" 4252d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 4262d61bbb3SSatish Balay { 4272d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 4282d61bbb3SSatish Balay Mat C; 4292d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 430273d9f13SBarry Smith int *rows,*cols,bs2=a->bs2; 43187828ca2SBarry Smith PetscScalar *array; 4322d61bbb3SSatish Balay 4332d61bbb3SSatish Balay PetscFunctionBegin; 43429bbc08cSBarry Smith if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place"); 435b0a32e0cSBarry Smith ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr); 436549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 4372d61bbb3SSatish Balay 438375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 43987828ca2SBarry Smith ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr); 44087828ca2SBarry Smith for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i]; 441375fe846SBarry Smith #else 4423eda8832SBarry Smith array = a->a; 443375fe846SBarry Smith #endif 444375fe846SBarry Smith 4452d61bbb3SSatish Balay for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 446273d9f13SBarry Smith ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 447606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 448b0a32e0cSBarry Smith ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr); 4492d61bbb3SSatish Balay cols = rows + bs; 4502d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4512d61bbb3SSatish Balay cols[0] = i*bs; 4522d61bbb3SSatish Balay for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 4532d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 4542d61bbb3SSatish Balay for (j=0; j<len; j++) { 4552d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 4562d61bbb3SSatish Balay for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 4572d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 4582d61bbb3SSatish Balay array += bs2; 4592d61bbb3SSatish Balay } 4602d61bbb3SSatish Balay } 461606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 462375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 463375fe846SBarry Smith ierr = PetscFree(array); 464375fe846SBarry Smith #endif 4652d61bbb3SSatish Balay 4662d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4672d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4682d61bbb3SSatish Balay 469c4992f7dSBarry Smith if (B) { 4702d61bbb3SSatish Balay *B = C; 4712d61bbb3SSatish Balay } else { 472273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 4732d61bbb3SSatish Balay } 4742d61bbb3SSatish Balay PetscFunctionReturn(0); 4752d61bbb3SSatish Balay } 4762d61bbb3SSatish Balay 4774a2ae208SSatish Balay #undef __FUNCT__ 4784a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary" 479b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer) 4802593348eSBarry Smith { 481b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 4829df24120SSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 48387828ca2SBarry Smith PetscScalar *aa; 484ce6f0cecSBarry Smith FILE *file; 4852593348eSBarry Smith 4863a40ed3dSBarry Smith PetscFunctionBegin; 487b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 488b0a32e0cSBarry Smith ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 489552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 4903b2fbd54SBarry Smith 491273d9f13SBarry Smith col_lens[1] = A->m; 492273d9f13SBarry Smith col_lens[2] = A->n; 4937e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 4942593348eSBarry Smith 4952593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 496b6490206SBarry Smith count = 0; 497b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 498b6490206SBarry Smith for (j=0; j<bs; j++) { 499b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 500b6490206SBarry Smith } 5012593348eSBarry Smith } 502273d9f13SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 503606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 5042593348eSBarry Smith 5052593348eSBarry Smith /* store column indices (zero start index) */ 506b0a32e0cSBarry Smith ierr = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 507b6490206SBarry Smith count = 0; 508b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 509b6490206SBarry Smith for (j=0; j<bs; j++) { 510b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 511b6490206SBarry Smith for (l=0; l<bs; l++) { 512b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 5132593348eSBarry Smith } 5142593348eSBarry Smith } 515b6490206SBarry Smith } 516b6490206SBarry Smith } 5170752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 518606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 5192593348eSBarry Smith 5202593348eSBarry Smith /* store nonzero values */ 52187828ca2SBarry Smith ierr = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 522b6490206SBarry Smith count = 0; 523b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 524b6490206SBarry Smith for (j=0; j<bs; j++) { 525b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 526b6490206SBarry Smith for (l=0; l<bs; l++) { 5277e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 528b6490206SBarry Smith } 529b6490206SBarry Smith } 530b6490206SBarry Smith } 531b6490206SBarry Smith } 5320752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 533606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 534ce6f0cecSBarry Smith 535b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 536ce6f0cecSBarry Smith if (file) { 537ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 538ce6f0cecSBarry Smith } 5393a40ed3dSBarry Smith PetscFunctionReturn(0); 5402593348eSBarry Smith } 5412593348eSBarry Smith 5424a2ae208SSatish Balay #undef __FUNCT__ 5434a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII" 544b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer) 5452593348eSBarry Smith { 546b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 547fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 548f3ef73ceSBarry Smith PetscViewerFormat format; 5492593348eSBarry Smith 5503a40ed3dSBarry Smith PetscFunctionBegin; 551b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 552456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 553b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 554fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 555bcd9e38bSBarry Smith Mat aij; 556bcd9e38bSBarry Smith ierr = MatConvert(A,MATSEQAIJ,&aij);CHKERRQ(ierr); 557bcd9e38bSBarry Smith ierr = MatView(aij,viewer);CHKERRQ(ierr); 558bcd9e38bSBarry Smith ierr = MatDestroy(aij);CHKERRQ(ierr); 55904929863SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 56004929863SHong Zhang PetscFunctionReturn(0); 561fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 562b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 56344cd7ae7SLois Curfman McInnes for (i=0; i<a->mbs; i++) { 56444cd7ae7SLois Curfman McInnes for (j=0; j<bs; j++) { 565b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 56644cd7ae7SLois Curfman McInnes for (k=a->i[i]; k<a->i[i+1]; k++) { 56744cd7ae7SLois Curfman McInnes for (l=0; l<bs; l++) { 568aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5690e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 57062b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %gi) ",bs*a->j[k]+l, 5710e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5720e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 57362b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %gi) ",bs*a->j[k]+l, 5740e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5750e6d2581SBarry Smith } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 57662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5770ef38995SBarry Smith } 57844cd7ae7SLois Curfman McInnes #else 5790ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 58062b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 5810ef38995SBarry Smith } 58244cd7ae7SLois Curfman McInnes #endif 58344cd7ae7SLois Curfman McInnes } 58444cd7ae7SLois Curfman McInnes } 585b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 58644cd7ae7SLois Curfman McInnes } 58744cd7ae7SLois Curfman McInnes } 588b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 5890ef38995SBarry Smith } else { 590b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 591b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 592b6490206SBarry Smith for (j=0; j<bs; j++) { 593b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 594b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 595b6490206SBarry Smith for (l=0; l<bs; l++) { 596aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5970e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 59862b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 5990e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 6000e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 60162b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 6020e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 6030ef38995SBarry Smith } else { 60462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 60588685aaeSLois Curfman McInnes } 60688685aaeSLois Curfman McInnes #else 60762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 60888685aaeSLois Curfman McInnes #endif 6092593348eSBarry Smith } 6102593348eSBarry Smith } 611b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6122593348eSBarry Smith } 6132593348eSBarry Smith } 614b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 615b6490206SBarry Smith } 616b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 6182593348eSBarry Smith } 6192593348eSBarry Smith 6204a2ae208SSatish Balay #undef __FUNCT__ 6214a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom" 622b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 6233270192aSSatish Balay { 62477ed5343SBarry Smith Mat A = (Mat) Aa; 6253270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 626b65db4caSBarry Smith int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2; 6270e6d2581SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 6283f1db9ecSBarry Smith MatScalar *aa; 629b0a32e0cSBarry Smith PetscViewer viewer; 6303270192aSSatish Balay 6313a40ed3dSBarry Smith PetscFunctionBegin; 6323270192aSSatish Balay 633b65db4caSBarry Smith /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/ 63477ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 63577ed5343SBarry Smith 636b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 63777ed5343SBarry Smith 6383270192aSSatish Balay /* loop over matrix elements drawing boxes */ 639b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 6403270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 6413270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 642273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 6433270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 6443270192aSSatish Balay aa = a->a + j*bs2; 6453270192aSSatish Balay for (k=0; k<bs; k++) { 6463270192aSSatish Balay for (l=0; l<bs; l++) { 6470e6d2581SBarry Smith if (PetscRealPart(*aa++) >= 0.) continue; 648b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 6493270192aSSatish Balay } 6503270192aSSatish Balay } 6513270192aSSatish Balay } 6523270192aSSatish Balay } 653b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 6543270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 6553270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 656273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 6573270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 6583270192aSSatish Balay aa = a->a + j*bs2; 6593270192aSSatish Balay for (k=0; k<bs; k++) { 6603270192aSSatish Balay for (l=0; l<bs; l++) { 6610e6d2581SBarry Smith if (PetscRealPart(*aa++) != 0.) continue; 662b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 6633270192aSSatish Balay } 6643270192aSSatish Balay } 6653270192aSSatish Balay } 6663270192aSSatish Balay } 6673270192aSSatish Balay 668b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 6693270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 6703270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 671273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 6723270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 6733270192aSSatish Balay aa = a->a + j*bs2; 6743270192aSSatish Balay for (k=0; k<bs; k++) { 6753270192aSSatish Balay for (l=0; l<bs; l++) { 6760e6d2581SBarry Smith if (PetscRealPart(*aa++) <= 0.) continue; 677b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 6783270192aSSatish Balay } 6793270192aSSatish Balay } 6803270192aSSatish Balay } 6813270192aSSatish Balay } 68277ed5343SBarry Smith PetscFunctionReturn(0); 68377ed5343SBarry Smith } 6843270192aSSatish Balay 6854a2ae208SSatish Balay #undef __FUNCT__ 6864a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw" 687b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer) 68877ed5343SBarry Smith { 6897dce120fSSatish Balay int ierr; 6900e6d2581SBarry Smith PetscReal xl,yl,xr,yr,w,h; 691b0a32e0cSBarry Smith PetscDraw draw; 69277ed5343SBarry Smith PetscTruth isnull; 6933270192aSSatish Balay 69477ed5343SBarry Smith PetscFunctionBegin; 69577ed5343SBarry Smith 696b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 697b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 69877ed5343SBarry Smith 69977ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 700273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 70177ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 702b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 703b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 70477ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 7053a40ed3dSBarry Smith PetscFunctionReturn(0); 7063270192aSSatish Balay } 7073270192aSSatish Balay 7084a2ae208SSatish Balay #undef __FUNCT__ 7094a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ" 710b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer) 7112593348eSBarry Smith { 71219bcc07fSBarry Smith int ierr; 7136831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 7142593348eSBarry Smith 7153a40ed3dSBarry Smith PetscFunctionBegin; 716b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 717b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 718fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 719fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 7200f5bd95cSBarry Smith if (issocket) { 72129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 7220f5bd95cSBarry Smith } else if (isascii){ 7233a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 7240f5bd95cSBarry Smith } else if (isbinary) { 7253a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 7260f5bd95cSBarry Smith } else if (isdraw) { 7273a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 7285cd90555SBarry Smith } else { 72929bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 7302593348eSBarry Smith } 7313a40ed3dSBarry Smith PetscFunctionReturn(0); 7322593348eSBarry Smith } 733b6490206SBarry Smith 734cd0e1443SSatish Balay 7354a2ae208SSatish Balay #undef __FUNCT__ 7364a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ" 737f15d580aSBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[]) 738cd0e1443SSatish Balay { 739cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 7402d61bbb3SSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 7412d61bbb3SSatish Balay int *ai = a->i,*ailen = a->ilen; 7422d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 7433f1db9ecSBarry Smith MatScalar *ap,*aa = a->a,zero = 0.0; 744cd0e1443SSatish Balay 7453a40ed3dSBarry Smith PetscFunctionBegin; 7462d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 747cd0e1443SSatish Balay row = im[k]; brow = row/bs; 74829bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 749a45adfd6SMatthew Knepley if (row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d too large", row); 7502d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 7512c3acbe9SBarry Smith nrow = ailen[brow]; 7522d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 75329bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 754a45adfd6SMatthew Knepley if (in[l] >= A->n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column %d too large", in[l]); 7552d61bbb3SSatish Balay col = in[l] ; 7562d61bbb3SSatish Balay bcol = col/bs; 7572d61bbb3SSatish Balay cidx = col%bs; 7582d61bbb3SSatish Balay ridx = row%bs; 7592d61bbb3SSatish Balay high = nrow; 7602d61bbb3SSatish Balay low = 0; /* assume unsorted */ 7612d61bbb3SSatish Balay while (high-low > 5) { 762cd0e1443SSatish Balay t = (low+high)/2; 763cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 764cd0e1443SSatish Balay else low = t; 765cd0e1443SSatish Balay } 766cd0e1443SSatish Balay for (i=low; i<high; i++) { 767cd0e1443SSatish Balay if (rp[i] > bcol) break; 768cd0e1443SSatish Balay if (rp[i] == bcol) { 7692d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 7702d61bbb3SSatish Balay goto finished; 771cd0e1443SSatish Balay } 772cd0e1443SSatish Balay } 7732d61bbb3SSatish Balay *v++ = zero; 7742d61bbb3SSatish Balay finished:; 775cd0e1443SSatish Balay } 776cd0e1443SSatish Balay } 7773a40ed3dSBarry Smith PetscFunctionReturn(0); 778cd0e1443SSatish Balay } 779cd0e1443SSatish Balay 78095d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 7814a2ae208SSatish Balay #undef __FUNCT__ 7824a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 783f15d580aSBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode addv) 78495d5f7c2SBarry Smith { 785563b5814SBarry Smith Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data; 786563b5814SBarry Smith int ierr,i,N = m*n*b->bs2; 787563b5814SBarry Smith MatScalar *vsingle; 78895d5f7c2SBarry Smith 78995d5f7c2SBarry Smith PetscFunctionBegin; 790563b5814SBarry Smith if (N > b->setvalueslen) { 791563b5814SBarry Smith if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);} 792b0a32e0cSBarry Smith ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr); 793563b5814SBarry Smith b->setvalueslen = N; 794563b5814SBarry Smith } 795563b5814SBarry Smith vsingle = b->setvaluescopy; 79695d5f7c2SBarry Smith for (i=0; i<N; i++) { 79795d5f7c2SBarry Smith vsingle[i] = v[i]; 79895d5f7c2SBarry Smith } 79995d5f7c2SBarry Smith ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr); 80095d5f7c2SBarry Smith PetscFunctionReturn(0); 80195d5f7c2SBarry Smith } 80295d5f7c2SBarry Smith #endif 80395d5f7c2SBarry Smith 8042d61bbb3SSatish Balay 8054a2ae208SSatish Balay #undef __FUNCT__ 8064a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 807f15d580aSBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,const int im[],int n,const int in[],const MatScalar v[],InsertMode is) 80892c4ed94SBarry Smith { 80992c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 8108a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 811273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 812549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 813273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 814f15d580aSBarry Smith const MatScalar *value = v; 815f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 81692c4ed94SBarry Smith 8173a40ed3dSBarry Smith PetscFunctionBegin; 8180e324ae4SSatish Balay if (roworiented) { 8190e324ae4SSatish Balay stepval = (n-1)*bs; 8200e324ae4SSatish Balay } else { 8210e324ae4SSatish Balay stepval = (m-1)*bs; 8220e324ae4SSatish Balay } 82392c4ed94SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 82492c4ed94SBarry Smith row = im[k]; 8255ef9f2a5SBarry Smith if (row < 0) continue; 826aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 827590ac198SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1); 82892c4ed94SBarry Smith #endif 82992c4ed94SBarry Smith rp = aj + ai[row]; 83092c4ed94SBarry Smith ap = aa + bs2*ai[row]; 83192c4ed94SBarry Smith rmax = imax[row]; 83292c4ed94SBarry Smith nrow = ailen[row]; 83392c4ed94SBarry Smith low = 0; 83492c4ed94SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 8355ef9f2a5SBarry Smith if (in[l] < 0) continue; 836aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 837590ac198SBarry Smith if (in[l] >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],a->nbs-1); 83892c4ed94SBarry Smith #endif 83992c4ed94SBarry Smith col = in[l]; 84092c4ed94SBarry Smith if (roworiented) { 8410e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 8420e324ae4SSatish Balay } else { 8430e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 84492c4ed94SBarry Smith } 84592c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 84692c4ed94SBarry Smith while (high-low > 7) { 84792c4ed94SBarry Smith t = (low+high)/2; 84892c4ed94SBarry Smith if (rp[t] > col) high = t; 84992c4ed94SBarry Smith else low = t; 85092c4ed94SBarry Smith } 85192c4ed94SBarry Smith for (i=low; i<high; i++) { 85292c4ed94SBarry Smith if (rp[i] > col) break; 85392c4ed94SBarry Smith if (rp[i] == col) { 8548a84c255SSatish Balay bap = ap + bs2*i; 8550e324ae4SSatish Balay if (roworiented) { 8568a84c255SSatish Balay if (is == ADD_VALUES) { 857dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 858dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 8598a84c255SSatish Balay bap[jj] += *value++; 860dd9472c6SBarry Smith } 861dd9472c6SBarry Smith } 8620e324ae4SSatish Balay } else { 863dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 864dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 8650e324ae4SSatish Balay bap[jj] = *value++; 8668a84c255SSatish Balay } 867dd9472c6SBarry Smith } 868dd9472c6SBarry Smith } 8690e324ae4SSatish Balay } else { 8700e324ae4SSatish Balay if (is == ADD_VALUES) { 871dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 872dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 8730e324ae4SSatish Balay *bap++ += *value++; 874dd9472c6SBarry Smith } 875dd9472c6SBarry Smith } 8760e324ae4SSatish Balay } else { 877dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 878dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 8790e324ae4SSatish Balay *bap++ = *value++; 8800e324ae4SSatish Balay } 8818a84c255SSatish Balay } 882dd9472c6SBarry Smith } 883dd9472c6SBarry Smith } 884f1241b54SBarry Smith goto noinsert2; 88592c4ed94SBarry Smith } 88692c4ed94SBarry Smith } 88789280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 888a45adfd6SMatthew Knepley else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 88992c4ed94SBarry Smith if (nrow >= rmax) { 89092c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 89192c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 8923f1db9ecSBarry Smith MatScalar *new_a; 89392c4ed94SBarry Smith 894a45adfd6SMatthew Knepley if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 89589280ab3SLois Curfman McInnes 89692c4ed94SBarry Smith /* malloc new storage space */ 8973f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 898b0a32e0cSBarry Smith ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 89992c4ed94SBarry Smith new_j = (int*)(new_a + bs2*new_nz); 90092c4ed94SBarry Smith new_i = new_j + new_nz; 90192c4ed94SBarry Smith 90292c4ed94SBarry Smith /* copy over old data into new slots */ 90392c4ed94SBarry Smith for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 90492c4ed94SBarry Smith for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 905549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 90692c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 907549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 908549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 909549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 910549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 91192c4ed94SBarry Smith /* free up old matrix storage */ 912606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 913606d414cSSatish Balay if (!a->singlemalloc) { 914606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 915606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 916606d414cSSatish Balay } 91792c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 918c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 91992c4ed94SBarry Smith 92092c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 92192c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 922b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 92392c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 92492c4ed94SBarry Smith a->reallocs++; 92592c4ed94SBarry Smith a->nz++; 92692c4ed94SBarry Smith } 92792c4ed94SBarry Smith N = nrow++ - 1; 92892c4ed94SBarry Smith /* shift up all the later entries in this row */ 92992c4ed94SBarry Smith for (ii=N; ii>=i; ii--) { 93092c4ed94SBarry Smith rp[ii+1] = rp[ii]; 931549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 93292c4ed94SBarry Smith } 933549d3d68SSatish Balay if (N >= i) { 934549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 935549d3d68SSatish Balay } 93692c4ed94SBarry Smith rp[i] = col; 9378a84c255SSatish Balay bap = ap + bs2*i; 9380e324ae4SSatish Balay if (roworiented) { 939dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 940dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 9410e324ae4SSatish Balay bap[jj] = *value++; 942dd9472c6SBarry Smith } 943dd9472c6SBarry Smith } 9440e324ae4SSatish Balay } else { 945dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 946dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 9470e324ae4SSatish Balay *bap++ = *value++; 9480e324ae4SSatish Balay } 949dd9472c6SBarry Smith } 950dd9472c6SBarry Smith } 951f1241b54SBarry Smith noinsert2:; 95292c4ed94SBarry Smith low = i; 95392c4ed94SBarry Smith } 95492c4ed94SBarry Smith ailen[row] = nrow; 95592c4ed94SBarry Smith } 9563a40ed3dSBarry Smith PetscFunctionReturn(0); 95792c4ed94SBarry Smith } 95892c4ed94SBarry Smith 9594a2ae208SSatish Balay #undef __FUNCT__ 9604a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ" 9618f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 962584200bdSSatish Balay { 963584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 964584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 965273d9f13SBarry Smith int m = A->m,*ip,N,*ailen = a->ilen; 966549d3d68SSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 9673f1db9ecSBarry Smith MatScalar *aa = a->a,*ap; 968584200bdSSatish Balay 9693a40ed3dSBarry Smith PetscFunctionBegin; 9703a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 971584200bdSSatish Balay 97243ee02c3SBarry Smith if (m) rmax = ailen[0]; 973584200bdSSatish Balay for (i=1; i<mbs; i++) { 974584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 975584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 976d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 977584200bdSSatish Balay if (fshift) { 978a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 979584200bdSSatish Balay N = ailen[i]; 980584200bdSSatish Balay for (j=0; j<N; j++) { 981584200bdSSatish Balay ip[j-fshift] = ip[j]; 982549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 983584200bdSSatish Balay } 984584200bdSSatish Balay } 985584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 986584200bdSSatish Balay } 987584200bdSSatish Balay if (mbs) { 988584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 989584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 990584200bdSSatish Balay } 991584200bdSSatish Balay /* reset ilen and imax for each row */ 992584200bdSSatish Balay for (i=0; i<mbs; i++) { 993584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 994584200bdSSatish Balay } 995a7c10996SSatish Balay a->nz = ai[mbs]; 996584200bdSSatish Balay 997584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 998584200bdSSatish Balay if (fshift && a->diag) { 999606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 1000b424e231SHong Zhang PetscLogObjectMemory(A,-(mbs+1)*sizeof(int)); 1001584200bdSSatish Balay a->diag = 0; 1002584200bdSSatish Balay } 1003bba1ac68SSatish Balay PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",m,A->n,a->bs,fshift*bs2,a->nz*bs2); 1004bba1ac68SSatish Balay PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs); 1005b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 1006e2f3b5e9SSatish Balay a->reallocs = 0; 10070e6d2581SBarry Smith A->info.nz_unneeded = (PetscReal)fshift*bs2; 1008cf4441caSHong Zhang 10093a40ed3dSBarry Smith PetscFunctionReturn(0); 1010584200bdSSatish Balay } 1011584200bdSSatish Balay 10125a838052SSatish Balay 1013bea157c4SSatish Balay 1014bea157c4SSatish Balay /* 1015bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 1016bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 1017bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 1018bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 1019bea157c4SSatish Balay */ 10204a2ae208SSatish Balay #undef __FUNCT__ 10214a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks" 1022bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 1023d9b7c43dSSatish Balay { 1024bea157c4SSatish Balay int i,j,k,row; 1025bea157c4SSatish Balay PetscTruth flg; 10263a40ed3dSBarry Smith 1027433994e6SBarry Smith PetscFunctionBegin; 1028bea157c4SSatish Balay for (i=0,j=0; i<n; j++) { 1029bea157c4SSatish Balay row = idx[i]; 1030bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 1031bea157c4SSatish Balay sizes[j] = 1; 1032bea157c4SSatish Balay i++; 1033e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 1034bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 1035bea157c4SSatish Balay i++; 1036bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 1037bea157c4SSatish Balay flg = PETSC_TRUE; 1038bea157c4SSatish Balay for (k=1; k<bs; k++) { 1039bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 1040bea157c4SSatish Balay flg = PETSC_FALSE; 1041bea157c4SSatish Balay break; 1042d9b7c43dSSatish Balay } 1043bea157c4SSatish Balay } 1044bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 1045bea157c4SSatish Balay sizes[j] = bs; 1046bea157c4SSatish Balay i+= bs; 1047bea157c4SSatish Balay } else { 1048bea157c4SSatish Balay sizes[j] = 1; 1049bea157c4SSatish Balay i++; 1050bea157c4SSatish Balay } 1051bea157c4SSatish Balay } 1052bea157c4SSatish Balay } 1053bea157c4SSatish Balay *bs_max = j; 10543a40ed3dSBarry Smith PetscFunctionReturn(0); 1055d9b7c43dSSatish Balay } 1056d9b7c43dSSatish Balay 10574a2ae208SSatish Balay #undef __FUNCT__ 10584a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ" 1059268466fbSBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,const PetscScalar *diag) 1060d9b7c43dSSatish Balay { 1061d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 1062b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 1063bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 106487828ca2SBarry Smith PetscScalar zero = 0.0; 10653f1db9ecSBarry Smith MatScalar *aa; 1066d9b7c43dSSatish Balay 10673a40ed3dSBarry Smith PetscFunctionBegin; 1068d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 1069b9b97703SBarry Smith ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr); 1070d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 1071d9b7c43dSSatish Balay 1072bea157c4SSatish Balay /* allocate memory for rows,sizes */ 1073b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 1074bea157c4SSatish Balay sizes = rows + is_n; 1075bea157c4SSatish Balay 1076563b5814SBarry Smith /* copy IS values to rows, and sort them */ 1077bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 1078bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 1079dffd3267SBarry Smith if (baij->keepzeroedrows) { 1080dffd3267SBarry Smith for (i=0; i<is_n; i++) { sizes[i] = 1; } 1081dffd3267SBarry Smith bs_max = is_n; 1082dffd3267SBarry Smith } else { 1083bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 1084dffd3267SBarry Smith } 1085b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 1086bea157c4SSatish Balay 1087bea157c4SSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 1088bea157c4SSatish Balay row = rows[j]; 1089273d9f13SBarry Smith if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 1090bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 1091bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 1092dffd3267SBarry Smith if (sizes[i] == bs && !baij->keepzeroedrows) { 1093bea157c4SSatish Balay if (diag) { 1094bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 1095bea157c4SSatish Balay baij->ilen[row/bs] = 1; 1096bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 1097bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 1098a07cd24cSSatish Balay } 1099563b5814SBarry Smith /* Now insert all the diagonal values for this bs */ 1100bea157c4SSatish Balay for (k=0; k<bs; k++) { 1101bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 1102bea157c4SSatish Balay } 1103bea157c4SSatish Balay } else { /* (!diag) */ 1104bea157c4SSatish Balay baij->ilen[row/bs] = 0; 1105bea157c4SSatish Balay } /* end (!diag) */ 1106bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 1107aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 110829bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 1109bea157c4SSatish Balay #endif 1110bea157c4SSatish Balay for (k=0; k<count; k++) { 1111d9b7c43dSSatish Balay aa[0] = zero; 1112d9b7c43dSSatish Balay aa += bs; 1113d9b7c43dSSatish Balay } 1114d9b7c43dSSatish Balay if (diag) { 1115f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 1116d9b7c43dSSatish Balay } 1117d9b7c43dSSatish Balay } 1118bea157c4SSatish Balay } 1119bea157c4SSatish Balay 1120606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 11219a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11223a40ed3dSBarry Smith PetscFunctionReturn(0); 1123d9b7c43dSSatish Balay } 11241c351548SSatish Balay 11254a2ae208SSatish Balay #undef __FUNCT__ 11264a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ" 1127f15d580aSBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is) 11282d61bbb3SSatish Balay { 11292d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 11302d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 1131273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 11322d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 1133549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 1134273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 11353f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 11362d61bbb3SSatish Balay 11372d61bbb3SSatish Balay PetscFunctionBegin; 11382d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 11392d61bbb3SSatish Balay row = im[k]; brow = row/bs; 11405ef9f2a5SBarry Smith if (row < 0) continue; 1141aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1142590ac198SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1); 11432d61bbb3SSatish Balay #endif 11442d61bbb3SSatish Balay rp = aj + ai[brow]; 11452d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 11462d61bbb3SSatish Balay rmax = imax[brow]; 11472d61bbb3SSatish Balay nrow = ailen[brow]; 11482d61bbb3SSatish Balay low = 0; 11492d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 11505ef9f2a5SBarry Smith if (in[l] < 0) continue; 1151aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1152590ac198SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1); 11532d61bbb3SSatish Balay #endif 11542d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 11552d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 11562d61bbb3SSatish Balay if (roworiented) { 11575ef9f2a5SBarry Smith value = v[l + k*n]; 11582d61bbb3SSatish Balay } else { 11592d61bbb3SSatish Balay value = v[k + l*m]; 11602d61bbb3SSatish Balay } 11612d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 11622d61bbb3SSatish Balay while (high-low > 7) { 11632d61bbb3SSatish Balay t = (low+high)/2; 11642d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 11652d61bbb3SSatish Balay else low = t; 11662d61bbb3SSatish Balay } 11672d61bbb3SSatish Balay for (i=low; i<high; i++) { 11682d61bbb3SSatish Balay if (rp[i] > bcol) break; 11692d61bbb3SSatish Balay if (rp[i] == bcol) { 11702d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 11712d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 11722d61bbb3SSatish Balay else *bap = value; 11732d61bbb3SSatish Balay goto noinsert1; 11742d61bbb3SSatish Balay } 11752d61bbb3SSatish Balay } 11762d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 1177a45adfd6SMatthew Knepley else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 11782d61bbb3SSatish Balay if (nrow >= rmax) { 11792d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 11802d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 11813f1db9ecSBarry Smith MatScalar *new_a; 11822d61bbb3SSatish Balay 1183a45adfd6SMatthew Knepley if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 11842d61bbb3SSatish Balay 11852d61bbb3SSatish Balay /* Malloc new storage space */ 11863f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 1187b0a32e0cSBarry Smith ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 11882d61bbb3SSatish Balay new_j = (int*)(new_a + bs2*new_nz); 11892d61bbb3SSatish Balay new_i = new_j + new_nz; 11902d61bbb3SSatish Balay 11912d61bbb3SSatish Balay /* copy over old data into new slots */ 11922d61bbb3SSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 11932d61bbb3SSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1194549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 11952d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1196549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1197549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1198549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1199549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 12002d61bbb3SSatish Balay /* free up old matrix storage */ 1201606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1202606d414cSSatish Balay if (!a->singlemalloc) { 1203606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1204606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1205606d414cSSatish Balay } 12062d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1207c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 12082d61bbb3SSatish Balay 12092d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 12102d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 1211b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 12122d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 12132d61bbb3SSatish Balay a->reallocs++; 12142d61bbb3SSatish Balay a->nz++; 12152d61bbb3SSatish Balay } 12162d61bbb3SSatish Balay N = nrow++ - 1; 12172d61bbb3SSatish Balay /* shift up all the later entries in this row */ 12182d61bbb3SSatish Balay for (ii=N; ii>=i; ii--) { 12192d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1220549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 12212d61bbb3SSatish Balay } 1222549d3d68SSatish Balay if (N>=i) { 1223549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1224549d3d68SSatish Balay } 12252d61bbb3SSatish Balay rp[i] = bcol; 12262d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 12272d61bbb3SSatish Balay noinsert1:; 12282d61bbb3SSatish Balay low = i; 12292d61bbb3SSatish Balay } 12302d61bbb3SSatish Balay ailen[brow] = nrow; 12312d61bbb3SSatish Balay } 12322d61bbb3SSatish Balay PetscFunctionReturn(0); 12332d61bbb3SSatish Balay } 12342d61bbb3SSatish Balay 12352d61bbb3SSatish Balay 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ" 1238b380c88cSHong Zhang int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatFactorInfo *info) 12392d61bbb3SSatish Balay { 12402d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 12412d61bbb3SSatish Balay Mat outA; 12422d61bbb3SSatish Balay int ierr; 1243667159a5SBarry Smith PetscTruth row_identity,col_identity; 12442d61bbb3SSatish Balay 12452d61bbb3SSatish Balay PetscFunctionBegin; 1246d3d32019SBarry Smith if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 1247667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1248667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1249667159a5SBarry Smith if (!row_identity || !col_identity) { 125029bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ILU"); 1251667159a5SBarry Smith } 12522d61bbb3SSatish Balay 12532d61bbb3SSatish Balay outA = inA; 12542d61bbb3SSatish Balay inA->factor = FACTOR_LU; 12552d61bbb3SSatish Balay 12562d61bbb3SSatish Balay if (!a->diag) { 1257c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 12582d61bbb3SSatish Balay } 1259cf242676SKris Buschelman 1260c38d4ed2SBarry Smith a->row = row; 1261c38d4ed2SBarry Smith a->col = col; 1262c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1263c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1264c38d4ed2SBarry Smith 1265c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 12664c49b128SBarry Smith ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); 1267b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 1268c38d4ed2SBarry Smith 1269cf242676SKris Buschelman /* 1270cf242676SKris Buschelman Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 1271cf242676SKris Buschelman for ILU(0) factorization with natural ordering 1272cf242676SKris Buschelman */ 1273cf242676SKris Buschelman if (a->bs < 8) { 1274cf242676SKris Buschelman ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr); 1275cf242676SKris Buschelman } else { 1276c38d4ed2SBarry Smith if (!a->solve_work) { 127787828ca2SBarry Smith ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 127887828ca2SBarry Smith PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar)); 1279c38d4ed2SBarry Smith } 12802d61bbb3SSatish Balay } 12812d61bbb3SSatish Balay 1282667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1283667159a5SBarry Smith 12842d61bbb3SSatish Balay PetscFunctionReturn(0); 12852d61bbb3SSatish Balay } 12864a2ae208SSatish Balay #undef __FUNCT__ 12874a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ" 1288ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1289ba4ca20aSSatish Balay { 1290c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1291ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1292d132466eSBarry Smith int ierr; 1293ba4ca20aSSatish Balay 12943a40ed3dSBarry Smith PetscFunctionBegin; 1295c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1296d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1297d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12983a40ed3dSBarry Smith PetscFunctionReturn(0); 1299ba4ca20aSSatish Balay } 1300d9b7c43dSSatish Balay 1301fb2e594dSBarry Smith EXTERN_C_BEGIN 13024a2ae208SSatish Balay #undef __FUNCT__ 13034a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 130427a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 130527a8da17SBarry Smith { 130627a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 130714db4f74SSatish Balay int i,nz,nbs; 130827a8da17SBarry Smith 130927a8da17SBarry Smith PetscFunctionBegin; 131014db4f74SSatish Balay nz = baij->maxnz/baij->bs2; 131114db4f74SSatish Balay nbs = baij->nbs; 131227a8da17SBarry Smith for (i=0; i<nz; i++) { 131327a8da17SBarry Smith baij->j[i] = indices[i]; 131427a8da17SBarry Smith } 131527a8da17SBarry Smith baij->nz = nz; 131614db4f74SSatish Balay for (i=0; i<nbs; i++) { 131727a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 131827a8da17SBarry Smith } 131927a8da17SBarry Smith 132027a8da17SBarry Smith PetscFunctionReturn(0); 132127a8da17SBarry Smith } 1322fb2e594dSBarry Smith EXTERN_C_END 132327a8da17SBarry Smith 13244a2ae208SSatish Balay #undef __FUNCT__ 13254a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices" 132627a8da17SBarry Smith /*@ 132727a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 132827a8da17SBarry Smith in the matrix. 132927a8da17SBarry Smith 133027a8da17SBarry Smith Input Parameters: 133127a8da17SBarry Smith + mat - the SeqBAIJ matrix 133227a8da17SBarry Smith - indices - the column indices 133327a8da17SBarry Smith 133415091d37SBarry Smith Level: advanced 133515091d37SBarry Smith 133627a8da17SBarry Smith Notes: 133727a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 133827a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 133927a8da17SBarry Smith of the MatSetValues() operation. 134027a8da17SBarry Smith 134127a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 134227a8da17SBarry Smith MatCreateSeqBAIJ(). 134327a8da17SBarry Smith 134427a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 134527a8da17SBarry Smith 134627a8da17SBarry Smith @*/ 134727a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 134827a8da17SBarry Smith { 134927a8da17SBarry Smith int ierr,(*f)(Mat,int *); 135027a8da17SBarry Smith 135127a8da17SBarry Smith PetscFunctionBegin; 135227a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1353c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 135427a8da17SBarry Smith if (f) { 135527a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 135627a8da17SBarry Smith } else { 135729bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 135827a8da17SBarry Smith } 135927a8da17SBarry Smith PetscFunctionReturn(0); 136027a8da17SBarry Smith } 136127a8da17SBarry Smith 13624a2ae208SSatish Balay #undef __FUNCT__ 13634a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ" 1364273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v) 1365273d9f13SBarry Smith { 1366273d9f13SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1367273d9f13SBarry Smith int ierr,i,j,n,row,bs,*ai,*aj,mbs; 1368273d9f13SBarry Smith PetscReal atmp; 136987828ca2SBarry Smith PetscScalar *x,zero = 0.0; 1370273d9f13SBarry Smith MatScalar *aa; 1371273d9f13SBarry Smith int ncols,brow,krow,kcol; 1372273d9f13SBarry Smith 1373273d9f13SBarry Smith PetscFunctionBegin; 1374273d9f13SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1375273d9f13SBarry Smith bs = a->bs; 1376273d9f13SBarry Smith aa = a->a; 1377273d9f13SBarry Smith ai = a->i; 1378273d9f13SBarry Smith aj = a->j; 1379273d9f13SBarry Smith mbs = a->mbs; 1380273d9f13SBarry Smith 1381273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 1382273d9f13SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1383273d9f13SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1384273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1385273d9f13SBarry Smith for (i=0; i<mbs; i++) { 1386273d9f13SBarry Smith ncols = ai[1] - ai[0]; ai++; 1387273d9f13SBarry Smith brow = bs*i; 1388273d9f13SBarry Smith for (j=0; j<ncols; j++){ 1389273d9f13SBarry Smith /* bcol = bs*(*aj); */ 1390273d9f13SBarry Smith for (kcol=0; kcol<bs; kcol++){ 1391273d9f13SBarry Smith for (krow=0; krow<bs; krow++){ 1392273d9f13SBarry Smith atmp = PetscAbsScalar(*aa); aa++; 1393273d9f13SBarry Smith row = brow + krow; /* row index */ 1394273d9f13SBarry Smith /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */ 1395273d9f13SBarry Smith if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp; 1396273d9f13SBarry Smith } 1397273d9f13SBarry Smith } 1398273d9f13SBarry Smith aj++; 1399273d9f13SBarry Smith } 1400273d9f13SBarry Smith } 1401273d9f13SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1402273d9f13SBarry Smith PetscFunctionReturn(0); 1403273d9f13SBarry Smith } 1404273d9f13SBarry Smith 14054a2ae208SSatish Balay #undef __FUNCT__ 14064a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ" 1407273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A) 1408273d9f13SBarry Smith { 1409273d9f13SBarry Smith int ierr; 1410273d9f13SBarry Smith 1411273d9f13SBarry Smith PetscFunctionBegin; 1412273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1413273d9f13SBarry Smith PetscFunctionReturn(0); 1414273d9f13SBarry Smith } 1415273d9f13SBarry Smith 14164a2ae208SSatish Balay #undef __FUNCT__ 14174a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ" 14184e7234bfSBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar *array[]) 1419f2a5309cSSatish Balay { 1420f2a5309cSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1421f2a5309cSSatish Balay PetscFunctionBegin; 1422f2a5309cSSatish Balay *array = a->a; 1423f2a5309cSSatish Balay PetscFunctionReturn(0); 1424f2a5309cSSatish Balay } 1425f2a5309cSSatish Balay 14264a2ae208SSatish Balay #undef __FUNCT__ 14274a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ" 14284e7234bfSBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar *array[]) 1429f2a5309cSSatish Balay { 1430f2a5309cSSatish Balay PetscFunctionBegin; 1431f2a5309cSSatish Balay PetscFunctionReturn(0); 1432f2a5309cSSatish Balay } 1433f2a5309cSSatish Balay 143442ee4b1aSHong Zhang #include "petscblaslapack.h" 143542ee4b1aSHong Zhang #undef __FUNCT__ 143642ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqBAIJ" 1437268466fbSBarry Smith int MatAXPY_SeqBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 143842ee4b1aSHong Zhang { 143942ee4b1aSHong Zhang Mat_SeqBAIJ *x = (Mat_SeqBAIJ *)X->data,*y = (Mat_SeqBAIJ *)Y->data; 14406550863bSHong Zhang int ierr,one=1,i,bs=y->bs,j,bs2; 144142ee4b1aSHong Zhang 144242ee4b1aSHong Zhang PetscFunctionBegin; 144342ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1444268466fbSBarry Smith BLaxpy_(&x->nz,(PetscScalar*)a,x->a,&one,y->a,&one); 1445c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1446c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1447c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1448c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1449c537a176SHong Zhang } 1450c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1451c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1452c4319e64SHong Zhang y->XtoY = X; 1453c537a176SHong Zhang } 1454c4319e64SHong Zhang bs2 = bs*bs; 1455c537a176SHong Zhang for (i=0; i<x->nz; i++) { 1456c4319e64SHong Zhang j = 0; 1457c4319e64SHong Zhang while (j < bs2){ 14586550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1459c4319e64SHong Zhang j++; 1460c537a176SHong Zhang } 1461c4319e64SHong Zhang } 1462c4319e64SHong Zhang PetscLogInfo(0,"MatAXPY_SeqBAIJ: ratio of nnz(X)/nnz(Y): %d/%d = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz)); 146342ee4b1aSHong Zhang } else { 146442ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 146542ee4b1aSHong Zhang } 146642ee4b1aSHong Zhang PetscFunctionReturn(0); 146742ee4b1aSHong Zhang } 146842ee4b1aSHong Zhang 14692593348eSBarry Smith /* -------------------------------------------------------------------*/ 1470cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1471cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1472cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1473cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 147497304618SKris Buschelman /* 4*/ MatMultAdd_SeqBAIJ_N, 14757c922b88SBarry Smith MatMultTranspose_SeqBAIJ, 14767c922b88SBarry Smith MatMultTransposeAdd_SeqBAIJ, 1477cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1478cc2dc46cSBarry Smith 0, 1479cc2dc46cSBarry Smith 0, 148097304618SKris Buschelman /*10*/ 0, 1481cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1482cc2dc46cSBarry Smith 0, 1483b6490206SBarry Smith 0, 1484f2501298SSatish Balay MatTranspose_SeqBAIJ, 148597304618SKris Buschelman /*15*/ MatGetInfo_SeqBAIJ, 1486cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1487cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1488cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1489cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 149097304618SKris Buschelman /*20*/ 0, 1491cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1492cc2dc46cSBarry Smith 0, 1493cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1494cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 149597304618SKris Buschelman /*25*/ MatZeroRows_SeqBAIJ, 1496cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1497cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1498cc2dc46cSBarry Smith 0, 1499cc2dc46cSBarry Smith 0, 150097304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqBAIJ, 1501cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1502cc2dc46cSBarry Smith 0, 1503f2a5309cSSatish Balay MatGetArray_SeqBAIJ, 1504f2a5309cSSatish Balay MatRestoreArray_SeqBAIJ, 150597304618SKris Buschelman /*35*/ MatDuplicate_SeqBAIJ, 1506cc2dc46cSBarry Smith 0, 1507cc2dc46cSBarry Smith 0, 1508cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1509cc2dc46cSBarry Smith 0, 151097304618SKris Buschelman /*40*/ MatAXPY_SeqBAIJ, 1511cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1512cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1513cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1514cc2dc46cSBarry Smith 0, 151597304618SKris Buschelman /*45*/ MatPrintHelp_SeqBAIJ, 1516cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1517cc2dc46cSBarry Smith 0, 1518cc2dc46cSBarry Smith 0, 1519cc2dc46cSBarry Smith 0, 152097304618SKris Buschelman /*50*/ MatGetBlockSize_SeqBAIJ, 15213b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 152292c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1523cc2dc46cSBarry Smith 0, 1524cc2dc46cSBarry Smith 0, 152597304618SKris Buschelman /*55*/ 0, 1526cc2dc46cSBarry Smith 0, 1527cc2dc46cSBarry Smith 0, 1528cc2dc46cSBarry Smith 0, 1529d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 153097304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqBAIJ, 1531b9b97703SBarry Smith MatDestroy_SeqBAIJ, 1532b9b97703SBarry Smith MatView_SeqBAIJ, 15333a64cc32SBarry Smith MatGetPetscMaps_Petsc, 1534273d9f13SBarry Smith 0, 153597304618SKris Buschelman /*65*/ 0, 1536273d9f13SBarry Smith 0, 1537273d9f13SBarry Smith 0, 1538273d9f13SBarry Smith 0, 1539273d9f13SBarry Smith 0, 154097304618SKris Buschelman /*70*/ MatGetRowMax_SeqBAIJ, 154197304618SKris Buschelman MatConvert_Basic, 1542273d9f13SBarry Smith 0, 154397304618SKris Buschelman 0, 154497304618SKris Buschelman 0, 154597304618SKris Buschelman /*75*/ 0, 154697304618SKris Buschelman 0, 154797304618SKris Buschelman 0, 154897304618SKris Buschelman 0, 154997304618SKris Buschelman 0, 155097304618SKris Buschelman /*80*/ 0, 155197304618SKris Buschelman 0, 155297304618SKris Buschelman 0, 155397304618SKris Buschelman 0, 155497304618SKris Buschelman 0, 155597304618SKris Buschelman /*85*/ MatLoad_SeqBAIJ 155697304618SKris Buschelman }; 15572593348eSBarry Smith 15583e90b805SBarry Smith EXTERN_C_BEGIN 15594a2ae208SSatish Balay #undef __FUNCT__ 15604a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ" 15613e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 15623e90b805SBarry Smith { 15633e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1564273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2; 1565d132466eSBarry Smith int ierr; 15663e90b805SBarry Smith 15673e90b805SBarry Smith PetscFunctionBegin; 15683e90b805SBarry Smith if (aij->nonew != 1) { 156929bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 15703e90b805SBarry Smith } 15713e90b805SBarry Smith 15723e90b805SBarry Smith /* allocate space for values if not already there */ 15733e90b805SBarry Smith if (!aij->saved_values) { 157487828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 15753e90b805SBarry Smith } 15763e90b805SBarry Smith 15773e90b805SBarry Smith /* copy values over */ 157887828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 15793e90b805SBarry Smith PetscFunctionReturn(0); 15803e90b805SBarry Smith } 15813e90b805SBarry Smith EXTERN_C_END 15823e90b805SBarry Smith 15833e90b805SBarry Smith EXTERN_C_BEGIN 15844a2ae208SSatish Balay #undef __FUNCT__ 15854a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ" 158632e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 15873e90b805SBarry Smith { 15883e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1589273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 15903e90b805SBarry Smith 15913e90b805SBarry Smith PetscFunctionBegin; 15923e90b805SBarry Smith if (aij->nonew != 1) { 159329bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 15943e90b805SBarry Smith } 15953e90b805SBarry Smith if (!aij->saved_values) { 159629bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 15973e90b805SBarry Smith } 15983e90b805SBarry Smith 15993e90b805SBarry Smith /* copy values over */ 160087828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 16013e90b805SBarry Smith PetscFunctionReturn(0); 16023e90b805SBarry Smith } 16033e90b805SBarry Smith EXTERN_C_END 16043e90b805SBarry Smith 1605273d9f13SBarry Smith EXTERN_C_BEGIN 1606273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*); 1607273d9f13SBarry Smith EXTERN_C_END 1608273d9f13SBarry Smith 1609273d9f13SBarry Smith EXTERN_C_BEGIN 16104a2ae208SSatish Balay #undef __FUNCT__ 1611a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqBAIJSetPreallocation_SeqBAIJ" 1612a23d5eceSKris Buschelman int MatSeqBAIJSetPreallocation_SeqBAIJ(Mat B,int bs,int nz,int *nnz) 1613a23d5eceSKris Buschelman { 1614a23d5eceSKris Buschelman Mat_SeqBAIJ *b; 1615a23d5eceSKris Buschelman int i,len,ierr,mbs,nbs,bs2,newbs = bs; 1616a23d5eceSKris Buschelman PetscTruth flg; 1617a23d5eceSKris Buschelman 1618a23d5eceSKris Buschelman PetscFunctionBegin; 1619a23d5eceSKris Buschelman 1620a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 1621a23d5eceSKris Buschelman ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr); 1622a23d5eceSKris Buschelman if (nnz && newbs != bs) { 1623a23d5eceSKris Buschelman SETERRQ(1,"Cannot change blocksize from command line if setting nnz"); 1624a23d5eceSKris Buschelman } 1625a23d5eceSKris Buschelman bs = newbs; 1626a23d5eceSKris Buschelman 1627a23d5eceSKris Buschelman mbs = B->m/bs; 1628a23d5eceSKris Buschelman nbs = B->n/bs; 1629a23d5eceSKris Buschelman bs2 = bs*bs; 1630a23d5eceSKris Buschelman 1631a23d5eceSKris Buschelman if (mbs*bs!=B->m || nbs*bs!=B->n) { 1632a23d5eceSKris Buschelman SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs); 1633a23d5eceSKris Buschelman } 1634a23d5eceSKris Buschelman 1635a23d5eceSKris Buschelman if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1636a23d5eceSKris Buschelman if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 1637a23d5eceSKris Buschelman if (nnz) { 1638a23d5eceSKris Buschelman for (i=0; i<mbs; i++) { 1639a23d5eceSKris Buschelman if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1640a23d5eceSKris Buschelman if (nnz[i] > nbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %d value %d rowlength %d",i,nnz[i],nbs); 1641a23d5eceSKris Buschelman } 1642a23d5eceSKris Buschelman } 1643a23d5eceSKris Buschelman 1644a23d5eceSKris Buschelman b = (Mat_SeqBAIJ*)B->data; 1645a23d5eceSKris Buschelman ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 1646a23d5eceSKris Buschelman B->ops->solve = MatSolve_SeqBAIJ_Update; 1647a23d5eceSKris Buschelman B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_Update; 1648a23d5eceSKris Buschelman if (!flg) { 1649a23d5eceSKris Buschelman switch (bs) { 1650a23d5eceSKris Buschelman case 1: 1651a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1652a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_1; 1653a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_1; 1654a23d5eceSKris Buschelman break; 1655a23d5eceSKris Buschelman case 2: 1656a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1657a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_2; 1658a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_2; 1659a23d5eceSKris Buschelman break; 1660a23d5eceSKris Buschelman case 3: 1661a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1662a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_3; 1663a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_3; 1664a23d5eceSKris Buschelman break; 1665a23d5eceSKris Buschelman case 4: 1666a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1667a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_4; 1668a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_4; 1669a23d5eceSKris Buschelman break; 1670a23d5eceSKris Buschelman case 5: 1671a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1672a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_5; 1673a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_5; 1674a23d5eceSKris Buschelman break; 1675a23d5eceSKris Buschelman case 6: 1676a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 1677a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_6; 1678a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_6; 1679a23d5eceSKris Buschelman break; 1680a23d5eceSKris Buschelman case 7: 1681a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7; 1682a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_7; 1683a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_7; 1684a23d5eceSKris Buschelman break; 1685a23d5eceSKris Buschelman default: 1686a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 1687a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_N; 1688a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_N; 1689a23d5eceSKris Buschelman break; 1690a23d5eceSKris Buschelman } 1691a23d5eceSKris Buschelman } 1692a23d5eceSKris Buschelman b->bs = bs; 1693a23d5eceSKris Buschelman b->mbs = mbs; 1694a23d5eceSKris Buschelman b->nbs = nbs; 1695a23d5eceSKris Buschelman ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 1696a23d5eceSKris Buschelman if (!nnz) { 1697a23d5eceSKris Buschelman if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1698a23d5eceSKris Buschelman else if (nz <= 0) nz = 1; 1699a23d5eceSKris Buschelman for (i=0; i<mbs; i++) b->imax[i] = nz; 1700a23d5eceSKris Buschelman nz = nz*mbs; 1701a23d5eceSKris Buschelman } else { 1702a23d5eceSKris Buschelman nz = 0; 1703a23d5eceSKris Buschelman for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 1704a23d5eceSKris Buschelman } 1705a23d5eceSKris Buschelman 1706a23d5eceSKris Buschelman /* allocate the matrix space */ 1707a23d5eceSKris Buschelman len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 1708a23d5eceSKris Buschelman ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 1709a23d5eceSKris Buschelman ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1710a23d5eceSKris Buschelman b->j = (int*)(b->a + nz*bs2); 1711a23d5eceSKris Buschelman ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 1712a23d5eceSKris Buschelman b->i = b->j + nz; 1713a23d5eceSKris Buschelman b->singlemalloc = PETSC_TRUE; 1714a23d5eceSKris Buschelman 1715a23d5eceSKris Buschelman b->i[0] = 0; 1716a23d5eceSKris Buschelman for (i=1; i<mbs+1; i++) { 1717a23d5eceSKris Buschelman b->i[i] = b->i[i-1] + b->imax[i-1]; 1718a23d5eceSKris Buschelman } 1719a23d5eceSKris Buschelman 1720a23d5eceSKris Buschelman /* b->ilen will count nonzeros in each block row so far. */ 1721a23d5eceSKris Buschelman ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1722a23d5eceSKris Buschelman PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1723a23d5eceSKris Buschelman for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 1724a23d5eceSKris Buschelman 1725a23d5eceSKris Buschelman b->bs = bs; 1726a23d5eceSKris Buschelman b->bs2 = bs2; 1727a23d5eceSKris Buschelman b->mbs = mbs; 1728a23d5eceSKris Buschelman b->nz = 0; 1729a23d5eceSKris Buschelman b->maxnz = nz*bs2; 1730a23d5eceSKris Buschelman B->info.nz_unneeded = (PetscReal)b->maxnz; 1731a23d5eceSKris Buschelman PetscFunctionReturn(0); 1732a23d5eceSKris Buschelman } 1733a23d5eceSKris Buschelman EXTERN_C_END 1734a23d5eceSKris Buschelman 17350bad9183SKris Buschelman /*MC 1736fafad747SKris Buschelman MATSEQBAIJ - MATSEQBAIJ = "seqbaij" - A matrix type to be used for sequential block sparse matrices, based on 17370bad9183SKris Buschelman block sparse compressed row format. 17380bad9183SKris Buschelman 17390bad9183SKris Buschelman Options Database Keys: 17400bad9183SKris Buschelman . -mat_type seqbaij - sets the matrix type to "seqbaij" during a call to MatSetFromOptions() 17410bad9183SKris Buschelman 17420bad9183SKris Buschelman Level: beginner 17430bad9183SKris Buschelman 17440bad9183SKris Buschelman .seealso: MatCreateSeqBAIJ 17450bad9183SKris Buschelman M*/ 17460bad9183SKris Buschelman 1747a23d5eceSKris Buschelman EXTERN_C_BEGIN 1748a23d5eceSKris Buschelman #undef __FUNCT__ 17494a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ" 1750273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B) 17512593348eSBarry Smith { 1752273d9f13SBarry Smith int ierr,size; 1753b6490206SBarry Smith Mat_SeqBAIJ *b; 17543b2fbd54SBarry Smith 17553a40ed3dSBarry Smith PetscFunctionBegin; 1756273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 175729bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1758b6490206SBarry Smith 1759273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1760273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1761b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr); 1762b0a32e0cSBarry Smith B->data = (void*)b; 1763549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1764549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 17652593348eSBarry Smith B->factor = 0; 17662593348eSBarry Smith B->lupivotthreshold = 1.0; 176790f02eecSBarry Smith B->mapping = 0; 17682593348eSBarry Smith b->row = 0; 17692593348eSBarry Smith b->col = 0; 1770e51c0b9cSSatish Balay b->icol = 0; 17712593348eSBarry Smith b->reallocs = 0; 17723e90b805SBarry Smith b->saved_values = 0; 1773cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 1774563b5814SBarry Smith b->setvalueslen = 0; 1775563b5814SBarry Smith b->setvaluescopy = PETSC_NULL; 1776563b5814SBarry Smith #endif 17772593348eSBarry Smith 17783a64cc32SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 17793a64cc32SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1780a5ae1ecdSBarry Smith 1781c4992f7dSBarry Smith b->sorted = PETSC_FALSE; 1782c4992f7dSBarry Smith b->roworiented = PETSC_TRUE; 17832593348eSBarry Smith b->nonew = 0; 17842593348eSBarry Smith b->diag = 0; 17852593348eSBarry Smith b->solve_work = 0; 1786de6a44a3SBarry Smith b->mult_work = 0; 17872a1b7f2aSHong Zhang B->spptr = 0; 17880e6d2581SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 1789883fce79SBarry Smith b->keepzeroedrows = PETSC_FALSE; 1790c4319e64SHong Zhang b->xtoy = 0; 1791c4319e64SHong Zhang b->XtoY = 0; 17924e220ebcSLois Curfman McInnes 1793f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 17943e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 1795bc4b532fSSatish Balay MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1796f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 17973e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 1798bc4b532fSSatish Balay MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1799f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 180027a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1801bc4b532fSSatish Balay MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 1802a6175056SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqaij_C", 1803273d9f13SBarry Smith "MatConvert_SeqBAIJ_SeqAIJ", 1804273d9f13SBarry Smith MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr); 1805a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetPreallocation_C", 1806a23d5eceSKris Buschelman "MatSeqBAIJSetPreallocation_SeqBAIJ", 1807a23d5eceSKris Buschelman MatSeqBAIJSetPreallocation_SeqBAIJ);CHKERRQ(ierr); 18083a40ed3dSBarry Smith PetscFunctionReturn(0); 18092593348eSBarry Smith } 1810273d9f13SBarry Smith EXTERN_C_END 18112593348eSBarry Smith 18124a2ae208SSatish Balay #undef __FUNCT__ 18134a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ" 18142e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 18152593348eSBarry Smith { 18162593348eSBarry Smith Mat C; 1817b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data; 181827a8da17SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1819de6a44a3SBarry Smith 18203a40ed3dSBarry Smith PetscFunctionBegin; 182129bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 18222593348eSBarry Smith 18232593348eSBarry Smith *B = 0; 1824273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1825273d9f13SBarry Smith ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr); 1826273d9f13SBarry Smith c = (Mat_SeqBAIJ*)C->data; 182744cd7ae7SLois Curfman McInnes 1828b6490206SBarry Smith c->bs = a->bs; 18299df24120SSatish Balay c->bs2 = a->bs2; 1830b6490206SBarry Smith c->mbs = a->mbs; 1831b6490206SBarry Smith c->nbs = a->nbs; 183235d8aa7fSBarry Smith ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 18332593348eSBarry Smith 1834b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1835b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 1836b6490206SBarry Smith for (i=0; i<mbs; i++) { 18372593348eSBarry Smith c->imax[i] = a->imax[i]; 18382593348eSBarry Smith c->ilen[i] = a->ilen[i]; 18392593348eSBarry Smith } 18402593348eSBarry Smith 18412593348eSBarry Smith /* allocate the matrix space */ 1842c4992f7dSBarry Smith c->singlemalloc = PETSC_TRUE; 18433f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 1844b0a32e0cSBarry Smith ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 18457e67e3f9SSatish Balay c->j = (int*)(c->a + nz*bs2); 1846de6a44a3SBarry Smith c->i = c->j + nz; 1847549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1848b6490206SBarry Smith if (mbs > 0) { 1849549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 18502e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1851549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 18522e8a6d31SBarry Smith } else { 1853549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 18542593348eSBarry Smith } 18552593348eSBarry Smith } 18562593348eSBarry Smith 1857b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 18582593348eSBarry Smith c->sorted = a->sorted; 18592593348eSBarry Smith c->roworiented = a->roworiented; 18602593348eSBarry Smith c->nonew = a->nonew; 18612593348eSBarry Smith 18622593348eSBarry Smith if (a->diag) { 1863b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1864b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 1865b6490206SBarry Smith for (i=0; i<mbs; i++) { 18662593348eSBarry Smith c->diag[i] = a->diag[i]; 18672593348eSBarry Smith } 186898305bb5SBarry Smith } else c->diag = 0; 18692593348eSBarry Smith c->nz = a->nz; 18702593348eSBarry Smith c->maxnz = a->maxnz; 18712593348eSBarry Smith c->solve_work = 0; 18722a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 18737fc0212eSBarry Smith c->mult_work = 0; 1874273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 1875273d9f13SBarry Smith C->assembled = PETSC_TRUE; 18762593348eSBarry Smith *B = C; 1877b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 18783a40ed3dSBarry Smith PetscFunctionReturn(0); 18792593348eSBarry Smith } 18802593348eSBarry Smith 18814a2ae208SSatish Balay #undef __FUNCT__ 18824a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ" 1883b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A) 18842593348eSBarry Smith { 1885b6490206SBarry Smith Mat_SeqBAIJ *a; 18862593348eSBarry Smith Mat B; 1887f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1888b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 188935aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1890a7c10996SSatish Balay int *masked,nmask,tmp,bs2,ishift; 189187828ca2SBarry Smith PetscScalar *aa; 189219bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 18932593348eSBarry Smith 18943a40ed3dSBarry Smith PetscFunctionBegin; 1895b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1896de6a44a3SBarry Smith bs2 = bs*bs; 1897b6490206SBarry Smith 1898d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 189929bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1900b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 19010752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1902552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 19032593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 19042593348eSBarry Smith 1905d64ed03dSBarry Smith if (header[3] < 0) { 190629bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ"); 1907d64ed03dSBarry Smith } 1908d64ed03dSBarry Smith 190929bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 191035aab85fSBarry Smith 191135aab85fSBarry Smith /* 191235aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 191335aab85fSBarry Smith divisible by the blocksize 191435aab85fSBarry Smith */ 1915b6490206SBarry Smith mbs = M/bs; 191635aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 191735aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 191835aab85fSBarry Smith else mbs++; 191935aab85fSBarry Smith if (extra_rows) { 1920b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 192135aab85fSBarry Smith } 1922b6490206SBarry Smith 19232593348eSBarry Smith /* read in row lengths */ 1924b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 19250752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 192635aab85fSBarry Smith for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 19272593348eSBarry Smith 1928b6490206SBarry Smith /* read in column indices */ 1929b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 19300752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 193135aab85fSBarry Smith for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 1932b6490206SBarry Smith 1933b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1934b0a32e0cSBarry Smith ierr = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 1935549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 1936b0a32e0cSBarry Smith ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 1937549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 193835aab85fSBarry Smith masked = mask + mbs; 1939b6490206SBarry Smith rowcount = 0; nzcount = 0; 1940b6490206SBarry Smith for (i=0; i<mbs; i++) { 194135aab85fSBarry Smith nmask = 0; 1942b6490206SBarry Smith for (j=0; j<bs; j++) { 1943b6490206SBarry Smith kmax = rowlengths[rowcount]; 1944b6490206SBarry Smith for (k=0; k<kmax; k++) { 194535aab85fSBarry Smith tmp = jj[nzcount++]/bs; 194635aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1947b6490206SBarry Smith } 1948b6490206SBarry Smith rowcount++; 1949b6490206SBarry Smith } 195035aab85fSBarry Smith browlengths[i] += nmask; 195135aab85fSBarry Smith /* zero out the mask elements we set */ 195235aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1953b6490206SBarry Smith } 1954b6490206SBarry Smith 19552593348eSBarry Smith /* create our matrix */ 1956dd23797bSKris Buschelman ierr = MatCreate(comm,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows,&B); 195778ae41b4SKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 195878ae41b4SKris Buschelman ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 1959b6490206SBarry Smith a = (Mat_SeqBAIJ*)B->data; 19602593348eSBarry Smith 19612593348eSBarry Smith /* set matrix "i" values */ 1962de6a44a3SBarry Smith a->i[0] = 0; 1963b6490206SBarry Smith for (i=1; i<= mbs; i++) { 1964b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1965b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 19662593348eSBarry Smith } 1967b6490206SBarry Smith a->nz = 0; 1968b6490206SBarry Smith for (i=0; i<mbs; i++) a->nz += browlengths[i]; 19692593348eSBarry Smith 1970b6490206SBarry Smith /* read in nonzero values */ 197187828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 19720752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 197335aab85fSBarry Smith for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 1974b6490206SBarry Smith 1975b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1976b6490206SBarry Smith nzcount = 0; jcount = 0; 1977b6490206SBarry Smith for (i=0; i<mbs; i++) { 1978b6490206SBarry Smith nzcountb = nzcount; 197935aab85fSBarry Smith nmask = 0; 1980b6490206SBarry Smith for (j=0; j<bs; j++) { 1981b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1982b6490206SBarry Smith for (k=0; k<kmax; k++) { 198335aab85fSBarry Smith tmp = jj[nzcount++]/bs; 198435aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1985b6490206SBarry Smith } 1986b6490206SBarry Smith } 1987de6a44a3SBarry Smith /* sort the masked values */ 1988433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1989de6a44a3SBarry Smith 1990b6490206SBarry Smith /* set "j" values into matrix */ 1991b6490206SBarry Smith maskcount = 1; 199235aab85fSBarry Smith for (j=0; j<nmask; j++) { 199335aab85fSBarry Smith a->j[jcount++] = masked[j]; 1994de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1995b6490206SBarry Smith } 1996b6490206SBarry Smith /* set "a" values into matrix */ 1997de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1998b6490206SBarry Smith for (j=0; j<bs; j++) { 1999b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 2000b6490206SBarry Smith for (k=0; k<kmax; k++) { 2001de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 2002de6a44a3SBarry Smith block = mask[tmp] - 1; 2003de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 2004de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 2005375fe846SBarry Smith a->a[idx] = (MatScalar)aa[nzcountb++]; 2006b6490206SBarry Smith } 2007b6490206SBarry Smith } 200835aab85fSBarry Smith /* zero out the mask elements we set */ 200935aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 2010b6490206SBarry Smith } 201129bbc08cSBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 2012b6490206SBarry Smith 2013606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2014606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 2015606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 2016606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 2017606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 2018b6490206SBarry Smith 201978ae41b4SKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 202078ae41b4SKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20219c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 202278ae41b4SKris Buschelman 202378ae41b4SKris Buschelman *A = B; 20243a40ed3dSBarry Smith PetscFunctionReturn(0); 20252593348eSBarry Smith } 20262593348eSBarry Smith 20274a2ae208SSatish Balay #undef __FUNCT__ 20284a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ" 2029273d9f13SBarry Smith /*@C 2030273d9f13SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 2031273d9f13SBarry Smith compressed row) format. For good matrix assembly performance the 2032273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 2033273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 2034273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 20352593348eSBarry Smith 2036273d9f13SBarry Smith Collective on MPI_Comm 2037273d9f13SBarry Smith 2038273d9f13SBarry Smith Input Parameters: 2039273d9f13SBarry Smith + comm - MPI communicator, set to PETSC_COMM_SELF 2040273d9f13SBarry Smith . bs - size of block 2041273d9f13SBarry Smith . m - number of rows 2042273d9f13SBarry Smith . n - number of columns 204335d8aa7fSBarry Smith . nz - number of nonzero blocks per block row (same for all rows) 204435d8aa7fSBarry Smith - nnz - array containing the number of nonzero blocks in the various block rows 2045273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 2046273d9f13SBarry Smith 2047273d9f13SBarry Smith Output Parameter: 2048273d9f13SBarry Smith . A - the matrix 2049273d9f13SBarry Smith 2050273d9f13SBarry Smith Options Database Keys: 2051273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 2052273d9f13SBarry Smith block calculations (much slower) 2053273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 2054273d9f13SBarry Smith 2055273d9f13SBarry Smith Level: intermediate 2056273d9f13SBarry Smith 2057273d9f13SBarry Smith Notes: 205835d8aa7fSBarry Smith A nonzero block is any block that as 1 or more nonzeros in it 205935d8aa7fSBarry Smith 2060273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 2061273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 2062273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 2063273d9f13SBarry Smith 2064273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2065273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 2066273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 2067273d9f13SBarry Smith matrices. 2068273d9f13SBarry Smith 2069273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 2070273d9f13SBarry Smith @*/ 2071ca01db9bSBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A) 2072273d9f13SBarry Smith { 2073273d9f13SBarry Smith int ierr; 2074273d9f13SBarry Smith 2075273d9f13SBarry Smith PetscFunctionBegin; 2076273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 2077273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr); 2078273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 2079273d9f13SBarry Smith PetscFunctionReturn(0); 2080273d9f13SBarry Smith } 2081273d9f13SBarry Smith 20824a2ae208SSatish Balay #undef __FUNCT__ 20834a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation" 2084273d9f13SBarry Smith /*@C 2085273d9f13SBarry Smith MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros 2086273d9f13SBarry Smith per row in the matrix. For good matrix assembly performance the 2087273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 2088273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 2089273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 2090273d9f13SBarry Smith 2091273d9f13SBarry Smith Collective on MPI_Comm 2092273d9f13SBarry Smith 2093273d9f13SBarry Smith Input Parameters: 2094273d9f13SBarry Smith + A - the matrix 2095273d9f13SBarry Smith . bs - size of block 2096273d9f13SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 2097273d9f13SBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 2098273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 2099273d9f13SBarry Smith 2100273d9f13SBarry Smith Options Database Keys: 2101273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 2102273d9f13SBarry Smith block calculations (much slower) 2103273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 2104273d9f13SBarry Smith 2105273d9f13SBarry Smith Level: intermediate 2106273d9f13SBarry Smith 2107273d9f13SBarry Smith Notes: 2108273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 2109273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 2110273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 2111273d9f13SBarry Smith 2112273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2113273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 2114273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 2115273d9f13SBarry Smith matrices. 2116273d9f13SBarry Smith 2117273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 2118273d9f13SBarry Smith @*/ 2119ca01db9bSBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) 2120273d9f13SBarry Smith { 2121ca01db9bSBarry Smith int ierr,(*f)(Mat,int,int,const int[]); 2122273d9f13SBarry Smith 2123273d9f13SBarry Smith PetscFunctionBegin; 2124a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 2125a23d5eceSKris Buschelman if (f) { 2126a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 2127273d9f13SBarry Smith } 2128273d9f13SBarry Smith PetscFunctionReturn(0); 2129273d9f13SBarry Smith } 2130