1*8661488fSKris Buschelman /*$Id: baij.c,v 1.233 2001/07/11 19:38:40 buschelm Exp $*/ 22593348eSBarry Smith 32593348eSBarry Smith /* 4b6490206SBarry Smith Defines the basic matrix operations for the BAIJ (compressed row) 52593348eSBarry Smith matrix storage format. 62593348eSBarry Smith */ 7a2ccead7SSatish Balay #include "petscsys.h" /*I "petscmat.h" I*/ 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 1295d5f7c2SBarry Smith /* UGLY, ugly, ugly 1395d5f7c2SBarry Smith When MatScalar == Scalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does 143477af42SKris Buschelman not exist. Otherwise ..._MatScalar() takes matrix dlements in single precision and 1595d5f7c2SBarry Smith inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ() 1695d5f7c2SBarry Smith converts the entries into single precision and then calls ..._MatScalar() to put them 1795d5f7c2SBarry Smith into the single precision data structures. 1895d5f7c2SBarry Smith */ 1995d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 20ca44d042SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode); 2195d5f7c2SBarry Smith #else 2295d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ 2395d5f7c2SBarry Smith #endif 2495d5f7c2SBarry Smith 252d61bbb3SSatish Balay #define CHUNKSIZE 10 26de6a44a3SBarry Smith 27be5855fcSBarry Smith /* 28be5855fcSBarry Smith Checks for missing diagonals 29be5855fcSBarry Smith */ 304a2ae208SSatish Balay #undef __FUNCT__ 314a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ" 32c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A) 33be5855fcSBarry Smith { 34be5855fcSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 35883fce79SBarry Smith int *diag,*jj = a->j,i,ierr; 36be5855fcSBarry Smith 37be5855fcSBarry Smith PetscFunctionBegin; 38c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 39883fce79SBarry Smith diag = a->diag; 400e8e8aceSBarry Smith for (i=0; i<a->mbs; i++) { 41be5855fcSBarry Smith if (jj[diag[i]] != i) { 4229bbc08cSBarry Smith SETERRQ1(1,"Matrix is missing diagonal number %d",i); 43be5855fcSBarry Smith } 44be5855fcSBarry Smith } 45be5855fcSBarry Smith PetscFunctionReturn(0); 46be5855fcSBarry Smith } 47be5855fcSBarry Smith 484a2ae208SSatish Balay #undef __FUNCT__ 494a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ" 50c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A) 51de6a44a3SBarry Smith { 52de6a44a3SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 5382502324SSatish Balay int i,j,*diag,m = a->mbs,ierr; 54de6a44a3SBarry Smith 553a40ed3dSBarry Smith PetscFunctionBegin; 56883fce79SBarry Smith if (a->diag) PetscFunctionReturn(0); 57883fce79SBarry Smith 58b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 59b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 607fc0212eSBarry Smith for (i=0; i<m; i++) { 61ecc615b2SBarry Smith diag[i] = a->i[i+1]; 62de6a44a3SBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 63de6a44a3SBarry Smith if (a->j[j] == i) { 64de6a44a3SBarry Smith diag[i] = j; 65de6a44a3SBarry Smith break; 66de6a44a3SBarry Smith } 67de6a44a3SBarry Smith } 68de6a44a3SBarry Smith } 69de6a44a3SBarry Smith a->diag = diag; 703a40ed3dSBarry Smith PetscFunctionReturn(0); 71de6a44a3SBarry Smith } 722593348eSBarry Smith 732593348eSBarry Smith 74ca44d042SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 753b2fbd54SBarry Smith 764a2ae208SSatish Balay #undef __FUNCT__ 774a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqBAIJ" 780e6d2581SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 793b2fbd54SBarry Smith { 803b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 813b2fbd54SBarry Smith int ierr,n = a->mbs,i; 823b2fbd54SBarry Smith 833a40ed3dSBarry Smith PetscFunctionBegin; 843b2fbd54SBarry Smith *nn = n; 853a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 863b2fbd54SBarry Smith if (symmetric) { 873b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 883b2fbd54SBarry Smith } else if (oshift == 1) { 893b2fbd54SBarry Smith /* temporarily add 1 to i and j indices */ 90435da068SBarry Smith int nz = a->i[n]; 913b2fbd54SBarry Smith for (i=0; i<nz; i++) a->j[i]++; 923b2fbd54SBarry Smith for (i=0; i<n+1; i++) a->i[i]++; 933b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 943b2fbd54SBarry Smith } else { 953b2fbd54SBarry Smith *ia = a->i; *ja = a->j; 963b2fbd54SBarry Smith } 973b2fbd54SBarry Smith 983a40ed3dSBarry Smith PetscFunctionReturn(0); 993b2fbd54SBarry Smith } 1003b2fbd54SBarry Smith 1014a2ae208SSatish Balay #undef __FUNCT__ 1024a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ" 103435da068SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done) 1043b2fbd54SBarry Smith { 1053b2fbd54SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 106606d414cSSatish Balay int i,n = a->mbs,ierr; 1073b2fbd54SBarry Smith 1083a40ed3dSBarry Smith PetscFunctionBegin; 1093a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 1103b2fbd54SBarry Smith if (symmetric) { 111606d414cSSatish Balay ierr = PetscFree(*ia);CHKERRQ(ierr); 112606d414cSSatish Balay ierr = PetscFree(*ja);CHKERRQ(ierr); 113af5da2bfSBarry Smith } else if (oshift == 1) { 114435da068SBarry Smith int nz = a->i[n]-1; 1153b2fbd54SBarry Smith for (i=0; i<nz; i++) a->j[i]--; 1163b2fbd54SBarry Smith for (i=0; i<n+1; i++) a->i[i]--; 1173b2fbd54SBarry Smith } 1183a40ed3dSBarry Smith PetscFunctionReturn(0); 1193b2fbd54SBarry Smith } 1203b2fbd54SBarry Smith 1214a2ae208SSatish Balay #undef __FUNCT__ 1224a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqBAIJ" 1232d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs) 1242d61bbb3SSatish Balay { 1252d61bbb3SSatish Balay Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data; 1262d61bbb3SSatish Balay 1272d61bbb3SSatish Balay PetscFunctionBegin; 1282d61bbb3SSatish Balay *bs = baij->bs; 1292d61bbb3SSatish Balay PetscFunctionReturn(0); 1302d61bbb3SSatish Balay } 1312d61bbb3SSatish Balay 1322d61bbb3SSatish Balay 1334a2ae208SSatish Balay #undef __FUNCT__ 1344a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqBAIJ" 135e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A) 1362d61bbb3SSatish Balay { 1372d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 138e51c0b9cSSatish Balay int ierr; 1392d61bbb3SSatish Balay 140433994e6SBarry Smith PetscFunctionBegin; 141aa482453SBarry Smith #if defined(PETSC_USE_LOG) 142b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz); 1432d61bbb3SSatish Balay #endif 144606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 145606d414cSSatish Balay if (!a->singlemalloc) { 146606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 147606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 148606d414cSSatish Balay } 149c38d4ed2SBarry Smith if (a->row) { 150c38d4ed2SBarry Smith ierr = ISDestroy(a->row);CHKERRQ(ierr); 151c38d4ed2SBarry Smith } 152c38d4ed2SBarry Smith if (a->col) { 153c38d4ed2SBarry Smith ierr = ISDestroy(a->col);CHKERRQ(ierr); 154c38d4ed2SBarry Smith } 155606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 156606d414cSSatish Balay if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);} 157606d414cSSatish Balay if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);} 158606d414cSSatish Balay if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);} 159606d414cSSatish Balay if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);} 160e51c0b9cSSatish Balay if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 161606d414cSSatish Balay if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);} 162563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 163563b5814SBarry Smith if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);} 164563b5814SBarry Smith #endif 165606d414cSSatish Balay ierr = PetscFree(a);CHKERRQ(ierr); 1662d61bbb3SSatish Balay PetscFunctionReturn(0); 1672d61bbb3SSatish Balay } 1682d61bbb3SSatish Balay 1694a2ae208SSatish Balay #undef __FUNCT__ 1704a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqBAIJ" 1712d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op) 1722d61bbb3SSatish Balay { 1732d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1742d61bbb3SSatish Balay 1752d61bbb3SSatish Balay PetscFunctionBegin; 176aa275fccSKris Buschelman switch (op) { 177aa275fccSKris Buschelman case MAT_ROW_ORIENTED: 178aa275fccSKris Buschelman a->roworiented = PETSC_TRUE; 179aa275fccSKris Buschelman break; 180aa275fccSKris Buschelman case MAT_COLUMN_ORIENTED: 181aa275fccSKris Buschelman a->roworiented = PETSC_FALSE; 182aa275fccSKris Buschelman break; 183aa275fccSKris Buschelman case MAT_COLUMNS_SORTED: 184aa275fccSKris Buschelman a->sorted = PETSC_TRUE; 185aa275fccSKris Buschelman break; 186aa275fccSKris Buschelman case MAT_COLUMNS_UNSORTED: 187aa275fccSKris Buschelman a->sorted = PETSC_FALSE; 188aa275fccSKris Buschelman break; 189aa275fccSKris Buschelman case MAT_KEEP_ZEROED_ROWS: 190aa275fccSKris Buschelman a->keepzeroedrows = PETSC_TRUE; 191aa275fccSKris Buschelman break; 192aa275fccSKris Buschelman case MAT_NO_NEW_NONZERO_LOCATIONS: 193aa275fccSKris Buschelman a->nonew = 1; 194aa275fccSKris Buschelman break; 195aa275fccSKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 196aa275fccSKris Buschelman a->nonew = -1; 197aa275fccSKris Buschelman break; 198aa275fccSKris Buschelman case MAT_NEW_NONZERO_ALLOCATION_ERR: 199aa275fccSKris Buschelman a->nonew = -2; 200aa275fccSKris Buschelman break; 201aa275fccSKris Buschelman case MAT_YES_NEW_NONZERO_LOCATIONS: 202aa275fccSKris Buschelman a->nonew = 0; 203aa275fccSKris Buschelman break; 204aa275fccSKris Buschelman case MAT_ROWS_SORTED: 205aa275fccSKris Buschelman case MAT_ROWS_UNSORTED: 206aa275fccSKris Buschelman case MAT_YES_NEW_DIAGONALS: 207aa275fccSKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 208aa275fccSKris Buschelman case MAT_USE_HASH_TABLE: 209b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n"); 210aa275fccSKris Buschelman break; 211aa275fccSKris Buschelman case MAT_NO_NEW_DIAGONALS: 21229bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); 213bd648c37SKris Buschelman break; 214bd648c37SKris Buschelman case MAT_USE_SINGLE_PRECISION_SOLVES: 215bd648c37SKris Buschelman #if defined(PETSC_USE_MAT_SINGLE) 216bd648c37SKris Buschelman if (a->bs==4) { 217*8661488fSKris Buschelman PetscTruth sse_enabled,single_prec,flg; 218bd648c37SKris Buschelman int ierr; 219bd648c37SKris Buschelman ierr = PetscSSEIsEnabled(&sse_enabled);CHKERRQ(ierr); 220bd648c37SKris Buschelman if (sse_enabled) { 221*8661488fSKris Buschelman a->single_precision_solves = PETSC_TRUE; 222*8661488fSKris Buschelman ierr = PetscOptionsGetLogical(PETSC_NULL,"-mat_single_precision_solves",&single_prec,&flg);CHKERRQ(ierr); 223*8661488fSKris Buschelman if (flg) { 224*8661488fSKris Buschelman if (!single_prec) { 225*8661488fSKris Buschelman a->single_precision_solves = PETSC_FALSE; 226*8661488fSKris Buschelman } 227*8661488fSKris Buschelman } 228bd648c37SKris Buschelman } else { 229bd648c37SKris Buschelman PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n"); 230bd648c37SKris Buschelman } 231bd648c37SKris Buschelman } else { 232bd648c37SKris Buschelman PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n"); 233bd648c37SKris Buschelman } 234bd648c37SKris Buschelman #else 235bd648c37SKris Buschelman PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n"); 236bd648c37SKris Buschelman #endif 237bd648c37SKris Buschelman break; 238aa275fccSKris Buschelman default: 23929bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"unknown option"); 2402d61bbb3SSatish Balay } 2412d61bbb3SSatish Balay PetscFunctionReturn(0); 2422d61bbb3SSatish Balay } 2432d61bbb3SSatish Balay 2444a2ae208SSatish Balay #undef __FUNCT__ 2454a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_SeqBAIJ" 2462d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 2472d61bbb3SSatish Balay { 2482d61bbb3SSatish Balay PetscFunctionBegin; 2494c49b128SBarry Smith if (m) *m = 0; 250273d9f13SBarry Smith if (n) *n = A->m; 2512d61bbb3SSatish Balay PetscFunctionReturn(0); 2522d61bbb3SSatish Balay } 2532d61bbb3SSatish Balay 2544a2ae208SSatish Balay #undef __FUNCT__ 2554a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ" 2562d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 2572d61bbb3SSatish Balay { 2582d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 25982502324SSatish Balay int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr; 2603f1db9ecSBarry Smith MatScalar *aa,*aa_i; 2613f1db9ecSBarry Smith Scalar *v_i; 2622d61bbb3SSatish Balay 2632d61bbb3SSatish Balay PetscFunctionBegin; 2642d61bbb3SSatish Balay bs = a->bs; 2652d61bbb3SSatish Balay ai = a->i; 2662d61bbb3SSatish Balay aj = a->j; 2672d61bbb3SSatish Balay aa = a->a; 2682d61bbb3SSatish Balay bs2 = a->bs2; 2692d61bbb3SSatish Balay 270273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 2712d61bbb3SSatish Balay 2722d61bbb3SSatish Balay bn = row/bs; /* Block number */ 2732d61bbb3SSatish Balay bp = row % bs; /* Block Position */ 2742d61bbb3SSatish Balay M = ai[bn+1] - ai[bn]; 2752d61bbb3SSatish Balay *nz = bs*M; 2762d61bbb3SSatish Balay 2772d61bbb3SSatish Balay if (v) { 2782d61bbb3SSatish Balay *v = 0; 2792d61bbb3SSatish Balay if (*nz) { 280b0a32e0cSBarry Smith ierr = PetscMalloc((*nz)*sizeof(Scalar),v);CHKERRQ(ierr); 2812d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 2822d61bbb3SSatish Balay v_i = *v + i*bs; 2832d61bbb3SSatish Balay aa_i = aa + bs2*(ai[bn] + i); 2842d61bbb3SSatish Balay for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];} 2852d61bbb3SSatish Balay } 2862d61bbb3SSatish Balay } 2872d61bbb3SSatish Balay } 2882d61bbb3SSatish Balay 2892d61bbb3SSatish Balay if (idx) { 2902d61bbb3SSatish Balay *idx = 0; 2912d61bbb3SSatish Balay if (*nz) { 292b0a32e0cSBarry Smith ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr); 2932d61bbb3SSatish Balay for (i=0; i<M; i++) { /* for each block in the block row */ 2942d61bbb3SSatish Balay idx_i = *idx + i*bs; 2952d61bbb3SSatish Balay itmp = bs*aj[ai[bn] + i]; 2962d61bbb3SSatish Balay for (j=0; j<bs; j++) {idx_i[j] = itmp++;} 2972d61bbb3SSatish Balay } 2982d61bbb3SSatish Balay } 2992d61bbb3SSatish Balay } 3002d61bbb3SSatish Balay PetscFunctionReturn(0); 3012d61bbb3SSatish Balay } 3022d61bbb3SSatish Balay 3034a2ae208SSatish Balay #undef __FUNCT__ 3044a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ" 3052d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 3062d61bbb3SSatish Balay { 307606d414cSSatish Balay int ierr; 308606d414cSSatish Balay 3092d61bbb3SSatish Balay PetscFunctionBegin; 310606d414cSSatish Balay if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}} 311606d414cSSatish Balay if (v) {if (*v) {ierr = PetscFree(*v);CHKERRQ(ierr);}} 3122d61bbb3SSatish Balay PetscFunctionReturn(0); 3132d61bbb3SSatish Balay } 3142d61bbb3SSatish Balay 3154a2ae208SSatish Balay #undef __FUNCT__ 3164a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ" 3172d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B) 3182d61bbb3SSatish Balay { 3192d61bbb3SSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 3202d61bbb3SSatish Balay Mat C; 3212d61bbb3SSatish Balay int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 322273d9f13SBarry Smith int *rows,*cols,bs2=a->bs2; 323375fe846SBarry Smith Scalar *array; 3242d61bbb3SSatish Balay 3252d61bbb3SSatish Balay PetscFunctionBegin; 32629bbc08cSBarry Smith if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place"); 327b0a32e0cSBarry Smith ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr); 328549d3d68SSatish Balay ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr); 3292d61bbb3SSatish Balay 330375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 331b0a32e0cSBarry Smith ierr = PetscMalloc(a->bs2*a->nz*sizeof(Scalar),&array);CHKERRQ(ierr); 332375fe846SBarry Smith for (i=0; i<a->bs2*a->nz; i++) array[i] = (Scalar)a->a[i]; 333375fe846SBarry Smith #else 3343eda8832SBarry Smith array = a->a; 335375fe846SBarry Smith #endif 336375fe846SBarry Smith 3372d61bbb3SSatish Balay for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1; 338273d9f13SBarry Smith ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr); 339606d414cSSatish Balay ierr = PetscFree(col);CHKERRQ(ierr); 340b0a32e0cSBarry Smith ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr); 3412d61bbb3SSatish Balay cols = rows + bs; 3422d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3432d61bbb3SSatish Balay cols[0] = i*bs; 3442d61bbb3SSatish Balay for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1; 3452d61bbb3SSatish Balay len = ai[i+1] - ai[i]; 3462d61bbb3SSatish Balay for (j=0; j<len; j++) { 3472d61bbb3SSatish Balay rows[0] = (*aj++)*bs; 3482d61bbb3SSatish Balay for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1; 3492d61bbb3SSatish Balay ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr); 3502d61bbb3SSatish Balay array += bs2; 3512d61bbb3SSatish Balay } 3522d61bbb3SSatish Balay } 353606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 354375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 355375fe846SBarry Smith ierr = PetscFree(array); 356375fe846SBarry Smith #endif 3572d61bbb3SSatish Balay 3582d61bbb3SSatish Balay ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3592d61bbb3SSatish Balay ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3602d61bbb3SSatish Balay 361c4992f7dSBarry Smith if (B) { 3622d61bbb3SSatish Balay *B = C; 3632d61bbb3SSatish Balay } else { 364273d9f13SBarry Smith ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 3652d61bbb3SSatish Balay } 3662d61bbb3SSatish Balay PetscFunctionReturn(0); 3672d61bbb3SSatish Balay } 3682d61bbb3SSatish Balay 3694a2ae208SSatish Balay #undef __FUNCT__ 3704a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary" 371b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer) 3722593348eSBarry Smith { 373b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 3749df24120SSatish Balay int i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 375b6490206SBarry Smith Scalar *aa; 376ce6f0cecSBarry Smith FILE *file; 3772593348eSBarry Smith 3783a40ed3dSBarry Smith PetscFunctionBegin; 379b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 380b0a32e0cSBarry Smith ierr = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr); 3812593348eSBarry Smith col_lens[0] = MAT_COOKIE; 3823b2fbd54SBarry Smith 383273d9f13SBarry Smith col_lens[1] = A->m; 384273d9f13SBarry Smith col_lens[2] = A->n; 3857e67e3f9SSatish Balay col_lens[3] = a->nz*bs2; 3862593348eSBarry Smith 3872593348eSBarry Smith /* store lengths of each row and write (including header) to file */ 388b6490206SBarry Smith count = 0; 389b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 390b6490206SBarry Smith for (j=0; j<bs; j++) { 391b6490206SBarry Smith col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 392b6490206SBarry Smith } 3932593348eSBarry Smith } 394273d9f13SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr); 395606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 3962593348eSBarry Smith 3972593348eSBarry Smith /* store column indices (zero start index) */ 398b0a32e0cSBarry Smith ierr = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr); 399b6490206SBarry Smith count = 0; 400b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 401b6490206SBarry Smith for (j=0; j<bs; j++) { 402b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 403b6490206SBarry Smith for (l=0; l<bs; l++) { 404b6490206SBarry Smith jj[count++] = bs*a->j[k] + l; 4052593348eSBarry Smith } 4062593348eSBarry Smith } 407b6490206SBarry Smith } 408b6490206SBarry Smith } 4090752156aSBarry Smith ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr); 410606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 4112593348eSBarry Smith 4122593348eSBarry Smith /* store nonzero values */ 413b0a32e0cSBarry Smith ierr = PetscMalloc((a->nz+1)*bs2*sizeof(Scalar),&aa);CHKERRQ(ierr); 414b6490206SBarry Smith count = 0; 415b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 416b6490206SBarry Smith for (j=0; j<bs; j++) { 417b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 418b6490206SBarry Smith for (l=0; l<bs; l++) { 4197e67e3f9SSatish Balay aa[count++] = a->a[bs2*k + l*bs + j]; 420b6490206SBarry Smith } 421b6490206SBarry Smith } 422b6490206SBarry Smith } 423b6490206SBarry Smith } 4240752156aSBarry Smith ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr); 425606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 426ce6f0cecSBarry Smith 427b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr); 428ce6f0cecSBarry Smith if (file) { 429ce6f0cecSBarry Smith fprintf(file,"-matload_block_size %d\n",a->bs); 430ce6f0cecSBarry Smith } 4313a40ed3dSBarry Smith PetscFunctionReturn(0); 4322593348eSBarry Smith } 4332593348eSBarry Smith 4344a2ae208SSatish Balay #undef __FUNCT__ 4354a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII" 436b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer) 4372593348eSBarry Smith { 438b6490206SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 439fb9695e5SSatish Balay int ierr,i,j,bs = a->bs,k,l,bs2=a->bs2; 440f3ef73ceSBarry Smith PetscViewerFormat format; 4412593348eSBarry Smith 4423a40ed3dSBarry Smith PetscFunctionBegin; 443b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 444fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) { 445b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block size is %d\n",bs);CHKERRQ(ierr); 446fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 44729bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 448fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 449b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 45044cd7ae7SLois Curfman McInnes for (i=0; i<a->mbs; i++) { 45144cd7ae7SLois Curfman McInnes for (j=0; j<bs; j++) { 452b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 45344cd7ae7SLois Curfman McInnes for (k=a->i[i]; k<a->i[i+1]; k++) { 45444cd7ae7SLois Curfman McInnes for (l=0; l<bs; l++) { 455aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4560e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 457b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l, 4580e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4590e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 460b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l, 4610e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4620e6d2581SBarry Smith } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) { 463b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4640ef38995SBarry Smith } 46544cd7ae7SLois Curfman McInnes #else 4660ef38995SBarry Smith if (a->a[bs2*k + l*bs + j] != 0.0) { 467b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 4680ef38995SBarry Smith } 46944cd7ae7SLois Curfman McInnes #endif 47044cd7ae7SLois Curfman McInnes } 47144cd7ae7SLois Curfman McInnes } 472b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 47344cd7ae7SLois Curfman McInnes } 47444cd7ae7SLois Curfman McInnes } 475b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 4760ef38995SBarry Smith } else { 477b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 478b6490206SBarry Smith for (i=0; i<a->mbs; i++) { 479b6490206SBarry Smith for (j=0; j<bs; j++) { 480b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr); 481b6490206SBarry Smith for (k=a->i[i]; k<a->i[i+1]; k++) { 482b6490206SBarry Smith for (l=0; l<bs; l++) { 483aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 4840e6d2581SBarry Smith if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) { 485b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l, 4860e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4870e6d2581SBarry Smith } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) { 488b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l, 4890e6d2581SBarry Smith PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 4900ef38995SBarry Smith } else { 491b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr); 49288685aaeSLois Curfman McInnes } 49388685aaeSLois Curfman McInnes #else 494b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr); 49588685aaeSLois Curfman McInnes #endif 4962593348eSBarry Smith } 4972593348eSBarry Smith } 498b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 4992593348eSBarry Smith } 5002593348eSBarry Smith } 501b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 502b6490206SBarry Smith } 503b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 5043a40ed3dSBarry Smith PetscFunctionReturn(0); 5052593348eSBarry Smith } 5062593348eSBarry Smith 5074a2ae208SSatish Balay #undef __FUNCT__ 5084a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom" 509b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa) 5103270192aSSatish Balay { 51177ed5343SBarry Smith Mat A = (Mat) Aa; 5123270192aSSatish Balay Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 513b65db4caSBarry Smith int row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2; 5140e6d2581SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 5153f1db9ecSBarry Smith MatScalar *aa; 516b0a32e0cSBarry Smith PetscViewer viewer; 5173270192aSSatish Balay 5183a40ed3dSBarry Smith PetscFunctionBegin; 5193270192aSSatish Balay 520b65db4caSBarry Smith /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/ 52177ed5343SBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 52277ed5343SBarry Smith 523b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 52477ed5343SBarry Smith 5253270192aSSatish Balay /* loop over matrix elements drawing boxes */ 526b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 5273270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5283270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 529273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 5303270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5313270192aSSatish Balay aa = a->a + j*bs2; 5323270192aSSatish Balay for (k=0; k<bs; k++) { 5333270192aSSatish Balay for (l=0; l<bs; l++) { 5340e6d2581SBarry Smith if (PetscRealPart(*aa++) >= 0.) continue; 535b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5363270192aSSatish Balay } 5373270192aSSatish Balay } 5383270192aSSatish Balay } 5393270192aSSatish Balay } 540b0a32e0cSBarry Smith color = PETSC_DRAW_CYAN; 5413270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5423270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 543273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 5443270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5453270192aSSatish Balay aa = a->a + j*bs2; 5463270192aSSatish Balay for (k=0; k<bs; k++) { 5473270192aSSatish Balay for (l=0; l<bs; l++) { 5480e6d2581SBarry Smith if (PetscRealPart(*aa++) != 0.) continue; 549b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5503270192aSSatish Balay } 5513270192aSSatish Balay } 5523270192aSSatish Balay } 5533270192aSSatish Balay } 5543270192aSSatish Balay 555b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 5563270192aSSatish Balay for (i=0,row=0; i<mbs; i++,row+=bs) { 5573270192aSSatish Balay for (j=a->i[i]; j<a->i[i+1]; j++) { 558273d9f13SBarry Smith y_l = A->m - row - 1.0; y_r = y_l + 1.0; 5593270192aSSatish Balay x_l = a->j[j]*bs; x_r = x_l + 1.0; 5603270192aSSatish Balay aa = a->a + j*bs2; 5613270192aSSatish Balay for (k=0; k<bs; k++) { 5623270192aSSatish Balay for (l=0; l<bs; l++) { 5630e6d2581SBarry Smith if (PetscRealPart(*aa++) <= 0.) continue; 564b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr); 5653270192aSSatish Balay } 5663270192aSSatish Balay } 5673270192aSSatish Balay } 5683270192aSSatish Balay } 56977ed5343SBarry Smith PetscFunctionReturn(0); 57077ed5343SBarry Smith } 5713270192aSSatish Balay 5724a2ae208SSatish Balay #undef __FUNCT__ 5734a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw" 574b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer) 57577ed5343SBarry Smith { 5767dce120fSSatish Balay int ierr; 5770e6d2581SBarry Smith PetscReal xl,yl,xr,yr,w,h; 578b0a32e0cSBarry Smith PetscDraw draw; 57977ed5343SBarry Smith PetscTruth isnull; 5803270192aSSatish Balay 58177ed5343SBarry Smith PetscFunctionBegin; 58277ed5343SBarry Smith 583b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 584b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 58577ed5343SBarry Smith 58677ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 587273d9f13SBarry Smith xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 58877ed5343SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 589b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 590b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr); 59177ed5343SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5923a40ed3dSBarry Smith PetscFunctionReturn(0); 5933270192aSSatish Balay } 5943270192aSSatish Balay 5954a2ae208SSatish Balay #undef __FUNCT__ 5964a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ" 597b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer) 5982593348eSBarry Smith { 59919bcc07fSBarry Smith int ierr; 6006831982aSBarry Smith PetscTruth issocket,isascii,isbinary,isdraw; 6012593348eSBarry Smith 6023a40ed3dSBarry Smith PetscFunctionBegin; 603b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 604b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 605fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 606fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 6070f5bd95cSBarry Smith if (issocket) { 60829bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported"); 6090f5bd95cSBarry Smith } else if (isascii){ 6103a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr); 6110f5bd95cSBarry Smith } else if (isbinary) { 6123a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr); 6130f5bd95cSBarry Smith } else if (isdraw) { 6143a40ed3dSBarry Smith ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr); 6155cd90555SBarry Smith } else { 61629bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name); 6172593348eSBarry Smith } 6183a40ed3dSBarry Smith PetscFunctionReturn(0); 6192593348eSBarry Smith } 620b6490206SBarry Smith 621cd0e1443SSatish Balay 6224a2ae208SSatish Balay #undef __FUNCT__ 6234a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ" 6242d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 625cd0e1443SSatish Balay { 626cd0e1443SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6272d61bbb3SSatish Balay int *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; 6282d61bbb3SSatish Balay int *ai = a->i,*ailen = a->ilen; 6292d61bbb3SSatish Balay int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 6303f1db9ecSBarry Smith MatScalar *ap,*aa = a->a,zero = 0.0; 631cd0e1443SSatish Balay 6323a40ed3dSBarry Smith PetscFunctionBegin; 6332d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over rows */ 634cd0e1443SSatish Balay row = im[k]; brow = row/bs; 63529bbc08cSBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row"); 636273d9f13SBarry Smith if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 6372d61bbb3SSatish Balay rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 6382c3acbe9SBarry Smith nrow = ailen[brow]; 6392d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over columns */ 64029bbc08cSBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column"); 641273d9f13SBarry Smith if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 6422d61bbb3SSatish Balay col = in[l] ; 6432d61bbb3SSatish Balay bcol = col/bs; 6442d61bbb3SSatish Balay cidx = col%bs; 6452d61bbb3SSatish Balay ridx = row%bs; 6462d61bbb3SSatish Balay high = nrow; 6472d61bbb3SSatish Balay low = 0; /* assume unsorted */ 6482d61bbb3SSatish Balay while (high-low > 5) { 649cd0e1443SSatish Balay t = (low+high)/2; 650cd0e1443SSatish Balay if (rp[t] > bcol) high = t; 651cd0e1443SSatish Balay else low = t; 652cd0e1443SSatish Balay } 653cd0e1443SSatish Balay for (i=low; i<high; i++) { 654cd0e1443SSatish Balay if (rp[i] > bcol) break; 655cd0e1443SSatish Balay if (rp[i] == bcol) { 6562d61bbb3SSatish Balay *v++ = ap[bs2*i+bs*cidx+ridx]; 6572d61bbb3SSatish Balay goto finished; 658cd0e1443SSatish Balay } 659cd0e1443SSatish Balay } 6602d61bbb3SSatish Balay *v++ = zero; 6612d61bbb3SSatish Balay finished:; 662cd0e1443SSatish Balay } 663cd0e1443SSatish Balay } 6643a40ed3dSBarry Smith PetscFunctionReturn(0); 665cd0e1443SSatish Balay } 666cd0e1443SSatish Balay 66795d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 6684a2ae208SSatish Balay #undef __FUNCT__ 6694a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 67095d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,Scalar *v,InsertMode addv) 67195d5f7c2SBarry Smith { 672563b5814SBarry Smith Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data; 673563b5814SBarry Smith int ierr,i,N = m*n*b->bs2; 674563b5814SBarry Smith MatScalar *vsingle; 67595d5f7c2SBarry Smith 67695d5f7c2SBarry Smith PetscFunctionBegin; 677563b5814SBarry Smith if (N > b->setvalueslen) { 678563b5814SBarry Smith if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);} 679b0a32e0cSBarry Smith ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr); 680563b5814SBarry Smith b->setvalueslen = N; 681563b5814SBarry Smith } 682563b5814SBarry Smith vsingle = b->setvaluescopy; 68395d5f7c2SBarry Smith for (i=0; i<N; i++) { 68495d5f7c2SBarry Smith vsingle[i] = v[i]; 68595d5f7c2SBarry Smith } 68695d5f7c2SBarry Smith ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr); 68795d5f7c2SBarry Smith PetscFunctionReturn(0); 68895d5f7c2SBarry Smith } 68995d5f7c2SBarry Smith #endif 69095d5f7c2SBarry Smith 6912d61bbb3SSatish Balay 6924a2ae208SSatish Balay #undef __FUNCT__ 6934a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ" 69495d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is) 69592c4ed94SBarry Smith { 69692c4ed94SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6978a84c255SSatish Balay int *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 698273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 699549d3d68SSatish Balay int *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr; 700273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 70195d5f7c2SBarry Smith MatScalar *value = v,*ap,*aa = a->a,*bap; 70292c4ed94SBarry Smith 7033a40ed3dSBarry Smith PetscFunctionBegin; 7040e324ae4SSatish Balay if (roworiented) { 7050e324ae4SSatish Balay stepval = (n-1)*bs; 7060e324ae4SSatish Balay } else { 7070e324ae4SSatish Balay stepval = (m-1)*bs; 7080e324ae4SSatish Balay } 70992c4ed94SBarry Smith for (k=0; k<m; k++) { /* loop over added rows */ 71092c4ed94SBarry Smith row = im[k]; 7115ef9f2a5SBarry Smith if (row < 0) continue; 712aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 71329bbc08cSBarry Smith if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); 71492c4ed94SBarry Smith #endif 71592c4ed94SBarry Smith rp = aj + ai[row]; 71692c4ed94SBarry Smith ap = aa + bs2*ai[row]; 71792c4ed94SBarry Smith rmax = imax[row]; 71892c4ed94SBarry Smith nrow = ailen[row]; 71992c4ed94SBarry Smith low = 0; 72092c4ed94SBarry Smith for (l=0; l<n; l++) { /* loop over added columns */ 7215ef9f2a5SBarry Smith if (in[l] < 0) continue; 722aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 72329bbc08cSBarry Smith if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); 72492c4ed94SBarry Smith #endif 72592c4ed94SBarry Smith col = in[l]; 72692c4ed94SBarry Smith if (roworiented) { 7270e324ae4SSatish Balay value = v + k*(stepval+bs)*bs + l*bs; 7280e324ae4SSatish Balay } else { 7290e324ae4SSatish Balay value = v + l*(stepval+bs)*bs + k*bs; 73092c4ed94SBarry Smith } 73192c4ed94SBarry Smith if (!sorted) low = 0; high = nrow; 73292c4ed94SBarry Smith while (high-low > 7) { 73392c4ed94SBarry Smith t = (low+high)/2; 73492c4ed94SBarry Smith if (rp[t] > col) high = t; 73592c4ed94SBarry Smith else low = t; 73692c4ed94SBarry Smith } 73792c4ed94SBarry Smith for (i=low; i<high; i++) { 73892c4ed94SBarry Smith if (rp[i] > col) break; 73992c4ed94SBarry Smith if (rp[i] == col) { 7408a84c255SSatish Balay bap = ap + bs2*i; 7410e324ae4SSatish Balay if (roworiented) { 7428a84c255SSatish Balay if (is == ADD_VALUES) { 743dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 744dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7458a84c255SSatish Balay bap[jj] += *value++; 746dd9472c6SBarry Smith } 747dd9472c6SBarry Smith } 7480e324ae4SSatish Balay } else { 749dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 750dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 7510e324ae4SSatish Balay bap[jj] = *value++; 7528a84c255SSatish Balay } 753dd9472c6SBarry Smith } 754dd9472c6SBarry Smith } 7550e324ae4SSatish Balay } else { 7560e324ae4SSatish Balay if (is == ADD_VALUES) { 757dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 758dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7590e324ae4SSatish Balay *bap++ += *value++; 760dd9472c6SBarry Smith } 761dd9472c6SBarry Smith } 7620e324ae4SSatish Balay } else { 763dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 764dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 7650e324ae4SSatish Balay *bap++ = *value++; 7660e324ae4SSatish Balay } 7678a84c255SSatish Balay } 768dd9472c6SBarry Smith } 769dd9472c6SBarry Smith } 770f1241b54SBarry Smith goto noinsert2; 77192c4ed94SBarry Smith } 77292c4ed94SBarry Smith } 77389280ab3SLois Curfman McInnes if (nonew == 1) goto noinsert2; 77429bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 77592c4ed94SBarry Smith if (nrow >= rmax) { 77692c4ed94SBarry Smith /* there is no extra room in row, therefore enlarge */ 77792c4ed94SBarry Smith int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 7783f1db9ecSBarry Smith MatScalar *new_a; 77992c4ed94SBarry Smith 78029bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 78189280ab3SLois Curfman McInnes 78292c4ed94SBarry Smith /* malloc new storage space */ 7833f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 784b0a32e0cSBarry Smith ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 78592c4ed94SBarry Smith new_j = (int*)(new_a + bs2*new_nz); 78692c4ed94SBarry Smith new_i = new_j + new_nz; 78792c4ed94SBarry Smith 78892c4ed94SBarry Smith /* copy over old data into new slots */ 78992c4ed94SBarry Smith for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];} 79092c4ed94SBarry Smith for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 791549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr); 79292c4ed94SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow); 793549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr); 794549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 795549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 796549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 79792c4ed94SBarry Smith /* free up old matrix storage */ 798606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 799606d414cSSatish Balay if (!a->singlemalloc) { 800606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 801606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 802606d414cSSatish Balay } 80392c4ed94SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 804c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 80592c4ed94SBarry Smith 80692c4ed94SBarry Smith rp = aj + ai[row]; ap = aa + bs2*ai[row]; 80792c4ed94SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 808b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 80992c4ed94SBarry Smith a->maxnz += bs2*CHUNKSIZE; 81092c4ed94SBarry Smith a->reallocs++; 81192c4ed94SBarry Smith a->nz++; 81292c4ed94SBarry Smith } 81392c4ed94SBarry Smith N = nrow++ - 1; 81492c4ed94SBarry Smith /* shift up all the later entries in this row */ 81592c4ed94SBarry Smith for (ii=N; ii>=i; ii--) { 81692c4ed94SBarry Smith rp[ii+1] = rp[ii]; 817549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 81892c4ed94SBarry Smith } 819549d3d68SSatish Balay if (N >= i) { 820549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 821549d3d68SSatish Balay } 82292c4ed94SBarry Smith rp[i] = col; 8238a84c255SSatish Balay bap = ap + bs2*i; 8240e324ae4SSatish Balay if (roworiented) { 825dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 826dd9472c6SBarry Smith for (jj=ii; jj<bs2; jj+=bs) { 8270e324ae4SSatish Balay bap[jj] = *value++; 828dd9472c6SBarry Smith } 829dd9472c6SBarry Smith } 8300e324ae4SSatish Balay } else { 831dd9472c6SBarry Smith for (ii=0; ii<bs; ii++,value+=stepval) { 832dd9472c6SBarry Smith for (jj=0; jj<bs; jj++) { 8330e324ae4SSatish Balay *bap++ = *value++; 8340e324ae4SSatish Balay } 835dd9472c6SBarry Smith } 836dd9472c6SBarry Smith } 837f1241b54SBarry Smith noinsert2:; 83892c4ed94SBarry Smith low = i; 83992c4ed94SBarry Smith } 84092c4ed94SBarry Smith ailen[row] = nrow; 84192c4ed94SBarry Smith } 8423a40ed3dSBarry Smith PetscFunctionReturn(0); 84392c4ed94SBarry Smith } 84492c4ed94SBarry Smith 845f2501298SSatish Balay 8464a2ae208SSatish Balay #undef __FUNCT__ 8474a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ" 8488f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 849584200bdSSatish Balay { 850584200bdSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 851584200bdSSatish Balay int fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; 852273d9f13SBarry Smith int m = A->m,*ip,N,*ailen = a->ilen; 853549d3d68SSatish Balay int mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr; 8543f1db9ecSBarry Smith MatScalar *aa = a->a,*ap; 855584200bdSSatish Balay 8563a40ed3dSBarry Smith PetscFunctionBegin; 8573a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 858584200bdSSatish Balay 85943ee02c3SBarry Smith if (m) rmax = ailen[0]; 860584200bdSSatish Balay for (i=1; i<mbs; i++) { 861584200bdSSatish Balay /* move each row back by the amount of empty slots (fshift) before it*/ 862584200bdSSatish Balay fshift += imax[i-1] - ailen[i-1]; 863d402145bSBarry Smith rmax = PetscMax(rmax,ailen[i]); 864584200bdSSatish Balay if (fshift) { 865a7c10996SSatish Balay ip = aj + ai[i]; ap = aa + bs2*ai[i]; 866584200bdSSatish Balay N = ailen[i]; 867584200bdSSatish Balay for (j=0; j<N; j++) { 868584200bdSSatish Balay ip[j-fshift] = ip[j]; 869549d3d68SSatish Balay ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 870584200bdSSatish Balay } 871584200bdSSatish Balay } 872584200bdSSatish Balay ai[i] = ai[i-1] + ailen[i-1]; 873584200bdSSatish Balay } 874584200bdSSatish Balay if (mbs) { 875584200bdSSatish Balay fshift += imax[mbs-1] - ailen[mbs-1]; 876584200bdSSatish Balay ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 877584200bdSSatish Balay } 878584200bdSSatish Balay /* reset ilen and imax for each row */ 879584200bdSSatish Balay for (i=0; i<mbs; i++) { 880584200bdSSatish Balay ailen[i] = imax[i] = ai[i+1] - ai[i]; 881584200bdSSatish Balay } 882a7c10996SSatish Balay a->nz = ai[mbs]; 883584200bdSSatish Balay 884584200bdSSatish Balay /* diagonals may have moved, so kill the diagonal pointers */ 885584200bdSSatish Balay if (fshift && a->diag) { 886606d414cSSatish Balay ierr = PetscFree(a->diag);CHKERRQ(ierr); 887b0a32e0cSBarry Smith PetscLogObjectMemory(A,-(m+1)*sizeof(int)); 888584200bdSSatish Balay a->diag = 0; 889584200bdSSatish Balay } 890b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n", 891273d9f13SBarry Smith m,A->n,a->bs,fshift*bs2,a->nz*bs2); 892b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n", 893584200bdSSatish Balay a->reallocs); 894b0a32e0cSBarry Smith PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax); 895e2f3b5e9SSatish Balay a->reallocs = 0; 8960e6d2581SBarry Smith A->info.nz_unneeded = (PetscReal)fshift*bs2; 8974e220ebcSLois Curfman McInnes 8983a40ed3dSBarry Smith PetscFunctionReturn(0); 899584200bdSSatish Balay } 900584200bdSSatish Balay 9015a838052SSatish Balay 902bea157c4SSatish Balay 903bea157c4SSatish Balay /* 904bea157c4SSatish Balay This function returns an array of flags which indicate the locations of contiguous 905bea157c4SSatish Balay blocks that should be zeroed. for eg: if bs = 3 and is = [0,1,2,3,5,6,7,8,9] 906bea157c4SSatish Balay then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)] 907bea157c4SSatish Balay Assume: sizes should be long enough to hold all the values. 908bea157c4SSatish Balay */ 9094a2ae208SSatish Balay #undef __FUNCT__ 9104a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks" 911bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max) 912d9b7c43dSSatish Balay { 913bea157c4SSatish Balay int i,j,k,row; 914bea157c4SSatish Balay PetscTruth flg; 9153a40ed3dSBarry Smith 916433994e6SBarry Smith PetscFunctionBegin; 917bea157c4SSatish Balay for (i=0,j=0; i<n; j++) { 918bea157c4SSatish Balay row = idx[i]; 919bea157c4SSatish Balay if (row%bs!=0) { /* Not the begining of a block */ 920bea157c4SSatish Balay sizes[j] = 1; 921bea157c4SSatish Balay i++; 922e4fda26cSSatish Balay } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */ 923bea157c4SSatish Balay sizes[j] = 1; /* Also makes sure atleast 'bs' values exist for next else */ 924bea157c4SSatish Balay i++; 925bea157c4SSatish Balay } else { /* Begining of the block, so check if the complete block exists */ 926bea157c4SSatish Balay flg = PETSC_TRUE; 927bea157c4SSatish Balay for (k=1; k<bs; k++) { 928bea157c4SSatish Balay if (row+k != idx[i+k]) { /* break in the block */ 929bea157c4SSatish Balay flg = PETSC_FALSE; 930bea157c4SSatish Balay break; 931d9b7c43dSSatish Balay } 932bea157c4SSatish Balay } 933bea157c4SSatish Balay if (flg == PETSC_TRUE) { /* No break in the bs */ 934bea157c4SSatish Balay sizes[j] = bs; 935bea157c4SSatish Balay i+= bs; 936bea157c4SSatish Balay } else { 937bea157c4SSatish Balay sizes[j] = 1; 938bea157c4SSatish Balay i++; 939bea157c4SSatish Balay } 940bea157c4SSatish Balay } 941bea157c4SSatish Balay } 942bea157c4SSatish Balay *bs_max = j; 9433a40ed3dSBarry Smith PetscFunctionReturn(0); 944d9b7c43dSSatish Balay } 945d9b7c43dSSatish Balay 9464a2ae208SSatish Balay #undef __FUNCT__ 9474a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ" 9488f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is,Scalar *diag) 949d9b7c43dSSatish Balay { 950d9b7c43dSSatish Balay Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data; 951b6815fffSBarry Smith int ierr,i,j,k,count,is_n,*is_idx,*rows; 952bea157c4SSatish Balay int bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max; 9533f1db9ecSBarry Smith Scalar zero = 0.0; 9543f1db9ecSBarry Smith MatScalar *aa; 955d9b7c43dSSatish Balay 9563a40ed3dSBarry Smith PetscFunctionBegin; 957d9b7c43dSSatish Balay /* Make a copy of the IS and sort it */ 958b9b97703SBarry Smith ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr); 959d9b7c43dSSatish Balay ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr); 960d9b7c43dSSatish Balay 961bea157c4SSatish Balay /* allocate memory for rows,sizes */ 962b0a32e0cSBarry Smith ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr); 963bea157c4SSatish Balay sizes = rows + is_n; 964bea157c4SSatish Balay 965563b5814SBarry Smith /* copy IS values to rows, and sort them */ 966bea157c4SSatish Balay for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; } 967bea157c4SSatish Balay ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr); 968dffd3267SBarry Smith if (baij->keepzeroedrows) { 969dffd3267SBarry Smith for (i=0; i<is_n; i++) { sizes[i] = 1; } 970dffd3267SBarry Smith bs_max = is_n; 971dffd3267SBarry Smith } else { 972bea157c4SSatish Balay ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr); 973dffd3267SBarry Smith } 974b97a56abSSatish Balay ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr); 975bea157c4SSatish Balay 976bea157c4SSatish Balay for (i=0,j=0; i<bs_max; j+=sizes[i],i++) { 977bea157c4SSatish Balay row = rows[j]; 978273d9f13SBarry Smith if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row); 979bea157c4SSatish Balay count = (baij->i[row/bs +1] - baij->i[row/bs])*bs; 980bea157c4SSatish Balay aa = baij->a + baij->i[row/bs]*bs2 + (row%bs); 981dffd3267SBarry Smith if (sizes[i] == bs && !baij->keepzeroedrows) { 982bea157c4SSatish Balay if (diag) { 983bea157c4SSatish Balay if (baij->ilen[row/bs] > 0) { 984bea157c4SSatish Balay baij->ilen[row/bs] = 1; 985bea157c4SSatish Balay baij->j[baij->i[row/bs]] = row/bs; 986bea157c4SSatish Balay ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr); 987a07cd24cSSatish Balay } 988563b5814SBarry Smith /* Now insert all the diagonal values for this bs */ 989bea157c4SSatish Balay for (k=0; k<bs; k++) { 990bea157c4SSatish Balay ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr); 991bea157c4SSatish Balay } 992bea157c4SSatish Balay } else { /* (!diag) */ 993bea157c4SSatish Balay baij->ilen[row/bs] = 0; 994bea157c4SSatish Balay } /* end (!diag) */ 995bea157c4SSatish Balay } else { /* (sizes[i] != bs) */ 996aa482453SBarry Smith #if defined (PETSC_USE_DEBUG) 99729bbc08cSBarry Smith if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1"); 998bea157c4SSatish Balay #endif 999bea157c4SSatish Balay for (k=0; k<count; k++) { 1000d9b7c43dSSatish Balay aa[0] = zero; 1001d9b7c43dSSatish Balay aa += bs; 1002d9b7c43dSSatish Balay } 1003d9b7c43dSSatish Balay if (diag) { 1004f830108cSBarry Smith ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr); 1005d9b7c43dSSatish Balay } 1006d9b7c43dSSatish Balay } 1007bea157c4SSatish Balay } 1008bea157c4SSatish Balay 1009606d414cSSatish Balay ierr = PetscFree(rows);CHKERRQ(ierr); 10109a8dea36SBarry Smith ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10113a40ed3dSBarry Smith PetscFunctionReturn(0); 1012d9b7c43dSSatish Balay } 10131c351548SSatish Balay 10144a2ae208SSatish Balay #undef __FUNCT__ 10154a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ" 10162d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 10172d61bbb3SSatish Balay { 10182d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 10192d61bbb3SSatish Balay int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 1020273d9f13SBarry Smith int *imax=a->imax,*ai=a->i,*ailen=a->ilen; 10212d61bbb3SSatish Balay int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 1022549d3d68SSatish Balay int ridx,cidx,bs2=a->bs2,ierr; 1023273d9f13SBarry Smith PetscTruth roworiented=a->roworiented; 10243f1db9ecSBarry Smith MatScalar *ap,value,*aa=a->a,*bap; 10252d61bbb3SSatish Balay 10262d61bbb3SSatish Balay PetscFunctionBegin; 10272d61bbb3SSatish Balay for (k=0; k<m; k++) { /* loop over added rows */ 10282d61bbb3SSatish Balay row = im[k]; brow = row/bs; 10295ef9f2a5SBarry Smith if (row < 0) continue; 1030aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1031273d9f13SBarry Smith if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m); 10322d61bbb3SSatish Balay #endif 10332d61bbb3SSatish Balay rp = aj + ai[brow]; 10342d61bbb3SSatish Balay ap = aa + bs2*ai[brow]; 10352d61bbb3SSatish Balay rmax = imax[brow]; 10362d61bbb3SSatish Balay nrow = ailen[brow]; 10372d61bbb3SSatish Balay low = 0; 10382d61bbb3SSatish Balay for (l=0; l<n; l++) { /* loop over added columns */ 10395ef9f2a5SBarry Smith if (in[l] < 0) continue; 1040aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 1041273d9f13SBarry Smith if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n); 10422d61bbb3SSatish Balay #endif 10432d61bbb3SSatish Balay col = in[l]; bcol = col/bs; 10442d61bbb3SSatish Balay ridx = row % bs; cidx = col % bs; 10452d61bbb3SSatish Balay if (roworiented) { 10465ef9f2a5SBarry Smith value = v[l + k*n]; 10472d61bbb3SSatish Balay } else { 10482d61bbb3SSatish Balay value = v[k + l*m]; 10492d61bbb3SSatish Balay } 10502d61bbb3SSatish Balay if (!sorted) low = 0; high = nrow; 10512d61bbb3SSatish Balay while (high-low > 7) { 10522d61bbb3SSatish Balay t = (low+high)/2; 10532d61bbb3SSatish Balay if (rp[t] > bcol) high = t; 10542d61bbb3SSatish Balay else low = t; 10552d61bbb3SSatish Balay } 10562d61bbb3SSatish Balay for (i=low; i<high; i++) { 10572d61bbb3SSatish Balay if (rp[i] > bcol) break; 10582d61bbb3SSatish Balay if (rp[i] == bcol) { 10592d61bbb3SSatish Balay bap = ap + bs2*i + bs*cidx + ridx; 10602d61bbb3SSatish Balay if (is == ADD_VALUES) *bap += value; 10612d61bbb3SSatish Balay else *bap = value; 10622d61bbb3SSatish Balay goto noinsert1; 10632d61bbb3SSatish Balay } 10642d61bbb3SSatish Balay } 10652d61bbb3SSatish Balay if (nonew == 1) goto noinsert1; 106629bbc08cSBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 10672d61bbb3SSatish Balay if (nrow >= rmax) { 10682d61bbb3SSatish Balay /* there is no extra room in row, therefore enlarge */ 10692d61bbb3SSatish Balay int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 10703f1db9ecSBarry Smith MatScalar *new_a; 10712d61bbb3SSatish Balay 107229bbc08cSBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); 10732d61bbb3SSatish Balay 10742d61bbb3SSatish Balay /* Malloc new storage space */ 10753f1db9ecSBarry Smith len = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int); 1076b0a32e0cSBarry Smith ierr = PetscMalloc(len,&new_a);CHKERRQ(ierr); 10772d61bbb3SSatish Balay new_j = (int*)(new_a + bs2*new_nz); 10782d61bbb3SSatish Balay new_i = new_j + new_nz; 10792d61bbb3SSatish Balay 10802d61bbb3SSatish Balay /* copy over old data into new slots */ 10812d61bbb3SSatish Balay for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];} 10822d61bbb3SSatish Balay for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;} 1083549d3d68SSatish Balay ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr); 10842d61bbb3SSatish Balay len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 1085549d3d68SSatish Balay ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr); 1086549d3d68SSatish Balay ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1087549d3d68SSatish Balay ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr); 1088549d3d68SSatish Balay ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr); 10892d61bbb3SSatish Balay /* free up old matrix storage */ 1090606d414cSSatish Balay ierr = PetscFree(a->a);CHKERRQ(ierr); 1091606d414cSSatish Balay if (!a->singlemalloc) { 1092606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 1093606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 1094606d414cSSatish Balay } 10952d61bbb3SSatish Balay aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 1096c4992f7dSBarry Smith a->singlemalloc = PETSC_TRUE; 10972d61bbb3SSatish Balay 10982d61bbb3SSatish Balay rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 10992d61bbb3SSatish Balay rmax = imax[brow] = imax[brow] + CHUNKSIZE; 1100b0a32e0cSBarry Smith PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar))); 11012d61bbb3SSatish Balay a->maxnz += bs2*CHUNKSIZE; 11022d61bbb3SSatish Balay a->reallocs++; 11032d61bbb3SSatish Balay a->nz++; 11042d61bbb3SSatish Balay } 11052d61bbb3SSatish Balay N = nrow++ - 1; 11062d61bbb3SSatish Balay /* shift up all the later entries in this row */ 11072d61bbb3SSatish Balay for (ii=N; ii>=i; ii--) { 11082d61bbb3SSatish Balay rp[ii+1] = rp[ii]; 1109549d3d68SSatish Balay ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); 11102d61bbb3SSatish Balay } 1111549d3d68SSatish Balay if (N>=i) { 1112549d3d68SSatish Balay ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1113549d3d68SSatish Balay } 11142d61bbb3SSatish Balay rp[i] = bcol; 11152d61bbb3SSatish Balay ap[bs2*i + bs*cidx + ridx] = value; 11162d61bbb3SSatish Balay noinsert1:; 11172d61bbb3SSatish Balay low = i; 11182d61bbb3SSatish Balay } 11192d61bbb3SSatish Balay ailen[brow] = nrow; 11202d61bbb3SSatish Balay } 11212d61bbb3SSatish Balay PetscFunctionReturn(0); 11222d61bbb3SSatish Balay } 11232d61bbb3SSatish Balay 1124b9b97703SBarry Smith EXTERN int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,MatLUInfo*,Mat*); 1125b9b97703SBarry Smith EXTERN int MatLUFactor_SeqBAIJ(Mat,IS,IS,MatLUInfo*); 1126ca44d042SBarry Smith EXTERN int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1127ca44d042SBarry Smith EXTERN int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*); 1128ca44d042SBarry Smith EXTERN int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**); 1129ca44d042SBarry Smith EXTERN int MatMultTranspose_SeqBAIJ(Mat,Vec,Vec); 1130ca44d042SBarry Smith EXTERN int MatMultTransposeAdd_SeqBAIJ(Mat,Vec,Vec,Vec); 1131ca44d042SBarry Smith EXTERN int MatScale_SeqBAIJ(Scalar*,Mat); 1132ca44d042SBarry Smith EXTERN int MatNorm_SeqBAIJ(Mat,NormType,PetscReal *); 1133ca44d042SBarry Smith EXTERN int MatEqual_SeqBAIJ(Mat,Mat,PetscTruth*); 1134ca44d042SBarry Smith EXTERN int MatGetDiagonal_SeqBAIJ(Mat,Vec); 1135ca44d042SBarry Smith EXTERN int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec); 1136ca44d042SBarry Smith EXTERN int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *); 1137ca44d042SBarry Smith EXTERN int MatZeroEntries_SeqBAIJ(Mat); 11382d61bbb3SSatish Balay 1139ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 1140ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 1141ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 1142ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 1143ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 1144ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 1145ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_6(Mat,Vec,Vec); 1146ca44d042SBarry Smith EXTERN int MatSolve_SeqBAIJ_7(Mat,Vec,Vec); 1147ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_7(Mat,Vec,Vec); 1148ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_6(Mat,Vec,Vec); 1149ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_5(Mat,Vec,Vec); 1150ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_4(Mat,Vec,Vec); 1151ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_3(Mat,Vec,Vec); 1152ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_2(Mat,Vec,Vec); 1153ca44d042SBarry Smith EXTERN int MatSolveTranspose_SeqBAIJ_1(Mat,Vec,Vec); 11542d61bbb3SSatish Balay 1155ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 1156ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 1157ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 1158ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 1159ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 1160ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 1161ca44d042SBarry Smith EXTERN int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*); 11622d61bbb3SSatish Balay 1163ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_1(Mat,Vec,Vec); 1164ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_2(Mat,Vec,Vec); 1165ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_3(Mat,Vec,Vec); 1166ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_4(Mat,Vec,Vec); 1167ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_5(Mat,Vec,Vec); 1168ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_6(Mat,Vec,Vec); 1169ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_7(Mat,Vec,Vec); 1170ca44d042SBarry Smith EXTERN int MatMult_SeqBAIJ_N(Mat,Vec,Vec); 11712d61bbb3SSatish Balay 1172ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec); 1173ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec); 1174ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec); 1175ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec); 1176ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec); 1177ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec); 1178ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec); 1179ca44d042SBarry Smith EXTERN int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec); 11802d61bbb3SSatish Balay 11814a2ae208SSatish Balay #undef __FUNCT__ 11824a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ" 11835ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 11842d61bbb3SSatish Balay { 11852d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 11862d61bbb3SSatish Balay Mat outA; 11872d61bbb3SSatish Balay int ierr; 1188667159a5SBarry Smith PetscTruth row_identity,col_identity; 11892d61bbb3SSatish Balay 11902d61bbb3SSatish Balay PetscFunctionBegin; 119129bbc08cSBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU"); 1192667159a5SBarry Smith ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); 1193667159a5SBarry Smith ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); 1194667159a5SBarry Smith if (!row_identity || !col_identity) { 119529bbc08cSBarry Smith SETERRQ(1,"Row and column permutations must be identity for in-place ILU"); 1196667159a5SBarry Smith } 11972d61bbb3SSatish Balay 11982d61bbb3SSatish Balay outA = inA; 11992d61bbb3SSatish Balay inA->factor = FACTOR_LU; 12002d61bbb3SSatish Balay 12012d61bbb3SSatish Balay if (!a->diag) { 1202c4992f7dSBarry Smith ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr); 12032d61bbb3SSatish Balay } 12042d61bbb3SSatish Balay /* 120515091d37SBarry Smith Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver 120615091d37SBarry Smith for ILU(0) factorization with natural ordering 12072d61bbb3SSatish Balay */ 1208*8661488fSKris Buschelman if (a->bs < 8) { 1209*8661488fSKris Buschelman ierr = MatSeqBAIJ_UpdateFactorizerSolver_NaturalOrdering(inA);CHKERRQ(ierr); 12103477af42SKris Buschelman } else { 1211c38d4ed2SBarry Smith a->row = row; 1212c38d4ed2SBarry Smith a->col = col; 1213c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); 1214c38d4ed2SBarry Smith ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); 1215c38d4ed2SBarry Smith 1216c38d4ed2SBarry Smith /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 12174c49b128SBarry Smith ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); 1218b0a32e0cSBarry Smith PetscLogObjectParent(inA,a->icol); 1219c38d4ed2SBarry Smith 1220c38d4ed2SBarry Smith if (!a->solve_work) { 1221b0a32e0cSBarry Smith ierr = PetscMalloc((inA->m+a->bs)*sizeof(Scalar),&a->solve_work);CHKERRQ(ierr); 1222b0a32e0cSBarry Smith PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(Scalar)); 1223c38d4ed2SBarry Smith } 12242d61bbb3SSatish Balay } 12252d61bbb3SSatish Balay 1226667159a5SBarry Smith ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr); 1227667159a5SBarry Smith 12282d61bbb3SSatish Balay PetscFunctionReturn(0); 12292d61bbb3SSatish Balay } 12304a2ae208SSatish Balay #undef __FUNCT__ 12314a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ" 1232ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A) 1233ba4ca20aSSatish Balay { 1234c38d4ed2SBarry Smith static PetscTruth called = PETSC_FALSE; 1235ba4ca20aSSatish Balay MPI_Comm comm = A->comm; 1236d132466eSBarry Smith int ierr; 1237ba4ca20aSSatish Balay 12383a40ed3dSBarry Smith PetscFunctionBegin; 1239c38d4ed2SBarry Smith if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE; 1240d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr); 1241d132466eSBarry Smith ierr = (*PetscHelpPrintf)(comm," -mat_block_size <block_size>\n");CHKERRQ(ierr); 12423a40ed3dSBarry Smith PetscFunctionReturn(0); 1243ba4ca20aSSatish Balay } 1244d9b7c43dSSatish Balay 1245fb2e594dSBarry Smith EXTERN_C_BEGIN 12464a2ae208SSatish Balay #undef __FUNCT__ 12474a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ" 124827a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices) 124927a8da17SBarry Smith { 125027a8da17SBarry Smith Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data; 125127a8da17SBarry Smith int i,nz,n; 125227a8da17SBarry Smith 125327a8da17SBarry Smith PetscFunctionBegin; 125427a8da17SBarry Smith nz = baij->maxnz; 1255273d9f13SBarry Smith n = mat->n; 125627a8da17SBarry Smith for (i=0; i<nz; i++) { 125727a8da17SBarry Smith baij->j[i] = indices[i]; 125827a8da17SBarry Smith } 125927a8da17SBarry Smith baij->nz = nz; 126027a8da17SBarry Smith for (i=0; i<n; i++) { 126127a8da17SBarry Smith baij->ilen[i] = baij->imax[i]; 126227a8da17SBarry Smith } 126327a8da17SBarry Smith 126427a8da17SBarry Smith PetscFunctionReturn(0); 126527a8da17SBarry Smith } 1266fb2e594dSBarry Smith EXTERN_C_END 126727a8da17SBarry Smith 12684a2ae208SSatish Balay #undef __FUNCT__ 12694a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices" 127027a8da17SBarry Smith /*@ 127127a8da17SBarry Smith MatSeqBAIJSetColumnIndices - Set the column indices for all the rows 127227a8da17SBarry Smith in the matrix. 127327a8da17SBarry Smith 127427a8da17SBarry Smith Input Parameters: 127527a8da17SBarry Smith + mat - the SeqBAIJ matrix 127627a8da17SBarry Smith - indices - the column indices 127727a8da17SBarry Smith 127815091d37SBarry Smith Level: advanced 127915091d37SBarry Smith 128027a8da17SBarry Smith Notes: 128127a8da17SBarry Smith This can be called if you have precomputed the nonzero structure of the 128227a8da17SBarry Smith matrix and want to provide it to the matrix object to improve the performance 128327a8da17SBarry Smith of the MatSetValues() operation. 128427a8da17SBarry Smith 128527a8da17SBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 128627a8da17SBarry Smith MatCreateSeqBAIJ(). 128727a8da17SBarry Smith 128827a8da17SBarry Smith MUST be called before any calls to MatSetValues(); 128927a8da17SBarry Smith 129027a8da17SBarry Smith @*/ 129127a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices) 129227a8da17SBarry Smith { 129327a8da17SBarry Smith int ierr,(*f)(Mat,int *); 129427a8da17SBarry Smith 129527a8da17SBarry Smith PetscFunctionBegin; 129627a8da17SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1297b9617806SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)())&f);CHKERRQ(ierr); 129827a8da17SBarry Smith if (f) { 129927a8da17SBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 130027a8da17SBarry Smith } else { 130129bbc08cSBarry Smith SETERRQ(1,"Wrong type of matrix to set column indices"); 130227a8da17SBarry Smith } 130327a8da17SBarry Smith PetscFunctionReturn(0); 130427a8da17SBarry Smith } 130527a8da17SBarry Smith 13064a2ae208SSatish Balay #undef __FUNCT__ 13074a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ" 1308273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v) 1309273d9f13SBarry Smith { 1310273d9f13SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1311273d9f13SBarry Smith int ierr,i,j,n,row,bs,*ai,*aj,mbs; 1312273d9f13SBarry Smith PetscReal atmp; 1313273d9f13SBarry Smith Scalar *x,zero = 0.0; 1314273d9f13SBarry Smith MatScalar *aa; 1315273d9f13SBarry Smith int ncols,brow,krow,kcol; 1316273d9f13SBarry Smith 1317273d9f13SBarry Smith PetscFunctionBegin; 1318273d9f13SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1319273d9f13SBarry Smith bs = a->bs; 1320273d9f13SBarry Smith aa = a->a; 1321273d9f13SBarry Smith ai = a->i; 1322273d9f13SBarry Smith aj = a->j; 1323273d9f13SBarry Smith mbs = a->mbs; 1324273d9f13SBarry Smith 1325273d9f13SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 1326273d9f13SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1327273d9f13SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1328273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1329273d9f13SBarry Smith for (i=0; i<mbs; i++) { 1330273d9f13SBarry Smith ncols = ai[1] - ai[0]; ai++; 1331273d9f13SBarry Smith brow = bs*i; 1332273d9f13SBarry Smith for (j=0; j<ncols; j++){ 1333273d9f13SBarry Smith /* bcol = bs*(*aj); */ 1334273d9f13SBarry Smith for (kcol=0; kcol<bs; kcol++){ 1335273d9f13SBarry Smith for (krow=0; krow<bs; krow++){ 1336273d9f13SBarry Smith atmp = PetscAbsScalar(*aa); aa++; 1337273d9f13SBarry Smith row = brow + krow; /* row index */ 1338273d9f13SBarry Smith /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */ 1339273d9f13SBarry Smith if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp; 1340273d9f13SBarry Smith } 1341273d9f13SBarry Smith } 1342273d9f13SBarry Smith aj++; 1343273d9f13SBarry Smith } 1344273d9f13SBarry Smith } 1345273d9f13SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1346273d9f13SBarry Smith PetscFunctionReturn(0); 1347273d9f13SBarry Smith } 1348273d9f13SBarry Smith 13494a2ae208SSatish Balay #undef __FUNCT__ 13504a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ" 1351273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A) 1352273d9f13SBarry Smith { 1353273d9f13SBarry Smith int ierr; 1354273d9f13SBarry Smith 1355273d9f13SBarry Smith PetscFunctionBegin; 1356273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr); 1357273d9f13SBarry Smith PetscFunctionReturn(0); 1358273d9f13SBarry Smith } 1359273d9f13SBarry Smith 13604a2ae208SSatish Balay #undef __FUNCT__ 13614a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ" 1362f2a5309cSSatish Balay int MatGetArray_SeqBAIJ(Mat A,Scalar **array) 1363f2a5309cSSatish Balay { 1364f2a5309cSSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1365f2a5309cSSatish Balay PetscFunctionBegin; 1366f2a5309cSSatish Balay *array = a->a; 1367f2a5309cSSatish Balay PetscFunctionReturn(0); 1368f2a5309cSSatish Balay } 1369f2a5309cSSatish Balay 13704a2ae208SSatish Balay #undef __FUNCT__ 13714a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ" 1372f2a5309cSSatish Balay int MatRestoreArray_SeqBAIJ(Mat A,Scalar **array) 1373f2a5309cSSatish Balay { 1374f2a5309cSSatish Balay PetscFunctionBegin; 1375f2a5309cSSatish Balay PetscFunctionReturn(0); 1376f2a5309cSSatish Balay } 1377f2a5309cSSatish Balay 13782593348eSBarry Smith /* -------------------------------------------------------------------*/ 1379cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ, 1380cc2dc46cSBarry Smith MatGetRow_SeqBAIJ, 1381cc2dc46cSBarry Smith MatRestoreRow_SeqBAIJ, 1382cc2dc46cSBarry Smith MatMult_SeqBAIJ_N, 1383cc2dc46cSBarry Smith MatMultAdd_SeqBAIJ_N, 13847c922b88SBarry Smith MatMultTranspose_SeqBAIJ, 13857c922b88SBarry Smith MatMultTransposeAdd_SeqBAIJ, 1386cc2dc46cSBarry Smith MatSolve_SeqBAIJ_N, 1387cc2dc46cSBarry Smith 0, 1388cc2dc46cSBarry Smith 0, 1389cc2dc46cSBarry Smith 0, 1390cc2dc46cSBarry Smith MatLUFactor_SeqBAIJ, 1391cc2dc46cSBarry Smith 0, 1392b6490206SBarry Smith 0, 1393f2501298SSatish Balay MatTranspose_SeqBAIJ, 1394cc2dc46cSBarry Smith MatGetInfo_SeqBAIJ, 1395cc2dc46cSBarry Smith MatEqual_SeqBAIJ, 1396cc2dc46cSBarry Smith MatGetDiagonal_SeqBAIJ, 1397cc2dc46cSBarry Smith MatDiagonalScale_SeqBAIJ, 1398cc2dc46cSBarry Smith MatNorm_SeqBAIJ, 1399b6490206SBarry Smith 0, 1400cc2dc46cSBarry Smith MatAssemblyEnd_SeqBAIJ, 1401cc2dc46cSBarry Smith 0, 1402cc2dc46cSBarry Smith MatSetOption_SeqBAIJ, 1403cc2dc46cSBarry Smith MatZeroEntries_SeqBAIJ, 1404cc2dc46cSBarry Smith MatZeroRows_SeqBAIJ, 1405cc2dc46cSBarry Smith MatLUFactorSymbolic_SeqBAIJ, 1406cc2dc46cSBarry Smith MatLUFactorNumeric_SeqBAIJ_N, 1407cc2dc46cSBarry Smith 0, 1408cc2dc46cSBarry Smith 0, 1409273d9f13SBarry Smith MatSetUpPreallocation_SeqBAIJ, 1410273d9f13SBarry Smith 0, 1411cc2dc46cSBarry Smith MatGetOwnershipRange_SeqBAIJ, 1412cc2dc46cSBarry Smith MatILUFactorSymbolic_SeqBAIJ, 1413cc2dc46cSBarry Smith 0, 1414f2a5309cSSatish Balay MatGetArray_SeqBAIJ, 1415f2a5309cSSatish Balay MatRestoreArray_SeqBAIJ, 14162e8a6d31SBarry Smith MatDuplicate_SeqBAIJ, 1417cc2dc46cSBarry Smith 0, 1418cc2dc46cSBarry Smith 0, 1419cc2dc46cSBarry Smith MatILUFactor_SeqBAIJ, 1420cc2dc46cSBarry Smith 0, 1421cc2dc46cSBarry Smith 0, 1422cc2dc46cSBarry Smith MatGetSubMatrices_SeqBAIJ, 1423cc2dc46cSBarry Smith MatIncreaseOverlap_SeqBAIJ, 1424cc2dc46cSBarry Smith MatGetValues_SeqBAIJ, 1425cc2dc46cSBarry Smith 0, 1426cc2dc46cSBarry Smith MatPrintHelp_SeqBAIJ, 1427cc2dc46cSBarry Smith MatScale_SeqBAIJ, 1428cc2dc46cSBarry Smith 0, 1429cc2dc46cSBarry Smith 0, 1430cc2dc46cSBarry Smith 0, 1431cc2dc46cSBarry Smith MatGetBlockSize_SeqBAIJ, 14323b2fbd54SBarry Smith MatGetRowIJ_SeqBAIJ, 143392c4ed94SBarry Smith MatRestoreRowIJ_SeqBAIJ, 1434cc2dc46cSBarry Smith 0, 1435cc2dc46cSBarry Smith 0, 1436cc2dc46cSBarry Smith 0, 1437cc2dc46cSBarry Smith 0, 1438cc2dc46cSBarry Smith 0, 1439cc2dc46cSBarry Smith 0, 1440d3825aa8SBarry Smith MatSetValuesBlocked_SeqBAIJ, 1441cc2dc46cSBarry Smith MatGetSubMatrix_SeqBAIJ, 1442b9b97703SBarry Smith MatDestroy_SeqBAIJ, 1443b9b97703SBarry Smith MatView_SeqBAIJ, 1444273d9f13SBarry Smith MatGetMaps_Petsc, 1445273d9f13SBarry Smith 0, 1446273d9f13SBarry Smith 0, 1447273d9f13SBarry Smith 0, 1448273d9f13SBarry Smith 0, 1449273d9f13SBarry Smith 0, 1450273d9f13SBarry Smith 0, 1451273d9f13SBarry Smith MatGetRowMax_SeqBAIJ, 1452273d9f13SBarry Smith MatConvert_Basic}; 14532593348eSBarry Smith 14543e90b805SBarry Smith EXTERN_C_BEGIN 14554a2ae208SSatish Balay #undef __FUNCT__ 14564a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ" 14573e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat) 14583e90b805SBarry Smith { 14593e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1460273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2; 1461d132466eSBarry Smith int ierr; 14623e90b805SBarry Smith 14633e90b805SBarry Smith PetscFunctionBegin; 14643e90b805SBarry Smith if (aij->nonew != 1) { 146529bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14663e90b805SBarry Smith } 14673e90b805SBarry Smith 14683e90b805SBarry Smith /* allocate space for values if not already there */ 14693e90b805SBarry Smith if (!aij->saved_values) { 1470f2a5309cSSatish Balay ierr = PetscMalloc((nz+1)*sizeof(Scalar),&aij->saved_values);CHKERRQ(ierr); 14713e90b805SBarry Smith } 14723e90b805SBarry Smith 14733e90b805SBarry Smith /* copy values over */ 1474d132466eSBarry Smith ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr); 14753e90b805SBarry Smith PetscFunctionReturn(0); 14763e90b805SBarry Smith } 14773e90b805SBarry Smith EXTERN_C_END 14783e90b805SBarry Smith 14793e90b805SBarry Smith EXTERN_C_BEGIN 14804a2ae208SSatish Balay #undef __FUNCT__ 14814a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ" 148232e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat) 14833e90b805SBarry Smith { 14843e90b805SBarry Smith Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data; 1485273d9f13SBarry Smith int nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr; 14863e90b805SBarry Smith 14873e90b805SBarry Smith PetscFunctionBegin; 14883e90b805SBarry Smith if (aij->nonew != 1) { 148929bbc08cSBarry Smith SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 14903e90b805SBarry Smith } 14913e90b805SBarry Smith if (!aij->saved_values) { 149229bbc08cSBarry Smith SETERRQ(1,"Must call MatStoreValues(A);first"); 14933e90b805SBarry Smith } 14943e90b805SBarry Smith 14953e90b805SBarry Smith /* copy values over */ 1496549d3d68SSatish Balay ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr); 14973e90b805SBarry Smith PetscFunctionReturn(0); 14983e90b805SBarry Smith } 14993e90b805SBarry Smith EXTERN_C_END 15003e90b805SBarry Smith 1501273d9f13SBarry Smith EXTERN_C_BEGIN 1502273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*); 1503273d9f13SBarry Smith EXTERN_C_END 1504273d9f13SBarry Smith 1505273d9f13SBarry Smith EXTERN_C_BEGIN 15064a2ae208SSatish Balay #undef __FUNCT__ 15074a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ" 1508273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B) 15092593348eSBarry Smith { 1510273d9f13SBarry Smith int ierr,size; 1511b6490206SBarry Smith Mat_SeqBAIJ *b; 15123b2fbd54SBarry Smith 15133a40ed3dSBarry Smith PetscFunctionBegin; 1514273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 151529bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1516b6490206SBarry Smith 1517273d9f13SBarry Smith B->m = B->M = PetscMax(B->m,B->M); 1518273d9f13SBarry Smith B->n = B->N = PetscMax(B->n,B->N); 1519b0a32e0cSBarry Smith ierr = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr); 1520b0a32e0cSBarry Smith B->data = (void*)b; 1521549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr); 1522549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 15232593348eSBarry Smith B->factor = 0; 15242593348eSBarry Smith B->lupivotthreshold = 1.0; 152590f02eecSBarry Smith B->mapping = 0; 15262593348eSBarry Smith b->row = 0; 15272593348eSBarry Smith b->col = 0; 1528e51c0b9cSSatish Balay b->icol = 0; 15292593348eSBarry Smith b->reallocs = 0; 15303e90b805SBarry Smith b->saved_values = 0; 1531cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 1532563b5814SBarry Smith b->setvalueslen = 0; 1533563b5814SBarry Smith b->setvaluescopy = PETSC_NULL; 1534563b5814SBarry Smith #endif 1535*8661488fSKris Buschelman b->single_precision_solves = PETSC_FALSE; 15362593348eSBarry Smith 1537273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 1538273d9f13SBarry Smith ierr = MapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1539a5ae1ecdSBarry Smith 1540c4992f7dSBarry Smith b->sorted = PETSC_FALSE; 1541c4992f7dSBarry Smith b->roworiented = PETSC_TRUE; 15422593348eSBarry Smith b->nonew = 0; 15432593348eSBarry Smith b->diag = 0; 15442593348eSBarry Smith b->solve_work = 0; 1545de6a44a3SBarry Smith b->mult_work = 0; 15462593348eSBarry Smith b->spptr = 0; 15470e6d2581SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 1548883fce79SBarry Smith b->keepzeroedrows = PETSC_FALSE; 15494e220ebcSLois Curfman McInnes 1550f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", 15513e90b805SBarry Smith "MatStoreValues_SeqBAIJ", 1552bc4b532fSSatish Balay MatStoreValues_SeqBAIJ);CHKERRQ(ierr); 1553f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", 15543e90b805SBarry Smith "MatRetrieveValues_SeqBAIJ", 1555bc4b532fSSatish Balay MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr); 1556f1af5d2fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C", 155727a8da17SBarry Smith "MatSeqBAIJSetColumnIndices_SeqBAIJ", 1558bc4b532fSSatish Balay MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr); 1559273d9f13SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_SeqBAIJ_SeqAIJ_C", 1560273d9f13SBarry Smith "MatConvert_SeqBAIJ_SeqAIJ", 1561273d9f13SBarry Smith MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr); 15623a40ed3dSBarry Smith PetscFunctionReturn(0); 15632593348eSBarry Smith } 1564273d9f13SBarry Smith EXTERN_C_END 15652593348eSBarry Smith 15664a2ae208SSatish Balay #undef __FUNCT__ 15674a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ" 15682e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 15692593348eSBarry Smith { 15702593348eSBarry Smith Mat C; 1571b6490206SBarry Smith Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data; 157227a8da17SBarry Smith int i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr; 1573de6a44a3SBarry Smith 15743a40ed3dSBarry Smith PetscFunctionBegin; 157529bbc08cSBarry Smith if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix"); 15762593348eSBarry Smith 15772593348eSBarry Smith *B = 0; 1578273d9f13SBarry Smith ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr); 1579273d9f13SBarry Smith ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr); 1580273d9f13SBarry Smith c = (Mat_SeqBAIJ*)C->data; 158144cd7ae7SLois Curfman McInnes 1582b6490206SBarry Smith c->bs = a->bs; 15839df24120SSatish Balay c->bs2 = a->bs2; 1584b6490206SBarry Smith c->mbs = a->mbs; 1585b6490206SBarry Smith c->nbs = a->nbs; 158635d8aa7fSBarry Smith ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 15872593348eSBarry Smith 1588b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr); 1589b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr); 1590b6490206SBarry Smith for (i=0; i<mbs; i++) { 15912593348eSBarry Smith c->imax[i] = a->imax[i]; 15922593348eSBarry Smith c->ilen[i] = a->ilen[i]; 15932593348eSBarry Smith } 15942593348eSBarry Smith 15952593348eSBarry Smith /* allocate the matrix space */ 1596c4992f7dSBarry Smith c->singlemalloc = PETSC_TRUE; 15973f1db9ecSBarry Smith len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int)); 1598b0a32e0cSBarry Smith ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr); 15997e67e3f9SSatish Balay c->j = (int*)(c->a + nz*bs2); 1600de6a44a3SBarry Smith c->i = c->j + nz; 1601549d3d68SSatish Balay ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr); 1602b6490206SBarry Smith if (mbs > 0) { 1603549d3d68SSatish Balay ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr); 16042e8a6d31SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 1605549d3d68SSatish Balay ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16062e8a6d31SBarry Smith } else { 1607549d3d68SSatish Balay ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr); 16082593348eSBarry Smith } 16092593348eSBarry Smith } 16102593348eSBarry Smith 1611b0a32e0cSBarry Smith PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 16122593348eSBarry Smith c->sorted = a->sorted; 16132593348eSBarry Smith c->roworiented = a->roworiented; 16142593348eSBarry Smith c->nonew = a->nonew; 16152593348eSBarry Smith 16162593348eSBarry Smith if (a->diag) { 1617b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr); 1618b0a32e0cSBarry Smith PetscLogObjectMemory(C,(mbs+1)*sizeof(int)); 1619b6490206SBarry Smith for (i=0; i<mbs; i++) { 16202593348eSBarry Smith c->diag[i] = a->diag[i]; 16212593348eSBarry Smith } 162298305bb5SBarry Smith } else c->diag = 0; 16232593348eSBarry Smith c->nz = a->nz; 16242593348eSBarry Smith c->maxnz = a->maxnz; 16252593348eSBarry Smith c->solve_work = 0; 16262593348eSBarry Smith c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 16277fc0212eSBarry Smith c->mult_work = 0; 1628273d9f13SBarry Smith C->preallocated = PETSC_TRUE; 1629273d9f13SBarry Smith C->assembled = PETSC_TRUE; 16302593348eSBarry Smith *B = C; 1631b0a32e0cSBarry Smith ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 16323a40ed3dSBarry Smith PetscFunctionReturn(0); 16332593348eSBarry Smith } 16342593348eSBarry Smith 1635273d9f13SBarry Smith EXTERN_C_BEGIN 16364a2ae208SSatish Balay #undef __FUNCT__ 16374a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ" 1638b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A) 16392593348eSBarry Smith { 1640b6490206SBarry Smith Mat_SeqBAIJ *a; 16412593348eSBarry Smith Mat B; 1642f1af5d2fSBarry Smith int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1; 1643b6490206SBarry Smith int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 164435aab85fSBarry Smith int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1645a7c10996SSatish Balay int *masked,nmask,tmp,bs2,ishift; 1646b6490206SBarry Smith Scalar *aa; 164719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 16482593348eSBarry Smith 16493a40ed3dSBarry Smith PetscFunctionBegin; 1650b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr); 1651de6a44a3SBarry Smith bs2 = bs*bs; 1652b6490206SBarry Smith 1653d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 165429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 1655b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 16560752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 165729bbc08cSBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object"); 16582593348eSBarry Smith M = header[1]; N = header[2]; nz = header[3]; 16592593348eSBarry Smith 1660d64ed03dSBarry Smith if (header[3] < 0) { 166129bbc08cSBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ"); 1662d64ed03dSBarry Smith } 1663d64ed03dSBarry Smith 166429bbc08cSBarry Smith if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices"); 166535aab85fSBarry Smith 166635aab85fSBarry Smith /* 166735aab85fSBarry Smith This code adds extra rows to make sure the number of rows is 166835aab85fSBarry Smith divisible by the blocksize 166935aab85fSBarry Smith */ 1670b6490206SBarry Smith mbs = M/bs; 167135aab85fSBarry Smith extra_rows = bs - M + bs*(mbs); 167235aab85fSBarry Smith if (extra_rows == bs) extra_rows = 0; 167335aab85fSBarry Smith else mbs++; 167435aab85fSBarry Smith if (extra_rows) { 1675b0a32e0cSBarry Smith PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n"); 167635aab85fSBarry Smith } 1677b6490206SBarry Smith 16782593348eSBarry Smith /* read in row lengths */ 1679b0a32e0cSBarry Smith ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr); 16800752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 168135aab85fSBarry Smith for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1; 16822593348eSBarry Smith 1683b6490206SBarry Smith /* read in column indices */ 1684b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr); 16850752156aSBarry Smith ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr); 168635aab85fSBarry Smith for (i=0; i<extra_rows; i++) jj[nz+i] = M+i; 1687b6490206SBarry Smith 1688b6490206SBarry Smith /* loop over row lengths determining block row lengths */ 1689b0a32e0cSBarry Smith ierr = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 1690549d3d68SSatish Balay ierr = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr); 1691b0a32e0cSBarry Smith ierr = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr); 1692549d3d68SSatish Balay ierr = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr); 169335aab85fSBarry Smith masked = mask + mbs; 1694b6490206SBarry Smith rowcount = 0; nzcount = 0; 1695b6490206SBarry Smith for (i=0; i<mbs; i++) { 169635aab85fSBarry Smith nmask = 0; 1697b6490206SBarry Smith for (j=0; j<bs; j++) { 1698b6490206SBarry Smith kmax = rowlengths[rowcount]; 1699b6490206SBarry Smith for (k=0; k<kmax; k++) { 170035aab85fSBarry Smith tmp = jj[nzcount++]/bs; 170135aab85fSBarry Smith if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1702b6490206SBarry Smith } 1703b6490206SBarry Smith rowcount++; 1704b6490206SBarry Smith } 170535aab85fSBarry Smith browlengths[i] += nmask; 170635aab85fSBarry Smith /* zero out the mask elements we set */ 170735aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1708b6490206SBarry Smith } 1709b6490206SBarry Smith 17102593348eSBarry Smith /* create our matrix */ 17113eee16b0SBarry Smith ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr); 17122593348eSBarry Smith B = *A; 1713b6490206SBarry Smith a = (Mat_SeqBAIJ*)B->data; 17142593348eSBarry Smith 17152593348eSBarry Smith /* set matrix "i" values */ 1716de6a44a3SBarry Smith a->i[0] = 0; 1717b6490206SBarry Smith for (i=1; i<= mbs; i++) { 1718b6490206SBarry Smith a->i[i] = a->i[i-1] + browlengths[i-1]; 1719b6490206SBarry Smith a->ilen[i-1] = browlengths[i-1]; 17202593348eSBarry Smith } 1721b6490206SBarry Smith a->nz = 0; 1722b6490206SBarry Smith for (i=0; i<mbs; i++) a->nz += browlengths[i]; 17232593348eSBarry Smith 1724b6490206SBarry Smith /* read in nonzero values */ 1725b0a32e0cSBarry Smith ierr = PetscMalloc((nz+extra_rows)*sizeof(Scalar),&aa);CHKERRQ(ierr); 17260752156aSBarry Smith ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr); 172735aab85fSBarry Smith for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0; 1728b6490206SBarry Smith 1729b6490206SBarry Smith /* set "a" and "j" values into matrix */ 1730b6490206SBarry Smith nzcount = 0; jcount = 0; 1731b6490206SBarry Smith for (i=0; i<mbs; i++) { 1732b6490206SBarry Smith nzcountb = nzcount; 173335aab85fSBarry Smith nmask = 0; 1734b6490206SBarry Smith for (j=0; j<bs; j++) { 1735b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1736b6490206SBarry Smith for (k=0; k<kmax; k++) { 173735aab85fSBarry Smith tmp = jj[nzcount++]/bs; 173835aab85fSBarry Smith if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1739b6490206SBarry Smith } 1740b6490206SBarry Smith } 1741de6a44a3SBarry Smith /* sort the masked values */ 1742433994e6SBarry Smith ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr); 1743de6a44a3SBarry Smith 1744b6490206SBarry Smith /* set "j" values into matrix */ 1745b6490206SBarry Smith maskcount = 1; 174635aab85fSBarry Smith for (j=0; j<nmask; j++) { 174735aab85fSBarry Smith a->j[jcount++] = masked[j]; 1748de6a44a3SBarry Smith mask[masked[j]] = maskcount++; 1749b6490206SBarry Smith } 1750b6490206SBarry Smith /* set "a" values into matrix */ 1751de6a44a3SBarry Smith ishift = bs2*a->i[i]; 1752b6490206SBarry Smith for (j=0; j<bs; j++) { 1753b6490206SBarry Smith kmax = rowlengths[i*bs+j]; 1754b6490206SBarry Smith for (k=0; k<kmax; k++) { 1755de6a44a3SBarry Smith tmp = jj[nzcountb]/bs ; 1756de6a44a3SBarry Smith block = mask[tmp] - 1; 1757de6a44a3SBarry Smith point = jj[nzcountb] - bs*tmp; 1758de6a44a3SBarry Smith idx = ishift + bs2*block + j + bs*point; 1759375fe846SBarry Smith a->a[idx] = (MatScalar)aa[nzcountb++]; 1760b6490206SBarry Smith } 1761b6490206SBarry Smith } 176235aab85fSBarry Smith /* zero out the mask elements we set */ 176335aab85fSBarry Smith for (j=0; j<nmask; j++) mask[masked[j]] = 0; 1764b6490206SBarry Smith } 176529bbc08cSBarry Smith if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix"); 1766b6490206SBarry Smith 1767606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 1768606d414cSSatish Balay ierr = PetscFree(browlengths);CHKERRQ(ierr); 1769606d414cSSatish Balay ierr = PetscFree(aa);CHKERRQ(ierr); 1770606d414cSSatish Balay ierr = PetscFree(jj);CHKERRQ(ierr); 1771606d414cSSatish Balay ierr = PetscFree(mask);CHKERRQ(ierr); 1772b6490206SBarry Smith 1773b6490206SBarry Smith B->assembled = PETSC_TRUE; 1774de6a44a3SBarry Smith 17759c01be13SBarry Smith ierr = MatView_Private(B);CHKERRQ(ierr); 17763a40ed3dSBarry Smith PetscFunctionReturn(0); 17772593348eSBarry Smith } 1778273d9f13SBarry Smith EXTERN_C_END 17792593348eSBarry Smith 17804a2ae208SSatish Balay #undef __FUNCT__ 17814a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ" 1782273d9f13SBarry Smith /*@C 1783273d9f13SBarry Smith MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 1784273d9f13SBarry Smith compressed row) format. For good matrix assembly performance the 1785273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 1786273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 1787273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 17882593348eSBarry Smith 1789273d9f13SBarry Smith Collective on MPI_Comm 1790273d9f13SBarry Smith 1791273d9f13SBarry Smith Input Parameters: 1792273d9f13SBarry Smith + comm - MPI communicator, set to PETSC_COMM_SELF 1793273d9f13SBarry Smith . bs - size of block 1794273d9f13SBarry Smith . m - number of rows 1795273d9f13SBarry Smith . n - number of columns 179635d8aa7fSBarry Smith . nz - number of nonzero blocks per block row (same for all rows) 179735d8aa7fSBarry Smith - nnz - array containing the number of nonzero blocks in the various block rows 1798273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 1799273d9f13SBarry Smith 1800273d9f13SBarry Smith Output Parameter: 1801273d9f13SBarry Smith . A - the matrix 1802273d9f13SBarry Smith 1803273d9f13SBarry Smith Options Database Keys: 1804273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 1805273d9f13SBarry Smith block calculations (much slower) 1806273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 1807273d9f13SBarry Smith 1808273d9f13SBarry Smith Level: intermediate 1809273d9f13SBarry Smith 1810273d9f13SBarry Smith Notes: 181135d8aa7fSBarry Smith A nonzero block is any block that as 1 or more nonzeros in it 181235d8aa7fSBarry Smith 1813273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 1814273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 1815273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 1816273d9f13SBarry Smith 1817273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1818273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1819273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 1820273d9f13SBarry Smith matrices. 1821273d9f13SBarry Smith 1822273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 1823273d9f13SBarry Smith @*/ 1824273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A) 1825273d9f13SBarry Smith { 1826273d9f13SBarry Smith int ierr; 1827273d9f13SBarry Smith 1828273d9f13SBarry Smith PetscFunctionBegin; 1829273d9f13SBarry Smith ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1830273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr); 1831273d9f13SBarry Smith ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr); 1832273d9f13SBarry Smith PetscFunctionReturn(0); 1833273d9f13SBarry Smith } 1834273d9f13SBarry Smith 18354a2ae208SSatish Balay #undef __FUNCT__ 18364a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation" 1837273d9f13SBarry Smith /*@C 1838273d9f13SBarry Smith MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros 1839273d9f13SBarry Smith per row in the matrix. For good matrix assembly performance the 1840273d9f13SBarry Smith user should preallocate the matrix storage by setting the parameter nz 1841273d9f13SBarry Smith (or the array nnz). By setting these parameters accurately, performance 1842273d9f13SBarry Smith during matrix assembly can be increased by more than a factor of 50. 1843273d9f13SBarry Smith 1844273d9f13SBarry Smith Collective on MPI_Comm 1845273d9f13SBarry Smith 1846273d9f13SBarry Smith Input Parameters: 1847273d9f13SBarry Smith + A - the matrix 1848273d9f13SBarry Smith . bs - size of block 1849273d9f13SBarry Smith . nz - number of block nonzeros per block row (same for all rows) 1850273d9f13SBarry Smith - nnz - array containing the number of block nonzeros in the various block rows 1851273d9f13SBarry Smith (possibly different for each block row) or PETSC_NULL 1852273d9f13SBarry Smith 1853273d9f13SBarry Smith Options Database Keys: 1854273d9f13SBarry Smith . -mat_no_unroll - uses code that does not unroll the loops in the 1855273d9f13SBarry Smith block calculations (much slower) 1856273d9f13SBarry Smith . -mat_block_size - size of the blocks to use 1857273d9f13SBarry Smith 1858273d9f13SBarry Smith Level: intermediate 1859273d9f13SBarry Smith 1860273d9f13SBarry Smith Notes: 1861273d9f13SBarry Smith The block AIJ format is fully compatible with standard Fortran 77 1862273d9f13SBarry Smith storage. That is, the stored row and column indices can begin at 1863273d9f13SBarry Smith either one (as in Fortran) or zero. See the users' manual for details. 1864273d9f13SBarry Smith 1865273d9f13SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1866273d9f13SBarry Smith Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1867273d9f13SBarry Smith allocation. For additional details, see the users manual chapter on 1868273d9f13SBarry Smith matrices. 1869273d9f13SBarry Smith 1870273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ() 1871273d9f13SBarry Smith @*/ 1872273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) 1873273d9f13SBarry Smith { 1874273d9f13SBarry Smith Mat_SeqBAIJ *b; 1875273d9f13SBarry Smith int i,len,ierr,mbs,nbs,bs2,newbs = bs; 1876273d9f13SBarry Smith PetscTruth flg; 1877273d9f13SBarry Smith 1878273d9f13SBarry Smith PetscFunctionBegin; 1879273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr); 1880273d9f13SBarry Smith if (!flg) PetscFunctionReturn(0); 1881273d9f13SBarry Smith 1882273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1883b0a32e0cSBarry Smith ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr); 1884273d9f13SBarry Smith if (nnz && newbs != bs) { 1885273d9f13SBarry Smith SETERRQ(1,"Cannot change blocksize from command line if setting nnz"); 1886273d9f13SBarry Smith } 1887273d9f13SBarry Smith bs = newbs; 1888273d9f13SBarry Smith 1889273d9f13SBarry Smith mbs = B->m/bs; 1890273d9f13SBarry Smith nbs = B->n/bs; 1891273d9f13SBarry Smith bs2 = bs*bs; 1892273d9f13SBarry Smith 1893273d9f13SBarry Smith if (mbs*bs!=B->m || nbs*bs!=B->n) { 189435d8aa7fSBarry Smith SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs); 1895273d9f13SBarry Smith } 1896273d9f13SBarry Smith 1897435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1898435da068SBarry Smith if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); 1899273d9f13SBarry Smith if (nnz) { 1900273d9f13SBarry Smith for (i=0; i<mbs; i++) { 1901273d9f13SBarry Smith if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); 19023a7fca6bSBarry Smith if (nnz[i] > nbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %d value %d rowlength %d",i,nnz[i],nbs); 1903273d9f13SBarry Smith } 1904273d9f13SBarry Smith } 1905273d9f13SBarry Smith 1906273d9f13SBarry Smith b = (Mat_SeqBAIJ*)B->data; 1907b0a32e0cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr); 1908273d9f13SBarry Smith if (!flg) { 1909273d9f13SBarry Smith switch (bs) { 1910273d9f13SBarry Smith case 1: 1911273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1912273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_1; 1913273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 1914273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_1; 1915273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_1; 1916273d9f13SBarry Smith break; 1917273d9f13SBarry Smith case 2: 1918273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1919273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_2; 1920273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 1921273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_2; 1922273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_2; 1923273d9f13SBarry Smith break; 1924273d9f13SBarry Smith case 3: 1925273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1926273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_3; 1927273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 1928273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_3; 1929273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_3; 1930273d9f13SBarry Smith break; 1931273d9f13SBarry Smith case 4: 1932273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1933273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_4; 1934273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4; 1935273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_4; 1936273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_4; 1937273d9f13SBarry Smith break; 1938273d9f13SBarry Smith case 5: 1939273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1940273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_5; 1941273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5; 1942273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_5; 1943273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_5; 1944273d9f13SBarry Smith break; 1945273d9f13SBarry Smith case 6: 1946273d9f13SBarry Smith B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6; 1947273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_6; 1948273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6; 1949273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_6; 1950273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_6; 1951273d9f13SBarry Smith break; 1952273d9f13SBarry Smith case 7: 1953273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_7; 1954273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_7; 1955273d9f13SBarry Smith B->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7; 1956273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_7; 1957273d9f13SBarry Smith break; 1958273d9f13SBarry Smith default: 1959273d9f13SBarry Smith B->ops->mult = MatMult_SeqBAIJ_N; 1960273d9f13SBarry Smith B->ops->solve = MatSolve_SeqBAIJ_N; 1961273d9f13SBarry Smith B->ops->multadd = MatMultAdd_SeqBAIJ_N; 1962273d9f13SBarry Smith break; 1963273d9f13SBarry Smith } 1964273d9f13SBarry Smith } 1965273d9f13SBarry Smith b->bs = bs; 1966273d9f13SBarry Smith b->mbs = mbs; 1967273d9f13SBarry Smith b->nbs = nbs; 1968b0a32e0cSBarry Smith ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr); 1969273d9f13SBarry Smith if (!nnz) { 1970435da068SBarry Smith if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; 1971273d9f13SBarry Smith else if (nz <= 0) nz = 1; 1972273d9f13SBarry Smith for (i=0; i<mbs; i++) b->imax[i] = nz; 1973273d9f13SBarry Smith nz = nz*mbs; 1974273d9f13SBarry Smith } else { 1975273d9f13SBarry Smith nz = 0; 1976273d9f13SBarry Smith for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} 1977273d9f13SBarry Smith } 1978273d9f13SBarry Smith 1979273d9f13SBarry Smith /* allocate the matrix space */ 1980273d9f13SBarry Smith len = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int); 1981b0a32e0cSBarry Smith ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr); 1982273d9f13SBarry Smith ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr); 1983273d9f13SBarry Smith b->j = (int*)(b->a + nz*bs2); 1984273d9f13SBarry Smith ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr); 1985273d9f13SBarry Smith b->i = b->j + nz; 1986273d9f13SBarry Smith b->singlemalloc = PETSC_TRUE; 1987273d9f13SBarry Smith 1988273d9f13SBarry Smith b->i[0] = 0; 1989273d9f13SBarry Smith for (i=1; i<mbs+1; i++) { 1990273d9f13SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 1991273d9f13SBarry Smith } 1992273d9f13SBarry Smith 1993273d9f13SBarry Smith /* b->ilen will count nonzeros in each block row so far. */ 199416d6aa21SSatish Balay ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr); 1995b0a32e0cSBarry Smith PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ)); 1996273d9f13SBarry Smith for (i=0; i<mbs; i++) { b->ilen[i] = 0;} 1997273d9f13SBarry Smith 1998273d9f13SBarry Smith b->bs = bs; 1999273d9f13SBarry Smith b->bs2 = bs2; 2000273d9f13SBarry Smith b->mbs = mbs; 2001273d9f13SBarry Smith b->nz = 0; 2002273d9f13SBarry Smith b->maxnz = nz*bs2; 2003273d9f13SBarry Smith B->info.nz_unneeded = (PetscReal)b->maxnz; 2004273d9f13SBarry Smith PetscFunctionReturn(0); 2005273d9f13SBarry Smith } 20062593348eSBarry Smith 2007