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"); 35377e54ba9SKris Buschelman case MAT_SYMMETRIC: 35477e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3559a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 3569a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 3579a4540c5SBarry Smith case MAT_HERMITIAN: 3589a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 3599a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 3609a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 36177e54ba9SKris Buschelman break; 362aa275fccSKris Buschelman default: 36329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 3642d61bbb3SSatish Balay } 3652d61bbb3SSatish Balay PetscFunctionReturn(0); 3662d61bbb3SSatish Balay } 3672d61bbb3SSatish Balay 3684a2ae208SSatish Balay #undef __FUNCT__ 3694a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ" 37087828ca2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 3712d61bbb3SSatish Balay { 3722d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 37382502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr; 3743f1db9ecSBarry Smith MatScalar *aa,*aa_i; 37587828ca2SBarry Smith PetscScalar *v_i; 3762d61bbb3SSatish Balay 3772d61bbb3SSatish Balay PetscFunctionBegin; 3782d61bbb3SSatish Balay bs = a->bs; 3792d61bbb3SSatish Balay ai = a->i; 3802d61bbb3SSatish Balay aj = a->j; 3812d61bbb3SSatish Balay aa = a->a; 3822d61bbb3SSatish Balay bs2 = a->bs2; 3832d61bbb3SSatish Balay 384a45adfd6SMatthew Knepley if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d out of range", row); 3852d61bbb3SSatish Balay 3862d61bbb3SSatish Balay bn = row/bs; /* Block number */ 3872d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 3882d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 3892d61bbb3SSatish Balay *nz = bs*M; 3902d61bbb3SSatish Balay 3912d61bbb3SSatish Balay if (v) { 3922d61bbb3SSatish Balay *v = 0; 3932d61bbb3SSatish Balay if (*nz) { 39487828ca2SBarry Smith ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr); 3952d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 3962d61bbb3SSatish Balay v_i = *v + i*bs; 3972d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 3982d61bbb3SSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 3992d61bbb3SSatish Balay } 4002d61bbb3SSatish Balay } 4012d61bbb3SSatish Balay } 4022d61bbb3SSatish Balay 4032d61bbb3SSatish Balay if (idx) { 4042d61bbb3SSatish Balay *idx = 0; 4052d61bbb3SSatish Balay if (*nz) { 406b0a32e0cSBarry Smith ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr); 4072d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 4082d61bbb3SSatish Balay idx_i = *idx + i*bs; 4092d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 4102d61bbb3SSatish Balay for (j=0; j<bs; j++) {idx_i[j] = itmp++;} 4112d61bbb3SSatish Balay } 4122d61bbb3SSatish Balay } 4132d61bbb3SSatish Balay } 4142d61bbb3SSatish Balay PetscFunctionReturn(0); 4152d61bbb3SSatish Balay } 4162d61bbb3SSatish Balay 4174a2ae208SSatish Balay #undef __FUNCT__ 4184a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ" 41987828ca2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v) 4202d61bbb3SSatish Balay { 421606d414cSSatish Balay int ierr; 422606d414cSSatish Balay 4232d61bbb3SSatish Balay PetscFunctionBegin; 424606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 425606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 4262d61bbb3SSatish Balay PetscFunctionReturn(0); 4272d61bbb3SSatish Balay } 4282d61bbb3SSatish Balay 4294a2ae208SSatish Balay #undef __FUNCT__ 4304a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ" 4312d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 4322d61bbb3SSatish Balay { 4332d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 4342d61bbb3SSatish Balay Mat C; 4352d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 436273d9f13SBarry Smith int *rows,*cols,bs2=a->bs2; 43787828ca2SBarry Smith PetscScalar *array; 4382d61bbb3SSatish Balay 4392d61bbb3SSatish Balay PetscFunctionBegin; 44029bbc08cSBarry Smith if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place"); 441b0a32e0cSBarry Smith ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr); 442549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 4432d61bbb3SSatish Balay 444375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 44587828ca2SBarry Smith ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr); 44687828ca2SBarry Smith for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i]; 447375fe846SBarry Smith #else 4483eda8832SBarry Smith array = a->a; 449375fe846SBarry Smith #endif 450375fe846SBarry Smith 4512d61bbb3SSatish Balay for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 452273d9f13SBarry Smith ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 453606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 454b0a32e0cSBarry Smith ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr); 4552d61bbb3SSatish Balay cols = rows + bs; 4562d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4572d61bbb3SSatish Balay cols[0] = i*bs; 4582d61bbb3SSatish Balay for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 4592d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 4602d61bbb3SSatish Balay for (j=0; j<len; j++) { 4612d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 4622d61bbb3SSatish Balay for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 4632d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 4642d61bbb3SSatish Balay array += bs2; 4652d61bbb3SSatish Balay } 4662d61bbb3SSatish Balay } 467606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 468375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 469375fe846SBarry Smith ierr = PetscFree(array); 470375fe846SBarry Smith #endif 4712d61bbb3SSatish Balay 4722d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4732d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4742d61bbb3SSatish Balay 475c4992f7dSBarry Smith if (B) { 4762d61bbb3SSatish Balay *B = C; 4772d61bbb3SSatish Balay } else { 478273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 4792d61bbb3SSatish Balay } 4802d61bbb3SSatish Balay PetscFunctionReturn(0); 4812d61bbb3SSatish Balay } 4822d61bbb3SSatish Balay 4834a2ae208SSatish Balay #undef __FUNCT__ 4844a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary" 485b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer) 4862593348eSBarry Smith { 487b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 4889df24120SSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 48987828ca2SBarry Smith PetscScalar *aa; 490ce6f0cecSBarry Smith FILE *file; 4912593348eSBarry Smith 4923a40ed3dSBarry Smith PetscFunctionBegin; 493b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 494b0a32e0cSBarry Smith ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 495552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 4963b2fbd54SBarry Smith 497273d9f13SBarry Smith col_lens[1] = A->m; 498273d9f13SBarry Smith col_lens[2] = A->n; 4997e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 5002593348eSBarry Smith 5012593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 502b6490206SBarry Smith count = 0; 503b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 504b6490206SBarry Smith for (j=0; j<bs; j++) { 505b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 506b6490206SBarry Smith } 5072593348eSBarry Smith } 508273d9f13SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 509606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 5102593348eSBarry Smith 5112593348eSBarry Smith /* store column indices (zero start index) */ 512b0a32e0cSBarry Smith ierr = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 513b6490206SBarry Smith count = 0; 514b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 515b6490206SBarry Smith for (j=0; j<bs; j++) { 516b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 517b6490206SBarry Smith for (l=0; l<bs; l++) { 518b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 5192593348eSBarry Smith } 5202593348eSBarry Smith } 521b6490206SBarry Smith } 522b6490206SBarry Smith } 5230752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 524606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 5252593348eSBarry Smith 5262593348eSBarry Smith /* store nonzero values */ 52787828ca2SBarry Smith ierr = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 528b6490206SBarry Smith count = 0; 529b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 530b6490206SBarry Smith for (j=0; j<bs; j++) { 531b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 532b6490206SBarry Smith for (l=0; l<bs; l++) { 5337e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 534b6490206SBarry Smith } 535b6490206SBarry Smith } 536b6490206SBarry Smith } 537b6490206SBarry Smith } 5380752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 539606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 540ce6f0cecSBarry Smith 541b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 542ce6f0cecSBarry Smith if (file) { 543ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 544ce6f0cecSBarry Smith } 5453a40ed3dSBarry Smith PetscFunctionReturn(0); 5462593348eSBarry Smith } 5472593348eSBarry Smith 5484a2ae208SSatish Balay #undef __FUNCT__ 5494a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII" 550b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer) 5512593348eSBarry Smith { 552b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 553fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 554f3ef73ceSBarry Smith PetscViewerFormat format; 5552593348eSBarry Smith 5563a40ed3dSBarry Smith PetscFunctionBegin; 557b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 558456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 559b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 560fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 561bcd9e38bSBarry Smith Mat aij; 562bcd9e38bSBarry Smith ierr = MatConvert(A,MATSEQAIJ,&aij);CHKERRQ(ierr); 563bcd9e38bSBarry Smith ierr = MatView(aij,viewer);CHKERRQ(ierr); 564bcd9e38bSBarry Smith ierr = MatDestroy(aij);CHKERRQ(ierr); 56504929863SHong Zhang } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 56604929863SHong Zhang PetscFunctionReturn(0); 567fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 568b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 56944cd7ae7SLois Curfman McInnes for (i=0; i<a->mbs; i++) { 57044cd7ae7SLois Curfman McInnes for (j=0; j<bs; j++) { 571b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 57244cd7ae7SLois Curfman McInnes for (k=a->i[i]; k<a->i[i+1]; k++) { 57344cd7ae7SLois Curfman McInnes for (l=0; l<bs; l++) { 574aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 5750e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 57662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %gi) ",bs*a->j[k]+l, 5770e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5780e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 57962b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %gi) ",bs*a->j[k]+l, 5800e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5810e6d2581SBarry Smith } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 58262b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 5830ef38995SBarry Smith } 58444cd7ae7SLois Curfman McInnes #else 5850ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 58662b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 5870ef38995SBarry Smith } 58844cd7ae7SLois Curfman McInnes #endif 58944cd7ae7SLois Curfman McInnes } 59044cd7ae7SLois Curfman McInnes } 591b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 59244cd7ae7SLois Curfman McInnes } 59344cd7ae7SLois Curfman McInnes } 594b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 5950ef38995SBarry Smith } else { 596b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 597b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 598b6490206SBarry Smith for (j=0; j<bs; j++) { 599b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 600b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 601b6490206SBarry Smith for (l=0; l<bs; l++) { 602aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 6030e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 60462b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l, 6050e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 6060e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 60762b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l, 6080e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 6090ef38995SBarry Smith } else { 61062b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 61188685aaeSLois Curfman McInnes } 61288685aaeSLois Curfman McInnes #else 61362b941d6SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 61488685aaeSLois Curfman McInnes #endif 6152593348eSBarry Smith } 6162593348eSBarry Smith } 617b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6182593348eSBarry Smith } 6192593348eSBarry Smith } 620b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 621b6490206SBarry Smith } 622b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 6233a40ed3dSBarry Smith PetscFunctionReturn(0); 6242593348eSBarry Smith } 6252593348eSBarry Smith 6264a2ae208SSatish Balay #undef __FUNCT__ 6274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom" 628b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 6293270192aSSatish Balay { 63077ed5343SBarry Smith Mat A = (Mat) Aa; 6313270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 632b65db4caSBarry Smith int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2; 6330e6d2581SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 6343f1db9ecSBarry Smith MatScalar *aa; 635b0a32e0cSBarry Smith PetscViewer viewer; 6363270192aSSatish Balay 6373a40ed3dSBarry Smith PetscFunctionBegin; 6383270192aSSatish Balay 639b65db4caSBarry Smith /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/ 64077ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 64177ed5343SBarry Smith 642b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 64377ed5343SBarry Smith 6443270192aSSatish Balay /* loop over matrix elements drawing boxes */ 645b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 6463270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 6473270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 648273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 6493270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 6503270192aSSatish Balay aa = a->a + j*bs2; 6513270192aSSatish Balay for (k=0; k<bs; k++) { 6523270192aSSatish Balay for (l=0; l<bs; l++) { 6530e6d2581SBarry Smith if (PetscRealPart(*aa++) >= 0.) continue; 654b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 6553270192aSSatish Balay } 6563270192aSSatish Balay } 6573270192aSSatish Balay } 6583270192aSSatish Balay } 659b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 6603270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 6613270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 662273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 6633270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 6643270192aSSatish Balay aa = a->a + j*bs2; 6653270192aSSatish Balay for (k=0; k<bs; k++) { 6663270192aSSatish Balay for (l=0; l<bs; l++) { 6670e6d2581SBarry Smith if (PetscRealPart(*aa++) != 0.) continue; 668b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 6693270192aSSatish Balay } 6703270192aSSatish Balay } 6713270192aSSatish Balay } 6723270192aSSatish Balay } 6733270192aSSatish Balay 674b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 6753270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 6763270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 677273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 6783270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 6793270192aSSatish Balay aa = a->a + j*bs2; 6803270192aSSatish Balay for (k=0; k<bs; k++) { 6813270192aSSatish Balay for (l=0; l<bs; l++) { 6820e6d2581SBarry Smith if (PetscRealPart(*aa++) <= 0.) continue; 683b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 6843270192aSSatish Balay } 6853270192aSSatish Balay } 6863270192aSSatish Balay } 6873270192aSSatish Balay } 68877ed5343SBarry Smith PetscFunctionReturn(0); 68977ed5343SBarry Smith } 6903270192aSSatish Balay 6914a2ae208SSatish Balay #undef __FUNCT__ 6924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw" 693b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer) 69477ed5343SBarry Smith { 6957dce120fSSatish Balay int ierr; 6960e6d2581SBarry Smith PetscReal xl,yl,xr,yr,w,h; 697b0a32e0cSBarry Smith PetscDraw draw; 69877ed5343SBarry Smith PetscTruth isnull; 6993270192aSSatish Balay 70077ed5343SBarry Smith PetscFunctionBegin; 70177ed5343SBarry Smith 702b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 703b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 70477ed5343SBarry Smith 70577ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 706273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 70777ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 708b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 709b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 71077ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 7113a40ed3dSBarry Smith PetscFunctionReturn(0); 7123270192aSSatish Balay } 7133270192aSSatish Balay 7144a2ae208SSatish Balay #undef __FUNCT__ 7154a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ" 716b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer) 7172593348eSBarry Smith { 71819bcc07fSBarry Smith int ierr; 719*a5e6ed63SBarry Smith PetscTruth isascii,isbinary,isdraw; 7202593348eSBarry Smith 7213a40ed3dSBarry Smith PetscFunctionBegin; 722b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 723fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 724fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 725*a5e6ed63SBarry Smith if (isascii){ 7263a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 7270f5bd95cSBarry Smith } else if (isbinary) { 7283a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 7290f5bd95cSBarry Smith } else if (isdraw) { 7303a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 7315cd90555SBarry Smith } else { 732*a5e6ed63SBarry Smith Mat B; 733*a5e6ed63SBarry Smith ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr); 734*a5e6ed63SBarry Smith ierr = MatView(B,viewer);CHKERRQ(ierr); 735*a5e6ed63SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 7362593348eSBarry Smith } 7373a40ed3dSBarry Smith PetscFunctionReturn(0); 7382593348eSBarry Smith } 739b6490206SBarry Smith 740cd0e1443SSatish Balay 7414a2ae208SSatish Balay #undef __FUNCT__ 7424a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ" 743f15d580aSBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[]) 744cd0e1443SSatish Balay { 745cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 7462d61bbb3SSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 7472d61bbb3SSatish Balay int *ai = a->i,*ailen = a->ilen; 7482d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 7493f1db9ecSBarry Smith MatScalar *ap,*aa = a->a,zero = 0.0; 750cd0e1443SSatish Balay 7513a40ed3dSBarry Smith PetscFunctionBegin; 7522d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 753cd0e1443SSatish Balay row = im[k]; brow = row/bs; 75429bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 755a45adfd6SMatthew Knepley if (row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %d too large", row); 7562d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 7572c3acbe9SBarry Smith nrow = ailen[brow]; 7582d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 75929bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 760a45adfd6SMatthew Knepley if (in[l] >= A->n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column %d too large", in[l]); 7612d61bbb3SSatish Balay col = in[l] ; 7622d61bbb3SSatish Balay bcol = col/bs; 7632d61bbb3SSatish Balay cidx = col%bs; 7642d61bbb3SSatish Balay ridx = row%bs; 7652d61bbb3SSatish Balay high = nrow; 7662d61bbb3SSatish Balay low = 0; /* assume unsorted */ 7672d61bbb3SSatish Balay while (high-low > 5) { 768cd0e1443SSatish Balay t = (low+high)/2; 769cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 770cd0e1443SSatish Balay else low = t; 771cd0e1443SSatish Balay } 772cd0e1443SSatish Balay for (i=low; i<high; i++) { 773cd0e1443SSatish Balay if (rp[i] > bcol) break; 774cd0e1443SSatish Balay if (rp[i] == bcol) { 7752d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 7762d61bbb3SSatish Balay goto finished; 777cd0e1443SSatish Balay } 778cd0e1443SSatish Balay } 7792d61bbb3SSatish Balay *v++ = zero; 7802d61bbb3SSatish Balay finished:; 781cd0e1443SSatish Balay } 782cd0e1443SSatish Balay } 7833a40ed3dSBarry Smith PetscFunctionReturn(0); 784cd0e1443SSatish Balay } 785cd0e1443SSatish Balay 78695d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 7874a2ae208SSatish Balay #undef __FUNCT__ 7884a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 789f15d580aSBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode addv) 79095d5f7c2SBarry Smith { 791563b5814SBarry Smith Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data; 792563b5814SBarry Smith int ierr,i,N = m*n*b->bs2; 793563b5814SBarry Smith MatScalar *vsingle; 79495d5f7c2SBarry Smith 79595d5f7c2SBarry Smith PetscFunctionBegin; 796563b5814SBarry Smith if (N > b->setvalueslen) { 797563b5814SBarry Smith if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);} 798b0a32e0cSBarry Smith ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr); 799563b5814SBarry Smith b->setvalueslen = N; 800563b5814SBarry Smith } 801563b5814SBarry Smith vsingle = b->setvaluescopy; 80295d5f7c2SBarry Smith for (i=0; i<N; i++) { 80395d5f7c2SBarry Smith vsingle[i] = v[i]; 80495d5f7c2SBarry Smith } 80595d5f7c2SBarry Smith ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr); 80695d5f7c2SBarry Smith PetscFunctionReturn(0); 80795d5f7c2SBarry Smith } 80895d5f7c2SBarry Smith #endif 80995d5f7c2SBarry Smith 8102d61bbb3SSatish Balay 8114a2ae208SSatish Balay #undef __FUNCT__ 8124a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 813f15d580aSBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,const int im[],int n,const int in[],const MatScalar v[],InsertMode is) 81492c4ed94SBarry Smith { 81592c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 8168a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 817273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 818549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 819273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 820f15d580aSBarry Smith const MatScalar *value = v; 821f15d580aSBarry Smith MatScalar *ap,*aa = a->a,*bap; 82292c4ed94SBarry Smith 8233a40ed3dSBarry Smith PetscFunctionBegin; 8240e324ae4SSatish Balay if (roworiented) { 8250e324ae4SSatish Balay stepval = (n-1)*bs; 8260e324ae4SSatish Balay } else { 8270e324ae4SSatish Balay stepval = (m-1)*bs; 8280e324ae4SSatish Balay } 82992c4ed94SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 83092c4ed94SBarry Smith row = im[k]; 8315ef9f2a5SBarry Smith if (row < 0) continue; 832aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 833590ac198SBarry Smith if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1); 83492c4ed94SBarry Smith #endif 83592c4ed94SBarry Smith rp = aj + ai[row]; 83692c4ed94SBarry Smith ap = aa + bs2*ai[row]; 83792c4ed94SBarry Smith rmax = imax[row]; 83892c4ed94SBarry Smith nrow = ailen[row]; 83992c4ed94SBarry Smith low = 0; 84092c4ed94SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 8415ef9f2a5SBarry Smith if (in[l] < 0) continue; 842aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 843590ac198SBarry Smith if (in[l] >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],a->nbs-1); 84492c4ed94SBarry Smith #endif 84592c4ed94SBarry Smith col = in[l]; 84692c4ed94SBarry Smith if (roworiented) { 8470e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 8480e324ae4SSatish Balay } else { 8490e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 85092c4ed94SBarry Smith } 85192c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 85292c4ed94SBarry Smith while (high-low > 7) { 85392c4ed94SBarry Smith t = (low+high)/2; 85492c4ed94SBarry Smith if (rp[t] > col) high = t; 85592c4ed94SBarry Smith else low = t; 85692c4ed94SBarry Smith } 85792c4ed94SBarry Smith for (i=low; i<high; i++) { 85892c4ed94SBarry Smith if (rp[i] > col) break; 85992c4ed94SBarry Smith if (rp[i] == col) { 8608a84c255SSatish Balay bap = ap + bs2*i; 8610e324ae4SSatish Balay if (roworiented) { 8628a84c255SSatish Balay if (is == ADD_VALUES) { 863dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 864dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 8658a84c255SSatish Balay bap[jj] += *value++; 866dd9472c6SBarry Smith } 867dd9472c6SBarry Smith } 8680e324ae4SSatish Balay } else { 869dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 870dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 8710e324ae4SSatish Balay bap[jj] = *value++; 8728a84c255SSatish Balay } 873dd9472c6SBarry Smith } 874dd9472c6SBarry Smith } 8750e324ae4SSatish Balay } else { 8760e324ae4SSatish Balay if (is == ADD_VALUES) { 877dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 878dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 8790e324ae4SSatish Balay *bap++ += *value++; 880dd9472c6SBarry Smith } 881dd9472c6SBarry Smith } 8820e324ae4SSatish Balay } else { 883dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 884dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 8850e324ae4SSatish Balay *bap++ = *value++; 8860e324ae4SSatish Balay } 8878a84c255SSatish Balay } 888dd9472c6SBarry Smith } 889dd9472c6SBarry Smith } 890f1241b54SBarry Smith goto noinsert2; 89192c4ed94SBarry Smith } 89292c4ed94SBarry Smith } 89389280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 894a45adfd6SMatthew Knepley else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 89592c4ed94SBarry Smith if (nrow >= rmax) { 89692c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 89792c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 8983f1db9ecSBarry Smith MatScalar *new_a; 89992c4ed94SBarry Smith 900a45adfd6SMatthew Knepley if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 90189280ab3SLois Curfman McInnes 90292c4ed94SBarry Smith /* malloc new storage space */ 9033f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 904b0a32e0cSBarry Smith ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 90592c4ed94SBarry Smith new_j = (int*)(new_a + bs2*new_nz); 90692c4ed94SBarry Smith new_i = new_j + new_nz; 90792c4ed94SBarry Smith 90892c4ed94SBarry Smith /* copy over old data into new slots */ 90992c4ed94SBarry Smith for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 91092c4ed94SBarry Smith for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 911549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 91292c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 913549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 914549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 915549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 916549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 91792c4ed94SBarry Smith /* free up old matrix storage */ 918606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 919606d414cSSatish Balay if (!a->singlemalloc) { 920606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 921606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 922606d414cSSatish Balay } 92392c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 924c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 92592c4ed94SBarry Smith 92692c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 92792c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 928b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 92992c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 93092c4ed94SBarry Smith a->reallocs++; 93192c4ed94SBarry Smith a->nz++; 93292c4ed94SBarry Smith } 93392c4ed94SBarry Smith N = nrow++ - 1; 93492c4ed94SBarry Smith /* shift up all the later entries in this row */ 93592c4ed94SBarry Smith for (ii=N; ii>=i; ii--) { 93692c4ed94SBarry Smith rp[ii+1] = rp[ii]; 937549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 93892c4ed94SBarry Smith } 939549d3d68SSatish Balay if (N >= i) { 940549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 941549d3d68SSatish Balay } 94292c4ed94SBarry Smith rp[i] = col; 9438a84c255SSatish Balay bap = ap + bs2*i; 9440e324ae4SSatish Balay if (roworiented) { 945dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 946dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 9470e324ae4SSatish Balay bap[jj] = *value++; 948dd9472c6SBarry Smith } 949dd9472c6SBarry Smith } 9500e324ae4SSatish Balay } else { 951dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 952dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 9530e324ae4SSatish Balay *bap++ = *value++; 9540e324ae4SSatish Balay } 955dd9472c6SBarry Smith } 956dd9472c6SBarry Smith } 957f1241b54SBarry Smith noinsert2:; 95892c4ed94SBarry Smith low = i; 95992c4ed94SBarry Smith } 96092c4ed94SBarry Smith ailen[row] = nrow; 96192c4ed94SBarry Smith } 9623a40ed3dSBarry Smith PetscFunctionReturn(0); 96392c4ed94SBarry Smith } 96492c4ed94SBarry Smith 9654a2ae208SSatish Balay #undef __FUNCT__ 9664a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ" 9678f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 968584200bdSSatish Balay { 969584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 970584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 971273d9f13SBarry Smith int m = A->m,*ip,N,*ailen = a->ilen; 972549d3d68SSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 9733f1db9ecSBarry Smith MatScalar *aa = a->a,*ap; 974584200bdSSatish Balay 9753a40ed3dSBarry Smith PetscFunctionBegin; 9763a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 977584200bdSSatish Balay 97843ee02c3SBarry Smith if (m) rmax = ailen[0]; 979584200bdSSatish Balay for (i=1; i<mbs; i++) { 980584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 981584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 982d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 983584200bdSSatish Balay if (fshift) { 984a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 985584200bdSSatish Balay N = ailen[i]; 986584200bdSSatish Balay for (j=0; j<N; j++) { 987584200bdSSatish Balay ip[j-fshift] = ip[j]; 988549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 989584200bdSSatish Balay } 990584200bdSSatish Balay } 991584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 992584200bdSSatish Balay } 993584200bdSSatish Balay if (mbs) { 994584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 995584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 996584200bdSSatish Balay } 997584200bdSSatish Balay /* reset ilen and imax for each row */ 998584200bdSSatish Balay for (i=0; i<mbs; i++) { 999584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 1000584200bdSSatish Balay } 1001a7c10996SSatish Balay a->nz = ai[mbs]; 1002584200bdSSatish Balay 1003584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 1004584200bdSSatish Balay if (fshift && a->diag) { 1005606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 1006b424e231SHong Zhang PetscLogObjectMemory(A,-(mbs+1)*sizeof(int)); 1007584200bdSSatish Balay a->diag = 0; 1008584200bdSSatish Balay } 1009bba1ac68SSatish 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); 1010bba1ac68SSatish Balay PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs); 1011b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 1012e2f3b5e9SSatish Balay a->reallocs = 0; 10130e6d2581SBarry Smith A->info.nz_unneeded = (PetscReal)fshift*bs2; 1014cf4441caSHong Zhang 10153a40ed3dSBarry Smith PetscFunctionReturn(0); 1016584200bdSSatish Balay } 1017584200bdSSatish Balay 10185a838052SSatish Balay 1019bea157c4SSatish Balay 1020bea157c4SSatish Balay /* 1021bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 1022bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 1023bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 1024bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 1025bea157c4SSatish Balay */ 10264a2ae208SSatish Balay #undef __FUNCT__ 10274a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks" 1028bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 1029d9b7c43dSSatish Balay { 1030bea157c4SSatish Balay int i,j,k,row; 1031bea157c4SSatish Balay PetscTruth flg; 10323a40ed3dSBarry Smith 1033433994e6SBarry Smith PetscFunctionBegin; 1034bea157c4SSatish Balay for (i=0,j=0; i<n; j++) { 1035bea157c4SSatish Balay row = idx[i]; 1036bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 1037bea157c4SSatish Balay sizes[j] = 1; 1038bea157c4SSatish Balay i++; 1039e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 1040bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 1041bea157c4SSatish Balay i++; 1042bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 1043bea157c4SSatish Balay flg = PETSC_TRUE; 1044bea157c4SSatish Balay for (k=1; k<bs; k++) { 1045bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 1046bea157c4SSatish Balay flg = PETSC_FALSE; 1047bea157c4SSatish Balay break; 1048d9b7c43dSSatish Balay } 1049bea157c4SSatish Balay } 1050bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 1051bea157c4SSatish Balay sizes[j] = bs; 1052bea157c4SSatish Balay i+= bs; 1053bea157c4SSatish Balay } else { 1054bea157c4SSatish Balay sizes[j] = 1; 1055bea157c4SSatish Balay i++; 1056bea157c4SSatish Balay } 1057bea157c4SSatish Balay } 1058bea157c4SSatish Balay } 1059bea157c4SSatish Balay *bs_max = j; 10603a40ed3dSBarry Smith PetscFunctionReturn(0); 1061d9b7c43dSSatish Balay } 1062d9b7c43dSSatish Balay 10634a2ae208SSatish Balay #undef __FUNCT__ 10644a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ" 1065268466fbSBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,const PetscScalar *diag) 1066d9b7c43dSSatish Balay { 1067d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 1068b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 1069bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 107087828ca2SBarry Smith PetscScalar zero = 0.0; 10713f1db9ecSBarry Smith MatScalar *aa; 1072d9b7c43dSSatish Balay 10733a40ed3dSBarry Smith PetscFunctionBegin; 1074d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 1075b9b97703SBarry Smith ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr); 1076d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 1077d9b7c43dSSatish Balay 1078bea157c4SSatish Balay /* allocate memory for rows,sizes */ 1079b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 1080bea157c4SSatish Balay sizes = rows + is_n; 1081bea157c4SSatish Balay 1082563b5814SBarry Smith /* copy IS values to rows, and sort them */ 1083bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 1084bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 1085dffd3267SBarry Smith if (baij->keepzeroedrows) { 1086dffd3267SBarry Smith for (i=0; i<is_n; i++) { sizes[i] = 1; } 1087dffd3267SBarry Smith bs_max = is_n; 1088dffd3267SBarry Smith } else { 1089bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 1090dffd3267SBarry Smith } 1091b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 1092bea157c4SSatish Balay 1093bea157c4SSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 1094bea157c4SSatish Balay row = rows[j]; 1095273d9f13SBarry Smith if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 1096bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 1097bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 1098dffd3267SBarry Smith if (sizes[i] == bs && !baij->keepzeroedrows) { 1099bea157c4SSatish Balay if (diag) { 1100bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 1101bea157c4SSatish Balay baij->ilen[row/bs] = 1; 1102bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 1103bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 1104a07cd24cSSatish Balay } 1105563b5814SBarry Smith /* Now insert all the diagonal values for this bs */ 1106bea157c4SSatish Balay for (k=0; k<bs; k++) { 1107bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 1108bea157c4SSatish Balay } 1109bea157c4SSatish Balay } else { /* (!diag) */ 1110bea157c4SSatish Balay baij->ilen[row/bs] = 0; 1111bea157c4SSatish Balay } /* end (!diag) */ 1112bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 1113aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 111429bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 1115bea157c4SSatish Balay #endif 1116bea157c4SSatish Balay for (k=0; k<count; k++) { 1117d9b7c43dSSatish Balay aa[0] = zero; 1118d9b7c43dSSatish Balay aa += bs; 1119d9b7c43dSSatish Balay } 1120d9b7c43dSSatish Balay if (diag) { 1121f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 1122d9b7c43dSSatish Balay } 1123d9b7c43dSSatish Balay } 1124bea157c4SSatish Balay } 1125bea157c4SSatish Balay 1126606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 11279a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11283a40ed3dSBarry Smith PetscFunctionReturn(0); 1129d9b7c43dSSatish Balay } 11301c351548SSatish Balay 11314a2ae208SSatish Balay #undef __FUNCT__ 11324a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ" 1133f15d580aSBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is) 11342d61bbb3SSatish Balay { 11352d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 11362d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 1137273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 11382d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 1139549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 1140273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 11413f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 11422d61bbb3SSatish Balay 11432d61bbb3SSatish Balay PetscFunctionBegin; 11442d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 11452d61bbb3SSatish Balay row = im[k]; brow = row/bs; 11465ef9f2a5SBarry Smith if (row < 0) continue; 1147aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1148590ac198SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1); 11492d61bbb3SSatish Balay #endif 11502d61bbb3SSatish Balay rp = aj + ai[brow]; 11512d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 11522d61bbb3SSatish Balay rmax = imax[brow]; 11532d61bbb3SSatish Balay nrow = ailen[brow]; 11542d61bbb3SSatish Balay low = 0; 11552d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 11565ef9f2a5SBarry Smith if (in[l] < 0) continue; 1157aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1158590ac198SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1); 11592d61bbb3SSatish Balay #endif 11602d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 11612d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 11622d61bbb3SSatish Balay if (roworiented) { 11635ef9f2a5SBarry Smith value = v[l + k*n]; 11642d61bbb3SSatish Balay } else { 11652d61bbb3SSatish Balay value = v[k + l*m]; 11662d61bbb3SSatish Balay } 11672d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 11682d61bbb3SSatish Balay while (high-low > 7) { 11692d61bbb3SSatish Balay t = (low+high)/2; 11702d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 11712d61bbb3SSatish Balay else low = t; 11722d61bbb3SSatish Balay } 11732d61bbb3SSatish Balay for (i=low; i<high; i++) { 11742d61bbb3SSatish Balay if (rp[i] > bcol) break; 11752d61bbb3SSatish Balay if (rp[i] == bcol) { 11762d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 11772d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 11782d61bbb3SSatish Balay else *bap = value; 11792d61bbb3SSatish Balay goto noinsert1; 11802d61bbb3SSatish Balay } 11812d61bbb3SSatish Balay } 11822d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 1183a45adfd6SMatthew Knepley else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 11842d61bbb3SSatish Balay if (nrow >= rmax) { 11852d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 11862d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 11873f1db9ecSBarry Smith MatScalar *new_a; 11882d61bbb3SSatish Balay 1189a45adfd6SMatthew Knepley if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col); 11902d61bbb3SSatish Balay 11912d61bbb3SSatish Balay /* Malloc new storage space */ 11923f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 1193b0a32e0cSBarry Smith ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 11942d61bbb3SSatish Balay new_j = (int*)(new_a + bs2*new_nz); 11952d61bbb3SSatish Balay new_i = new_j + new_nz; 11962d61bbb3SSatish Balay 11972d61bbb3SSatish Balay /* copy over old data into new slots */ 11982d61bbb3SSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 11992d61bbb3SSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1200549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 12012d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1202549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1203549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1204549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1205549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 12062d61bbb3SSatish Balay /* free up old matrix storage */ 1207606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1208606d414cSSatish Balay if (!a->singlemalloc) { 1209606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1210606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1211606d414cSSatish Balay } 12122d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1213c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 12142d61bbb3SSatish Balay 12152d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 12162d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 1217b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 12182d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 12192d61bbb3SSatish Balay a->reallocs++; 12202d61bbb3SSatish Balay a->nz++; 12212d61bbb3SSatish Balay } 12222d61bbb3SSatish Balay N = nrow++ - 1; 12232d61bbb3SSatish Balay /* shift up all the later entries in this row */ 12242d61bbb3SSatish Balay for (ii=N; ii>=i; ii--) { 12252d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1226549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 12272d61bbb3SSatish Balay } 1228549d3d68SSatish Balay if (N>=i) { 1229549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1230549d3d68SSatish Balay } 12312d61bbb3SSatish Balay rp[i] = bcol; 12322d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 12332d61bbb3SSatish Balay noinsert1:; 12342d61bbb3SSatish Balay low = i; 12352d61bbb3SSatish Balay } 12362d61bbb3SSatish Balay ailen[brow] = nrow; 12372d61bbb3SSatish Balay } 12382d61bbb3SSatish Balay PetscFunctionReturn(0); 12392d61bbb3SSatish Balay } 12402d61bbb3SSatish Balay 12412d61bbb3SSatish Balay 12424a2ae208SSatish Balay #undef __FUNCT__ 12434a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ" 1244b380c88cSHong Zhang int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatFactorInfo *info) 12452d61bbb3SSatish Balay { 12462d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 12472d61bbb3SSatish Balay Mat outA; 12482d61bbb3SSatish Balay int ierr; 1249667159a5SBarry Smith PetscTruth row_identity,col_identity; 12502d61bbb3SSatish Balay 12512d61bbb3SSatish Balay PetscFunctionBegin; 1252d3d32019SBarry Smith if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 1253667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1254667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1255667159a5SBarry Smith if (!row_identity || !col_identity) { 125629bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ILU"); 1257667159a5SBarry Smith } 12582d61bbb3SSatish Balay 12592d61bbb3SSatish Balay outA = inA; 12602d61bbb3SSatish Balay inA->factor = FACTOR_LU; 12612d61bbb3SSatish Balay 12622d61bbb3SSatish Balay if (!a->diag) { 1263c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 12642d61bbb3SSatish Balay } 1265cf242676SKris Buschelman 1266c38d4ed2SBarry Smith a->row = row; 1267c38d4ed2SBarry Smith a->col = col; 1268c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1269c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1270c38d4ed2SBarry Smith 1271c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 12724c49b128SBarry Smith ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); 1273b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 1274c38d4ed2SBarry Smith 1275cf242676SKris Buschelman /* 1276cf242676SKris Buschelman Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 1277cf242676SKris Buschelman for ILU(0) factorization with natural ordering 1278cf242676SKris Buschelman */ 1279cf242676SKris Buschelman if (a->bs < 8) { 1280cf242676SKris Buschelman ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr); 1281cf242676SKris Buschelman } else { 1282c38d4ed2SBarry Smith if (!a->solve_work) { 128387828ca2SBarry Smith ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); 128487828ca2SBarry Smith PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar)); 1285c38d4ed2SBarry Smith } 12862d61bbb3SSatish Balay } 12872d61bbb3SSatish Balay 1288667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1289667159a5SBarry Smith 12902d61bbb3SSatish Balay PetscFunctionReturn(0); 12912d61bbb3SSatish Balay } 12924a2ae208SSatish Balay #undef __FUNCT__ 12934a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ" 1294ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1295ba4ca20aSSatish Balay { 1296c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1297ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1298d132466eSBarry Smith int ierr; 1299ba4ca20aSSatish Balay 13003a40ed3dSBarry Smith PetscFunctionBegin; 1301c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1302d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1303d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 13043a40ed3dSBarry Smith PetscFunctionReturn(0); 1305ba4ca20aSSatish Balay } 1306d9b7c43dSSatish Balay 1307fb2e594dSBarry Smith EXTERN_C_BEGIN 13084a2ae208SSatish Balay #undef __FUNCT__ 13094a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 131027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 131127a8da17SBarry Smith { 131227a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 131314db4f74SSatish Balay int i,nz,nbs; 131427a8da17SBarry Smith 131527a8da17SBarry Smith PetscFunctionBegin; 131614db4f74SSatish Balay nz = baij->maxnz/baij->bs2; 131714db4f74SSatish Balay nbs = baij->nbs; 131827a8da17SBarry Smith for (i=0; i<nz; i++) { 131927a8da17SBarry Smith baij->j[i] = indices[i]; 132027a8da17SBarry Smith } 132127a8da17SBarry Smith baij->nz = nz; 132214db4f74SSatish Balay for (i=0; i<nbs; i++) { 132327a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 132427a8da17SBarry Smith } 132527a8da17SBarry Smith 132627a8da17SBarry Smith PetscFunctionReturn(0); 132727a8da17SBarry Smith } 1328fb2e594dSBarry Smith EXTERN_C_END 132927a8da17SBarry Smith 13304a2ae208SSatish Balay #undef __FUNCT__ 13314a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices" 133227a8da17SBarry Smith /*@ 133327a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 133427a8da17SBarry Smith in the matrix. 133527a8da17SBarry Smith 133627a8da17SBarry Smith Input Parameters: 133727a8da17SBarry Smith + mat - the SeqBAIJ matrix 133827a8da17SBarry Smith - indices - the column indices 133927a8da17SBarry Smith 134015091d37SBarry Smith Level: advanced 134115091d37SBarry Smith 134227a8da17SBarry Smith Notes: 134327a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 134427a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 134527a8da17SBarry Smith of the MatSetValues() operation. 134627a8da17SBarry Smith 134727a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 134827a8da17SBarry Smith MatCreateSeqBAIJ(). 134927a8da17SBarry Smith 135027a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 135127a8da17SBarry Smith 135227a8da17SBarry Smith @*/ 135327a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 135427a8da17SBarry Smith { 135527a8da17SBarry Smith int ierr,(*f)(Mat,int *); 135627a8da17SBarry Smith 135727a8da17SBarry Smith PetscFunctionBegin; 135827a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1359c134de8dSSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); 136027a8da17SBarry Smith if (f) { 136127a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 136227a8da17SBarry Smith } else { 136329bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 136427a8da17SBarry Smith } 136527a8da17SBarry Smith PetscFunctionReturn(0); 136627a8da17SBarry Smith } 136727a8da17SBarry Smith 13684a2ae208SSatish Balay #undef __FUNCT__ 13694a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ" 1370273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v) 1371273d9f13SBarry Smith { 1372273d9f13SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1373273d9f13SBarry Smith int ierr,i,j,n,row,bs,*ai,*aj,mbs; 1374273d9f13SBarry Smith PetscReal atmp; 137587828ca2SBarry Smith PetscScalar *x,zero = 0.0; 1376273d9f13SBarry Smith MatScalar *aa; 1377273d9f13SBarry Smith int ncols,brow,krow,kcol; 1378273d9f13SBarry Smith 1379273d9f13SBarry Smith PetscFunctionBegin; 1380273d9f13SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1381273d9f13SBarry Smith bs = a->bs; 1382273d9f13SBarry Smith aa = a->a; 1383273d9f13SBarry Smith ai = a->i; 1384273d9f13SBarry Smith aj = a->j; 1385273d9f13SBarry Smith mbs = a->mbs; 1386273d9f13SBarry Smith 1387273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 1388273d9f13SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1389273d9f13SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1390273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1391273d9f13SBarry Smith for (i=0; i<mbs; i++) { 1392273d9f13SBarry Smith ncols = ai[1] - ai[0]; ai++; 1393273d9f13SBarry Smith brow = bs*i; 1394273d9f13SBarry Smith for (j=0; j<ncols; j++){ 1395273d9f13SBarry Smith /* bcol = bs*(*aj); */ 1396273d9f13SBarry Smith for (kcol=0; kcol<bs; kcol++){ 1397273d9f13SBarry Smith for (krow=0; krow<bs; krow++){ 1398273d9f13SBarry Smith atmp = PetscAbsScalar(*aa); aa++; 1399273d9f13SBarry Smith row = brow + krow; /* row index */ 1400273d9f13SBarry Smith /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */ 1401273d9f13SBarry Smith if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp; 1402273d9f13SBarry Smith } 1403273d9f13SBarry Smith } 1404273d9f13SBarry Smith aj++; 1405273d9f13SBarry Smith } 1406273d9f13SBarry Smith } 1407273d9f13SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1408273d9f13SBarry Smith PetscFunctionReturn(0); 1409273d9f13SBarry Smith } 1410273d9f13SBarry Smith 14114a2ae208SSatish Balay #undef __FUNCT__ 14124a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ" 1413273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A) 1414273d9f13SBarry Smith { 1415273d9f13SBarry Smith int ierr; 1416273d9f13SBarry Smith 1417273d9f13SBarry Smith PetscFunctionBegin; 1418273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1419273d9f13SBarry Smith PetscFunctionReturn(0); 1420273d9f13SBarry Smith } 1421273d9f13SBarry Smith 14224a2ae208SSatish Balay #undef __FUNCT__ 14234a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ" 14244e7234bfSBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar *array[]) 1425f2a5309cSSatish Balay { 1426f2a5309cSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1427f2a5309cSSatish Balay PetscFunctionBegin; 1428f2a5309cSSatish Balay *array = a->a; 1429f2a5309cSSatish Balay PetscFunctionReturn(0); 1430f2a5309cSSatish Balay } 1431f2a5309cSSatish Balay 14324a2ae208SSatish Balay #undef __FUNCT__ 14334a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ" 14344e7234bfSBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar *array[]) 1435f2a5309cSSatish Balay { 1436f2a5309cSSatish Balay PetscFunctionBegin; 1437f2a5309cSSatish Balay PetscFunctionReturn(0); 1438f2a5309cSSatish Balay } 1439f2a5309cSSatish Balay 144042ee4b1aSHong Zhang #include "petscblaslapack.h" 144142ee4b1aSHong Zhang #undef __FUNCT__ 144242ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqBAIJ" 1443268466fbSBarry Smith int MatAXPY_SeqBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 144442ee4b1aSHong Zhang { 144542ee4b1aSHong Zhang Mat_SeqBAIJ *x = (Mat_SeqBAIJ *)X->data,*y = (Mat_SeqBAIJ *)Y->data; 14466550863bSHong Zhang int ierr,one=1,i,bs=y->bs,j,bs2; 144742ee4b1aSHong Zhang 144842ee4b1aSHong Zhang PetscFunctionBegin; 144942ee4b1aSHong Zhang if (str == SAME_NONZERO_PATTERN) { 1450268466fbSBarry Smith BLaxpy_(&x->nz,(PetscScalar*)a,x->a,&one,y->a,&one); 1451c537a176SHong Zhang } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ 1452c4319e64SHong Zhang if (y->xtoy && y->XtoY != X) { 1453c4319e64SHong Zhang ierr = PetscFree(y->xtoy);CHKERRQ(ierr); 1454c4319e64SHong Zhang ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); 1455c537a176SHong Zhang } 1456c4319e64SHong Zhang if (!y->xtoy) { /* get xtoy */ 1457c4319e64SHong Zhang ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); 1458c4319e64SHong Zhang y->XtoY = X; 1459c537a176SHong Zhang } 1460c4319e64SHong Zhang bs2 = bs*bs; 1461c537a176SHong Zhang for (i=0; i<x->nz; i++) { 1462c4319e64SHong Zhang j = 0; 1463c4319e64SHong Zhang while (j < bs2){ 14646550863bSHong Zhang y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]); 1465c4319e64SHong Zhang j++; 1466c537a176SHong Zhang } 1467c4319e64SHong Zhang } 1468c4319e64SHong 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)); 146942ee4b1aSHong Zhang } else { 147042ee4b1aSHong Zhang ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 147142ee4b1aSHong Zhang } 147242ee4b1aSHong Zhang PetscFunctionReturn(0); 147342ee4b1aSHong Zhang } 147442ee4b1aSHong Zhang 14752593348eSBarry Smith /* -------------------------------------------------------------------*/ 1476cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1477cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1478cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1479cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 148097304618SKris Buschelman /* 4*/ MatMultAdd_SeqBAIJ_N, 14817c922b88SBarry Smith MatMultTranspose_SeqBAIJ, 14827c922b88SBarry Smith MatMultTransposeAdd_SeqBAIJ, 1483cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1484cc2dc46cSBarry Smith 0, 1485cc2dc46cSBarry Smith 0, 148697304618SKris Buschelman /*10*/ 0, 1487cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1488cc2dc46cSBarry Smith 0, 1489b6490206SBarry Smith 0, 1490f2501298SSatish Balay MatTranspose_SeqBAIJ, 149197304618SKris Buschelman /*15*/ MatGetInfo_SeqBAIJ, 1492cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1493cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1494cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1495cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 149697304618SKris Buschelman /*20*/ 0, 1497cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1498cc2dc46cSBarry Smith 0, 1499cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1500cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 150197304618SKris Buschelman /*25*/ MatZeroRows_SeqBAIJ, 1502cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1503cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1504cc2dc46cSBarry Smith 0, 1505cc2dc46cSBarry Smith 0, 150697304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqBAIJ, 1507cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1508cc2dc46cSBarry Smith 0, 1509f2a5309cSSatish Balay MatGetArray_SeqBAIJ, 1510f2a5309cSSatish Balay MatRestoreArray_SeqBAIJ, 151197304618SKris Buschelman /*35*/ MatDuplicate_SeqBAIJ, 1512cc2dc46cSBarry Smith 0, 1513cc2dc46cSBarry Smith 0, 1514cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1515cc2dc46cSBarry Smith 0, 151697304618SKris Buschelman /*40*/ MatAXPY_SeqBAIJ, 1517cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1518cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1519cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1520cc2dc46cSBarry Smith 0, 152197304618SKris Buschelman /*45*/ MatPrintHelp_SeqBAIJ, 1522cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1523cc2dc46cSBarry Smith 0, 1524cc2dc46cSBarry Smith 0, 1525cc2dc46cSBarry Smith 0, 152697304618SKris Buschelman /*50*/ MatGetBlockSize_SeqBAIJ, 15273b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 152892c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1529cc2dc46cSBarry Smith 0, 1530cc2dc46cSBarry Smith 0, 153197304618SKris Buschelman /*55*/ 0, 1532cc2dc46cSBarry Smith 0, 1533cc2dc46cSBarry Smith 0, 1534cc2dc46cSBarry Smith 0, 1535d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 153697304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqBAIJ, 1537b9b97703SBarry Smith MatDestroy_SeqBAIJ, 1538b9b97703SBarry Smith MatView_SeqBAIJ, 15393a64cc32SBarry Smith MatGetPetscMaps_Petsc, 1540273d9f13SBarry Smith 0, 154197304618SKris Buschelman /*65*/ 0, 1542273d9f13SBarry Smith 0, 1543273d9f13SBarry Smith 0, 1544273d9f13SBarry Smith 0, 1545273d9f13SBarry Smith 0, 154697304618SKris Buschelman /*70*/ MatGetRowMax_SeqBAIJ, 154797304618SKris Buschelman MatConvert_Basic, 1548273d9f13SBarry Smith 0, 154997304618SKris Buschelman 0, 155097304618SKris Buschelman 0, 155197304618SKris Buschelman /*75*/ 0, 155297304618SKris Buschelman 0, 155397304618SKris Buschelman 0, 155497304618SKris Buschelman 0, 155597304618SKris Buschelman 0, 155697304618SKris Buschelman /*80*/ 0, 155797304618SKris Buschelman 0, 155897304618SKris Buschelman 0, 155997304618SKris Buschelman 0, 156097304618SKris Buschelman /*85*/ MatLoad_SeqBAIJ 156197304618SKris Buschelman }; 15622593348eSBarry Smith 15633e90b805SBarry Smith EXTERN_C_BEGIN 15644a2ae208SSatish Balay #undef __FUNCT__ 15654a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ" 15663e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 15673e90b805SBarry Smith { 15683e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1569273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2; 1570d132466eSBarry Smith int ierr; 15713e90b805SBarry Smith 15723e90b805SBarry Smith PetscFunctionBegin; 15733e90b805SBarry Smith if (aij->nonew != 1) { 157429bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 15753e90b805SBarry Smith } 15763e90b805SBarry Smith 15773e90b805SBarry Smith /* allocate space for values if not already there */ 15783e90b805SBarry Smith if (!aij->saved_values) { 157987828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); 15803e90b805SBarry Smith } 15813e90b805SBarry Smith 15823e90b805SBarry Smith /* copy values over */ 158387828ca2SBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); 15843e90b805SBarry Smith PetscFunctionReturn(0); 15853e90b805SBarry Smith } 15863e90b805SBarry Smith EXTERN_C_END 15873e90b805SBarry Smith 15883e90b805SBarry Smith EXTERN_C_BEGIN 15894a2ae208SSatish Balay #undef __FUNCT__ 15904a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ" 159132e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 15923e90b805SBarry Smith { 15933e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1594273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 15953e90b805SBarry Smith 15963e90b805SBarry Smith PetscFunctionBegin; 15973e90b805SBarry Smith if (aij->nonew != 1) { 159829bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 15993e90b805SBarry Smith } 16003e90b805SBarry Smith if (!aij->saved_values) { 160129bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 16023e90b805SBarry Smith } 16033e90b805SBarry Smith 16043e90b805SBarry Smith /* copy values over */ 160587828ca2SBarry Smith ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); 16063e90b805SBarry Smith PetscFunctionReturn(0); 16073e90b805SBarry Smith } 16083e90b805SBarry Smith EXTERN_C_END 16093e90b805SBarry Smith 1610273d9f13SBarry Smith EXTERN_C_BEGIN 1611f248c16bSBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,const MatType,Mat*); 1612273d9f13SBarry Smith EXTERN_C_END 1613273d9f13SBarry Smith 1614273d9f13SBarry Smith EXTERN_C_BEGIN 16154a2ae208SSatish Balay #undef __FUNCT__ 1616a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqBAIJSetPreallocation_SeqBAIJ" 1617a23d5eceSKris Buschelman int MatSeqBAIJSetPreallocation_SeqBAIJ(Mat B,int bs,int nz,int *nnz) 1618a23d5eceSKris Buschelman { 1619a23d5eceSKris Buschelman Mat_SeqBAIJ *b; 1620a23d5eceSKris Buschelman int i,len,ierr,mbs,nbs,bs2,newbs = bs; 1621a23d5eceSKris Buschelman PetscTruth flg; 1622a23d5eceSKris Buschelman 1623a23d5eceSKris Buschelman PetscFunctionBegin; 1624a23d5eceSKris Buschelman 1625a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 1626a23d5eceSKris Buschelman ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr); 1627a23d5eceSKris Buschelman if (nnz && newbs != bs) { 1628a23d5eceSKris Buschelman SETERRQ(1,"Cannot change blocksize from command line if setting nnz"); 1629a23d5eceSKris Buschelman } 1630a23d5eceSKris Buschelman bs = newbs; 1631a23d5eceSKris Buschelman 1632a23d5eceSKris Buschelman mbs = B->m/bs; 1633a23d5eceSKris Buschelman nbs = B->n/bs; 1634a23d5eceSKris Buschelman bs2 = bs*bs; 1635a23d5eceSKris Buschelman 1636a23d5eceSKris Buschelman if (mbs*bs!=B->m || nbs*bs!=B->n) { 1637a23d5eceSKris Buschelman SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs); 1638a23d5eceSKris Buschelman } 1639a23d5eceSKris Buschelman 1640a23d5eceSKris Buschelman if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1641a23d5eceSKris Buschelman if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 1642a23d5eceSKris Buschelman if (nnz) { 1643a23d5eceSKris Buschelman for (i=0; i<mbs; i++) { 1644a23d5eceSKris Buschelman if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1645a23d5eceSKris 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); 1646a23d5eceSKris Buschelman } 1647a23d5eceSKris Buschelman } 1648a23d5eceSKris Buschelman 1649a23d5eceSKris Buschelman b = (Mat_SeqBAIJ*)B->data; 1650a23d5eceSKris Buschelman ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 1651a23d5eceSKris Buschelman B->ops->solve = MatSolve_SeqBAIJ_Update; 1652a23d5eceSKris Buschelman B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_Update; 1653a23d5eceSKris Buschelman if (!flg) { 1654a23d5eceSKris Buschelman switch (bs) { 1655a23d5eceSKris Buschelman case 1: 1656a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1657a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_1; 1658a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_1; 1659a23d5eceSKris Buschelman break; 1660a23d5eceSKris Buschelman case 2: 1661a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1662a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_2; 1663a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_2; 1664a23d5eceSKris Buschelman break; 1665a23d5eceSKris Buschelman case 3: 1666a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1667a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_3; 1668a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_3; 1669a23d5eceSKris Buschelman break; 1670a23d5eceSKris Buschelman case 4: 1671a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1672a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_4; 1673a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_4; 1674a23d5eceSKris Buschelman break; 1675a23d5eceSKris Buschelman case 5: 1676a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1677a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_5; 1678a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_5; 1679a23d5eceSKris Buschelman break; 1680a23d5eceSKris Buschelman case 6: 1681a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 1682a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_6; 1683a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_6; 1684a23d5eceSKris Buschelman break; 1685a23d5eceSKris Buschelman case 7: 1686a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7; 1687a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_7; 1688a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_7; 1689a23d5eceSKris Buschelman break; 1690a23d5eceSKris Buschelman default: 1691a23d5eceSKris Buschelman B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N; 1692a23d5eceSKris Buschelman B->ops->mult = MatMult_SeqBAIJ_N; 1693a23d5eceSKris Buschelman B->ops->multadd = MatMultAdd_SeqBAIJ_N; 1694a23d5eceSKris Buschelman break; 1695a23d5eceSKris Buschelman } 1696a23d5eceSKris Buschelman } 1697a23d5eceSKris Buschelman b->bs = bs; 1698a23d5eceSKris Buschelman b->mbs = mbs; 1699a23d5eceSKris Buschelman b->nbs = nbs; 1700a23d5eceSKris Buschelman ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 1701a23d5eceSKris Buschelman if (!nnz) { 1702a23d5eceSKris Buschelman if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1703a23d5eceSKris Buschelman else if (nz <= 0) nz = 1; 1704a23d5eceSKris Buschelman for (i=0; i<mbs; i++) b->imax[i] = nz; 1705a23d5eceSKris Buschelman nz = nz*mbs; 1706a23d5eceSKris Buschelman } else { 1707a23d5eceSKris Buschelman nz = 0; 1708a23d5eceSKris Buschelman for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 1709a23d5eceSKris Buschelman } 1710a23d5eceSKris Buschelman 1711a23d5eceSKris Buschelman /* allocate the matrix space */ 1712a23d5eceSKris Buschelman len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 1713a23d5eceSKris Buschelman ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 1714a23d5eceSKris Buschelman ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1715a23d5eceSKris Buschelman b->j = (int*)(b->a + nz*bs2); 1716a23d5eceSKris Buschelman ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 1717a23d5eceSKris Buschelman b->i = b->j + nz; 1718a23d5eceSKris Buschelman b->singlemalloc = PETSC_TRUE; 1719a23d5eceSKris Buschelman 1720a23d5eceSKris Buschelman b->i[0] = 0; 1721a23d5eceSKris Buschelman for (i=1; i<mbs+1; i++) { 1722a23d5eceSKris Buschelman b->i[i] = b->i[i-1] + b->imax[i-1]; 1723a23d5eceSKris Buschelman } 1724a23d5eceSKris Buschelman 1725a23d5eceSKris Buschelman /* b->ilen will count nonzeros in each block row so far. */ 1726a23d5eceSKris Buschelman ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1727a23d5eceSKris Buschelman PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1728a23d5eceSKris Buschelman for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 1729a23d5eceSKris Buschelman 1730a23d5eceSKris Buschelman b->bs = bs; 1731a23d5eceSKris Buschelman b->bs2 = bs2; 1732a23d5eceSKris Buschelman b->mbs = mbs; 1733a23d5eceSKris Buschelman b->nz = 0; 1734a23d5eceSKris Buschelman b->maxnz = nz*bs2; 1735a23d5eceSKris Buschelman B->info.nz_unneeded = (PetscReal)b->maxnz; 1736a23d5eceSKris Buschelman PetscFunctionReturn(0); 1737a23d5eceSKris Buschelman } 1738a23d5eceSKris Buschelman EXTERN_C_END 1739a23d5eceSKris Buschelman 17400bad9183SKris Buschelman /*MC 1741fafad747SKris Buschelman MATSEQBAIJ - MATSEQBAIJ = "seqbaij" - A matrix type to be used for sequential block sparse matrices, based on 17420bad9183SKris Buschelman block sparse compressed row format. 17430bad9183SKris Buschelman 17440bad9183SKris Buschelman Options Database Keys: 17450bad9183SKris Buschelman . -mat_type seqbaij - sets the matrix type to "seqbaij" during a call to MatSetFromOptions() 17460bad9183SKris Buschelman 17470bad9183SKris Buschelman Level: beginner 17480bad9183SKris Buschelman 17490bad9183SKris Buschelman .seealso: MatCreateSeqBAIJ 17500bad9183SKris Buschelman M*/ 17510bad9183SKris Buschelman 1752a23d5eceSKris Buschelman EXTERN_C_BEGIN 1753a23d5eceSKris Buschelman #undef __FUNCT__ 17544a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ" 1755273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B) 17562593348eSBarry Smith { 1757273d9f13SBarry Smith int ierr,size; 1758b6490206SBarry Smith Mat_SeqBAIJ *b; 17593b2fbd54SBarry Smith 17603a40ed3dSBarry Smith PetscFunctionBegin; 1761273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 176229bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1763b6490206SBarry Smith 1764273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1765273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1766b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr); 1767b0a32e0cSBarry Smith B->data = (void*)b; 1768549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1769549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 17702593348eSBarry Smith B->factor = 0; 17712593348eSBarry Smith B->lupivotthreshold = 1.0; 177290f02eecSBarry Smith B->mapping = 0; 17732593348eSBarry Smith b->row = 0; 17742593348eSBarry Smith b->col = 0; 1775e51c0b9cSSatish Balay b->icol = 0; 17762593348eSBarry Smith b->reallocs = 0; 17773e90b805SBarry Smith b->saved_values = 0; 1778cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 1779563b5814SBarry Smith b->setvalueslen = 0; 1780563b5814SBarry Smith b->setvaluescopy = PETSC_NULL; 1781563b5814SBarry Smith #endif 17822593348eSBarry Smith 17833a64cc32SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 17843a64cc32SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1785a5ae1ecdSBarry Smith 1786c4992f7dSBarry Smith b->sorted = PETSC_FALSE; 1787c4992f7dSBarry Smith b->roworiented = PETSC_TRUE; 17882593348eSBarry Smith b->nonew = 0; 17892593348eSBarry Smith b->diag = 0; 17902593348eSBarry Smith b->solve_work = 0; 1791de6a44a3SBarry Smith b->mult_work = 0; 17922a1b7f2aSHong Zhang B->spptr = 0; 17930e6d2581SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 1794883fce79SBarry Smith b->keepzeroedrows = PETSC_FALSE; 1795c4319e64SHong Zhang b->xtoy = 0; 1796c4319e64SHong Zhang b->XtoY = 0; 17974e220ebcSLois Curfman McInnes 1798f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 17993e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 1800bc4b532fSSatish Balay MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1801f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 18023e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 1803bc4b532fSSatish Balay MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1804f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 180527a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1806bc4b532fSSatish Balay MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 1807a6175056SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqaij_C", 1808273d9f13SBarry Smith "MatConvert_SeqBAIJ_SeqAIJ", 1809273d9f13SBarry Smith MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr); 1810a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetPreallocation_C", 1811a23d5eceSKris Buschelman "MatSeqBAIJSetPreallocation_SeqBAIJ", 1812a23d5eceSKris Buschelman MatSeqBAIJSetPreallocation_SeqBAIJ);CHKERRQ(ierr); 18133a40ed3dSBarry Smith PetscFunctionReturn(0); 18142593348eSBarry Smith } 1815273d9f13SBarry Smith EXTERN_C_END 18162593348eSBarry Smith 18174a2ae208SSatish Balay #undef __FUNCT__ 18184a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ" 18192e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 18202593348eSBarry Smith { 18212593348eSBarry Smith Mat C; 1822b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data; 182327a8da17SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1824de6a44a3SBarry Smith 18253a40ed3dSBarry Smith PetscFunctionBegin; 182629bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 18272593348eSBarry Smith 18282593348eSBarry Smith *B = 0; 1829273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1830273d9f13SBarry Smith ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr); 1831273d9f13SBarry Smith c = (Mat_SeqBAIJ*)C->data; 183244cd7ae7SLois Curfman McInnes 1833b6490206SBarry Smith c->bs = a->bs; 18349df24120SSatish Balay c->bs2 = a->bs2; 1835b6490206SBarry Smith c->mbs = a->mbs; 1836b6490206SBarry Smith c->nbs = a->nbs; 183735d8aa7fSBarry Smith ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 18382593348eSBarry Smith 1839b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1840b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 1841b6490206SBarry Smith for (i=0; i<mbs; i++) { 18422593348eSBarry Smith c->imax[i] = a->imax[i]; 18432593348eSBarry Smith c->ilen[i] = a->ilen[i]; 18442593348eSBarry Smith } 18452593348eSBarry Smith 18462593348eSBarry Smith /* allocate the matrix space */ 1847c4992f7dSBarry Smith c->singlemalloc = PETSC_TRUE; 18483f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 1849b0a32e0cSBarry Smith ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 18507e67e3f9SSatish Balay c->j = (int*)(c->a + nz*bs2); 1851de6a44a3SBarry Smith c->i = c->j + nz; 1852549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1853b6490206SBarry Smith if (mbs > 0) { 1854549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 18552e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1856549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 18572e8a6d31SBarry Smith } else { 1858549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 18592593348eSBarry Smith } 18602593348eSBarry Smith } 18612593348eSBarry Smith 1862b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 18632593348eSBarry Smith c->sorted = a->sorted; 18642593348eSBarry Smith c->roworiented = a->roworiented; 18652593348eSBarry Smith c->nonew = a->nonew; 18662593348eSBarry Smith 18672593348eSBarry Smith if (a->diag) { 1868b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1869b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 1870b6490206SBarry Smith for (i=0; i<mbs; i++) { 18712593348eSBarry Smith c->diag[i] = a->diag[i]; 18722593348eSBarry Smith } 187398305bb5SBarry Smith } else c->diag = 0; 18742593348eSBarry Smith c->nz = a->nz; 18752593348eSBarry Smith c->maxnz = a->maxnz; 18762593348eSBarry Smith c->solve_work = 0; 18772a1b7f2aSHong Zhang C->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 18787fc0212eSBarry Smith c->mult_work = 0; 1879273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 1880273d9f13SBarry Smith C->assembled = PETSC_TRUE; 18812593348eSBarry Smith *B = C; 1882b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 18833a40ed3dSBarry Smith PetscFunctionReturn(0); 18842593348eSBarry Smith } 18852593348eSBarry Smith 18864a2ae208SSatish Balay #undef __FUNCT__ 18874a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ" 18888e9aea5cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,const MatType type,Mat *A) 18892593348eSBarry Smith { 1890b6490206SBarry Smith Mat_SeqBAIJ *a; 18912593348eSBarry Smith Mat B; 1892f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1893b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 189435aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1895a7c10996SSatish Balay int *masked,nmask,tmp,bs2,ishift; 189687828ca2SBarry Smith PetscScalar *aa; 189719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 18982593348eSBarry Smith 18993a40ed3dSBarry Smith PetscFunctionBegin; 1900b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1901de6a44a3SBarry Smith bs2 = bs*bs; 1902b6490206SBarry Smith 1903d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 190429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1905b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 19060752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1907552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 19082593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 19092593348eSBarry Smith 1910d64ed03dSBarry Smith if (header[3] < 0) { 191129bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ"); 1912d64ed03dSBarry Smith } 1913d64ed03dSBarry Smith 191429bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 191535aab85fSBarry Smith 191635aab85fSBarry Smith /* 191735aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 191835aab85fSBarry Smith divisible by the blocksize 191935aab85fSBarry Smith */ 1920b6490206SBarry Smith mbs = M/bs; 192135aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 192235aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 192335aab85fSBarry Smith else mbs++; 192435aab85fSBarry Smith if (extra_rows) { 1925b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 192635aab85fSBarry Smith } 1927b6490206SBarry Smith 19282593348eSBarry Smith /* read in row lengths */ 1929b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 19300752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 193135aab85fSBarry Smith for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 19322593348eSBarry Smith 1933b6490206SBarry Smith /* read in column indices */ 1934b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 19350752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 193635aab85fSBarry Smith for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 1937b6490206SBarry Smith 1938b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1939b0a32e0cSBarry Smith ierr = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 1940549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 1941b0a32e0cSBarry Smith ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 1942549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 194335aab85fSBarry Smith masked = mask + mbs; 1944b6490206SBarry Smith rowcount = 0; nzcount = 0; 1945b6490206SBarry Smith for (i=0; i<mbs; i++) { 194635aab85fSBarry Smith nmask = 0; 1947b6490206SBarry Smith for (j=0; j<bs; j++) { 1948b6490206SBarry Smith kmax = rowlengths[rowcount]; 1949b6490206SBarry Smith for (k=0; k<kmax; k++) { 195035aab85fSBarry Smith tmp = jj[nzcount++]/bs; 195135aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1952b6490206SBarry Smith } 1953b6490206SBarry Smith rowcount++; 1954b6490206SBarry Smith } 195535aab85fSBarry Smith browlengths[i] += nmask; 195635aab85fSBarry Smith /* zero out the mask elements we set */ 195735aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1958b6490206SBarry Smith } 1959b6490206SBarry Smith 19602593348eSBarry Smith /* create our matrix */ 1961dd23797bSKris Buschelman ierr = MatCreate(comm,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows,&B); 196278ae41b4SKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 196378ae41b4SKris Buschelman ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 1964b6490206SBarry Smith a = (Mat_SeqBAIJ*)B->data; 19652593348eSBarry Smith 19662593348eSBarry Smith /* set matrix "i" values */ 1967de6a44a3SBarry Smith a->i[0] = 0; 1968b6490206SBarry Smith for (i=1; i<= mbs; i++) { 1969b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1970b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 19712593348eSBarry Smith } 1972b6490206SBarry Smith a->nz = 0; 1973b6490206SBarry Smith for (i=0; i<mbs; i++) a->nz += browlengths[i]; 19742593348eSBarry Smith 1975b6490206SBarry Smith /* read in nonzero values */ 197687828ca2SBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 19770752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 197835aab85fSBarry Smith for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 1979b6490206SBarry Smith 1980b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1981b6490206SBarry Smith nzcount = 0; jcount = 0; 1982b6490206SBarry Smith for (i=0; i<mbs; i++) { 1983b6490206SBarry Smith nzcountb = nzcount; 198435aab85fSBarry Smith nmask = 0; 1985b6490206SBarry Smith for (j=0; j<bs; j++) { 1986b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1987b6490206SBarry Smith for (k=0; k<kmax; k++) { 198835aab85fSBarry Smith tmp = jj[nzcount++]/bs; 198935aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1990b6490206SBarry Smith } 1991b6490206SBarry Smith } 1992de6a44a3SBarry Smith /* sort the masked values */ 1993433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1994de6a44a3SBarry Smith 1995b6490206SBarry Smith /* set "j" values into matrix */ 1996b6490206SBarry Smith maskcount = 1; 199735aab85fSBarry Smith for (j=0; j<nmask; j++) { 199835aab85fSBarry Smith a->j[jcount++] = masked[j]; 1999de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 2000b6490206SBarry Smith } 2001b6490206SBarry Smith /* set "a" values into matrix */ 2002de6a44a3SBarry Smith ishift = bs2*a->i[i]; 2003b6490206SBarry Smith for (j=0; j<bs; j++) { 2004b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 2005b6490206SBarry Smith for (k=0; k<kmax; k++) { 2006de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 2007de6a44a3SBarry Smith block = mask[tmp] - 1; 2008de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 2009de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 2010375fe846SBarry Smith a->a[idx] = (MatScalar)aa[nzcountb++]; 2011b6490206SBarry Smith } 2012b6490206SBarry Smith } 201335aab85fSBarry Smith /* zero out the mask elements we set */ 201435aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 2015b6490206SBarry Smith } 201629bbc08cSBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 2017b6490206SBarry Smith 2018606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 2019606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 2020606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 2021606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 2022606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 2023b6490206SBarry Smith 202478ae41b4SKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 202578ae41b4SKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20269c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 202778ae41b4SKris Buschelman 202878ae41b4SKris Buschelman *A = B; 20293a40ed3dSBarry Smith PetscFunctionReturn(0); 20302593348eSBarry Smith } 20312593348eSBarry Smith 20324a2ae208SSatish Balay #undef __FUNCT__ 20334a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ" 2034273d9f13SBarry Smith /*@C 2035273d9f13SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 2036273d9f13SBarry Smith compressed row) format. For good matrix assembly performance the 2037273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 2038273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 2039273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 20402593348eSBarry Smith 2041273d9f13SBarry Smith Collective on MPI_Comm 2042273d9f13SBarry Smith 2043273d9f13SBarry Smith Input Parameters: 2044273d9f13SBarry Smith + comm - MPI communicator, set to PETSC_COMM_SELF 2045273d9f13SBarry Smith . bs - size of block 2046273d9f13SBarry Smith . m - number of rows 2047273d9f13SBarry Smith . n - number of columns 204835d8aa7fSBarry Smith . nz - number of nonzero blocks per block row (same for all rows) 204935d8aa7fSBarry Smith - nnz - array containing the number of nonzero blocks in the various block rows 2050273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 2051273d9f13SBarry Smith 2052273d9f13SBarry Smith Output Parameter: 2053273d9f13SBarry Smith . A - the matrix 2054273d9f13SBarry Smith 2055273d9f13SBarry Smith Options Database Keys: 2056273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 2057273d9f13SBarry Smith block calculations (much slower) 2058273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 2059273d9f13SBarry Smith 2060273d9f13SBarry Smith Level: intermediate 2061273d9f13SBarry Smith 2062273d9f13SBarry Smith Notes: 206335d8aa7fSBarry Smith A nonzero block is any block that as 1 or more nonzeros in it 206435d8aa7fSBarry Smith 2065273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 2066273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 2067273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 2068273d9f13SBarry Smith 2069273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2070273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 2071273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 2072273d9f13SBarry Smith matrices. 2073273d9f13SBarry Smith 2074273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 2075273d9f13SBarry Smith @*/ 2076ca01db9bSBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A) 2077273d9f13SBarry Smith { 2078273d9f13SBarry Smith int ierr; 2079273d9f13SBarry Smith 2080273d9f13SBarry Smith PetscFunctionBegin; 2081273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 2082273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr); 2083273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 2084273d9f13SBarry Smith PetscFunctionReturn(0); 2085273d9f13SBarry Smith } 2086273d9f13SBarry Smith 20874a2ae208SSatish Balay #undef __FUNCT__ 20884a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation" 2089273d9f13SBarry Smith /*@C 2090273d9f13SBarry Smith MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros 2091273d9f13SBarry Smith per row in the matrix. For good matrix assembly performance the 2092273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 2093273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 2094273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 2095273d9f13SBarry Smith 2096273d9f13SBarry Smith Collective on MPI_Comm 2097273d9f13SBarry Smith 2098273d9f13SBarry Smith Input Parameters: 2099273d9f13SBarry Smith + A - the matrix 2100273d9f13SBarry Smith . bs - size of block 2101273d9f13SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 2102273d9f13SBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 2103273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 2104273d9f13SBarry Smith 2105273d9f13SBarry Smith Options Database Keys: 2106273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 2107273d9f13SBarry Smith block calculations (much slower) 2108273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 2109273d9f13SBarry Smith 2110273d9f13SBarry Smith Level: intermediate 2111273d9f13SBarry Smith 2112273d9f13SBarry Smith Notes: 2113273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 2114273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 2115273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 2116273d9f13SBarry Smith 2117273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2118273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 2119273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 2120273d9f13SBarry Smith matrices. 2121273d9f13SBarry Smith 2122273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 2123273d9f13SBarry Smith @*/ 2124ca01db9bSBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) 2125273d9f13SBarry Smith { 2126ca01db9bSBarry Smith int ierr,(*f)(Mat,int,int,const int[]); 2127273d9f13SBarry Smith 2128273d9f13SBarry Smith PetscFunctionBegin; 2129a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 2130a23d5eceSKris Buschelman if (f) { 2131a23d5eceSKris Buschelman ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr); 2132273d9f13SBarry Smith } 2133273d9f13SBarry Smith PetscFunctionReturn(0); 2134273d9f13SBarry Smith } 2135