1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*6831982aSBarry Smith static char vcid[] = "$Id: baij.c,v 1.185 1999/10/13 20:37:28 bsmith Exp bsmith $"; 32593348eSBarry Smith #endif 42593348eSBarry Smith 52593348eSBarry Smith /* 6b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 72593348eSBarry Smith matrix storage format. 82593348eSBarry Smith */ 93369ce9aSBarry Smith #include "sys.h" 1070f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 111a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 121a6a6d98SBarry Smith #include "src/inline/spops.h" 133b547af2SSatish Balay 142d61bbb3SSatish Balay #define CHUNKSIZE 10 15de6a44a3SBarry Smith 16be5855fcSBarry Smith /* 17be5855fcSBarry Smith Checks for missing diagonals 18be5855fcSBarry Smith */ 19be5855fcSBarry Smith #undef __FUNC__ 20be5855fcSBarry Smith #define __FUNC__ "MatMissingDiag_SeqBAIJ" 21be5855fcSBarry Smith int MatMissingDiag_SeqBAIJ(Mat A) 22be5855fcSBarry Smith { 23be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 24be5855fcSBarry Smith int *diag = a->diag, *jj = a->j,i; 25be5855fcSBarry Smith 26be5855fcSBarry Smith PetscFunctionBegin; 270e8e8aceSBarry Smith for ( i=0; i<a->mbs; i++ ) { 28be5855fcSBarry Smith if (jj[diag[i]] != i) { 29be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 30be5855fcSBarry Smith } 31be5855fcSBarry Smith } 32be5855fcSBarry Smith PetscFunctionReturn(0); 33be5855fcSBarry Smith } 34be5855fcSBarry Smith 355615d1e5SSatish Balay #undef __FUNC__ 365615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ" 37de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A) 38de6a44a3SBarry Smith { 39de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 407fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 41de6a44a3SBarry Smith 423a40ed3dSBarry Smith PetscFunctionBegin; 43de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 44de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 457fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 46ecc615b2SBarry Smith diag[i] = a->i[i+1]; 47de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 48de6a44a3SBarry Smith if (a->j[j] == i) { 49de6a44a3SBarry Smith diag[i] = j; 50de6a44a3SBarry Smith break; 51de6a44a3SBarry Smith } 52de6a44a3SBarry Smith } 53de6a44a3SBarry Smith } 54de6a44a3SBarry Smith a->diag = diag; 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56de6a44a3SBarry Smith } 572593348eSBarry Smith 582593348eSBarry Smith 593b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 603b2fbd54SBarry Smith 615615d1e5SSatish Balay #undef __FUNC__ 625615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 633b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 643b2fbd54SBarry Smith PetscTruth *done) 653b2fbd54SBarry Smith { 663b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 673b2fbd54SBarry Smith int ierr, n = a->mbs,i; 683b2fbd54SBarry Smith 693a40ed3dSBarry Smith PetscFunctionBegin; 703b2fbd54SBarry Smith *nn = n; 713a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 723b2fbd54SBarry Smith if (symmetric) { 733b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 743b2fbd54SBarry Smith } else if (oshift == 1) { 753b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 763b2fbd54SBarry Smith int nz = a->i[n] + 1; 773b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 783b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 793b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 803b2fbd54SBarry Smith } else { 813b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 823b2fbd54SBarry Smith } 833b2fbd54SBarry Smith 843a40ed3dSBarry Smith PetscFunctionReturn(0); 853b2fbd54SBarry Smith } 863b2fbd54SBarry Smith 875615d1e5SSatish Balay #undef __FUNC__ 88d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 893b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 903b2fbd54SBarry Smith PetscTruth *done) 913b2fbd54SBarry Smith { 923b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 93606d414cSSatish Balay int i,n = a->mbs,ierr; 943b2fbd54SBarry Smith 953a40ed3dSBarry Smith PetscFunctionBegin; 963a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 973b2fbd54SBarry Smith if (symmetric) { 98606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 99606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 100af5da2bfSBarry Smith } else if (oshift == 1) { 1013b2fbd54SBarry Smith int nz = a->i[n]; 1023b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 1033b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 1043b2fbd54SBarry Smith } 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 1063b2fbd54SBarry Smith } 1073b2fbd54SBarry Smith 1082d61bbb3SSatish Balay #undef __FUNC__ 1092d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 1102d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 1112d61bbb3SSatish Balay { 1122d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 1132d61bbb3SSatish Balay 1142d61bbb3SSatish Balay PetscFunctionBegin; 1152d61bbb3SSatish Balay *bs = baij->bs; 1162d61bbb3SSatish Balay PetscFunctionReturn(0); 1172d61bbb3SSatish Balay } 1182d61bbb3SSatish Balay 1192d61bbb3SSatish Balay 1202d61bbb3SSatish Balay #undef __FUNC__ 1212d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 122e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1232d61bbb3SSatish Balay { 1242d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 125e51c0b9cSSatish Balay int ierr; 1262d61bbb3SSatish Balay 127433994e6SBarry Smith PetscFunctionBegin; 12894d884c6SBarry Smith if (A->mapping) { 12994d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 13094d884c6SBarry Smith } 13194d884c6SBarry Smith if (A->bmapping) { 13294d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 13394d884c6SBarry Smith } 13461b13de0SBarry Smith if (A->rmap) { 13561b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 13661b13de0SBarry Smith } 13761b13de0SBarry Smith if (A->cmap) { 13861b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 13961b13de0SBarry Smith } 140aa482453SBarry Smith #if defined(PETSC_USE_LOG) 141e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1422d61bbb3SSatish Balay #endif 143606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 144606d414cSSatish Balay if (!a->singlemalloc) { 145606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 146606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 147606d414cSSatish Balay } 148606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 149606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 150606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 151606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 152606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 153e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 154606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 155606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1562d61bbb3SSatish Balay PLogObjectDestroy(A); 1572d61bbb3SSatish Balay PetscHeaderDestroy(A); 1582d61bbb3SSatish Balay PetscFunctionReturn(0); 1592d61bbb3SSatish Balay } 1602d61bbb3SSatish Balay 1612d61bbb3SSatish Balay #undef __FUNC__ 1622d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1632d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1642d61bbb3SSatish Balay { 1652d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1662d61bbb3SSatish Balay 1672d61bbb3SSatish Balay PetscFunctionBegin; 1682d61bbb3SSatish Balay if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 1692d61bbb3SSatish Balay else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 1702d61bbb3SSatish Balay else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 1712d61bbb3SSatish Balay else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 1722d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1734787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1744787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1752d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1762d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1772d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1782d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1792d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1802d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 181b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 182b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 183981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1842d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1852d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1862d61bbb3SSatish Balay } else { 1872d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1882d61bbb3SSatish Balay } 1892d61bbb3SSatish Balay PetscFunctionReturn(0); 1902d61bbb3SSatish Balay } 1912d61bbb3SSatish Balay 1922d61bbb3SSatish Balay 1932d61bbb3SSatish Balay #undef __FUNC__ 1942d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 1952d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 1962d61bbb3SSatish Balay { 1972d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1982d61bbb3SSatish Balay 1992d61bbb3SSatish Balay PetscFunctionBegin; 2002d61bbb3SSatish Balay if (m) *m = a->m; 2012d61bbb3SSatish Balay if (n) *n = a->n; 2022d61bbb3SSatish Balay PetscFunctionReturn(0); 2032d61bbb3SSatish Balay } 2042d61bbb3SSatish Balay 2052d61bbb3SSatish Balay #undef __FUNC__ 2062d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 2072d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2082d61bbb3SSatish Balay { 2092d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2102d61bbb3SSatish Balay 2112d61bbb3SSatish Balay PetscFunctionBegin; 2122d61bbb3SSatish Balay *m = 0; *n = a->m; 2132d61bbb3SSatish Balay PetscFunctionReturn(0); 2142d61bbb3SSatish Balay } 2152d61bbb3SSatish Balay 2162d61bbb3SSatish Balay #undef __FUNC__ 2172d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 2182d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2192d61bbb3SSatish Balay { 2202d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2212d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2223f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2233f1db9ecSBarry Smith Scalar *v_i; 2242d61bbb3SSatish Balay 2252d61bbb3SSatish Balay PetscFunctionBegin; 2262d61bbb3SSatish Balay bs = a->bs; 2272d61bbb3SSatish Balay ai = a->i; 2282d61bbb3SSatish Balay aj = a->j; 2292d61bbb3SSatish Balay aa = a->a; 2302d61bbb3SSatish Balay bs2 = a->bs2; 2312d61bbb3SSatish Balay 2322d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2332d61bbb3SSatish Balay 2342d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2352d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2362d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2372d61bbb3SSatish Balay *nz = bs*M; 2382d61bbb3SSatish Balay 2392d61bbb3SSatish Balay if (v) { 2402d61bbb3SSatish Balay *v = 0; 2412d61bbb3SSatish Balay if (*nz) { 2422d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v); 2432d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2442d61bbb3SSatish Balay v_i = *v + i*bs; 2452d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2462d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2472d61bbb3SSatish Balay } 2482d61bbb3SSatish Balay } 2492d61bbb3SSatish Balay } 2502d61bbb3SSatish Balay 2512d61bbb3SSatish Balay if (idx) { 2522d61bbb3SSatish Balay *idx = 0; 2532d61bbb3SSatish Balay if (*nz) { 2542d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx); 2552d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2562d61bbb3SSatish Balay idx_i = *idx + i*bs; 2572d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2582d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2592d61bbb3SSatish Balay } 2602d61bbb3SSatish Balay } 2612d61bbb3SSatish Balay } 2622d61bbb3SSatish Balay PetscFunctionReturn(0); 2632d61bbb3SSatish Balay } 2642d61bbb3SSatish Balay 2652d61bbb3SSatish Balay #undef __FUNC__ 2662d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2672d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2682d61bbb3SSatish Balay { 269606d414cSSatish Balay int ierr; 270606d414cSSatish Balay 2712d61bbb3SSatish Balay PetscFunctionBegin; 272606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 273606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2742d61bbb3SSatish Balay PetscFunctionReturn(0); 2752d61bbb3SSatish Balay } 2762d61bbb3SSatish Balay 2772d61bbb3SSatish Balay #undef __FUNC__ 2782d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2792d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2802d61bbb3SSatish Balay { 2812d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2822d61bbb3SSatish Balay Mat C; 2832d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2842d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2853f1db9ecSBarry Smith MatScalar *array = a->a; 2862d61bbb3SSatish Balay 2872d61bbb3SSatish Balay PetscFunctionBegin; 2882d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2892d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 290549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2912d61bbb3SSatish Balay 2922d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2932d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 294606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 2952d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 2962d61bbb3SSatish Balay cols = rows + bs; 2972d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 2982d61bbb3SSatish Balay cols[0] = i*bs; 2992d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 3002d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3012d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 3022d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3032d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 3042d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3052d61bbb3SSatish Balay array += bs2; 3062d61bbb3SSatish Balay } 3072d61bbb3SSatish Balay } 308606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 3092d61bbb3SSatish Balay 3102d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3112d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3122d61bbb3SSatish Balay 3132d61bbb3SSatish Balay if (B != PETSC_NULL) { 3142d61bbb3SSatish Balay *B = C; 3152d61bbb3SSatish Balay } else { 316f830108cSBarry Smith PetscOps *Abops; 317cc2dc46cSBarry Smith MatOps Aops; 318f830108cSBarry Smith 3192d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 320606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 321606d414cSSatish Balay if (!a->singlemalloc) { 322606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 323606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 324606d414cSSatish Balay } 325606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 326606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 327606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 328606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 329606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 330f830108cSBarry Smith 3317413bad6SBarry Smith 3327413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3337413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3347413bad6SBarry Smith 335f830108cSBarry Smith /* 336f830108cSBarry Smith This is horrible, horrible code. We need to keep the 337f830108cSBarry Smith A pointers for the bops and ops but copy everything 338f830108cSBarry Smith else from C. 339f830108cSBarry Smith */ 340f830108cSBarry Smith Abops = A->bops; 341f830108cSBarry Smith Aops = A->ops; 342549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 343f830108cSBarry Smith A->bops = Abops; 344f830108cSBarry Smith A->ops = Aops; 34527a8da17SBarry Smith A->qlist = 0; 346f830108cSBarry Smith 3472d61bbb3SSatish Balay PetscHeaderDestroy(C); 3482d61bbb3SSatish Balay } 3492d61bbb3SSatish Balay PetscFunctionReturn(0); 3502d61bbb3SSatish Balay } 3512d61bbb3SSatish Balay 3525615d1e5SSatish Balay #undef __FUNC__ 353d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 354b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3552593348eSBarry Smith { 356b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3579df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 358b6490206SBarry Smith Scalar *aa; 359ce6f0cecSBarry Smith FILE *file; 3602593348eSBarry Smith 3613a40ed3dSBarry Smith PetscFunctionBegin; 36290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3632593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3642593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3653b2fbd54SBarry Smith 3662593348eSBarry Smith col_lens[1] = a->m; 3672593348eSBarry Smith col_lens[2] = a->n; 3687e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3692593348eSBarry Smith 3702593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 371b6490206SBarry Smith count = 0; 372b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 373b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 374b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 375b6490206SBarry Smith } 3762593348eSBarry Smith } 3770752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 378606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3792593348eSBarry Smith 3802593348eSBarry Smith /* store column indices (zero start index) */ 38166963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj); 382b6490206SBarry Smith count = 0; 383b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 384b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 385b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 386b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 387b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3882593348eSBarry Smith } 3892593348eSBarry Smith } 390b6490206SBarry Smith } 391b6490206SBarry Smith } 3920752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 393606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 3942593348eSBarry Smith 3952593348eSBarry Smith /* store nonzero values */ 39666963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 397b6490206SBarry Smith count = 0; 398b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 399b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 400b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 401b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4027e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 403b6490206SBarry Smith } 404b6490206SBarry Smith } 405b6490206SBarry Smith } 406b6490206SBarry Smith } 4070752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 408606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 409ce6f0cecSBarry Smith 410ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 411ce6f0cecSBarry Smith if (file) { 412ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 413ce6f0cecSBarry Smith } 4143a40ed3dSBarry Smith PetscFunctionReturn(0); 4152593348eSBarry Smith } 4162593348eSBarry Smith 4175615d1e5SSatish Balay #undef __FUNC__ 418d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 419b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4202593348eSBarry Smith { 421b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4229df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4232593348eSBarry Smith char *outputname; 4242593348eSBarry Smith 4253a40ed3dSBarry Smith PetscFunctionBegin; 42677ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 427888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 428639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 429d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4300ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 431*6831982aSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 4320ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 433*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 43444cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 43544cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 436*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 43744cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 43844cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 439aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4400ef38995SBarry Smith if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 441*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 442*6831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4430ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 444*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 445*6831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4460ef38995SBarry Smith } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 447*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4480ef38995SBarry Smith } 44944cd7ae7SLois Curfman McInnes #else 4500ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 451*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 4520ef38995SBarry Smith } 45344cd7ae7SLois Curfman McInnes #endif 45444cd7ae7SLois Curfman McInnes } 45544cd7ae7SLois Curfman McInnes } 456*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 45744cd7ae7SLois Curfman McInnes } 45844cd7ae7SLois Curfman McInnes } 459*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 4600ef38995SBarry Smith } else { 461*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 462b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 463b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 464*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 465b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 466b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 467aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 468e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 469*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 470*6831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4710ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 472*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 473*6831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4740ef38995SBarry Smith } else { 475*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 47688685aaeSLois Curfman McInnes } 47788685aaeSLois Curfman McInnes #else 478*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 47988685aaeSLois Curfman McInnes #endif 4802593348eSBarry Smith } 4812593348eSBarry Smith } 482*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 4832593348eSBarry Smith } 4842593348eSBarry Smith } 485*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 486b6490206SBarry Smith } 487*6831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 4883a40ed3dSBarry Smith PetscFunctionReturn(0); 4892593348eSBarry Smith } 4902593348eSBarry Smith 4915615d1e5SSatish Balay #undef __FUNC__ 49277ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 49377ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 4943270192aSSatish Balay { 49577ed5343SBarry Smith Mat A = (Mat) Aa; 4963270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 4977dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 498fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 4993f1db9ecSBarry Smith MatScalar *aa; 50077ed5343SBarry Smith MPI_Comm comm; 50177ed5343SBarry Smith Viewer viewer; 5023270192aSSatish Balay 5033a40ed3dSBarry Smith PetscFunctionBegin; 50477ed5343SBarry Smith /* 50577ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 50677ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 50777ed5343SBarry Smith rest should return immediately. 50877ed5343SBarry Smith */ 50977ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 510d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 51177ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 5123270192aSSatish Balay 51377ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 51477ed5343SBarry Smith 51577ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 51677ed5343SBarry Smith 5173270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5183270192aSSatish Balay color = DRAW_BLUE; 5193270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5203270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5213270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5223270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5233270192aSSatish Balay aa = a->a + j*bs2; 5243270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5253270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5263270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 527433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5283270192aSSatish Balay } 5293270192aSSatish Balay } 5303270192aSSatish Balay } 5313270192aSSatish Balay } 5323270192aSSatish Balay color = DRAW_CYAN; 5333270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5343270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5353270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5363270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5373270192aSSatish Balay aa = a->a + j*bs2; 5383270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5393270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5403270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 541433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5423270192aSSatish Balay } 5433270192aSSatish Balay } 5443270192aSSatish Balay } 5453270192aSSatish Balay } 5463270192aSSatish Balay 5473270192aSSatish Balay color = DRAW_RED; 5483270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5493270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5503270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5513270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5523270192aSSatish Balay aa = a->a + j*bs2; 5533270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5543270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5553270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 556433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5573270192aSSatish Balay } 5583270192aSSatish Balay } 5593270192aSSatish Balay } 5603270192aSSatish Balay } 56177ed5343SBarry Smith PetscFunctionReturn(0); 56277ed5343SBarry Smith } 5633270192aSSatish Balay 56477ed5343SBarry Smith #undef __FUNC__ 56577ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 56677ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 56777ed5343SBarry Smith { 56877ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5697dce120fSSatish Balay int ierr; 5707dce120fSSatish Balay double xl,yl,xr,yr,w,h; 57177ed5343SBarry Smith Draw draw; 57277ed5343SBarry Smith PetscTruth isnull; 5733270192aSSatish Balay 57477ed5343SBarry Smith PetscFunctionBegin; 57577ed5343SBarry Smith 57677ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 57777ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 57877ed5343SBarry Smith 57977ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 58077ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 58177ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5823270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 58377ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 58477ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5853a40ed3dSBarry Smith PetscFunctionReturn(0); 5863270192aSSatish Balay } 5873270192aSSatish Balay 5885615d1e5SSatish Balay #undef __FUNC__ 589d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 590e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5912593348eSBarry Smith { 59219bcc07fSBarry Smith int ierr; 593*6831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 5942593348eSBarry Smith 5953a40ed3dSBarry Smith PetscFunctionBegin; 596*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 597*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 598*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 599*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 6000f5bd95cSBarry Smith if (issocket) { 6017b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 6020f5bd95cSBarry Smith } else if (isascii){ 6033a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6040f5bd95cSBarry Smith } else if (isbinary) { 6053a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6060f5bd95cSBarry Smith } else if (isdraw) { 6073a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6085cd90555SBarry Smith } else { 6090f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 6102593348eSBarry Smith } 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 6122593348eSBarry Smith } 613b6490206SBarry Smith 614cd0e1443SSatish Balay 6155615d1e5SSatish Balay #undef __FUNC__ 6162d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6172d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 618cd0e1443SSatish Balay { 619cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6202d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6212d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6222d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6233f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 624cd0e1443SSatish Balay 6253a40ed3dSBarry Smith PetscFunctionBegin; 6262d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 627cd0e1443SSatish Balay row = im[k]; brow = row/bs; 628a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 629a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6302d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6312c3acbe9SBarry Smith nrow = ailen[brow]; 6322d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 633a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 634a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6352d61bbb3SSatish Balay col = in[l] ; 6362d61bbb3SSatish Balay bcol = col/bs; 6372d61bbb3SSatish Balay cidx = col%bs; 6382d61bbb3SSatish Balay ridx = row%bs; 6392d61bbb3SSatish Balay high = nrow; 6402d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6412d61bbb3SSatish Balay while (high-low > 5) { 642cd0e1443SSatish Balay t = (low+high)/2; 643cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 644cd0e1443SSatish Balay else low = t; 645cd0e1443SSatish Balay } 646cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 647cd0e1443SSatish Balay if (rp[i] > bcol) break; 648cd0e1443SSatish Balay if (rp[i] == bcol) { 6492d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6502d61bbb3SSatish Balay goto finished; 651cd0e1443SSatish Balay } 652cd0e1443SSatish Balay } 6532d61bbb3SSatish Balay *v++ = zero; 6542d61bbb3SSatish Balay finished:; 655cd0e1443SSatish Balay } 656cd0e1443SSatish Balay } 6573a40ed3dSBarry Smith PetscFunctionReturn(0); 658cd0e1443SSatish Balay } 659cd0e1443SSatish Balay 6602d61bbb3SSatish Balay 6615615d1e5SSatish Balay #undef __FUNC__ 66205a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 66392c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 66492c4ed94SBarry Smith { 66592c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6668a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 667dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 668549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6693f1db9ecSBarry Smith Scalar *value = v; 6703f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 67192c4ed94SBarry Smith 6723a40ed3dSBarry Smith PetscFunctionBegin; 6730e324ae4SSatish Balay if (roworiented) { 6740e324ae4SSatish Balay stepval = (n-1)*bs; 6750e324ae4SSatish Balay } else { 6760e324ae4SSatish Balay stepval = (m-1)*bs; 6770e324ae4SSatish Balay } 67892c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 67992c4ed94SBarry Smith row = im[k]; 6805ef9f2a5SBarry Smith if (row < 0) continue; 681aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 682a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 68392c4ed94SBarry Smith #endif 68492c4ed94SBarry Smith rp = aj + ai[row]; 68592c4ed94SBarry Smith ap = aa + bs2*ai[row]; 68692c4ed94SBarry Smith rmax = imax[row]; 68792c4ed94SBarry Smith nrow = ailen[row]; 68892c4ed94SBarry Smith low = 0; 68992c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6905ef9f2a5SBarry Smith if (in[l] < 0) continue; 691aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 692a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 69392c4ed94SBarry Smith #endif 69492c4ed94SBarry Smith col = in[l]; 69592c4ed94SBarry Smith if (roworiented) { 6960e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 6970e324ae4SSatish Balay } else { 6980e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 69992c4ed94SBarry Smith } 70092c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 70192c4ed94SBarry Smith while (high-low > 7) { 70292c4ed94SBarry Smith t = (low+high)/2; 70392c4ed94SBarry Smith if (rp[t] > col) high = t; 70492c4ed94SBarry Smith else low = t; 70592c4ed94SBarry Smith } 70692c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 70792c4ed94SBarry Smith if (rp[i] > col) break; 70892c4ed94SBarry Smith if (rp[i] == col) { 7098a84c255SSatish Balay bap = ap + bs2*i; 7100e324ae4SSatish Balay if (roworiented) { 7118a84c255SSatish Balay if (is == ADD_VALUES) { 712dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 713dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7148a84c255SSatish Balay bap[jj] += *value++; 715dd9472c6SBarry Smith } 716dd9472c6SBarry Smith } 7170e324ae4SSatish Balay } else { 718dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 719dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7200e324ae4SSatish Balay bap[jj] = *value++; 7218a84c255SSatish Balay } 722dd9472c6SBarry Smith } 723dd9472c6SBarry Smith } 7240e324ae4SSatish Balay } else { 7250e324ae4SSatish Balay if (is == ADD_VALUES) { 726dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 727dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7280e324ae4SSatish Balay *bap++ += *value++; 729dd9472c6SBarry Smith } 730dd9472c6SBarry Smith } 7310e324ae4SSatish Balay } else { 732dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 733dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7340e324ae4SSatish Balay *bap++ = *value++; 7350e324ae4SSatish Balay } 7368a84c255SSatish Balay } 737dd9472c6SBarry Smith } 738dd9472c6SBarry Smith } 739f1241b54SBarry Smith goto noinsert2; 74092c4ed94SBarry Smith } 74192c4ed94SBarry Smith } 74289280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 743a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 74492c4ed94SBarry Smith if (nrow >= rmax) { 74592c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 74692c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7473f1db9ecSBarry Smith MatScalar *new_a; 74892c4ed94SBarry Smith 749a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 75089280ab3SLois Curfman McInnes 75192c4ed94SBarry Smith /* malloc new storage space */ 7523f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7533f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 75492c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 75592c4ed94SBarry Smith new_i = new_j + new_nz; 75692c4ed94SBarry Smith 75792c4ed94SBarry Smith /* copy over old data into new slots */ 75892c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 75992c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 760549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 76192c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 762549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 763549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 764549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 765549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 76692c4ed94SBarry Smith /* free up old matrix storage */ 767606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 768606d414cSSatish Balay if (!a->singlemalloc) { 769606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 770606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 771606d414cSSatish Balay } 77292c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 77392c4ed94SBarry Smith a->singlemalloc = 1; 77492c4ed94SBarry Smith 77592c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 77692c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7773f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 77892c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 77992c4ed94SBarry Smith a->reallocs++; 78092c4ed94SBarry Smith a->nz++; 78192c4ed94SBarry Smith } 78292c4ed94SBarry Smith N = nrow++ - 1; 78392c4ed94SBarry Smith /* shift up all the later entries in this row */ 78492c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 78592c4ed94SBarry Smith rp[ii+1] = rp[ii]; 786549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 78792c4ed94SBarry Smith } 788549d3d68SSatish Balay if (N >= i) { 789549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 790549d3d68SSatish Balay } 79192c4ed94SBarry Smith rp[i] = col; 7928a84c255SSatish Balay bap = ap + bs2*i; 7930e324ae4SSatish Balay if (roworiented) { 794dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 795dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7960e324ae4SSatish Balay bap[jj] = *value++; 797dd9472c6SBarry Smith } 798dd9472c6SBarry Smith } 7990e324ae4SSatish Balay } else { 800dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 801dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 8020e324ae4SSatish Balay *bap++ = *value++; 8030e324ae4SSatish Balay } 804dd9472c6SBarry Smith } 805dd9472c6SBarry Smith } 806f1241b54SBarry Smith noinsert2:; 80792c4ed94SBarry Smith low = i; 80892c4ed94SBarry Smith } 80992c4ed94SBarry Smith ailen[row] = nrow; 81092c4ed94SBarry Smith } 8113a40ed3dSBarry Smith PetscFunctionReturn(0); 81292c4ed94SBarry Smith } 81392c4ed94SBarry Smith 814f2501298SSatish Balay 8155615d1e5SSatish Balay #undef __FUNC__ 8165615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8178f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 818584200bdSSatish Balay { 819584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 820584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 821a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 822549d3d68SSatish Balay int mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr; 8233f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 824584200bdSSatish Balay 8253a40ed3dSBarry Smith PetscFunctionBegin; 8263a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 827584200bdSSatish Balay 82843ee02c3SBarry Smith if (m) rmax = ailen[0]; 829584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 830584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 831584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 832d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 833584200bdSSatish Balay if (fshift) { 834a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 835584200bdSSatish Balay N = ailen[i]; 836584200bdSSatish Balay for ( j=0; j<N; j++ ) { 837584200bdSSatish Balay ip[j-fshift] = ip[j]; 838549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 839584200bdSSatish Balay } 840584200bdSSatish Balay } 841584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 842584200bdSSatish Balay } 843584200bdSSatish Balay if (mbs) { 844584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 845584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 846584200bdSSatish Balay } 847584200bdSSatish Balay /* reset ilen and imax for each row */ 848584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 849584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 850584200bdSSatish Balay } 851a7c10996SSatish Balay a->nz = ai[mbs]; 852584200bdSSatish Balay 853584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 854584200bdSSatish Balay if (fshift && a->diag) { 855606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 856584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 857584200bdSSatish Balay a->diag = 0; 858584200bdSSatish Balay } 8593d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 860219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8613d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 862584200bdSSatish Balay a->reallocs); 863d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 864e2f3b5e9SSatish Balay a->reallocs = 0; 8654e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8664e220ebcSLois Curfman McInnes 8673a40ed3dSBarry Smith PetscFunctionReturn(0); 868584200bdSSatish Balay } 869584200bdSSatish Balay 8705a838052SSatish Balay 871bea157c4SSatish Balay 872bea157c4SSatish Balay /* 873bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 874bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 875bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 876bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 877bea157c4SSatish Balay */ 8785615d1e5SSatish Balay #undef __FUNC__ 879bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks" 880bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 881d9b7c43dSSatish Balay { 882bea157c4SSatish Balay int i,j,k,row; 883bea157c4SSatish Balay PetscTruth flg; 8843a40ed3dSBarry Smith 885433994e6SBarry Smith PetscFunctionBegin; 886bea157c4SSatish Balay for ( i=0,j=0; i<n; j++ ) { 887bea157c4SSatish Balay row = idx[i]; 888bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 889bea157c4SSatish Balay sizes[j] = 1; 890bea157c4SSatish Balay i++; 891e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 892bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 893bea157c4SSatish Balay i++; 894bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 895bea157c4SSatish Balay flg = PETSC_TRUE; 896bea157c4SSatish Balay for ( k=1; k<bs; k++ ) { 897bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 898bea157c4SSatish Balay flg = PETSC_FALSE; 899bea157c4SSatish Balay break; 900d9b7c43dSSatish Balay } 901bea157c4SSatish Balay } 902bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 903bea157c4SSatish Balay sizes[j] = bs; 904bea157c4SSatish Balay i+= bs; 905bea157c4SSatish Balay } else { 906bea157c4SSatish Balay sizes[j] = 1; 907bea157c4SSatish Balay i++; 908bea157c4SSatish Balay } 909bea157c4SSatish Balay } 910bea157c4SSatish Balay } 911bea157c4SSatish Balay *bs_max = j; 9123a40ed3dSBarry Smith PetscFunctionReturn(0); 913d9b7c43dSSatish Balay } 914d9b7c43dSSatish Balay 9155615d1e5SSatish Balay #undef __FUNC__ 9165615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 9178f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 918d9b7c43dSSatish Balay { 919d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 920b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 921bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9223f1db9ecSBarry Smith Scalar zero = 0.0; 9233f1db9ecSBarry Smith MatScalar *aa; 924d9b7c43dSSatish Balay 9253a40ed3dSBarry Smith PetscFunctionBegin; 926d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 927d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 928d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 929d9b7c43dSSatish Balay 930bea157c4SSatish Balay /* allocate memory for rows,sizes */ 931bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 932bea157c4SSatish Balay sizes = rows + is_n; 933bea157c4SSatish Balay 934bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 935bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 936bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 937bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 938b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 939bea157c4SSatish Balay 940bea157c4SSatish Balay for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) { 941bea157c4SSatish Balay row = rows[j]; 942b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 943bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 944bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 945bea157c4SSatish Balay if (sizes[i] == bs) { 946bea157c4SSatish Balay if (diag) { 947bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 948bea157c4SSatish Balay baij->ilen[row/bs] = 1; 949bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 950bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 951a07cd24cSSatish Balay } 952bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 953bea157c4SSatish Balay for ( k=0; k<bs; k++ ) { 954bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 955bea157c4SSatish Balay } 956bea157c4SSatish Balay } else { /* (!diag) */ 957bea157c4SSatish Balay baij->ilen[row/bs] = 0; 958bea157c4SSatish Balay } /* end (!diag) */ 959bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 960aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 961bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 962bea157c4SSatish Balay #endif 963bea157c4SSatish Balay for ( k=0; k<count; k++ ) { 964d9b7c43dSSatish Balay aa[0] = zero; 965d9b7c43dSSatish Balay aa+=bs; 966d9b7c43dSSatish Balay } 967d9b7c43dSSatish Balay if (diag) { 968f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 969d9b7c43dSSatish Balay } 970d9b7c43dSSatish Balay } 971bea157c4SSatish Balay } 972bea157c4SSatish Balay 973606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 9749a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9753a40ed3dSBarry Smith PetscFunctionReturn(0); 976d9b7c43dSSatish Balay } 9771c351548SSatish Balay 9785615d1e5SSatish Balay #undef __FUNC__ 9792d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9802d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9812d61bbb3SSatish Balay { 9822d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9832d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9842d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9852d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 986549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 9873f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9882d61bbb3SSatish Balay 9892d61bbb3SSatish Balay PetscFunctionBegin; 9902d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9912d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9925ef9f2a5SBarry Smith if (row < 0) continue; 993aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 99432e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 9952d61bbb3SSatish Balay #endif 9962d61bbb3SSatish Balay rp = aj + ai[brow]; 9972d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 9982d61bbb3SSatish Balay rmax = imax[brow]; 9992d61bbb3SSatish Balay nrow = ailen[brow]; 10002d61bbb3SSatish Balay low = 0; 10012d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 10025ef9f2a5SBarry Smith if (in[l] < 0) continue; 1003aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 100432e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 10052d61bbb3SSatish Balay #endif 10062d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10072d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10082d61bbb3SSatish Balay if (roworiented) { 10095ef9f2a5SBarry Smith value = v[l + k*n]; 10102d61bbb3SSatish Balay } else { 10112d61bbb3SSatish Balay value = v[k + l*m]; 10122d61bbb3SSatish Balay } 10132d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10142d61bbb3SSatish Balay while (high-low > 7) { 10152d61bbb3SSatish Balay t = (low+high)/2; 10162d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10172d61bbb3SSatish Balay else low = t; 10182d61bbb3SSatish Balay } 10192d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 10202d61bbb3SSatish Balay if (rp[i] > bcol) break; 10212d61bbb3SSatish Balay if (rp[i] == bcol) { 10222d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10232d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10242d61bbb3SSatish Balay else *bap = value; 10252d61bbb3SSatish Balay goto noinsert1; 10262d61bbb3SSatish Balay } 10272d61bbb3SSatish Balay } 10282d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10292d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10302d61bbb3SSatish Balay if (nrow >= rmax) { 10312d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10322d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10333f1db9ecSBarry Smith MatScalar *new_a; 10342d61bbb3SSatish Balay 10352d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10362d61bbb3SSatish Balay 10372d61bbb3SSatish Balay /* Malloc new storage space */ 10383f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10393f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 10402d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 10412d61bbb3SSatish Balay new_i = new_j + new_nz; 10422d61bbb3SSatish Balay 10432d61bbb3SSatish Balay /* copy over old data into new slots */ 10442d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10452d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1046549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10472d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1048549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1049549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1050549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1051549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10522d61bbb3SSatish Balay /* free up old matrix storage */ 1053606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1054606d414cSSatish Balay if (!a->singlemalloc) { 1055606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1056606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1057606d414cSSatish Balay } 10582d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 10592d61bbb3SSatish Balay a->singlemalloc = 1; 10602d61bbb3SSatish Balay 10612d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10622d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10633f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10642d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10652d61bbb3SSatish Balay a->reallocs++; 10662d61bbb3SSatish Balay a->nz++; 10672d61bbb3SSatish Balay } 10682d61bbb3SSatish Balay N = nrow++ - 1; 10692d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10702d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10712d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1072549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10732d61bbb3SSatish Balay } 1074549d3d68SSatish Balay if (N>=i) { 1075549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1076549d3d68SSatish Balay } 10772d61bbb3SSatish Balay rp[i] = bcol; 10782d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10792d61bbb3SSatish Balay noinsert1:; 10802d61bbb3SSatish Balay low = i; 10812d61bbb3SSatish Balay } 10822d61bbb3SSatish Balay ailen[brow] = nrow; 10832d61bbb3SSatish Balay } 10842d61bbb3SSatish Balay PetscFunctionReturn(0); 10852d61bbb3SSatish Balay } 10862d61bbb3SSatish Balay 10872d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10882d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10892d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 10907b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 10917b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 10922d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10932d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 10942d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 10952d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 10962d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 10972d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 10982d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 10992d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 11002d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 11012d61bbb3SSatish Balay 11022d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 11032d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 11042d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 11052d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 11062d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 11072d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 110815091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 11092d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 11102d61bbb3SSatish Balay 11112d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 11122d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 11132d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 11142d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 11152d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 11162d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 111715091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11182d61bbb3SSatish Balay 11192d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11202d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11212d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11222d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11232d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 112415091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11252d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11262d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11272d61bbb3SSatish Balay 11282d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11292d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11302d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11312d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11322d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 113315091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11342d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11352d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11362d61bbb3SSatish Balay 11372d61bbb3SSatish Balay #undef __FUNC__ 11382d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 11395ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11402d61bbb3SSatish Balay { 11412d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 11422d61bbb3SSatish Balay Mat outA; 11432d61bbb3SSatish Balay int ierr; 1144667159a5SBarry Smith PetscTruth row_identity, col_identity; 11452d61bbb3SSatish Balay 11462d61bbb3SSatish Balay PetscFunctionBegin; 1147b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1148667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1149667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1150667159a5SBarry Smith if (!row_identity || !col_identity) { 1151b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1152667159a5SBarry Smith } 11532d61bbb3SSatish Balay 11542d61bbb3SSatish Balay outA = inA; 11552d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11562d61bbb3SSatish Balay a->row = row; 11572d61bbb3SSatish Balay a->col = col; 11582d61bbb3SSatish Balay 1159e51c0b9cSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1160e51c0b9cSSatish Balay ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr); 11611e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1162e51c0b9cSSatish Balay 11632d61bbb3SSatish Balay if (!a->solve_work) { 11642d61bbb3SSatish Balay a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 11652d61bbb3SSatish Balay PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 11662d61bbb3SSatish Balay } 11672d61bbb3SSatish Balay 11682d61bbb3SSatish Balay if (!a->diag) { 11692d61bbb3SSatish Balay ierr = MatMarkDiag_SeqBAIJ(inA);CHKERRQ(ierr); 11702d61bbb3SSatish Balay } 11712d61bbb3SSatish Balay /* 117215091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 117315091d37SBarry Smith for ILU(0) factorization with natural ordering 11742d61bbb3SSatish Balay */ 117515091d37SBarry Smith switch (a->bs) { 117615091d37SBarry Smith case 2: 117715091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 117815091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 117915091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 118015091d37SBarry Smith break; 118115091d37SBarry Smith case 3: 118215091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 118315091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 118415091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 118515091d37SBarry Smith break; 118615091d37SBarry Smith case 4: 1187667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1188f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 118915091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 119015091d37SBarry Smith break; 119115091d37SBarry Smith case 5: 1192667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1193667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 119415091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 119515091d37SBarry Smith break; 119615091d37SBarry Smith case 6: 119715091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 119815091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 119915091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 120015091d37SBarry Smith break; 120115091d37SBarry Smith case 7: 120215091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 120315091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 120415091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 120515091d37SBarry Smith break; 12062d61bbb3SSatish Balay } 12072d61bbb3SSatish Balay 1208667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1209667159a5SBarry Smith 12102d61bbb3SSatish Balay PetscFunctionReturn(0); 12112d61bbb3SSatish Balay } 12122d61bbb3SSatish Balay #undef __FUNC__ 1213d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1214ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1215ba4ca20aSSatish Balay { 1216ba4ca20aSSatish Balay static int called = 0; 1217ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1218d132466eSBarry Smith int ierr; 1219ba4ca20aSSatish Balay 12203a40ed3dSBarry Smith PetscFunctionBegin; 12213a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 1222d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1223d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12243a40ed3dSBarry Smith PetscFunctionReturn(0); 1225ba4ca20aSSatish Balay } 1226d9b7c43dSSatish Balay 1227fb2e594dSBarry Smith EXTERN_C_BEGIN 122827a8da17SBarry Smith #undef __FUNC__ 122927a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 123027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 123127a8da17SBarry Smith { 123227a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 123327a8da17SBarry Smith int i,nz,n; 123427a8da17SBarry Smith 123527a8da17SBarry Smith PetscFunctionBegin; 123627a8da17SBarry Smith nz = baij->maxnz; 123727a8da17SBarry Smith n = baij->n; 123827a8da17SBarry Smith for (i=0; i<nz; i++) { 123927a8da17SBarry Smith baij->j[i] = indices[i]; 124027a8da17SBarry Smith } 124127a8da17SBarry Smith baij->nz = nz; 124227a8da17SBarry Smith for ( i=0; i<n; i++ ) { 124327a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 124427a8da17SBarry Smith } 124527a8da17SBarry Smith 124627a8da17SBarry Smith PetscFunctionReturn(0); 124727a8da17SBarry Smith } 1248fb2e594dSBarry Smith EXTERN_C_END 124927a8da17SBarry Smith 125027a8da17SBarry Smith #undef __FUNC__ 125127a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 125227a8da17SBarry Smith /*@ 125327a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 125427a8da17SBarry Smith in the matrix. 125527a8da17SBarry Smith 125627a8da17SBarry Smith Input Parameters: 125727a8da17SBarry Smith + mat - the SeqBAIJ matrix 125827a8da17SBarry Smith - indices - the column indices 125927a8da17SBarry Smith 126015091d37SBarry Smith Level: advanced 126115091d37SBarry Smith 126227a8da17SBarry Smith Notes: 126327a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 126427a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 126527a8da17SBarry Smith of the MatSetValues() operation. 126627a8da17SBarry Smith 126727a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 126827a8da17SBarry Smith MatCreateSeqBAIJ(). 126927a8da17SBarry Smith 127027a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 127127a8da17SBarry Smith 127227a8da17SBarry Smith @*/ 127327a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 127427a8da17SBarry Smith { 127527a8da17SBarry Smith int ierr,(*f)(Mat,int *); 127627a8da17SBarry Smith 127727a8da17SBarry Smith PetscFunctionBegin; 127827a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1279549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 128027a8da17SBarry Smith if (f) { 128127a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 128227a8da17SBarry Smith } else { 128327a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 128427a8da17SBarry Smith } 128527a8da17SBarry Smith PetscFunctionReturn(0); 128627a8da17SBarry Smith } 128727a8da17SBarry Smith 12882593348eSBarry Smith /* -------------------------------------------------------------------*/ 1289cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1290cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1291cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1292cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1293cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1294cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1295cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1296cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1297cc2dc46cSBarry Smith 0, 1298cc2dc46cSBarry Smith 0, 1299cc2dc46cSBarry Smith 0, 1300cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1301cc2dc46cSBarry Smith 0, 1302b6490206SBarry Smith 0, 1303f2501298SSatish Balay MatTranspose_SeqBAIJ, 1304cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1305cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1306cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1307cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1308cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1309b6490206SBarry Smith 0, 1310cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1311cc2dc46cSBarry Smith 0, 1312cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1313cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1314cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1315cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1316cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1317cc2dc46cSBarry Smith 0, 1318cc2dc46cSBarry Smith 0, 1319cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1320cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1321cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1322cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1323cc2dc46cSBarry Smith 0, 1324cc2dc46cSBarry Smith 0, 1325cc2dc46cSBarry Smith 0, 13262e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1327cc2dc46cSBarry Smith 0, 1328cc2dc46cSBarry Smith 0, 1329cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1330cc2dc46cSBarry Smith 0, 1331cc2dc46cSBarry Smith 0, 1332cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1333cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1334cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1335cc2dc46cSBarry Smith 0, 1336cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1337cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1338cc2dc46cSBarry Smith 0, 1339cc2dc46cSBarry Smith 0, 1340cc2dc46cSBarry Smith 0, 1341cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 13423b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 134392c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1344cc2dc46cSBarry Smith 0, 1345cc2dc46cSBarry Smith 0, 1346cc2dc46cSBarry Smith 0, 1347cc2dc46cSBarry Smith 0, 1348cc2dc46cSBarry Smith 0, 1349cc2dc46cSBarry Smith 0, 1350d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1351cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1352cc2dc46cSBarry Smith 0, 1353cc2dc46cSBarry Smith 0, 1354cc2dc46cSBarry Smith MatGetMaps_Petsc}; 13552593348eSBarry Smith 13563e90b805SBarry Smith EXTERN_C_BEGIN 13573e90b805SBarry Smith #undef __FUNC__ 13583e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ" 13593e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 13603e90b805SBarry Smith { 13613e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 13623e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1363d132466eSBarry Smith int ierr; 13643e90b805SBarry Smith 13653e90b805SBarry Smith PetscFunctionBegin; 13663e90b805SBarry Smith if (aij->nonew != 1) { 13673e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13683e90b805SBarry Smith } 13693e90b805SBarry Smith 13703e90b805SBarry Smith /* allocate space for values if not already there */ 13713e90b805SBarry Smith if (!aij->saved_values) { 13723e90b805SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 13733e90b805SBarry Smith } 13743e90b805SBarry Smith 13753e90b805SBarry Smith /* copy values over */ 1376d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 13773e90b805SBarry Smith PetscFunctionReturn(0); 13783e90b805SBarry Smith } 13793e90b805SBarry Smith EXTERN_C_END 13803e90b805SBarry Smith 13813e90b805SBarry Smith EXTERN_C_BEGIN 13823e90b805SBarry Smith #undef __FUNC__ 138332e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ" 138432e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 13853e90b805SBarry Smith { 13863e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1387549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 13883e90b805SBarry Smith 13893e90b805SBarry Smith PetscFunctionBegin; 13903e90b805SBarry Smith if (aij->nonew != 1) { 13913e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13923e90b805SBarry Smith } 13933e90b805SBarry Smith if (!aij->saved_values) { 13943e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 13953e90b805SBarry Smith } 13963e90b805SBarry Smith 13973e90b805SBarry Smith /* copy values over */ 1398549d3d68SSatish Balay ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 13993e90b805SBarry Smith PetscFunctionReturn(0); 14003e90b805SBarry Smith } 14013e90b805SBarry Smith EXTERN_C_END 14023e90b805SBarry Smith 14035615d1e5SSatish Balay #undef __FUNC__ 14045615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 14052593348eSBarry Smith /*@C 140644cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 140744cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 140844cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 14097fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 14102bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 14112593348eSBarry Smith 1412db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1413db81eaa0SLois Curfman McInnes 14142593348eSBarry Smith Input Parameters: 1415db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1416b6490206SBarry Smith . bs - size of block 14172593348eSBarry Smith . m - number of rows 14182593348eSBarry Smith . n - number of columns 1419b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14207fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14212bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14222593348eSBarry Smith 14232593348eSBarry Smith Output Parameter: 14242593348eSBarry Smith . A - the matrix 14252593348eSBarry Smith 14260513a670SBarry Smith Options Database Keys: 1427db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1428db81eaa0SLois Curfman McInnes block calculations (much slower) 1429db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14300513a670SBarry Smith 143115091d37SBarry Smith Level: intermediate 143215091d37SBarry Smith 14332593348eSBarry Smith Notes: 143444cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 14352593348eSBarry Smith storage. That is, the stored row and column indices can begin at 143644cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 14372593348eSBarry Smith 14382593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 14392593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14402593348eSBarry Smith allocation. For additional details, see the users manual chapter on 14416da5968aSLois Curfman McInnes matrices. 14422593348eSBarry Smith 1443db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 14442593348eSBarry Smith @*/ 1445b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 14462593348eSBarry Smith { 14472593348eSBarry Smith Mat B; 1448b6490206SBarry Smith Mat_SeqBAIJ *b; 1449302e33e3SBarry Smith int i,len,ierr,flg,mbs,nbs,bs2,size; 14503b2fbd54SBarry Smith 14513a40ed3dSBarry Smith PetscFunctionBegin; 1452d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1453a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1454b6490206SBarry Smith 1455962fb4a0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1456302e33e3SBarry Smith mbs = m/bs; 1457302e33e3SBarry Smith nbs = n/bs; 1458302e33e3SBarry Smith bs2 = bs*bs; 14590513a670SBarry Smith 14603a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1461a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 14623a40ed3dSBarry Smith } 14632593348eSBarry Smith 1464b73539f3SBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz); 1465b73539f3SBarry Smith if (nnz) { 1466faf3f1fcSBarry Smith for (i=0; i<mbs; i++) { 1467b73539f3SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 1468b73539f3SBarry Smith } 1469b73539f3SBarry Smith } 1470b73539f3SBarry Smith 14712593348eSBarry Smith *A = 0; 14723f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 14732593348eSBarry Smith PLogObjectCreate(B); 1474b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1475549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1476549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 14771a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 14781a6a6d98SBarry Smith if (!flg) { 14797fc0212eSBarry Smith switch (bs) { 14807fc0212eSBarry Smith case 1: 1481f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1482f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1483f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1484f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 14857fc0212eSBarry Smith break; 14864eeb42bcSBarry Smith case 2: 1487f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1488f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1489f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1490f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 14916d84be18SBarry Smith break; 1492f327dd97SBarry Smith case 3: 1493f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1494f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1495f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1496f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 14974eeb42bcSBarry Smith break; 14986d84be18SBarry Smith case 4: 1499f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1500f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1501f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1502f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 15036d84be18SBarry Smith break; 15046d84be18SBarry Smith case 5: 1505f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1506f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1507f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1508f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 15096d84be18SBarry Smith break; 151015091d37SBarry Smith case 6: 151115091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 151215091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 151315091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 151415091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 151515091d37SBarry Smith break; 151648e9ddb2SSatish Balay case 7: 1517f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1518f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1519f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 152048e9ddb2SSatish Balay break; 15217fc0212eSBarry Smith } 15221a6a6d98SBarry Smith } 1523e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1524e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15252593348eSBarry Smith B->factor = 0; 15262593348eSBarry Smith B->lupivotthreshold = 1.0; 152790f02eecSBarry Smith B->mapping = 0; 15282593348eSBarry Smith b->row = 0; 15292593348eSBarry Smith b->col = 0; 1530e51c0b9cSSatish Balay b->icol = 0; 15312593348eSBarry Smith b->reallocs = 0; 15323e90b805SBarry Smith b->saved_values = 0; 15332593348eSBarry Smith 153444cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 153544cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1536a5ae1ecdSBarry Smith 15377413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 15387413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1539a5ae1ecdSBarry Smith 1540b6490206SBarry Smith b->mbs = mbs; 1541f2501298SSatish Balay b->nbs = nbs; 1542b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax); 15432593348eSBarry Smith if (nnz == PETSC_NULL) { 1544119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 15452593348eSBarry Smith else if (nz <= 0) nz = 1; 1546b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1547b6490206SBarry Smith nz = nz*mbs; 15483a40ed3dSBarry Smith } else { 15492593348eSBarry Smith nz = 0; 1550b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 15512593348eSBarry Smith } 15522593348eSBarry Smith 15532593348eSBarry Smith /* allocate the matrix space */ 15543f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 15553f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a); 1556549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15577e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 1558549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15592593348eSBarry Smith b->i = b->j + nz; 15602593348eSBarry Smith b->singlemalloc = 1; 15612593348eSBarry Smith 1562de6a44a3SBarry Smith b->i[0] = 0; 1563b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 15642593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 15652593348eSBarry Smith } 15662593348eSBarry Smith 1567b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1568b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1569f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1570b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 15712593348eSBarry Smith 1572b6490206SBarry Smith b->bs = bs; 15739df24120SSatish Balay b->bs2 = bs2; 1574b6490206SBarry Smith b->mbs = mbs; 15752593348eSBarry Smith b->nz = 0; 157618e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 15772593348eSBarry Smith b->sorted = 0; 15782593348eSBarry Smith b->roworiented = 1; 15792593348eSBarry Smith b->nonew = 0; 15802593348eSBarry Smith b->diag = 0; 15812593348eSBarry Smith b->solve_work = 0; 1582de6a44a3SBarry Smith b->mult_work = 0; 15832593348eSBarry Smith b->spptr = 0; 15844e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 15854e220ebcSLois Curfman McInnes 15862593348eSBarry Smith *A = B; 15872593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 15882593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 158927a8da17SBarry Smith 15903e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 15913e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 15923e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 15933e90b805SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 15943e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 15953e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 159627a8da17SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 159727a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 159827a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 15993a40ed3dSBarry Smith PetscFunctionReturn(0); 16002593348eSBarry Smith } 16012593348eSBarry Smith 16025615d1e5SSatish Balay #undef __FUNC__ 16032e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 16042e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 16052593348eSBarry Smith { 16062593348eSBarry Smith Mat C; 1607b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 160827a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1609de6a44a3SBarry Smith 16103a40ed3dSBarry Smith PetscFunctionBegin; 1611a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 16122593348eSBarry Smith 16132593348eSBarry Smith *B = 0; 16143f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 16152593348eSBarry Smith PLogObjectCreate(C); 1616b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1617549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1618e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1619e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 16202593348eSBarry Smith C->factor = A->factor; 16212593348eSBarry Smith c->row = 0; 16222593348eSBarry Smith c->col = 0; 1623cac97260SSatish Balay c->icol = 0; 162432e87ba3SBarry Smith c->saved_values = 0; 16252593348eSBarry Smith C->assembled = PETSC_TRUE; 16262593348eSBarry Smith 162744cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 162844cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 162944cd7ae7SLois Curfman McInnes C->M = a->m; 163044cd7ae7SLois Curfman McInnes C->N = a->n; 163144cd7ae7SLois Curfman McInnes 1632b6490206SBarry Smith c->bs = a->bs; 16339df24120SSatish Balay c->bs2 = a->bs2; 1634b6490206SBarry Smith c->mbs = a->mbs; 1635b6490206SBarry Smith c->nbs = a->nbs; 16362593348eSBarry Smith 1637b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1638b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1639b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16402593348eSBarry Smith c->imax[i] = a->imax[i]; 16412593348eSBarry Smith c->ilen[i] = a->ilen[i]; 16422593348eSBarry Smith } 16432593348eSBarry Smith 16442593348eSBarry Smith /* allocate the matrix space */ 16452593348eSBarry Smith c->singlemalloc = 1; 16463f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 16473f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a); 16487e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1649de6a44a3SBarry Smith c->i = c->j + nz; 1650549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1651b6490206SBarry Smith if (mbs > 0) { 1652549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16532e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1654549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16552e8a6d31SBarry Smith } else { 1656549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16572593348eSBarry Smith } 16582593348eSBarry Smith } 16592593348eSBarry Smith 1660f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16612593348eSBarry Smith c->sorted = a->sorted; 16622593348eSBarry Smith c->roworiented = a->roworiented; 16632593348eSBarry Smith c->nonew = a->nonew; 16642593348eSBarry Smith 16652593348eSBarry Smith if (a->diag) { 1666b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag); 1667b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1668b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16692593348eSBarry Smith c->diag[i] = a->diag[i]; 16702593348eSBarry Smith } 167198305bb5SBarry Smith } else c->diag = 0; 16722593348eSBarry Smith c->nz = a->nz; 16732593348eSBarry Smith c->maxnz = a->maxnz; 16742593348eSBarry Smith c->solve_work = 0; 16752593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 16767fc0212eSBarry Smith c->mult_work = 0; 16772593348eSBarry Smith *B = C; 16787b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 16793a40ed3dSBarry Smith PetscFunctionReturn(0); 16802593348eSBarry Smith } 16812593348eSBarry Smith 16825615d1e5SSatish Balay #undef __FUNC__ 16835615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 168419bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 16852593348eSBarry Smith { 1686b6490206SBarry Smith Mat_SeqBAIJ *a; 16872593348eSBarry Smith Mat B; 1688de6a44a3SBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1689b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 169035aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1691a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1692b6490206SBarry Smith Scalar *aa; 169319bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 16942593348eSBarry Smith 16953a40ed3dSBarry Smith PetscFunctionBegin; 16961a6a6d98SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1697de6a44a3SBarry Smith bs2 = bs*bs; 1698b6490206SBarry Smith 1699d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1700cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 170190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 17020752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1703a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 17042593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 17052593348eSBarry Smith 1706d64ed03dSBarry Smith if (header[3] < 0) { 1707a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1708d64ed03dSBarry Smith } 1709d64ed03dSBarry Smith 1710a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 171135aab85fSBarry Smith 171235aab85fSBarry Smith /* 171335aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 171435aab85fSBarry Smith divisible by the blocksize 171535aab85fSBarry Smith */ 1716b6490206SBarry Smith mbs = M/bs; 171735aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 171835aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 171935aab85fSBarry Smith else mbs++; 172035aab85fSBarry Smith if (extra_rows) { 1721537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 172235aab85fSBarry Smith } 1723b6490206SBarry Smith 17242593348eSBarry Smith /* read in row lengths */ 172535aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 17260752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 172735aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 17282593348eSBarry Smith 1729b6490206SBarry Smith /* read in column indices */ 173035aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj); 17310752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 173235aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1733b6490206SBarry Smith 1734b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1735b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1736549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 173735aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask); 1738549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 173935aab85fSBarry Smith masked = mask + mbs; 1740b6490206SBarry Smith rowcount = 0; nzcount = 0; 1741b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 174235aab85fSBarry Smith nmask = 0; 1743b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1744b6490206SBarry Smith kmax = rowlengths[rowcount]; 1745b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 174635aab85fSBarry Smith tmp = jj[nzcount++]/bs; 174735aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1748b6490206SBarry Smith } 1749b6490206SBarry Smith rowcount++; 1750b6490206SBarry Smith } 175135aab85fSBarry Smith browlengths[i] += nmask; 175235aab85fSBarry Smith /* zero out the mask elements we set */ 175335aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1754b6490206SBarry Smith } 1755b6490206SBarry Smith 17562593348eSBarry Smith /* create our matrix */ 17573eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17582593348eSBarry Smith B = *A; 1759b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 17602593348eSBarry Smith 17612593348eSBarry Smith /* set matrix "i" values */ 1762de6a44a3SBarry Smith a->i[0] = 0; 1763b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1764b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1765b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17662593348eSBarry Smith } 1767b6490206SBarry Smith a->nz = 0; 1768b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 17692593348eSBarry Smith 1770b6490206SBarry Smith /* read in nonzero values */ 177135aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 17720752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 177335aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1774b6490206SBarry Smith 1775b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1776b6490206SBarry Smith nzcount = 0; jcount = 0; 1777b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1778b6490206SBarry Smith nzcountb = nzcount; 177935aab85fSBarry Smith nmask = 0; 1780b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1781b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1782b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 178335aab85fSBarry Smith tmp = jj[nzcount++]/bs; 178435aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1785b6490206SBarry Smith } 1786b6490206SBarry Smith rowcount++; 1787b6490206SBarry Smith } 1788de6a44a3SBarry Smith /* sort the masked values */ 1789433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1790de6a44a3SBarry Smith 1791b6490206SBarry Smith /* set "j" values into matrix */ 1792b6490206SBarry Smith maskcount = 1; 179335aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 179435aab85fSBarry Smith a->j[jcount++] = masked[j]; 1795de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1796b6490206SBarry Smith } 1797b6490206SBarry Smith /* set "a" values into matrix */ 1798de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1799b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1800b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1801b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1802de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1803de6a44a3SBarry Smith block = mask[tmp] - 1; 1804de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1805de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1806b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1807b6490206SBarry Smith } 1808b6490206SBarry Smith } 180935aab85fSBarry Smith /* zero out the mask elements we set */ 181035aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1811b6490206SBarry Smith } 1812a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1813b6490206SBarry Smith 1814606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1815606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1816606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1817606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1818606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1819b6490206SBarry Smith 1820b6490206SBarry Smith B->assembled = PETSC_TRUE; 1821de6a44a3SBarry Smith 18229c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 18233a40ed3dSBarry Smith PetscFunctionReturn(0); 18242593348eSBarry Smith } 18252593348eSBarry Smith 18262593348eSBarry Smith 18272593348eSBarry Smith 1828