1*c4992f7dSBarry Smith /*$Id: baij.c,v 1.189 1999/11/10 03:19:27 bsmith Exp bsmith $*/ 22593348eSBarry Smith 32593348eSBarry Smith /* 4b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 52593348eSBarry Smith matrix storage format. 62593348eSBarry Smith */ 73369ce9aSBarry Smith #include "sys.h" 870f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 91a6a6d98SBarry Smith #include "src/vec/vecimpl.h" 101a6a6d98SBarry Smith #include "src/inline/spops.h" 113b547af2SSatish Balay 122d61bbb3SSatish Balay #define CHUNKSIZE 10 13de6a44a3SBarry Smith 14be5855fcSBarry Smith /* 15be5855fcSBarry Smith Checks for missing diagonals 16be5855fcSBarry Smith */ 17be5855fcSBarry Smith #undef __FUNC__ 18*c4992f7dSBarry Smith #define __FUNC__ "MatMissingDiagonal_SeqBAIJ" 19*c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A) 20be5855fcSBarry Smith { 21be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 22*c4992f7dSBarry Smith int *diag = a->diag, *jj = a->j,i,ierr; 23be5855fcSBarry Smith 24be5855fcSBarry Smith PetscFunctionBegin; 25*c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 260e8e8aceSBarry Smith for ( i=0; i<a->mbs; i++ ) { 27be5855fcSBarry Smith if (jj[diag[i]] != i) { 28be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 29be5855fcSBarry Smith } 30be5855fcSBarry Smith } 31be5855fcSBarry Smith PetscFunctionReturn(0); 32be5855fcSBarry Smith } 33be5855fcSBarry Smith 345615d1e5SSatish Balay #undef __FUNC__ 35*c4992f7dSBarry Smith #define __FUNC__ "MatMarkDiagonal_SeqBAIJ" 36*c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A) 37de6a44a3SBarry Smith { 38de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 397fc0212eSBarry Smith int i,j, *diag, m = a->mbs; 40de6a44a3SBarry Smith 413a40ed3dSBarry Smith PetscFunctionBegin; 42de6a44a3SBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int));CHKPTRQ(diag); 43de6a44a3SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 447fc0212eSBarry Smith for ( i=0; i<m; i++ ) { 45ecc615b2SBarry Smith diag[i] = a->i[i+1]; 46de6a44a3SBarry Smith for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 47de6a44a3SBarry Smith if (a->j[j] == i) { 48de6a44a3SBarry Smith diag[i] = j; 49de6a44a3SBarry Smith break; 50de6a44a3SBarry Smith } 51de6a44a3SBarry Smith } 52de6a44a3SBarry Smith } 53de6a44a3SBarry Smith a->diag = diag; 543a40ed3dSBarry Smith PetscFunctionReturn(0); 55de6a44a3SBarry Smith } 562593348eSBarry Smith 572593348eSBarry Smith 583b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 593b2fbd54SBarry Smith 605615d1e5SSatish Balay #undef __FUNC__ 615615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ" 623b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 633b2fbd54SBarry Smith PetscTruth *done) 643b2fbd54SBarry Smith { 653b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 663b2fbd54SBarry Smith int ierr, n = a->mbs,i; 673b2fbd54SBarry Smith 683a40ed3dSBarry Smith PetscFunctionBegin; 693b2fbd54SBarry Smith *nn = n; 703a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 713b2fbd54SBarry Smith if (symmetric) { 723b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 733b2fbd54SBarry Smith } else if (oshift == 1) { 743b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 753b2fbd54SBarry Smith int nz = a->i[n] + 1; 763b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 773b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]++; 783b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 793b2fbd54SBarry Smith } else { 803b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 813b2fbd54SBarry Smith } 823b2fbd54SBarry Smith 833a40ed3dSBarry Smith PetscFunctionReturn(0); 843b2fbd54SBarry Smith } 853b2fbd54SBarry Smith 865615d1e5SSatish Balay #undef __FUNC__ 87d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ" 883b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 893b2fbd54SBarry Smith PetscTruth *done) 903b2fbd54SBarry Smith { 913b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 92606d414cSSatish Balay int i,n = a->mbs,ierr; 933b2fbd54SBarry Smith 943a40ed3dSBarry Smith PetscFunctionBegin; 953a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 963b2fbd54SBarry Smith if (symmetric) { 97606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 98606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 99af5da2bfSBarry Smith } else if (oshift == 1) { 1003b2fbd54SBarry Smith int nz = a->i[n]; 1013b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 1023b2fbd54SBarry Smith for ( i=0; i<n+1; i++ ) a->i[i]--; 1033b2fbd54SBarry Smith } 1043a40ed3dSBarry Smith PetscFunctionReturn(0); 1053b2fbd54SBarry Smith } 1063b2fbd54SBarry Smith 1072d61bbb3SSatish Balay #undef __FUNC__ 1082d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ" 1092d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs) 1102d61bbb3SSatish Balay { 1112d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data; 1122d61bbb3SSatish Balay 1132d61bbb3SSatish Balay PetscFunctionBegin; 1142d61bbb3SSatish Balay *bs = baij->bs; 1152d61bbb3SSatish Balay PetscFunctionReturn(0); 1162d61bbb3SSatish Balay } 1172d61bbb3SSatish Balay 1182d61bbb3SSatish Balay 1192d61bbb3SSatish Balay #undef __FUNC__ 1202d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ" 121e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1222d61bbb3SSatish Balay { 1232d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 124e51c0b9cSSatish Balay int ierr; 1252d61bbb3SSatish Balay 126433994e6SBarry Smith PetscFunctionBegin; 12794d884c6SBarry Smith if (A->mapping) { 12894d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr); 12994d884c6SBarry Smith } 13094d884c6SBarry Smith if (A->bmapping) { 13194d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr); 13294d884c6SBarry Smith } 13361b13de0SBarry Smith if (A->rmap) { 13461b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 13561b13de0SBarry Smith } 13661b13de0SBarry Smith if (A->cmap) { 13761b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 13861b13de0SBarry Smith } 139aa482453SBarry Smith #if defined(PETSC_USE_LOG) 140e1311b90SBarry Smith PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 1412d61bbb3SSatish Balay #endif 142606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 143606d414cSSatish Balay if (!a->singlemalloc) { 144606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 145606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 146606d414cSSatish Balay } 147c38d4ed2SBarry Smith if (a->row) { 148c38d4ed2SBarry Smith ierr = ISDestroy(a->row);CHKERRQ(ierr); 149c38d4ed2SBarry Smith } 150c38d4ed2SBarry Smith if (a->col) { 151c38d4ed2SBarry Smith ierr = ISDestroy(a->col);CHKERRQ(ierr); 152c38d4ed2SBarry Smith } 153606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 154606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 155606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 156606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 157606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 158e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 159606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 160606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1612d61bbb3SSatish Balay PLogObjectDestroy(A); 1622d61bbb3SSatish Balay PetscHeaderDestroy(A); 1632d61bbb3SSatish Balay PetscFunctionReturn(0); 1642d61bbb3SSatish Balay } 1652d61bbb3SSatish Balay 1662d61bbb3SSatish Balay #undef __FUNC__ 1672d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ" 1682d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1692d61bbb3SSatish Balay { 1702d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1712d61bbb3SSatish Balay 1722d61bbb3SSatish Balay PetscFunctionBegin; 173*c4992f7dSBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = PETSC_TRUE; 174*c4992f7dSBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = PETSC_FALSE; 175*c4992f7dSBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = PETSC_TRUE; 176*c4992f7dSBarry Smith else if (op == MAT_COLUMNS_UNSORTED) a->sorted = PETSC_FALSE; 177*c4992f7dSBarry Smith else if (op == MAT_KEEP_ZEROED_ROWS) a->keepzeroedrows = PETSC_TRUE; 1782d61bbb3SSatish Balay else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 1794787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 1804787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 1812d61bbb3SSatish Balay else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 1822d61bbb3SSatish Balay else if (op == MAT_ROWS_SORTED || 1832d61bbb3SSatish Balay op == MAT_ROWS_UNSORTED || 1842d61bbb3SSatish Balay op == MAT_SYMMETRIC || 1852d61bbb3SSatish Balay op == MAT_STRUCTURALLY_SYMMETRIC || 1862d61bbb3SSatish Balay op == MAT_YES_NEW_DIAGONALS || 187b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 188b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) { 189981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 1902d61bbb3SSatish Balay } else if (op == MAT_NO_NEW_DIAGONALS) { 1912d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 1922d61bbb3SSatish Balay } else { 1932d61bbb3SSatish Balay SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1942d61bbb3SSatish Balay } 1952d61bbb3SSatish Balay PetscFunctionReturn(0); 1962d61bbb3SSatish Balay } 1972d61bbb3SSatish Balay 1982d61bbb3SSatish Balay 1992d61bbb3SSatish Balay #undef __FUNC__ 2002d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ" 2012d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 2022d61bbb3SSatish Balay { 2032d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2042d61bbb3SSatish Balay 2052d61bbb3SSatish Balay PetscFunctionBegin; 2062d61bbb3SSatish Balay if (m) *m = a->m; 2072d61bbb3SSatish Balay if (n) *n = a->n; 2082d61bbb3SSatish Balay PetscFunctionReturn(0); 2092d61bbb3SSatish Balay } 2102d61bbb3SSatish Balay 2112d61bbb3SSatish Balay #undef __FUNC__ 2122d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ" 2132d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2142d61bbb3SSatish Balay { 2152d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2162d61bbb3SSatish Balay 2172d61bbb3SSatish Balay PetscFunctionBegin; 2182d61bbb3SSatish Balay *m = 0; *n = a->m; 2192d61bbb3SSatish Balay PetscFunctionReturn(0); 2202d61bbb3SSatish Balay } 2212d61bbb3SSatish Balay 2222d61bbb3SSatish Balay #undef __FUNC__ 2232d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ" 2242d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2252d61bbb3SSatish Balay { 2262d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 2272d61bbb3SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 2283f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2293f1db9ecSBarry Smith Scalar *v_i; 2302d61bbb3SSatish Balay 2312d61bbb3SSatish Balay PetscFunctionBegin; 2322d61bbb3SSatish Balay bs = a->bs; 2332d61bbb3SSatish Balay ai = a->i; 2342d61bbb3SSatish Balay aj = a->j; 2352d61bbb3SSatish Balay aa = a->a; 2362d61bbb3SSatish Balay bs2 = a->bs2; 2372d61bbb3SSatish Balay 2382d61bbb3SSatish Balay if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 2392d61bbb3SSatish Balay 2402d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2412d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2422d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2432d61bbb3SSatish Balay *nz = bs*M; 2442d61bbb3SSatish Balay 2452d61bbb3SSatish Balay if (v) { 2462d61bbb3SSatish Balay *v = 0; 2472d61bbb3SSatish Balay if (*nz) { 2482d61bbb3SSatish Balay *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) );CHKPTRQ(*v); 2492d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2502d61bbb3SSatish Balay v_i = *v + i*bs; 2512d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2522d61bbb3SSatish Balay for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 2532d61bbb3SSatish Balay } 2542d61bbb3SSatish Balay } 2552d61bbb3SSatish Balay } 2562d61bbb3SSatish Balay 2572d61bbb3SSatish Balay if (idx) { 2582d61bbb3SSatish Balay *idx = 0; 2592d61bbb3SSatish Balay if (*nz) { 2602d61bbb3SSatish Balay *idx = (int *) PetscMalloc( (*nz)*sizeof(int) );CHKPTRQ(*idx); 2612d61bbb3SSatish Balay for ( i=0; i<M; i++ ) { /* for each block in the block row */ 2622d61bbb3SSatish Balay idx_i = *idx + i*bs; 2632d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2642d61bbb3SSatish Balay for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 2652d61bbb3SSatish Balay } 2662d61bbb3SSatish Balay } 2672d61bbb3SSatish Balay } 2682d61bbb3SSatish Balay PetscFunctionReturn(0); 2692d61bbb3SSatish Balay } 2702d61bbb3SSatish Balay 2712d61bbb3SSatish Balay #undef __FUNC__ 2722d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ" 2732d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2742d61bbb3SSatish Balay { 275606d414cSSatish Balay int ierr; 276606d414cSSatish Balay 2772d61bbb3SSatish Balay PetscFunctionBegin; 278606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 279606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 2802d61bbb3SSatish Balay PetscFunctionReturn(0); 2812d61bbb3SSatish Balay } 2822d61bbb3SSatish Balay 2832d61bbb3SSatish Balay #undef __FUNC__ 2842d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ" 2852d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 2862d61bbb3SSatish Balay { 2872d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 2882d61bbb3SSatish Balay Mat C; 2892d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 2902d61bbb3SSatish Balay int *rows,*cols,bs2=a->bs2; 2913f1db9ecSBarry Smith MatScalar *array = a->a; 2922d61bbb3SSatish Balay 2932d61bbb3SSatish Balay PetscFunctionBegin; 2942d61bbb3SSatish Balay if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place"); 2952d61bbb3SSatish Balay col = (int *) PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col); 296549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 2972d61bbb3SSatish Balay 2982d61bbb3SSatish Balay for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 2992d61bbb3SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 300606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 3012d61bbb3SSatish Balay rows = (int *) PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows); 3022d61bbb3SSatish Balay cols = rows + bs; 3032d61bbb3SSatish Balay for ( i=0; i<mbs; i++ ) { 3042d61bbb3SSatish Balay cols[0] = i*bs; 3052d61bbb3SSatish Balay for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 3062d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3072d61bbb3SSatish Balay for ( j=0; j<len; j++ ) { 3082d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3092d61bbb3SSatish Balay for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 3102d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3112d61bbb3SSatish Balay array += bs2; 3122d61bbb3SSatish Balay } 3132d61bbb3SSatish Balay } 314606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 3152d61bbb3SSatish Balay 3162d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3172d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3182d61bbb3SSatish Balay 319*c4992f7dSBarry Smith if (B) { 3202d61bbb3SSatish Balay *B = C; 3212d61bbb3SSatish Balay } else { 322f830108cSBarry Smith PetscOps *Abops; 323cc2dc46cSBarry Smith MatOps Aops; 324f830108cSBarry Smith 3252d61bbb3SSatish Balay /* This isn't really an in-place transpose */ 326606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 327606d414cSSatish Balay if (!a->singlemalloc) { 328606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 329606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 330606d414cSSatish Balay } 331606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 332606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 333606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 334606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 335606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 336f830108cSBarry Smith 3377413bad6SBarry Smith 3387413bad6SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 3397413bad6SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 3407413bad6SBarry Smith 341f830108cSBarry Smith /* 342f830108cSBarry Smith This is horrible, horrible code. We need to keep the 343f830108cSBarry Smith A pointers for the bops and ops but copy everything 344f830108cSBarry Smith else from C. 345f830108cSBarry Smith */ 346f830108cSBarry Smith Abops = A->bops; 347f830108cSBarry Smith Aops = A->ops; 348549d3d68SSatish Balay ierr = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr); 349f830108cSBarry Smith A->bops = Abops; 350f830108cSBarry Smith A->ops = Aops; 35127a8da17SBarry Smith A->qlist = 0; 352f830108cSBarry Smith 3532d61bbb3SSatish Balay PetscHeaderDestroy(C); 3542d61bbb3SSatish Balay } 3552d61bbb3SSatish Balay PetscFunctionReturn(0); 3562d61bbb3SSatish Balay } 3572d61bbb3SSatish Balay 3585615d1e5SSatish Balay #undef __FUNC__ 359d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary" 360b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 3612593348eSBarry Smith { 362b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 3639df24120SSatish Balay int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 364b6490206SBarry Smith Scalar *aa; 365ce6f0cecSBarry Smith FILE *file; 3662593348eSBarry Smith 3673a40ed3dSBarry Smith PetscFunctionBegin; 36890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 3692593348eSBarry Smith col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 3702593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3713b2fbd54SBarry Smith 3722593348eSBarry Smith col_lens[1] = a->m; 3732593348eSBarry Smith col_lens[2] = a->n; 3747e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3752593348eSBarry Smith 3762593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 377b6490206SBarry Smith count = 0; 378b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 379b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 380b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 381b6490206SBarry Smith } 3822593348eSBarry Smith } 3830752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr); 384606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3852593348eSBarry Smith 3862593348eSBarry Smith /* store column indices (zero start index) */ 38766963ce1SSatish Balay jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) );CHKPTRQ(jj); 388b6490206SBarry Smith count = 0; 389b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 390b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 391b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 392b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 393b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 3942593348eSBarry Smith } 3952593348eSBarry Smith } 396b6490206SBarry Smith } 397b6490206SBarry Smith } 3980752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 399606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 4002593348eSBarry Smith 4012593348eSBarry Smith /* store nonzero values */ 40266963ce1SSatish Balay aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa); 403b6490206SBarry Smith count = 0; 404b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 405b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 406b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 407b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 4087e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 409b6490206SBarry Smith } 410b6490206SBarry Smith } 411b6490206SBarry Smith } 412b6490206SBarry Smith } 4130752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 414606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 415ce6f0cecSBarry Smith 416ce6f0cecSBarry Smith ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 417ce6f0cecSBarry Smith if (file) { 418ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 419ce6f0cecSBarry Smith } 4203a40ed3dSBarry Smith PetscFunctionReturn(0); 4212593348eSBarry Smith } 4222593348eSBarry Smith 4235615d1e5SSatish Balay #undef __FUNC__ 424d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII" 425b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 4262593348eSBarry Smith { 427b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 4289df24120SSatish Balay int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 4292593348eSBarry Smith char *outputname; 4302593348eSBarry Smith 4313a40ed3dSBarry Smith PetscFunctionBegin; 43277ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 433888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 434639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 435d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 4360ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 4376831982aSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported"); 4380ef38995SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 4396831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 44044cd7ae7SLois Curfman McInnes for ( i=0; i<a->mbs; i++ ) { 44144cd7ae7SLois Curfman McInnes for ( j=0; j<bs; j++ ) { 4426831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 44344cd7ae7SLois Curfman McInnes for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 44444cd7ae7SLois Curfman McInnes for ( l=0; l<bs; l++ ) { 445aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4460ef38995SBarry Smith if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 4476831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 4486831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4490ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 4506831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 4516831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4520ef38995SBarry Smith } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) { 4536831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4540ef38995SBarry Smith } 45544cd7ae7SLois Curfman McInnes #else 4560ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 4576831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 4580ef38995SBarry Smith } 45944cd7ae7SLois Curfman McInnes #endif 46044cd7ae7SLois Curfman McInnes } 46144cd7ae7SLois Curfman McInnes } 4626831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 46344cd7ae7SLois Curfman McInnes } 46444cd7ae7SLois Curfman McInnes } 4656831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 4660ef38995SBarry Smith } else { 4676831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 468b6490206SBarry Smith for ( i=0; i<a->mbs; i++ ) { 469b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 4706831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 471b6490206SBarry Smith for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 472b6490206SBarry Smith for ( l=0; l<bs; l++ ) { 473aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 474e20fef11SSatish Balay if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) { 4756831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 4766831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4770ef38995SBarry Smith } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) { 4786831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 4796831982aSBarry Smith PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4800ef38995SBarry Smith } else { 4816831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 48288685aaeSLois Curfman McInnes } 48388685aaeSLois Curfman McInnes #else 4846831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 48588685aaeSLois Curfman McInnes #endif 4862593348eSBarry Smith } 4872593348eSBarry Smith } 4886831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 4892593348eSBarry Smith } 4902593348eSBarry Smith } 4916831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 492b6490206SBarry Smith } 4936831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 4943a40ed3dSBarry Smith PetscFunctionReturn(0); 4952593348eSBarry Smith } 4962593348eSBarry Smith 4975615d1e5SSatish Balay #undef __FUNC__ 49877ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom" 49977ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa) 5003270192aSSatish Balay { 50177ed5343SBarry Smith Mat A = (Mat) Aa; 5023270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5037dce120fSSatish Balay int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank; 504fae8c767SSatish Balay double xl,yl,xr,yr,x_l,x_r,y_l,y_r; 5053f1db9ecSBarry Smith MatScalar *aa; 50677ed5343SBarry Smith MPI_Comm comm; 50777ed5343SBarry Smith Viewer viewer; 5083270192aSSatish Balay 5093a40ed3dSBarry Smith PetscFunctionBegin; 51077ed5343SBarry Smith /* 51177ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 51277ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 51377ed5343SBarry Smith rest should return immediately. 51477ed5343SBarry Smith */ 51577ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 516d132466eSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 51777ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 5183270192aSSatish Balay 51977ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 52077ed5343SBarry Smith 52177ed5343SBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 52277ed5343SBarry Smith 5233270192aSSatish Balay /* loop over matrix elements drawing boxes */ 5243270192aSSatish Balay color = DRAW_BLUE; 5253270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5263270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5273270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5283270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5293270192aSSatish Balay aa = a->a + j*bs2; 5303270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5313270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5323270192aSSatish Balay if (PetscReal(*aa++) >= 0.) continue; 533433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5343270192aSSatish Balay } 5353270192aSSatish Balay } 5363270192aSSatish Balay } 5373270192aSSatish Balay } 5383270192aSSatish Balay color = DRAW_CYAN; 5393270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5403270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5413270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5423270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5433270192aSSatish Balay aa = a->a + j*bs2; 5443270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5453270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5463270192aSSatish Balay if (PetscReal(*aa++) != 0.) continue; 547433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5483270192aSSatish Balay } 5493270192aSSatish Balay } 5503270192aSSatish Balay } 5513270192aSSatish Balay } 5523270192aSSatish Balay 5533270192aSSatish Balay color = DRAW_RED; 5543270192aSSatish Balay for ( i=0,row=0; i<mbs; i++,row+=bs ) { 5553270192aSSatish Balay for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 5563270192aSSatish Balay y_l = a->m - row - 1.0; y_r = y_l + 1.0; 5573270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5583270192aSSatish Balay aa = a->a + j*bs2; 5593270192aSSatish Balay for ( k=0; k<bs; k++ ) { 5603270192aSSatish Balay for ( l=0; l<bs; l++ ) { 5613270192aSSatish Balay if (PetscReal(*aa++) <= 0.) continue; 562433994e6SBarry Smith ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5633270192aSSatish Balay } 5643270192aSSatish Balay } 5653270192aSSatish Balay } 5663270192aSSatish Balay } 56777ed5343SBarry Smith PetscFunctionReturn(0); 56877ed5343SBarry Smith } 5693270192aSSatish Balay 57077ed5343SBarry Smith #undef __FUNC__ 57177ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw" 57277ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer) 57377ed5343SBarry Smith { 57477ed5343SBarry Smith Mat_SeqBAIJ *a=(Mat_SeqBAIJ *) A->data; 5757dce120fSSatish Balay int ierr; 5767dce120fSSatish Balay double xl,yl,xr,yr,w,h; 57777ed5343SBarry Smith Draw draw; 57877ed5343SBarry Smith PetscTruth isnull; 5793270192aSSatish Balay 58077ed5343SBarry Smith PetscFunctionBegin; 58177ed5343SBarry Smith 58277ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 58377ed5343SBarry Smith ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 58477ed5343SBarry Smith 58577ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 58677ed5343SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 58777ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 5883270192aSSatish Balay ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 58977ed5343SBarry Smith ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 59077ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5913a40ed3dSBarry Smith PetscFunctionReturn(0); 5923270192aSSatish Balay } 5933270192aSSatish Balay 5945615d1e5SSatish Balay #undef __FUNC__ 595d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ" 596e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer) 5972593348eSBarry Smith { 59819bcc07fSBarry Smith int ierr; 5996831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 6002593348eSBarry Smith 6013a40ed3dSBarry Smith PetscFunctionBegin; 6026831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 6036831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 6046831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 6056831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr); 6060f5bd95cSBarry Smith if (issocket) { 6077b2a1423SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported"); 6080f5bd95cSBarry Smith } else if (isascii){ 6093a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6100f5bd95cSBarry Smith } else if (isbinary) { 6113a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6120f5bd95cSBarry Smith } else if (isdraw) { 6133a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6145cd90555SBarry Smith } else { 6150f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 6162593348eSBarry Smith } 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 6182593348eSBarry Smith } 619b6490206SBarry Smith 620cd0e1443SSatish Balay 6215615d1e5SSatish Balay #undef __FUNC__ 6222d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ" 6232d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 624cd0e1443SSatish Balay { 625cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6262d61bbb3SSatish Balay int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 6272d61bbb3SSatish Balay int *ai = a->i, *ailen = a->ilen; 6282d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6293f1db9ecSBarry Smith MatScalar *ap, *aa = a->a, zero = 0.0; 630cd0e1443SSatish Balay 6313a40ed3dSBarry Smith PetscFunctionBegin; 6322d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over rows */ 633cd0e1443SSatish Balay row = im[k]; brow = row/bs; 634a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 635a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 6362d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6372c3acbe9SBarry Smith nrow = ailen[brow]; 6382d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over columns */ 639a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 640a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 6412d61bbb3SSatish Balay col = in[l] ; 6422d61bbb3SSatish Balay bcol = col/bs; 6432d61bbb3SSatish Balay cidx = col%bs; 6442d61bbb3SSatish Balay ridx = row%bs; 6452d61bbb3SSatish Balay high = nrow; 6462d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6472d61bbb3SSatish Balay while (high-low > 5) { 648cd0e1443SSatish Balay t = (low+high)/2; 649cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 650cd0e1443SSatish Balay else low = t; 651cd0e1443SSatish Balay } 652cd0e1443SSatish Balay for ( i=low; i<high; i++ ) { 653cd0e1443SSatish Balay if (rp[i] > bcol) break; 654cd0e1443SSatish Balay if (rp[i] == bcol) { 6552d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6562d61bbb3SSatish Balay goto finished; 657cd0e1443SSatish Balay } 658cd0e1443SSatish Balay } 6592d61bbb3SSatish Balay *v++ = zero; 6602d61bbb3SSatish Balay finished:; 661cd0e1443SSatish Balay } 662cd0e1443SSatish Balay } 6633a40ed3dSBarry Smith PetscFunctionReturn(0); 664cd0e1443SSatish Balay } 665cd0e1443SSatish Balay 6662d61bbb3SSatish Balay 6675615d1e5SSatish Balay #undef __FUNC__ 66805a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ" 66992c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 67092c4ed94SBarry Smith { 67192c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 6728a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 673dd9472c6SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented; 674549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 6753f1db9ecSBarry Smith Scalar *value = v; 6763f1db9ecSBarry Smith MatScalar *ap,*aa=a->a,*bap; 67792c4ed94SBarry Smith 6783a40ed3dSBarry Smith PetscFunctionBegin; 6790e324ae4SSatish Balay if (roworiented) { 6800e324ae4SSatish Balay stepval = (n-1)*bs; 6810e324ae4SSatish Balay } else { 6820e324ae4SSatish Balay stepval = (m-1)*bs; 6830e324ae4SSatish Balay } 68492c4ed94SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 68592c4ed94SBarry Smith row = im[k]; 6865ef9f2a5SBarry Smith if (row < 0) continue; 687aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 688a8c6a408SBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 68992c4ed94SBarry Smith #endif 69092c4ed94SBarry Smith rp = aj + ai[row]; 69192c4ed94SBarry Smith ap = aa + bs2*ai[row]; 69292c4ed94SBarry Smith rmax = imax[row]; 69392c4ed94SBarry Smith nrow = ailen[row]; 69492c4ed94SBarry Smith low = 0; 69592c4ed94SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 6965ef9f2a5SBarry Smith if (in[l] < 0) continue; 697aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 698a8c6a408SBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 69992c4ed94SBarry Smith #endif 70092c4ed94SBarry Smith col = in[l]; 70192c4ed94SBarry Smith if (roworiented) { 7020e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 7030e324ae4SSatish Balay } else { 7040e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 70592c4ed94SBarry Smith } 70692c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 70792c4ed94SBarry Smith while (high-low > 7) { 70892c4ed94SBarry Smith t = (low+high)/2; 70992c4ed94SBarry Smith if (rp[t] > col) high = t; 71092c4ed94SBarry Smith else low = t; 71192c4ed94SBarry Smith } 71292c4ed94SBarry Smith for ( i=low; i<high; i++ ) { 71392c4ed94SBarry Smith if (rp[i] > col) break; 71492c4ed94SBarry Smith if (rp[i] == col) { 7158a84c255SSatish Balay bap = ap + bs2*i; 7160e324ae4SSatish Balay if (roworiented) { 7178a84c255SSatish Balay if (is == ADD_VALUES) { 718dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 719dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7208a84c255SSatish Balay bap[jj] += *value++; 721dd9472c6SBarry Smith } 722dd9472c6SBarry Smith } 7230e324ae4SSatish Balay } else { 724dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 725dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 7260e324ae4SSatish Balay bap[jj] = *value++; 7278a84c255SSatish Balay } 728dd9472c6SBarry Smith } 729dd9472c6SBarry Smith } 7300e324ae4SSatish Balay } else { 7310e324ae4SSatish Balay if (is == ADD_VALUES) { 732dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 733dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7340e324ae4SSatish Balay *bap++ += *value++; 735dd9472c6SBarry Smith } 736dd9472c6SBarry Smith } 7370e324ae4SSatish Balay } else { 738dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 739dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 7400e324ae4SSatish Balay *bap++ = *value++; 7410e324ae4SSatish Balay } 7428a84c255SSatish Balay } 743dd9472c6SBarry Smith } 744dd9472c6SBarry Smith } 745f1241b54SBarry Smith goto noinsert2; 74692c4ed94SBarry Smith } 74792c4ed94SBarry Smith } 74889280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 749a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 75092c4ed94SBarry Smith if (nrow >= rmax) { 75192c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 75292c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7533f1db9ecSBarry Smith MatScalar *new_a; 75492c4ed94SBarry Smith 755a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 75689280ab3SLois Curfman McInnes 75792c4ed94SBarry Smith /* malloc new storage space */ 7583f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 7593f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 76092c4ed94SBarry Smith new_j = (int *) (new_a + bs2*new_nz); 76192c4ed94SBarry Smith new_i = new_j + new_nz; 76292c4ed94SBarry Smith 76392c4ed94SBarry Smith /* copy over old data into new slots */ 76492c4ed94SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 76592c4ed94SBarry Smith for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 766549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 76792c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 768549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 769549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 770549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 771549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 77292c4ed94SBarry Smith /* free up old matrix storage */ 773606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 774606d414cSSatish Balay if (!a->singlemalloc) { 775606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 776606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 777606d414cSSatish Balay } 77892c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 779*c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 78092c4ed94SBarry Smith 78192c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 78292c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 7833f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 78492c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 78592c4ed94SBarry Smith a->reallocs++; 78692c4ed94SBarry Smith a->nz++; 78792c4ed94SBarry Smith } 78892c4ed94SBarry Smith N = nrow++ - 1; 78992c4ed94SBarry Smith /* shift up all the later entries in this row */ 79092c4ed94SBarry Smith for ( ii=N; ii>=i; ii-- ) { 79192c4ed94SBarry Smith rp[ii+1] = rp[ii]; 792549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 79392c4ed94SBarry Smith } 794549d3d68SSatish Balay if (N >= i) { 795549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 796549d3d68SSatish Balay } 79792c4ed94SBarry Smith rp[i] = col; 7988a84c255SSatish Balay bap = ap + bs2*i; 7990e324ae4SSatish Balay if (roworiented) { 800dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval) { 801dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs ) { 8020e324ae4SSatish Balay bap[jj] = *value++; 803dd9472c6SBarry Smith } 804dd9472c6SBarry Smith } 8050e324ae4SSatish Balay } else { 806dd9472c6SBarry Smith for ( ii=0; ii<bs; ii++,value+=stepval ) { 807dd9472c6SBarry Smith for (jj=0; jj<bs; jj++ ) { 8080e324ae4SSatish Balay *bap++ = *value++; 8090e324ae4SSatish Balay } 810dd9472c6SBarry Smith } 811dd9472c6SBarry Smith } 812f1241b54SBarry Smith noinsert2:; 81392c4ed94SBarry Smith low = i; 81492c4ed94SBarry Smith } 81592c4ed94SBarry Smith ailen[row] = nrow; 81692c4ed94SBarry Smith } 8173a40ed3dSBarry Smith PetscFunctionReturn(0); 81892c4ed94SBarry Smith } 81992c4ed94SBarry Smith 820f2501298SSatish Balay 8215615d1e5SSatish Balay #undef __FUNC__ 8225615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ" 8238f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 824584200bdSSatish Balay { 825584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 826584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 827a7c10996SSatish Balay int m = a->m,*ip, N, *ailen = a->ilen; 828549d3d68SSatish Balay int mbs = a->mbs, bs2 = a->bs2,rmax = 0,ierr; 8293f1db9ecSBarry Smith MatScalar *aa = a->a, *ap; 830584200bdSSatish Balay 8313a40ed3dSBarry Smith PetscFunctionBegin; 8323a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 833584200bdSSatish Balay 83443ee02c3SBarry Smith if (m) rmax = ailen[0]; 835584200bdSSatish Balay for ( i=1; i<mbs; i++ ) { 836584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 837584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 838d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 839584200bdSSatish Balay if (fshift) { 840a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 841584200bdSSatish Balay N = ailen[i]; 842584200bdSSatish Balay for ( j=0; j<N; j++ ) { 843584200bdSSatish Balay ip[j-fshift] = ip[j]; 844549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 845584200bdSSatish Balay } 846584200bdSSatish Balay } 847584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 848584200bdSSatish Balay } 849584200bdSSatish Balay if (mbs) { 850584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 851584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 852584200bdSSatish Balay } 853584200bdSSatish Balay /* reset ilen and imax for each row */ 854584200bdSSatish Balay for ( i=0; i<mbs; i++ ) { 855584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 856584200bdSSatish Balay } 857a7c10996SSatish Balay a->nz = ai[mbs]; 858584200bdSSatish Balay 859584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 860584200bdSSatish Balay if (fshift && a->diag) { 861606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 862584200bdSSatish Balay PLogObjectMemory(A,-(m+1)*sizeof(int)); 863584200bdSSatish Balay a->diag = 0; 864584200bdSSatish Balay } 8653d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 866219d9a1aSLois Curfman McInnes m,a->n,a->bs,fshift*bs2,a->nz*bs2); 8673d2013a6SLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 868584200bdSSatish Balay a->reallocs); 869d402145bSBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 870e2f3b5e9SSatish Balay a->reallocs = 0; 8714e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift*bs2; 8724e220ebcSLois Curfman McInnes 8733a40ed3dSBarry Smith PetscFunctionReturn(0); 874584200bdSSatish Balay } 875584200bdSSatish Balay 8765a838052SSatish Balay 877bea157c4SSatish Balay 878bea157c4SSatish Balay /* 879bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 880bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 881bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 882bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 883bea157c4SSatish Balay */ 8845615d1e5SSatish Balay #undef __FUNC__ 885bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks" 886bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 887d9b7c43dSSatish Balay { 888bea157c4SSatish Balay int i,j,k,row; 889bea157c4SSatish Balay PetscTruth flg; 8903a40ed3dSBarry Smith 891433994e6SBarry Smith PetscFunctionBegin; 892bea157c4SSatish Balay for ( i=0,j=0; i<n; j++ ) { 893bea157c4SSatish Balay row = idx[i]; 894bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 895bea157c4SSatish Balay sizes[j] = 1; 896bea157c4SSatish Balay i++; 897e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 898bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 899bea157c4SSatish Balay i++; 900bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 901bea157c4SSatish Balay flg = PETSC_TRUE; 902bea157c4SSatish Balay for ( k=1; k<bs; k++ ) { 903bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 904bea157c4SSatish Balay flg = PETSC_FALSE; 905bea157c4SSatish Balay break; 906d9b7c43dSSatish Balay } 907bea157c4SSatish Balay } 908bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 909bea157c4SSatish Balay sizes[j] = bs; 910bea157c4SSatish Balay i+= bs; 911bea157c4SSatish Balay } else { 912bea157c4SSatish Balay sizes[j] = 1; 913bea157c4SSatish Balay i++; 914bea157c4SSatish Balay } 915bea157c4SSatish Balay } 916bea157c4SSatish Balay } 917bea157c4SSatish Balay *bs_max = j; 9183a40ed3dSBarry Smith PetscFunctionReturn(0); 919d9b7c43dSSatish Balay } 920d9b7c43dSSatish Balay 9215615d1e5SSatish Balay #undef __FUNC__ 9225615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ" 9238f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag) 924d9b7c43dSSatish Balay { 925d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 926b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 927bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9283f1db9ecSBarry Smith Scalar zero = 0.0; 9293f1db9ecSBarry Smith MatScalar *aa; 930d9b7c43dSSatish Balay 9313a40ed3dSBarry Smith PetscFunctionBegin; 932d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 933d9b7c43dSSatish Balay ierr = ISGetSize(is,&is_n);CHKERRQ(ierr); 934d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 935d9b7c43dSSatish Balay 936bea157c4SSatish Balay /* allocate memory for rows,sizes */ 937bea157c4SSatish Balay rows = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows); 938bea157c4SSatish Balay sizes = rows + is_n; 939bea157c4SSatish Balay 940bea157c4SSatish Balay /* initialize copy IS valurs to rows, and sort them */ 941bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 942bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 943bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 944b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 945bea157c4SSatish Balay 946bea157c4SSatish Balay for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) { 947bea157c4SSatish Balay row = rows[j]; 948b6815fffSBarry Smith if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row); 949bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 950bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 951bea157c4SSatish Balay if (sizes[i] == bs) { 952bea157c4SSatish Balay if (diag) { 953bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 954bea157c4SSatish Balay baij->ilen[row/bs] = 1; 955bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 956bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 957a07cd24cSSatish Balay } 958bea157c4SSatish Balay /* Now insert all the diagoanl values for this bs */ 959bea157c4SSatish Balay for ( k=0; k<bs; k++ ) { 960bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 961bea157c4SSatish Balay } 962bea157c4SSatish Balay } else { /* (!diag) */ 963bea157c4SSatish Balay baij->ilen[row/bs] = 0; 964bea157c4SSatish Balay } /* end (!diag) */ 965bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 966aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 967bea157c4SSatish Balay if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1"); 968bea157c4SSatish Balay #endif 969bea157c4SSatish Balay for ( k=0; k<count; k++ ) { 970d9b7c43dSSatish Balay aa[0] = zero; 971d9b7c43dSSatish Balay aa+=bs; 972d9b7c43dSSatish Balay } 973d9b7c43dSSatish Balay if (diag) { 974f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 975d9b7c43dSSatish Balay } 976d9b7c43dSSatish Balay } 977bea157c4SSatish Balay } 978bea157c4SSatish Balay 979606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 9809a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9813a40ed3dSBarry Smith PetscFunctionReturn(0); 982d9b7c43dSSatish Balay } 9831c351548SSatish Balay 9845615d1e5SSatish Balay #undef __FUNC__ 9852d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ" 9862d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 9872d61bbb3SSatish Balay { 9882d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 9892d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 9902d61bbb3SSatish Balay int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 9912d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 992549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 9933f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 9942d61bbb3SSatish Balay 9952d61bbb3SSatish Balay PetscFunctionBegin; 9962d61bbb3SSatish Balay for ( k=0; k<m; k++ ) { /* loop over added rows */ 9972d61bbb3SSatish Balay row = im[k]; brow = row/bs; 9985ef9f2a5SBarry Smith if (row < 0) continue; 999aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 100032e87ba3SBarry Smith if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m); 10012d61bbb3SSatish Balay #endif 10022d61bbb3SSatish Balay rp = aj + ai[brow]; 10032d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 10042d61bbb3SSatish Balay rmax = imax[brow]; 10052d61bbb3SSatish Balay nrow = ailen[brow]; 10062d61bbb3SSatish Balay low = 0; 10072d61bbb3SSatish Balay for ( l=0; l<n; l++ ) { /* loop over added columns */ 10085ef9f2a5SBarry Smith if (in[l] < 0) continue; 1009aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 101032e87ba3SBarry Smith if (in[l] >= a->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->n); 10112d61bbb3SSatish Balay #endif 10122d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10132d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10142d61bbb3SSatish Balay if (roworiented) { 10155ef9f2a5SBarry Smith value = v[l + k*n]; 10162d61bbb3SSatish Balay } else { 10172d61bbb3SSatish Balay value = v[k + l*m]; 10182d61bbb3SSatish Balay } 10192d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10202d61bbb3SSatish Balay while (high-low > 7) { 10212d61bbb3SSatish Balay t = (low+high)/2; 10222d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10232d61bbb3SSatish Balay else low = t; 10242d61bbb3SSatish Balay } 10252d61bbb3SSatish Balay for ( i=low; i<high; i++ ) { 10262d61bbb3SSatish Balay if (rp[i] > bcol) break; 10272d61bbb3SSatish Balay if (rp[i] == bcol) { 10282d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10292d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10302d61bbb3SSatish Balay else *bap = value; 10312d61bbb3SSatish Balay goto noinsert1; 10322d61bbb3SSatish Balay } 10332d61bbb3SSatish Balay } 10342d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 10352d61bbb3SSatish Balay else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10362d61bbb3SSatish Balay if (nrow >= rmax) { 10372d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10382d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10393f1db9ecSBarry Smith MatScalar *new_a; 10402d61bbb3SSatish Balay 10412d61bbb3SSatish Balay if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 10422d61bbb3SSatish Balay 10432d61bbb3SSatish Balay /* Malloc new storage space */ 10443f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 10453f1db9ecSBarry Smith new_a = (MatScalar *) PetscMalloc( len );CHKPTRQ(new_a); 10462d61bbb3SSatish Balay new_j = (int *) (new_a + bs2*new_nz); 10472d61bbb3SSatish Balay new_i = new_j + new_nz; 10482d61bbb3SSatish Balay 10492d61bbb3SSatish Balay /* copy over old data into new slots */ 10502d61bbb3SSatish Balay for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 10512d61bbb3SSatish Balay for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1052549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10532d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1054549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1055549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1056549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1057549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10582d61bbb3SSatish Balay /* free up old matrix storage */ 1059606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1060606d414cSSatish Balay if (!a->singlemalloc) { 1061606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1062606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1063606d414cSSatish Balay } 10642d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1065*c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 10662d61bbb3SSatish Balay 10672d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10682d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 10693f1db9ecSBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 10702d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 10712d61bbb3SSatish Balay a->reallocs++; 10722d61bbb3SSatish Balay a->nz++; 10732d61bbb3SSatish Balay } 10742d61bbb3SSatish Balay N = nrow++ - 1; 10752d61bbb3SSatish Balay /* shift up all the later entries in this row */ 10762d61bbb3SSatish Balay for ( ii=N; ii>=i; ii-- ) { 10772d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1078549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 10792d61bbb3SSatish Balay } 1080549d3d68SSatish Balay if (N>=i) { 1081549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1082549d3d68SSatish Balay } 10832d61bbb3SSatish Balay rp[i] = bcol; 10842d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 10852d61bbb3SSatish Balay noinsert1:; 10862d61bbb3SSatish Balay low = i; 10872d61bbb3SSatish Balay } 10882d61bbb3SSatish Balay ailen[brow] = nrow; 10892d61bbb3SSatish Balay } 10902d61bbb3SSatish Balay PetscFunctionReturn(0); 10912d61bbb3SSatish Balay } 10922d61bbb3SSatish Balay 10932d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 10942d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 10952d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 10967b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 10977b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 10982d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec); 10992d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 11002d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat); 11012d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *); 11022d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*); 11032d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec); 11042d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 11052d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 11062d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat); 11072d61bbb3SSatish Balay 11082d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 11092d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 11102d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 11112d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 11122d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 11132d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 111415091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 11152d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 1116f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_7(Mat,Vec,Vec); 1117f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_6(Mat,Vec,Vec); 1118f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_5(Mat,Vec,Vec); 1119f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_4(Mat,Vec,Vec); 1120f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_3(Mat,Vec,Vec); 1121f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_2(Mat,Vec,Vec); 1122f1af5d2fSBarry Smith extern int MatSolveTrans_SeqBAIJ_1(Mat,Vec,Vec); 11232d61bbb3SSatish Balay 11242d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 11252d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 11262d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 11272d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 11282d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 11292d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 113015091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11312d61bbb3SSatish Balay 11322d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 11332d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 11342d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 11352d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 11362d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 113715091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 11382d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 11392d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11402d61bbb3SSatish Balay 11412d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 11422d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 11432d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 11442d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 11452d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 114615091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 11472d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 11482d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11492d61bbb3SSatish Balay 11502d61bbb3SSatish Balay #undef __FUNC__ 11512d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ" 11525ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11532d61bbb3SSatish Balay { 11542d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 11552d61bbb3SSatish Balay Mat outA; 11562d61bbb3SSatish Balay int ierr; 1157667159a5SBarry Smith PetscTruth row_identity, col_identity; 11582d61bbb3SSatish Balay 11592d61bbb3SSatish Balay PetscFunctionBegin; 1160b6d21a55SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU"); 1161667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1162667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1163667159a5SBarry Smith if (!row_identity || !col_identity) { 1164b008c796SBarry Smith SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU"); 1165667159a5SBarry Smith } 11662d61bbb3SSatish Balay 11672d61bbb3SSatish Balay outA = inA; 11682d61bbb3SSatish Balay inA->factor = FACTOR_LU; 11692d61bbb3SSatish Balay 11702d61bbb3SSatish Balay if (!a->diag) { 1171*c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 11722d61bbb3SSatish Balay } 11732d61bbb3SSatish Balay /* 117415091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 117515091d37SBarry Smith for ILU(0) factorization with natural ordering 11762d61bbb3SSatish Balay */ 117715091d37SBarry Smith switch (a->bs) { 1178f1af5d2fSBarry Smith case 1: 1179f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_2_NaturalOrdering; 1180f1af5d2fSBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n"); 118115091d37SBarry Smith case 2: 118215091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering; 118315091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 1184f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_2_NaturalOrdering; 118515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n"); 118615091d37SBarry Smith break; 118715091d37SBarry Smith case 3: 118815091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering; 118915091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 1190f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_3_NaturalOrdering; 119115091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n"); 119215091d37SBarry Smith break; 119315091d37SBarry Smith case 4: 1194667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering; 1195f830108cSBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_4_NaturalOrdering; 1196f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_4_NaturalOrdering; 119715091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n"); 119815091d37SBarry Smith break; 119915091d37SBarry Smith case 5: 1200667159a5SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering; 1201667159a5SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 1202f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_5_NaturalOrdering; 120315091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n"); 120415091d37SBarry Smith break; 120515091d37SBarry Smith case 6: 120615091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering; 120715091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 1208f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_6_NaturalOrdering; 120915091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n"); 121015091d37SBarry Smith break; 121115091d37SBarry Smith case 7: 121215091d37SBarry Smith inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering; 1213f1af5d2fSBarry Smith inA->ops->solvetrans = MatSolveTrans_SeqBAIJ_7_NaturalOrdering; 121415091d37SBarry Smith inA->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering; 121515091d37SBarry Smith PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n"); 121615091d37SBarry Smith break; 1217c38d4ed2SBarry Smith default: 1218c38d4ed2SBarry Smith a->row = row; 1219c38d4ed2SBarry Smith a->col = col; 1220c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1221c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1222c38d4ed2SBarry Smith 1223c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1224c38d4ed2SBarry Smith ierr = ISInvertPermutation(col,&(a->icol));CHKERRQ(ierr); 1225c38d4ed2SBarry Smith PLogObjectParent(inA,a->icol); 1226c38d4ed2SBarry Smith 1227c38d4ed2SBarry Smith if (!a->solve_work) { 1228c38d4ed2SBarry Smith a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1229c38d4ed2SBarry Smith PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar)); 1230c38d4ed2SBarry Smith } 12312d61bbb3SSatish Balay } 12322d61bbb3SSatish Balay 1233667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1234667159a5SBarry Smith 12352d61bbb3SSatish Balay PetscFunctionReturn(0); 12362d61bbb3SSatish Balay } 12372d61bbb3SSatish Balay #undef __FUNC__ 1238d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ" 1239ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1240ba4ca20aSSatish Balay { 1241c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1242ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1243d132466eSBarry Smith int ierr; 1244ba4ca20aSSatish Balay 12453a40ed3dSBarry Smith PetscFunctionBegin; 1246c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1247d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1248d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12493a40ed3dSBarry Smith PetscFunctionReturn(0); 1250ba4ca20aSSatish Balay } 1251d9b7c43dSSatish Balay 1252fb2e594dSBarry Smith EXTERN_C_BEGIN 125327a8da17SBarry Smith #undef __FUNC__ 125427a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 125527a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 125627a8da17SBarry Smith { 125727a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 125827a8da17SBarry Smith int i,nz,n; 125927a8da17SBarry Smith 126027a8da17SBarry Smith PetscFunctionBegin; 126127a8da17SBarry Smith nz = baij->maxnz; 126227a8da17SBarry Smith n = baij->n; 126327a8da17SBarry Smith for (i=0; i<nz; i++) { 126427a8da17SBarry Smith baij->j[i] = indices[i]; 126527a8da17SBarry Smith } 126627a8da17SBarry Smith baij->nz = nz; 126727a8da17SBarry Smith for ( i=0; i<n; i++ ) { 126827a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 126927a8da17SBarry Smith } 127027a8da17SBarry Smith 127127a8da17SBarry Smith PetscFunctionReturn(0); 127227a8da17SBarry Smith } 1273fb2e594dSBarry Smith EXTERN_C_END 127427a8da17SBarry Smith 127527a8da17SBarry Smith #undef __FUNC__ 127627a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices" 127727a8da17SBarry Smith /*@ 127827a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 127927a8da17SBarry Smith in the matrix. 128027a8da17SBarry Smith 128127a8da17SBarry Smith Input Parameters: 128227a8da17SBarry Smith + mat - the SeqBAIJ matrix 128327a8da17SBarry Smith - indices - the column indices 128427a8da17SBarry Smith 128515091d37SBarry Smith Level: advanced 128615091d37SBarry Smith 128727a8da17SBarry Smith Notes: 128827a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 128927a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 129027a8da17SBarry Smith of the MatSetValues() operation. 129127a8da17SBarry Smith 129227a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 129327a8da17SBarry Smith MatCreateSeqBAIJ(). 129427a8da17SBarry Smith 129527a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 129627a8da17SBarry Smith 129727a8da17SBarry Smith @*/ 129827a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 129927a8da17SBarry Smith { 130027a8da17SBarry Smith int ierr,(*f)(Mat,int *); 130127a8da17SBarry Smith 130227a8da17SBarry Smith PetscFunctionBegin; 130327a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1304549d3d68SSatish Balay ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr); 130527a8da17SBarry Smith if (f) { 130627a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 130727a8da17SBarry Smith } else { 130827a8da17SBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 130927a8da17SBarry Smith } 131027a8da17SBarry Smith PetscFunctionReturn(0); 131127a8da17SBarry Smith } 131227a8da17SBarry Smith 13132593348eSBarry Smith /* -------------------------------------------------------------------*/ 1314cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1315cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1316cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1317cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1318cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 1319cc2dc46cSBarry Smith MatMultTrans_SeqBAIJ, 1320cc2dc46cSBarry Smith MatMultTransAdd_SeqBAIJ, 1321cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1322cc2dc46cSBarry Smith 0, 1323cc2dc46cSBarry Smith 0, 1324cc2dc46cSBarry Smith 0, 1325cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1326cc2dc46cSBarry Smith 0, 1327b6490206SBarry Smith 0, 1328f2501298SSatish Balay MatTranspose_SeqBAIJ, 1329cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1330cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1331cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1332cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1333cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1334b6490206SBarry Smith 0, 1335cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1336cc2dc46cSBarry Smith 0, 1337cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1338cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1339cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1340cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1341cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1342cc2dc46cSBarry Smith 0, 1343cc2dc46cSBarry Smith 0, 1344cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1345cc2dc46cSBarry Smith MatGetSize_SeqBAIJ, 1346cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1347cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1348cc2dc46cSBarry Smith 0, 1349cc2dc46cSBarry Smith 0, 1350cc2dc46cSBarry Smith 0, 13512e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1352cc2dc46cSBarry Smith 0, 1353cc2dc46cSBarry Smith 0, 1354cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1355cc2dc46cSBarry Smith 0, 1356cc2dc46cSBarry Smith 0, 1357cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1358cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1359cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1360cc2dc46cSBarry Smith 0, 1361cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1362cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1363cc2dc46cSBarry Smith 0, 1364cc2dc46cSBarry Smith 0, 1365cc2dc46cSBarry Smith 0, 1366cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 13673b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 136892c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1369cc2dc46cSBarry Smith 0, 1370cc2dc46cSBarry Smith 0, 1371cc2dc46cSBarry Smith 0, 1372cc2dc46cSBarry Smith 0, 1373cc2dc46cSBarry Smith 0, 1374cc2dc46cSBarry Smith 0, 1375d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1376cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1377cc2dc46cSBarry Smith 0, 1378cc2dc46cSBarry Smith 0, 1379cc2dc46cSBarry Smith MatGetMaps_Petsc}; 13802593348eSBarry Smith 13813e90b805SBarry Smith EXTERN_C_BEGIN 13823e90b805SBarry Smith #undef __FUNC__ 13833e90b805SBarry Smith #define __FUNC__ "MatStoreValues_SeqBAIJ" 13843e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 13853e90b805SBarry Smith { 13863e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 13873e90b805SBarry Smith int nz = aij->i[aij->m]*aij->bs*aij->bs2; 1388d132466eSBarry Smith int ierr; 13893e90b805SBarry Smith 13903e90b805SBarry Smith PetscFunctionBegin; 13913e90b805SBarry Smith if (aij->nonew != 1) { 13923e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 13933e90b805SBarry Smith } 13943e90b805SBarry Smith 13953e90b805SBarry Smith /* allocate space for values if not already there */ 13963e90b805SBarry Smith if (!aij->saved_values) { 13973e90b805SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 13983e90b805SBarry Smith } 13993e90b805SBarry Smith 14003e90b805SBarry Smith /* copy values over */ 1401d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 14023e90b805SBarry Smith PetscFunctionReturn(0); 14033e90b805SBarry Smith } 14043e90b805SBarry Smith EXTERN_C_END 14053e90b805SBarry Smith 14063e90b805SBarry Smith EXTERN_C_BEGIN 14073e90b805SBarry Smith #undef __FUNC__ 140832e87ba3SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqBAIJ" 140932e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 14103e90b805SBarry Smith { 14113e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1412549d3d68SSatish Balay int nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr; 14133e90b805SBarry Smith 14143e90b805SBarry Smith PetscFunctionBegin; 14153e90b805SBarry Smith if (aij->nonew != 1) { 14163e90b805SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14173e90b805SBarry Smith } 14183e90b805SBarry Smith if (!aij->saved_values) { 14193e90b805SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 14203e90b805SBarry Smith } 14213e90b805SBarry Smith 14223e90b805SBarry Smith /* copy values over */ 1423549d3d68SSatish Balay ierr = PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 14243e90b805SBarry Smith PetscFunctionReturn(0); 14253e90b805SBarry Smith } 14263e90b805SBarry Smith EXTERN_C_END 14273e90b805SBarry Smith 14285615d1e5SSatish Balay #undef __FUNC__ 14295615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ" 14302593348eSBarry Smith /*@C 143144cd7ae7SLois Curfman McInnes MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 143244cd7ae7SLois Curfman McInnes compressed row) format. For good matrix assembly performance the 143344cd7ae7SLois Curfman McInnes user should preallocate the matrix storage by setting the parameter nz 14347fc3c18eSBarry Smith (or the array nnz). By setting these parameters accurately, performance 14352bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 14362593348eSBarry Smith 1437db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1438db81eaa0SLois Curfman McInnes 14392593348eSBarry Smith Input Parameters: 1440db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 1441b6490206SBarry Smith . bs - size of block 14422593348eSBarry Smith . m - number of rows 14432593348eSBarry Smith . n - number of columns 1444b6490206SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 14457fc3c18eSBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 14462bd5e0b2SLois Curfman McInnes (possibly different for each block row) or PETSC_NULL 14472593348eSBarry Smith 14482593348eSBarry Smith Output Parameter: 14492593348eSBarry Smith . A - the matrix 14502593348eSBarry Smith 14510513a670SBarry Smith Options Database Keys: 1452db81eaa0SLois Curfman McInnes . -mat_no_unroll - uses code that does not unroll the loops in the 1453db81eaa0SLois Curfman McInnes block calculations (much slower) 1454db81eaa0SLois Curfman McInnes . -mat_block_size - size of the blocks to use 14550513a670SBarry Smith 145615091d37SBarry Smith Level: intermediate 145715091d37SBarry Smith 14582593348eSBarry Smith Notes: 145944cd7ae7SLois Curfman McInnes The block AIJ format is fully compatible with standard Fortran 77 14602593348eSBarry Smith storage. That is, the stored row and column indices can begin at 146144cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 14622593348eSBarry Smith 14632593348eSBarry Smith Specify the preallocated storage with either nz or nnz (not both). 14642593348eSBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14652593348eSBarry Smith allocation. For additional details, see the users manual chapter on 14666da5968aSLois Curfman McInnes matrices. 14672593348eSBarry Smith 1468db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 14692593348eSBarry Smith @*/ 1470b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 14712593348eSBarry Smith { 14722593348eSBarry Smith Mat B; 1473b6490206SBarry Smith Mat_SeqBAIJ *b; 1474f1af5d2fSBarry Smith int i,len,ierr,mbs,nbs,bs2,size; 1475f1af5d2fSBarry Smith PetscTruth flg; 14763b2fbd54SBarry Smith 14773a40ed3dSBarry Smith PetscFunctionBegin; 1478d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1479a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1480b6490206SBarry Smith 1481962fb4a0SBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1482302e33e3SBarry Smith mbs = m/bs; 1483302e33e3SBarry Smith nbs = n/bs; 1484302e33e3SBarry Smith bs2 = bs*bs; 14850513a670SBarry Smith 14863a40ed3dSBarry Smith if (mbs*bs!=m || nbs*bs!=n) { 1487a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize"); 14883a40ed3dSBarry Smith } 14892593348eSBarry Smith 1490b73539f3SBarry Smith if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz); 1491b73539f3SBarry Smith if (nnz) { 1492faf3f1fcSBarry Smith for (i=0; i<mbs; i++) { 1493b73539f3SBarry 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]); 1494b73539f3SBarry Smith } 1495b73539f3SBarry Smith } 1496b73539f3SBarry Smith 14972593348eSBarry Smith *A = 0; 14983f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView); 14992593348eSBarry Smith PLogObjectCreate(B); 1500b6490206SBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqBAIJ));CHKPTRQ(b); 1501549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1502549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 15031a6a6d98SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 15041a6a6d98SBarry Smith if (!flg) { 15057fc0212eSBarry Smith switch (bs) { 15067fc0212eSBarry Smith case 1: 1507f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1508f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1509f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_1; 1510f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1511f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 15127fc0212eSBarry Smith break; 15134eeb42bcSBarry Smith case 2: 1514f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1515f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1516f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_2; 1517f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1518f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 15196d84be18SBarry Smith break; 1520f327dd97SBarry Smith case 3: 1521f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1522f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1523f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_3; 1524f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1525f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 15264eeb42bcSBarry Smith break; 15276d84be18SBarry Smith case 4: 1528f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1529f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1530f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_4; 1531f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1532f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 15336d84be18SBarry Smith break; 15346d84be18SBarry Smith case 5: 1535f830108cSBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1536f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1537f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_5; 1538f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1539f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 15406d84be18SBarry Smith break; 154115091d37SBarry Smith case 6: 154215091d37SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 154315091d37SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 1544f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_6; 154515091d37SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 154615091d37SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 154715091d37SBarry Smith break; 154848e9ddb2SSatish Balay case 7: 1549f830108cSBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1550f830108cSBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1551f1af5d2fSBarry Smith B->ops->solvetrans = MatSolveTrans_SeqBAIJ_7; 1552f830108cSBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 155348e9ddb2SSatish Balay break; 15547fc0212eSBarry Smith } 15551a6a6d98SBarry Smith } 1556e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqBAIJ; 1557e1311b90SBarry Smith B->ops->view = MatView_SeqBAIJ; 15582593348eSBarry Smith B->factor = 0; 15592593348eSBarry Smith B->lupivotthreshold = 1.0; 156090f02eecSBarry Smith B->mapping = 0; 15612593348eSBarry Smith b->row = 0; 15622593348eSBarry Smith b->col = 0; 1563e51c0b9cSSatish Balay b->icol = 0; 15642593348eSBarry Smith b->reallocs = 0; 15653e90b805SBarry Smith b->saved_values = 0; 15662593348eSBarry Smith 156744cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 156844cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1569a5ae1ecdSBarry Smith 15707413bad6SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 15717413bad6SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1572a5ae1ecdSBarry Smith 1573b6490206SBarry Smith b->mbs = mbs; 1574f2501298SSatish Balay b->nbs = nbs; 1575b6490206SBarry Smith b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(b->imax); 1576*c4992f7dSBarry Smith if (!nnz) { 1577119a934aSSatish Balay if (nz == PETSC_DEFAULT) nz = 5; 15782593348eSBarry Smith else if (nz <= 0) nz = 1; 1579b6490206SBarry Smith for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1580b6490206SBarry Smith nz = nz*mbs; 15813a40ed3dSBarry Smith } else { 15822593348eSBarry Smith nz = 0; 1583b6490206SBarry Smith for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 15842593348eSBarry Smith } 15852593348eSBarry Smith 15862593348eSBarry Smith /* allocate the matrix space */ 15873f1db9ecSBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int); 15883f1db9ecSBarry Smith b->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(b->a); 1589549d3d68SSatish Balay ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 15907e67e3f9SSatish Balay b->j = (int *) (b->a + nz*bs2); 1591549d3d68SSatish Balay ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 15922593348eSBarry Smith b->i = b->j + nz; 1593*c4992f7dSBarry Smith b->singlemalloc = PETSC_TRUE; 15942593348eSBarry Smith 1595de6a44a3SBarry Smith b->i[0] = 0; 1596b6490206SBarry Smith for (i=1; i<mbs+1; i++) { 15972593348eSBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 15982593348eSBarry Smith } 15992593348eSBarry Smith 1600b6490206SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 1601b6490206SBarry Smith b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1602f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1603b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 16042593348eSBarry Smith 1605b6490206SBarry Smith b->bs = bs; 16069df24120SSatish Balay b->bs2 = bs2; 1607b6490206SBarry Smith b->mbs = mbs; 16082593348eSBarry Smith b->nz = 0; 160918e7b8e4SLois Curfman McInnes b->maxnz = nz*bs2; 1610*c4992f7dSBarry Smith b->sorted = PETSC_FALSE; 1611*c4992f7dSBarry Smith b->roworiented = PETSC_TRUE; 16122593348eSBarry Smith b->nonew = 0; 16132593348eSBarry Smith b->diag = 0; 16142593348eSBarry Smith b->solve_work = 0; 1615de6a44a3SBarry Smith b->mult_work = 0; 16162593348eSBarry Smith b->spptr = 0; 16174e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 16184e220ebcSLois Curfman McInnes 16192593348eSBarry Smith *A = B; 16202593348eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 16212593348eSBarry Smith if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); } 162227a8da17SBarry Smith 1623f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 16243e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 16253e90b805SBarry Smith (void*)MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1626f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 16273e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 16283e90b805SBarry Smith (void*)MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1629f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 163027a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 163127a8da17SBarry Smith (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 16323a40ed3dSBarry Smith PetscFunctionReturn(0); 16332593348eSBarry Smith } 16342593348eSBarry Smith 16355615d1e5SSatish Balay #undef __FUNC__ 16362e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ" 16372e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 16382593348eSBarry Smith { 16392593348eSBarry Smith Mat C; 1640b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 164127a8da17SBarry Smith int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1642de6a44a3SBarry Smith 16433a40ed3dSBarry Smith PetscFunctionBegin; 1644a8c6a408SBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix"); 16452593348eSBarry Smith 16462593348eSBarry Smith *B = 0; 16473f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView); 16482593348eSBarry Smith PLogObjectCreate(C); 1649b6490206SBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqBAIJ));CHKPTRQ(c); 1650549d3d68SSatish Balay ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1651e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqBAIJ; 1652e1311b90SBarry Smith C->ops->view = MatView_SeqBAIJ; 16532593348eSBarry Smith C->factor = A->factor; 16542593348eSBarry Smith c->row = 0; 16552593348eSBarry Smith c->col = 0; 1656cac97260SSatish Balay c->icol = 0; 165732e87ba3SBarry Smith c->saved_values = 0; 16582593348eSBarry Smith C->assembled = PETSC_TRUE; 16592593348eSBarry Smith 166044cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 166144cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 166244cd7ae7SLois Curfman McInnes C->M = a->m; 166344cd7ae7SLois Curfman McInnes C->N = a->n; 166444cd7ae7SLois Curfman McInnes 1665b6490206SBarry Smith c->bs = a->bs; 16669df24120SSatish Balay c->bs2 = a->bs2; 1667b6490206SBarry Smith c->mbs = a->mbs; 1668b6490206SBarry Smith c->nbs = a->nbs; 16692593348eSBarry Smith 1670b6490206SBarry Smith c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax); 1671b6490206SBarry Smith c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen); 1672b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 16732593348eSBarry Smith c->imax[i] = a->imax[i]; 16742593348eSBarry Smith c->ilen[i] = a->ilen[i]; 16752593348eSBarry Smith } 16762593348eSBarry Smith 16772593348eSBarry Smith /* allocate the matrix space */ 1678*c4992f7dSBarry Smith c->singlemalloc = PETSC_TRUE; 16793f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 16803f1db9ecSBarry Smith c->a = (MatScalar *) PetscMalloc( len );CHKPTRQ(c->a); 16817e67e3f9SSatish Balay c->j = (int *) (c->a + nz*bs2); 1682de6a44a3SBarry Smith c->i = c->j + nz; 1683549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1684b6490206SBarry Smith if (mbs > 0) { 1685549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16862e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1687549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16882e8a6d31SBarry Smith } else { 1689549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16902593348eSBarry Smith } 16912593348eSBarry Smith } 16922593348eSBarry Smith 1693f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16942593348eSBarry Smith c->sorted = a->sorted; 16952593348eSBarry Smith c->roworiented = a->roworiented; 16962593348eSBarry Smith c->nonew = a->nonew; 16972593348eSBarry Smith 16982593348eSBarry Smith if (a->diag) { 1699b6490206SBarry Smith c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) );CHKPTRQ(c->diag); 1700b6490206SBarry Smith PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1701b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 17022593348eSBarry Smith c->diag[i] = a->diag[i]; 17032593348eSBarry Smith } 170498305bb5SBarry Smith } else c->diag = 0; 17052593348eSBarry Smith c->nz = a->nz; 17062593348eSBarry Smith c->maxnz = a->maxnz; 17072593348eSBarry Smith c->solve_work = 0; 17082593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 17097fc0212eSBarry Smith c->mult_work = 0; 17102593348eSBarry Smith *B = C; 17117b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 17123a40ed3dSBarry Smith PetscFunctionReturn(0); 17132593348eSBarry Smith } 17142593348eSBarry Smith 17155615d1e5SSatish Balay #undef __FUNC__ 17165615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ" 171719bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 17182593348eSBarry Smith { 1719b6490206SBarry Smith Mat_SeqBAIJ *a; 17202593348eSBarry Smith Mat B; 1721f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1722b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 172335aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1724a7c10996SSatish Balay int *masked, nmask,tmp,bs2,ishift; 1725b6490206SBarry Smith Scalar *aa; 172619bcc07fSBarry Smith MPI_Comm comm = ((PetscObject) viewer)->comm; 17272593348eSBarry Smith 17283a40ed3dSBarry Smith PetscFunctionBegin; 1729f1af5d2fSBarry Smith ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1730de6a44a3SBarry Smith bs2 = bs*bs; 1731b6490206SBarry Smith 1732d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1733cc0fae11SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 173490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 17350752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 1736a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object"); 17372593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 17382593348eSBarry Smith 1739d64ed03dSBarry Smith if (header[3] < 0) { 1740a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ"); 1741d64ed03dSBarry Smith } 1742d64ed03dSBarry Smith 1743a8c6a408SBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices"); 174435aab85fSBarry Smith 174535aab85fSBarry Smith /* 174635aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 174735aab85fSBarry Smith divisible by the blocksize 174835aab85fSBarry Smith */ 1749b6490206SBarry Smith mbs = M/bs; 175035aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 175135aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 175235aab85fSBarry Smith else mbs++; 175335aab85fSBarry Smith if (extra_rows) { 1754537820f0SBarry Smith PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 175535aab85fSBarry Smith } 1756b6490206SBarry Smith 17572593348eSBarry Smith /* read in row lengths */ 175835aab85fSBarry Smith rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 17590752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 176035aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 17612593348eSBarry Smith 1762b6490206SBarry Smith /* read in column indices */ 176335aab85fSBarry Smith jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) );CHKPTRQ(jj); 17640752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 176535aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1766b6490206SBarry Smith 1767b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1768b6490206SBarry Smith browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1769549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 177035aab85fSBarry Smith mask = (int *) PetscMalloc( 2*mbs*sizeof(int) );CHKPTRQ(mask); 1771549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 177235aab85fSBarry Smith masked = mask + mbs; 1773b6490206SBarry Smith rowcount = 0; nzcount = 0; 1774b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 177535aab85fSBarry Smith nmask = 0; 1776b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1777b6490206SBarry Smith kmax = rowlengths[rowcount]; 1778b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 177935aab85fSBarry Smith tmp = jj[nzcount++]/bs; 178035aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1781b6490206SBarry Smith } 1782b6490206SBarry Smith rowcount++; 1783b6490206SBarry Smith } 178435aab85fSBarry Smith browlengths[i] += nmask; 178535aab85fSBarry Smith /* zero out the mask elements we set */ 178635aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1787b6490206SBarry Smith } 1788b6490206SBarry Smith 17892593348eSBarry Smith /* create our matrix */ 17903eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17912593348eSBarry Smith B = *A; 1792b6490206SBarry Smith a = (Mat_SeqBAIJ *) B->data; 17932593348eSBarry Smith 17942593348eSBarry Smith /* set matrix "i" values */ 1795de6a44a3SBarry Smith a->i[0] = 0; 1796b6490206SBarry Smith for ( i=1; i<= mbs; i++ ) { 1797b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1798b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17992593348eSBarry Smith } 1800b6490206SBarry Smith a->nz = 0; 1801b6490206SBarry Smith for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 18022593348eSBarry Smith 1803b6490206SBarry Smith /* read in nonzero values */ 180435aab85fSBarry Smith aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 18050752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 180635aab85fSBarry Smith for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1807b6490206SBarry Smith 1808b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1809b6490206SBarry Smith nzcount = 0; jcount = 0; 1810b6490206SBarry Smith for ( i=0; i<mbs; i++ ) { 1811b6490206SBarry Smith nzcountb = nzcount; 181235aab85fSBarry Smith nmask = 0; 1813b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1814b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1815b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 181635aab85fSBarry Smith tmp = jj[nzcount++]/bs; 181735aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1818b6490206SBarry Smith } 1819b6490206SBarry Smith rowcount++; 1820b6490206SBarry Smith } 1821de6a44a3SBarry Smith /* sort the masked values */ 1822433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1823de6a44a3SBarry Smith 1824b6490206SBarry Smith /* set "j" values into matrix */ 1825b6490206SBarry Smith maskcount = 1; 182635aab85fSBarry Smith for ( j=0; j<nmask; j++ ) { 182735aab85fSBarry Smith a->j[jcount++] = masked[j]; 1828de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1829b6490206SBarry Smith } 1830b6490206SBarry Smith /* set "a" values into matrix */ 1831de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1832b6490206SBarry Smith for ( j=0; j<bs; j++ ) { 1833b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1834b6490206SBarry Smith for ( k=0; k<kmax; k++ ) { 1835de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1836de6a44a3SBarry Smith block = mask[tmp] - 1; 1837de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1838de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1839b6490206SBarry Smith a->a[idx] = aa[nzcountb++]; 1840b6490206SBarry Smith } 1841b6490206SBarry Smith } 184235aab85fSBarry Smith /* zero out the mask elements we set */ 184335aab85fSBarry Smith for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1844b6490206SBarry Smith } 1845a8c6a408SBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix"); 1846b6490206SBarry Smith 1847606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1848606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1849606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1850606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1851606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1852b6490206SBarry Smith 1853b6490206SBarry Smith B->assembled = PETSC_TRUE; 1854de6a44a3SBarry Smith 18559c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 18563a40ed3dSBarry Smith PetscFunctionReturn(0); 18572593348eSBarry Smith } 18582593348eSBarry Smith 18592593348eSBarry Smith 18602593348eSBarry Smith 1861