xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 15091d3721b14bd18f7efa625bb8738e103dca31)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*15091d37SBarry Smith static char vcid[] = "$Id: baij.c,v 1.166 1999/03/04 22:35:37 bsmith Exp bsmith $";
32593348eSBarry Smith #endif
42593348eSBarry Smith 
52593348eSBarry Smith /*
6b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
72593348eSBarry Smith   matrix storage format.
82593348eSBarry Smith */
93369ce9aSBarry Smith #include "sys.h"
1070f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
111a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
121a6a6d98SBarry Smith #include "src/inline/spops.h"
133b547af2SSatish Balay 
142d61bbb3SSatish Balay #define CHUNKSIZE  10
15de6a44a3SBarry Smith 
16be5855fcSBarry Smith /*
17be5855fcSBarry Smith      Checks for missing diagonals
18be5855fcSBarry Smith */
19be5855fcSBarry Smith #undef __FUNC__
20be5855fcSBarry Smith #define __FUNC__ "MatMissingDiag_SeqBAIJ"
21be5855fcSBarry Smith int MatMissingDiag_SeqBAIJ(Mat A)
22be5855fcSBarry Smith {
23be5855fcSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
24be5855fcSBarry Smith   int         *diag = a->diag, *jj = a->j,i;
25be5855fcSBarry Smith 
26be5855fcSBarry Smith   PetscFunctionBegin;
270e8e8aceSBarry Smith   for ( i=0; i<a->mbs; i++ ) {
28be5855fcSBarry Smith     if (jj[diag[i]] != i) {
29be5855fcSBarry Smith       SETERRQ1(1,1,"Matrix is missing diagonal number %d",i);
30be5855fcSBarry Smith     }
31be5855fcSBarry Smith   }
32be5855fcSBarry Smith   PetscFunctionReturn(0);
33be5855fcSBarry Smith }
34be5855fcSBarry Smith 
355615d1e5SSatish Balay #undef __FUNC__
365615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqBAIJ"
37de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A)
38de6a44a3SBarry Smith {
39de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
407fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
41de6a44a3SBarry Smith 
423a40ed3dSBarry Smith   PetscFunctionBegin;
43de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
44de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
457fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
46ecc615b2SBarry Smith     diag[i] = a->i[i+1];
47de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
48de6a44a3SBarry Smith       if (a->j[j] == i) {
49de6a44a3SBarry Smith         diag[i] = j;
50de6a44a3SBarry Smith         break;
51de6a44a3SBarry Smith       }
52de6a44a3SBarry Smith     }
53de6a44a3SBarry Smith   }
54de6a44a3SBarry Smith   a->diag = diag;
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56de6a44a3SBarry Smith }
572593348eSBarry Smith 
582593348eSBarry Smith 
593b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
603b2fbd54SBarry Smith 
615615d1e5SSatish Balay #undef __FUNC__
625615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqBAIJ"
633b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
643b2fbd54SBarry Smith                             PetscTruth *done)
653b2fbd54SBarry Smith {
663b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
673b2fbd54SBarry Smith   int         ierr, n = a->mbs,i;
683b2fbd54SBarry Smith 
693a40ed3dSBarry Smith   PetscFunctionBegin;
703b2fbd54SBarry Smith   *nn = n;
713a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
723b2fbd54SBarry Smith   if (symmetric) {
733b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
743b2fbd54SBarry Smith   } else if (oshift == 1) {
753b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
763b2fbd54SBarry Smith     int nz = a->i[n] + 1;
773b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]++;
783b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]++;
793b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
803b2fbd54SBarry Smith   } else {
813b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
823b2fbd54SBarry Smith   }
833b2fbd54SBarry Smith 
843a40ed3dSBarry Smith   PetscFunctionReturn(0);
853b2fbd54SBarry Smith }
863b2fbd54SBarry Smith 
875615d1e5SSatish Balay #undef __FUNC__
88d4bb536fSBarry Smith #define __FUNC__ "MatRestoreRowIJ_SeqBAIJ"
893b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
903b2fbd54SBarry Smith                                 PetscTruth *done)
913b2fbd54SBarry Smith {
923b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
933b2fbd54SBarry Smith   int         i,n = a->mbs;
943b2fbd54SBarry Smith 
953a40ed3dSBarry Smith   PetscFunctionBegin;
963a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
973b2fbd54SBarry Smith   if (symmetric) {
983b2fbd54SBarry Smith     PetscFree(*ia);
993b2fbd54SBarry Smith     PetscFree(*ja);
100af5da2bfSBarry Smith   } else if (oshift == 1) {
1013b2fbd54SBarry Smith     int nz = a->i[n];
1023b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]--;
1033b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]--;
1043b2fbd54SBarry Smith   }
1053a40ed3dSBarry Smith   PetscFunctionReturn(0);
1063b2fbd54SBarry Smith }
1073b2fbd54SBarry Smith 
1082d61bbb3SSatish Balay #undef __FUNC__
1092d61bbb3SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqBAIJ"
1102d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs)
1112d61bbb3SSatish Balay {
1122d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data;
1132d61bbb3SSatish Balay 
1142d61bbb3SSatish Balay   PetscFunctionBegin;
1152d61bbb3SSatish Balay   *bs = baij->bs;
1162d61bbb3SSatish Balay   PetscFunctionReturn(0);
1172d61bbb3SSatish Balay }
1182d61bbb3SSatish Balay 
1192d61bbb3SSatish Balay 
1202d61bbb3SSatish Balay #undef __FUNC__
1212d61bbb3SSatish Balay #define __FUNC__ "MatDestroy_SeqBAIJ"
122e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
1232d61bbb3SSatish Balay {
1242d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
125e51c0b9cSSatish Balay   int         ierr;
1262d61bbb3SSatish Balay 
12794d884c6SBarry Smith   if (--A->refct > 0) PetscFunctionReturn(0);
12894d884c6SBarry Smith 
12994d884c6SBarry Smith   if (A->mapping) {
13094d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr);
13194d884c6SBarry Smith   }
13294d884c6SBarry Smith   if (A->bmapping) {
13394d884c6SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr);
13494d884c6SBarry Smith   }
13561b13de0SBarry Smith   if (A->rmap) {
13661b13de0SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
13761b13de0SBarry Smith   }
13861b13de0SBarry Smith   if (A->cmap) {
13961b13de0SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
14061b13de0SBarry Smith   }
1412d61bbb3SSatish Balay #if defined(USE_PETSC_LOG)
142e1311b90SBarry Smith   PLogObjectState((PetscObject) A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
1432d61bbb3SSatish Balay #endif
1442d61bbb3SSatish Balay   PetscFree(a->a);
1452d61bbb3SSatish Balay   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
1462d61bbb3SSatish Balay   if (a->diag) PetscFree(a->diag);
1472d61bbb3SSatish Balay   if (a->ilen) PetscFree(a->ilen);
1482d61bbb3SSatish Balay   if (a->imax) PetscFree(a->imax);
1492d61bbb3SSatish Balay   if (a->solve_work) PetscFree(a->solve_work);
1502d61bbb3SSatish Balay   if (a->mult_work) PetscFree(a->mult_work);
151e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
1522d61bbb3SSatish Balay   PetscFree(a);
1532d61bbb3SSatish Balay   PLogObjectDestroy(A);
1542d61bbb3SSatish Balay   PetscHeaderDestroy(A);
1552d61bbb3SSatish Balay   PetscFunctionReturn(0);
1562d61bbb3SSatish Balay }
1572d61bbb3SSatish Balay 
1582d61bbb3SSatish Balay #undef __FUNC__
1592d61bbb3SSatish Balay #define __FUNC__ "MatSetOption_SeqBAIJ"
1602d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1612d61bbb3SSatish Balay {
1622d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1632d61bbb3SSatish Balay 
1642d61bbb3SSatish Balay   PetscFunctionBegin;
1652d61bbb3SSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented = 1;
1662d61bbb3SSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented = 0;
1672d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted      = 1;
1682d61bbb3SSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted      = 0;
1692d61bbb3SSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew       = 1;
1704787f768SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew       = -1;
1714787f768SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew       = -2;
1722d61bbb3SSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew       = 0;
1732d61bbb3SSatish Balay   else if (op == MAT_ROWS_SORTED ||
1742d61bbb3SSatish Balay            op == MAT_ROWS_UNSORTED ||
1752d61bbb3SSatish Balay            op == MAT_SYMMETRIC ||
1762d61bbb3SSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
1772d61bbb3SSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
178b51ba29fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
179b51ba29fSSatish Balay            op == MAT_USE_HASH_TABLE) {
180981c4779SBarry Smith     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
1812d61bbb3SSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
1822d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
1832d61bbb3SSatish Balay   } else {
1842d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
1852d61bbb3SSatish Balay   }
1862d61bbb3SSatish Balay   PetscFunctionReturn(0);
1872d61bbb3SSatish Balay }
1882d61bbb3SSatish Balay 
1892d61bbb3SSatish Balay 
1902d61bbb3SSatish Balay #undef __FUNC__
1912d61bbb3SSatish Balay #define __FUNC__ "MatGetSize_SeqBAIJ"
1922d61bbb3SSatish Balay int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
1932d61bbb3SSatish Balay {
1942d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1952d61bbb3SSatish Balay 
1962d61bbb3SSatish Balay   PetscFunctionBegin;
1972d61bbb3SSatish Balay   if (m) *m = a->m;
1982d61bbb3SSatish Balay   if (n) *n = a->n;
1992d61bbb3SSatish Balay   PetscFunctionReturn(0);
2002d61bbb3SSatish Balay }
2012d61bbb3SSatish Balay 
2022d61bbb3SSatish Balay #undef __FUNC__
2032d61bbb3SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqBAIJ"
2042d61bbb3SSatish Balay int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
2052d61bbb3SSatish Balay {
2062d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2072d61bbb3SSatish Balay 
2082d61bbb3SSatish Balay   PetscFunctionBegin;
2092d61bbb3SSatish Balay   *m = 0; *n = a->m;
2102d61bbb3SSatish Balay   PetscFunctionReturn(0);
2112d61bbb3SSatish Balay }
2122d61bbb3SSatish Balay 
2132d61bbb3SSatish Balay #undef __FUNC__
2142d61bbb3SSatish Balay #define __FUNC__ "MatGetRow_SeqBAIJ"
2152d61bbb3SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2162d61bbb3SSatish Balay {
2172d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data;
2182d61bbb3SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
2193f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
2203f1db9ecSBarry Smith   Scalar       *v_i;
2212d61bbb3SSatish Balay 
2222d61bbb3SSatish Balay   PetscFunctionBegin;
2232d61bbb3SSatish Balay   bs  = a->bs;
2242d61bbb3SSatish Balay   ai  = a->i;
2252d61bbb3SSatish Balay   aj  = a->j;
2262d61bbb3SSatish Balay   aa  = a->a;
2272d61bbb3SSatish Balay   bs2 = a->bs2;
2282d61bbb3SSatish Balay 
2292d61bbb3SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
2302d61bbb3SSatish Balay 
2312d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2322d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2332d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2342d61bbb3SSatish Balay   *nz = bs*M;
2352d61bbb3SSatish Balay 
2362d61bbb3SSatish Balay   if (v) {
2372d61bbb3SSatish Balay     *v = 0;
2382d61bbb3SSatish Balay     if (*nz) {
2392d61bbb3SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
2402d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2412d61bbb3SSatish Balay         v_i  = *v + i*bs;
2422d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2432d61bbb3SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2442d61bbb3SSatish Balay       }
2452d61bbb3SSatish Balay     }
2462d61bbb3SSatish Balay   }
2472d61bbb3SSatish Balay 
2482d61bbb3SSatish Balay   if (idx) {
2492d61bbb3SSatish Balay     *idx = 0;
2502d61bbb3SSatish Balay     if (*nz) {
2512d61bbb3SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
2522d61bbb3SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2532d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2542d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2552d61bbb3SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2562d61bbb3SSatish Balay       }
2572d61bbb3SSatish Balay     }
2582d61bbb3SSatish Balay   }
2592d61bbb3SSatish Balay   PetscFunctionReturn(0);
2602d61bbb3SSatish Balay }
2612d61bbb3SSatish Balay 
2622d61bbb3SSatish Balay #undef __FUNC__
2632d61bbb3SSatish Balay #define __FUNC__ "MatRestoreRow_SeqBAIJ"
2642d61bbb3SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2652d61bbb3SSatish Balay {
2662d61bbb3SSatish Balay   PetscFunctionBegin;
2672d61bbb3SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
2682d61bbb3SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
2692d61bbb3SSatish Balay   PetscFunctionReturn(0);
2702d61bbb3SSatish Balay }
2712d61bbb3SSatish Balay 
2722d61bbb3SSatish Balay #undef __FUNC__
2732d61bbb3SSatish Balay #define __FUNC__ "MatTranspose_SeqBAIJ"
2742d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
2752d61bbb3SSatish Balay {
2762d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
2772d61bbb3SSatish Balay   Mat         C;
2782d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
2792d61bbb3SSatish Balay   int         *rows,*cols,bs2=a->bs2;
2803f1db9ecSBarry Smith   MatScalar   *array = a->a;
2812d61bbb3SSatish Balay 
2822d61bbb3SSatish Balay   PetscFunctionBegin;
2832d61bbb3SSatish Balay   if (B==PETSC_NULL && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
2842d61bbb3SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
2852d61bbb3SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
2862d61bbb3SSatish Balay 
2872d61bbb3SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
2882d61bbb3SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
2892d61bbb3SSatish Balay   PetscFree(col);
2902d61bbb3SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
2912d61bbb3SSatish Balay   cols = rows + bs;
2922d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2932d61bbb3SSatish Balay     cols[0] = i*bs;
2942d61bbb3SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
2952d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
2962d61bbb3SSatish Balay     for ( j=0; j<len; j++ ) {
2972d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
2982d61bbb3SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
2992d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
3002d61bbb3SSatish Balay       array += bs2;
3012d61bbb3SSatish Balay     }
3022d61bbb3SSatish Balay   }
3032d61bbb3SSatish Balay   PetscFree(rows);
3042d61bbb3SSatish Balay 
3052d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
3062d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay   if (B != PETSC_NULL) {
3092d61bbb3SSatish Balay     *B = C;
3102d61bbb3SSatish Balay   } else {
311f830108cSBarry Smith     PetscOps *Abops;
312cc2dc46cSBarry Smith     MatOps   Aops;
313f830108cSBarry Smith 
3142d61bbb3SSatish Balay     /* This isn't really an in-place transpose */
3152d61bbb3SSatish Balay     PetscFree(a->a);
3162d61bbb3SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
3172d61bbb3SSatish Balay     if (a->diag) PetscFree(a->diag);
3182d61bbb3SSatish Balay     if (a->ilen) PetscFree(a->ilen);
3192d61bbb3SSatish Balay     if (a->imax) PetscFree(a->imax);
3202d61bbb3SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
3212d61bbb3SSatish Balay     PetscFree(a);
322f830108cSBarry Smith 
3237413bad6SBarry Smith 
3247413bad6SBarry Smith     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
3257413bad6SBarry Smith     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
3267413bad6SBarry Smith 
327f830108cSBarry Smith     /*
328f830108cSBarry Smith        This is horrible, horrible code. We need to keep the
329f830108cSBarry Smith       A pointers for the bops and ops but copy everything
330f830108cSBarry Smith       else from C.
331f830108cSBarry Smith     */
332f830108cSBarry Smith     Abops    = A->bops;
333f830108cSBarry Smith     Aops     = A->ops;
3342d61bbb3SSatish Balay     PetscMemcpy(A,C,sizeof(struct _p_Mat));
335f830108cSBarry Smith     A->bops  = Abops;
336f830108cSBarry Smith     A->ops   = Aops;
33727a8da17SBarry Smith     A->qlist = 0;
338f830108cSBarry Smith 
3392d61bbb3SSatish Balay     PetscHeaderDestroy(C);
3402d61bbb3SSatish Balay   }
3412d61bbb3SSatish Balay   PetscFunctionReturn(0);
3422d61bbb3SSatish Balay }
3432d61bbb3SSatish Balay 
3445615d1e5SSatish Balay #undef __FUNC__
345d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
346b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
3472593348eSBarry Smith {
348b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3499df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
350b6490206SBarry Smith   Scalar      *aa;
351ce6f0cecSBarry Smith   FILE        *file;
3522593348eSBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
35490ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
3552593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3562593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3573b2fbd54SBarry Smith 
3582593348eSBarry Smith   col_lens[1] = a->m;
3592593348eSBarry Smith   col_lens[2] = a->n;
3607e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3612593348eSBarry Smith 
3622593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
363b6490206SBarry Smith   count = 0;
364b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
365b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
366b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
367b6490206SBarry Smith     }
3682593348eSBarry Smith   }
3690752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3702593348eSBarry Smith   PetscFree(col_lens);
3712593348eSBarry Smith 
3722593348eSBarry Smith   /* store column indices (zero start index) */
37366963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj);
374b6490206SBarry Smith   count = 0;
375b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
376b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
377b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
378b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
379b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3802593348eSBarry Smith         }
3812593348eSBarry Smith       }
382b6490206SBarry Smith     }
383b6490206SBarry Smith   }
3840752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr);
385b6490206SBarry Smith   PetscFree(jj);
3862593348eSBarry Smith 
3872593348eSBarry Smith   /* store nonzero values */
38866963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa);
389b6490206SBarry Smith   count = 0;
390b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
391b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
392b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
393b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
3947e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
395b6490206SBarry Smith         }
396b6490206SBarry Smith       }
397b6490206SBarry Smith     }
398b6490206SBarry Smith   }
3990752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
400b6490206SBarry Smith   PetscFree(aa);
401ce6f0cecSBarry Smith 
402ce6f0cecSBarry Smith   ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
403ce6f0cecSBarry Smith   if (file) {
404ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
405ce6f0cecSBarry Smith   }
4063a40ed3dSBarry Smith   PetscFunctionReturn(0);
4072593348eSBarry Smith }
4082593348eSBarry Smith 
4095615d1e5SSatish Balay #undef __FUNC__
410d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
411b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
4122593348eSBarry Smith {
413b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4149df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
4152593348eSBarry Smith   FILE        *fd;
4162593348eSBarry Smith   char        *outputname;
4172593348eSBarry Smith 
4183a40ed3dSBarry Smith   PetscFunctionBegin;
41990ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
42077ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
42190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
422639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
4230ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);
4240ef38995SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
4257b2a1423SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Socket format not supported");
4260ef38995SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
42744cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
42844cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
42944cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
43044cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
43144cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
4323a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
4330ef38995SBarry Smith             if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
43444cd7ae7SLois Curfman McInnes               fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
435e20fef11SSatish Balay                       PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
4360ef38995SBarry Smith             } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
437766eeae4SLois Curfman McInnes               fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
438e20fef11SSatish Balay                       PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));
4390ef38995SBarry Smith             } else if (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
440e20fef11SSatish Balay               fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
4410ef38995SBarry Smith             }
44244cd7ae7SLois Curfman McInnes #else
4430ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
44444cd7ae7SLois Curfman McInnes               fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
4450ef38995SBarry Smith             }
44644cd7ae7SLois Curfman McInnes #endif
44744cd7ae7SLois Curfman McInnes           }
44844cd7ae7SLois Curfman McInnes         }
44944cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
45044cd7ae7SLois Curfman McInnes       }
45144cd7ae7SLois Curfman McInnes     }
4520ef38995SBarry Smith   } else {
453b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
454b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
455b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
456b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
457b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
4583a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
459e20fef11SSatish Balay             if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) {
46088685aaeSLois Curfman McInnes               fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
461e20fef11SSatish Balay                 PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
4620ef38995SBarry Smith             } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) {
463766eeae4SLois Curfman McInnes               fprintf(fd," %d %g - %g i",bs*a->j[k]+l,
464e20fef11SSatish Balay                 PetscReal(a->a[bs2*k + l*bs + j]),-PetscImaginary(a->a[bs2*k + l*bs + j]));
4650ef38995SBarry Smith             } else {
466e20fef11SSatish Balay               fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
46788685aaeSLois Curfman McInnes             }
46888685aaeSLois Curfman McInnes #else
4697e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
47088685aaeSLois Curfman McInnes #endif
4712593348eSBarry Smith           }
4722593348eSBarry Smith         }
4732593348eSBarry Smith         fprintf(fd,"\n");
4742593348eSBarry Smith       }
4752593348eSBarry Smith     }
476b6490206SBarry Smith   }
4772593348eSBarry Smith   fflush(fd);
4783a40ed3dSBarry Smith   PetscFunctionReturn(0);
4792593348eSBarry Smith }
4802593348eSBarry Smith 
4815615d1e5SSatish Balay #undef __FUNC__
48277ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom"
48377ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa)
4843270192aSSatish Balay {
48577ed5343SBarry Smith   Mat          A = (Mat) Aa;
4863270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
4877dce120fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
488fae8c767SSatish Balay   double       xl,yl,xr,yr,x_l,x_r,y_l,y_r;
4893f1db9ecSBarry Smith   MatScalar    *aa;
49077ed5343SBarry Smith   MPI_Comm     comm;
49177ed5343SBarry Smith   Viewer       viewer;
4923270192aSSatish Balay 
4933a40ed3dSBarry Smith   PetscFunctionBegin;
49477ed5343SBarry Smith   /*
49577ed5343SBarry Smith       This is nasty. If this is called from an originally parallel matrix
49677ed5343SBarry Smith    then all processes call this, but only the first has the matrix so the
49777ed5343SBarry Smith    rest should return immediately.
49877ed5343SBarry Smith   */
49977ed5343SBarry Smith   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
50077ed5343SBarry Smith   MPI_Comm_rank(comm,&rank);
50177ed5343SBarry Smith   if (rank) PetscFunctionReturn(0);
5023270192aSSatish Balay 
50377ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr);
50477ed5343SBarry Smith 
50577ed5343SBarry Smith   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr);
50677ed5343SBarry Smith 
5073270192aSSatish Balay   /* loop over matrix elements drawing boxes */
5083270192aSSatish Balay   color = DRAW_BLUE;
5093270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5103270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5113270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5123270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5133270192aSSatish Balay       aa = a->a + j*bs2;
5143270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5153270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5163270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
517b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5183270192aSSatish Balay         }
5193270192aSSatish Balay       }
5203270192aSSatish Balay     }
5213270192aSSatish Balay   }
5223270192aSSatish Balay   color = DRAW_CYAN;
5233270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5243270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5253270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5263270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5273270192aSSatish Balay       aa = a->a + j*bs2;
5283270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5293270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5303270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
531b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5323270192aSSatish Balay         }
5333270192aSSatish Balay       }
5343270192aSSatish Balay     }
5353270192aSSatish Balay   }
5363270192aSSatish Balay 
5373270192aSSatish Balay   color = DRAW_RED;
5383270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5393270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5403270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5413270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5423270192aSSatish Balay       aa = a->a + j*bs2;
5433270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5443270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5453270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
546b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5473270192aSSatish Balay         }
5483270192aSSatish Balay       }
5493270192aSSatish Balay     }
5503270192aSSatish Balay   }
55177ed5343SBarry Smith   PetscFunctionReturn(0);
55277ed5343SBarry Smith }
5533270192aSSatish Balay 
55477ed5343SBarry Smith #undef __FUNC__
55577ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
55677ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
55777ed5343SBarry Smith {
55877ed5343SBarry Smith   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
5597dce120fSSatish Balay   int          ierr;
5607dce120fSSatish Balay   double       xl,yl,xr,yr,w,h;
56177ed5343SBarry Smith   Draw         draw;
56277ed5343SBarry Smith   PetscTruth   isnull;
5633270192aSSatish Balay 
56477ed5343SBarry Smith   PetscFunctionBegin;
56577ed5343SBarry Smith 
56677ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
56777ed5343SBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
56877ed5343SBarry Smith 
56977ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
57077ed5343SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
57177ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
5723270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
57377ed5343SBarry Smith   ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A); CHKERRQ(ierr);
57477ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5753a40ed3dSBarry Smith   PetscFunctionReturn(0);
5763270192aSSatish Balay }
5773270192aSSatish Balay 
5785615d1e5SSatish Balay #undef __FUNC__
579d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
580e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer)
5812593348eSBarry Smith {
58219bcc07fSBarry Smith   ViewerType  vtype;
58319bcc07fSBarry Smith   int         ierr;
5842593348eSBarry Smith 
5853a40ed3dSBarry Smith   PetscFunctionBegin;
58619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5877b2a1423SBarry Smith   if (PetscTypeCompare(vtype,SOCKET_VIEWER)) {
5887b2a1423SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported");
5893f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){
5903a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
5913f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) {
5923a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
5933f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) {
5943a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
5955cd90555SBarry Smith   } else {
5965cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
5972593348eSBarry Smith   }
5983a40ed3dSBarry Smith   PetscFunctionReturn(0);
5992593348eSBarry Smith }
600b6490206SBarry Smith 
601cd0e1443SSatish Balay 
6025615d1e5SSatish Balay #undef __FUNC__
6032d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
6042d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
605cd0e1443SSatish Balay {
606cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6072d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
6082d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
6092d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6103f1db9ecSBarry Smith   MatScalar  *ap, *aa = a->a, zero = 0.0;
611cd0e1443SSatish Balay 
6123a40ed3dSBarry Smith   PetscFunctionBegin;
6132d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
614cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
615a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
616a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
6172d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6182c3acbe9SBarry Smith     nrow = ailen[brow];
6192d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
620a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
621a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6222d61bbb3SSatish Balay       col  = in[l] ;
6232d61bbb3SSatish Balay       bcol = col/bs;
6242d61bbb3SSatish Balay       cidx = col%bs;
6252d61bbb3SSatish Balay       ridx = row%bs;
6262d61bbb3SSatish Balay       high = nrow;
6272d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6282d61bbb3SSatish Balay       while (high-low > 5) {
629cd0e1443SSatish Balay         t = (low+high)/2;
630cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
631cd0e1443SSatish Balay         else             low  = t;
632cd0e1443SSatish Balay       }
633cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
634cd0e1443SSatish Balay         if (rp[i] > bcol) break;
635cd0e1443SSatish Balay         if (rp[i] == bcol) {
6362d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6372d61bbb3SSatish Balay           goto finished;
638cd0e1443SSatish Balay         }
639cd0e1443SSatish Balay       }
6402d61bbb3SSatish Balay       *v++ = zero;
6412d61bbb3SSatish Balay       finished:;
642cd0e1443SSatish Balay     }
643cd0e1443SSatish Balay   }
6443a40ed3dSBarry Smith   PetscFunctionReturn(0);
645cd0e1443SSatish Balay }
646cd0e1443SSatish Balay 
6472d61bbb3SSatish Balay 
6485615d1e5SSatish Balay #undef __FUNC__
64905a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
65092c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
65192c4ed94SBarry Smith {
65292c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6538a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
654dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
655dd9472c6SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
6563f1db9ecSBarry Smith   Scalar      *value = v;
6573f1db9ecSBarry Smith   MatScalar   *ap,*aa=a->a,*bap;
65892c4ed94SBarry Smith 
6593a40ed3dSBarry Smith   PetscFunctionBegin;
6600e324ae4SSatish Balay   if (roworiented) {
6610e324ae4SSatish Balay     stepval = (n-1)*bs;
6620e324ae4SSatish Balay   } else {
6630e324ae4SSatish Balay     stepval = (m-1)*bs;
6640e324ae4SSatish Balay   }
66592c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
66692c4ed94SBarry Smith     row  = im[k];
6675ef9f2a5SBarry Smith     if (row < 0) continue;
6683a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
669a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
67092c4ed94SBarry Smith #endif
67192c4ed94SBarry Smith     rp   = aj + ai[row];
67292c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
67392c4ed94SBarry Smith     rmax = imax[row];
67492c4ed94SBarry Smith     nrow = ailen[row];
67592c4ed94SBarry Smith     low  = 0;
67692c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6775ef9f2a5SBarry Smith       if (in[l] < 0) continue;
6783a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
679a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
68092c4ed94SBarry Smith #endif
68192c4ed94SBarry Smith       col = in[l];
68292c4ed94SBarry Smith       if (roworiented) {
6830e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
6840e324ae4SSatish Balay       } else {
6850e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
68692c4ed94SBarry Smith       }
68792c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
68892c4ed94SBarry Smith       while (high-low > 7) {
68992c4ed94SBarry Smith         t = (low+high)/2;
69092c4ed94SBarry Smith         if (rp[t] > col) high = t;
69192c4ed94SBarry Smith         else             low  = t;
69292c4ed94SBarry Smith       }
69392c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
69492c4ed94SBarry Smith         if (rp[i] > col) break;
69592c4ed94SBarry Smith         if (rp[i] == col) {
6968a84c255SSatish Balay           bap  = ap +  bs2*i;
6970e324ae4SSatish Balay           if (roworiented) {
6988a84c255SSatish Balay             if (is == ADD_VALUES) {
699dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
700dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7018a84c255SSatish Balay                   bap[jj] += *value++;
702dd9472c6SBarry Smith                 }
703dd9472c6SBarry Smith               }
7040e324ae4SSatish Balay             } else {
705dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
706dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7070e324ae4SSatish Balay                   bap[jj] = *value++;
7088a84c255SSatish Balay                 }
709dd9472c6SBarry Smith               }
710dd9472c6SBarry Smith             }
7110e324ae4SSatish Balay           } else {
7120e324ae4SSatish Balay             if (is == ADD_VALUES) {
713dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
714dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7150e324ae4SSatish Balay                   *bap++ += *value++;
716dd9472c6SBarry Smith                 }
717dd9472c6SBarry Smith               }
7180e324ae4SSatish Balay             } else {
719dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
720dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7210e324ae4SSatish Balay                   *bap++  = *value++;
7220e324ae4SSatish Balay                 }
7238a84c255SSatish Balay               }
724dd9472c6SBarry Smith             }
725dd9472c6SBarry Smith           }
726f1241b54SBarry Smith           goto noinsert2;
72792c4ed94SBarry Smith         }
72892c4ed94SBarry Smith       }
72989280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
730a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
73192c4ed94SBarry Smith       if (nrow >= rmax) {
73292c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
73392c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
7343f1db9ecSBarry Smith         MatScalar *new_a;
73592c4ed94SBarry Smith 
736a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
73789280ab3SLois Curfman McInnes 
73892c4ed94SBarry Smith         /* malloc new storage space */
7393f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
7403f1db9ecSBarry Smith         new_a   = (MatScalar *) PetscMalloc( len ); CHKPTRQ(new_a);
74192c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
74292c4ed94SBarry Smith         new_i   = new_j + new_nz;
74392c4ed94SBarry Smith 
74492c4ed94SBarry Smith         /* copy over old data into new slots */
74592c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
74692c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
74792c4ed94SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));
74892c4ed94SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow);
749dd9472c6SBarry Smith         PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));
7503f1db9ecSBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));
7513f1db9ecSBarry Smith         PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));
7523f1db9ecSBarry Smith         PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));
75392c4ed94SBarry Smith         /* free up old matrix storage */
75492c4ed94SBarry Smith         PetscFree(a->a);
75592c4ed94SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
75692c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
75792c4ed94SBarry Smith         a->singlemalloc = 1;
75892c4ed94SBarry Smith 
75992c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
76092c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
7613f1db9ecSBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
76292c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
76392c4ed94SBarry Smith         a->reallocs++;
76492c4ed94SBarry Smith         a->nz++;
76592c4ed94SBarry Smith       }
76692c4ed94SBarry Smith       N = nrow++ - 1;
76792c4ed94SBarry Smith       /* shift up all the later entries in this row */
76892c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
76992c4ed94SBarry Smith         rp[ii+1] = rp[ii];
7703f1db9ecSBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
77192c4ed94SBarry Smith       }
7723f1db9ecSBarry Smith       if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
77392c4ed94SBarry Smith       rp[i] = col;
7748a84c255SSatish Balay       bap   = ap +  bs2*i;
7750e324ae4SSatish Balay       if (roworiented) {
776dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
777dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
7780e324ae4SSatish Balay             bap[jj] = *value++;
779dd9472c6SBarry Smith           }
780dd9472c6SBarry Smith         }
7810e324ae4SSatish Balay       } else {
782dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
783dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
7840e324ae4SSatish Balay             *bap++  = *value++;
7850e324ae4SSatish Balay           }
786dd9472c6SBarry Smith         }
787dd9472c6SBarry Smith       }
788f1241b54SBarry Smith       noinsert2:;
78992c4ed94SBarry Smith       low = i;
79092c4ed94SBarry Smith     }
79192c4ed94SBarry Smith     ailen[row] = nrow;
79292c4ed94SBarry Smith   }
7933a40ed3dSBarry Smith   PetscFunctionReturn(0);
79492c4ed94SBarry Smith }
79592c4ed94SBarry Smith 
796f2501298SSatish Balay 
7975615d1e5SSatish Balay #undef __FUNC__
7985615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
7998f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
800584200bdSSatish Balay {
801584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
802584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
803a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
80443ee02c3SBarry Smith   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0;
8053f1db9ecSBarry Smith   MatScalar  *aa = a->a, *ap;
806584200bdSSatish Balay 
8073a40ed3dSBarry Smith   PetscFunctionBegin;
8083a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
809584200bdSSatish Balay 
81043ee02c3SBarry Smith   if (m) rmax = ailen[0];
811584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
812584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
813584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
814d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
815584200bdSSatish Balay     if (fshift) {
816a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
817584200bdSSatish Balay       N = ailen[i];
818584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
819584200bdSSatish Balay         ip[j-fshift] = ip[j];
8203f1db9ecSBarry Smith         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
821584200bdSSatish Balay       }
822584200bdSSatish Balay     }
823584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
824584200bdSSatish Balay   }
825584200bdSSatish Balay   if (mbs) {
826584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
827584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
828584200bdSSatish Balay   }
829584200bdSSatish Balay   /* reset ilen and imax for each row */
830584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
831584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
832584200bdSSatish Balay   }
833a7c10996SSatish Balay   a->nz = ai[mbs];
834584200bdSSatish Balay 
835584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
836584200bdSSatish Balay   if (fshift && a->diag) {
837584200bdSSatish Balay     PetscFree(a->diag);
838584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
839584200bdSSatish Balay     a->diag = 0;
840584200bdSSatish Balay   }
8413d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
842219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8433d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
844584200bdSSatish Balay            a->reallocs);
845d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
846e2f3b5e9SSatish Balay   a->reallocs          = 0;
8474e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8484e220ebcSLois Curfman McInnes 
8493a40ed3dSBarry Smith   PetscFunctionReturn(0);
850584200bdSSatish Balay }
851584200bdSSatish Balay 
8525a838052SSatish Balay 
853bea157c4SSatish Balay 
854bea157c4SSatish Balay /*
855bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
856bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
857bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
858bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
859bea157c4SSatish Balay */
8605615d1e5SSatish Balay #undef __FUNC__
861bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks"
862bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
863d9b7c43dSSatish Balay {
864bea157c4SSatish Balay   int i,j,k,row;
865bea157c4SSatish Balay   PetscTruth flg;
8663a40ed3dSBarry Smith 
867bea157c4SSatish Balay   /*   PetscFunctionBegin;*/
868bea157c4SSatish Balay   for ( i=0,j=0; i<n; j++ ) {
869bea157c4SSatish Balay     row = idx[i];
870bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
871bea157c4SSatish Balay       sizes[j] = 1;
872bea157c4SSatish Balay       i++;
873e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
874bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
875bea157c4SSatish Balay       i++;
876bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
877bea157c4SSatish Balay       flg = PETSC_TRUE;
878bea157c4SSatish Balay       for ( k=1; k<bs; k++ ) {
879bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
880bea157c4SSatish Balay           flg = PETSC_FALSE;
881bea157c4SSatish Balay           break;
882d9b7c43dSSatish Balay         }
883bea157c4SSatish Balay       }
884bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
885bea157c4SSatish Balay         sizes[j] = bs;
886bea157c4SSatish Balay         i+= bs;
887bea157c4SSatish Balay       } else {
888bea157c4SSatish Balay         sizes[j] = 1;
889bea157c4SSatish Balay         i++;
890bea157c4SSatish Balay       }
891bea157c4SSatish Balay     }
892bea157c4SSatish Balay   }
893bea157c4SSatish Balay   *bs_max = j;
8943a40ed3dSBarry Smith   PetscFunctionReturn(0);
895d9b7c43dSSatish Balay }
896d9b7c43dSSatish Balay 
8975615d1e5SSatish Balay #undef __FUNC__
8985615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
8998f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
900d9b7c43dSSatish Balay {
901d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
902b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
903bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
9043f1db9ecSBarry Smith   Scalar      zero = 0.0;
9053f1db9ecSBarry Smith   MatScalar   *aa;
906d9b7c43dSSatish Balay 
9073a40ed3dSBarry Smith   PetscFunctionBegin;
908d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
909d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
910d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
911d9b7c43dSSatish Balay 
912bea157c4SSatish Balay   /* allocate memory for rows,sizes */
913bea157c4SSatish Balay   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int)); CHKPTRQ(rows);
914bea157c4SSatish Balay   sizes = rows + is_n;
915bea157c4SSatish Balay 
916bea157c4SSatish Balay   /* initialize copy IS valurs to rows, and sort them */
917bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
918bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows); CHKERRQ(ierr);
919bea157c4SSatish Balay   ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max); CHKERRQ(ierr);
920b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
921bea157c4SSatish Balay 
922bea157c4SSatish Balay   for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) {
923bea157c4SSatish Balay     row   = rows[j];
924b6815fffSBarry Smith     if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row);
925bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
926bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
927bea157c4SSatish Balay     if (sizes[i] == bs) {
928bea157c4SSatish Balay       if (diag) {
929bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
930bea157c4SSatish Balay           baij->ilen[row/bs] = 1;
931bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
932bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar)); CHKERRQ(ierr);
933a07cd24cSSatish Balay         }
934bea157c4SSatish Balay         /* Now insert all the diagoanl values for this bs */
935bea157c4SSatish Balay         for ( k=0; k<bs; k++ ) {
936bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
937bea157c4SSatish Balay         }
938bea157c4SSatish Balay       } else { /* (!diag) */
939bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
940bea157c4SSatish Balay       } /* end (!diag) */
941bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
942bea157c4SSatish Balay #if defined (USE_PETSC_DEBUG)
943bea157c4SSatish Balay       if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1");
944bea157c4SSatish Balay #endif
945bea157c4SSatish Balay       for ( k=0; k<count; k++ ) {
946d9b7c43dSSatish Balay         aa[0] = zero;
947d9b7c43dSSatish Balay         aa+=bs;
948d9b7c43dSSatish Balay       }
949d9b7c43dSSatish Balay       if (diag) {
950f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
951d9b7c43dSSatish Balay       }
952d9b7c43dSSatish Balay     }
953bea157c4SSatish Balay   }
954bea157c4SSatish Balay 
955bea157c4SSatish Balay   PetscFree(rows);
9569a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9573a40ed3dSBarry Smith   PetscFunctionReturn(0);
958d9b7c43dSSatish Balay }
9591c351548SSatish Balay 
9605615d1e5SSatish Balay #undef __FUNC__
9612d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9622d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9632d61bbb3SSatish Balay {
9642d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9652d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9662d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9672d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
9682d61bbb3SSatish Balay   int         ridx,cidx,bs2=a->bs2;
9693f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
9702d61bbb3SSatish Balay 
9712d61bbb3SSatish Balay   PetscFunctionBegin;
9722d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
9732d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
9745ef9f2a5SBarry Smith     if (row < 0) continue;
9752d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9762d61bbb3SSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
9772d61bbb3SSatish Balay #endif
9782d61bbb3SSatish Balay     rp   = aj + ai[brow];
9792d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
9802d61bbb3SSatish Balay     rmax = imax[brow];
9812d61bbb3SSatish Balay     nrow = ailen[brow];
9822d61bbb3SSatish Balay     low  = 0;
9832d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
9845ef9f2a5SBarry Smith       if (in[l] < 0) continue;
9852d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9862d61bbb3SSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
9872d61bbb3SSatish Balay #endif
9882d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
9892d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
9902d61bbb3SSatish Balay       if (roworiented) {
9915ef9f2a5SBarry Smith         value = v[l + k*n];
9922d61bbb3SSatish Balay       } else {
9932d61bbb3SSatish Balay         value = v[k + l*m];
9942d61bbb3SSatish Balay       }
9952d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
9962d61bbb3SSatish Balay       while (high-low > 7) {
9972d61bbb3SSatish Balay         t = (low+high)/2;
9982d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
9992d61bbb3SSatish Balay         else              low  = t;
10002d61bbb3SSatish Balay       }
10012d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
10022d61bbb3SSatish Balay         if (rp[i] > bcol) break;
10032d61bbb3SSatish Balay         if (rp[i] == bcol) {
10042d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
10052d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
10062d61bbb3SSatish Balay           else                  *bap  = value;
10072d61bbb3SSatish Balay           goto noinsert1;
10082d61bbb3SSatish Balay         }
10092d61bbb3SSatish Balay       }
10102d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
10112d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
10122d61bbb3SSatish Balay       if (nrow >= rmax) {
10132d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
10142d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
10153f1db9ecSBarry Smith         MatScalar *new_a;
10162d61bbb3SSatish Balay 
10172d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
10182d61bbb3SSatish Balay 
10192d61bbb3SSatish Balay         /* Malloc new storage space */
10203f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
10213f1db9ecSBarry Smith         new_a   = (MatScalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10222d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
10232d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10242d61bbb3SSatish Balay 
10252d61bbb3SSatish Balay         /* copy over old data into new slots */
10262d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
10272d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
10282d61bbb3SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
10292d61bbb3SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
10302d61bbb3SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
10312d61bbb3SSatish Balay                                                            len*sizeof(int));
10323f1db9ecSBarry Smith         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));
10333f1db9ecSBarry Smith         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));
10342d61bbb3SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
10353f1db9ecSBarry Smith                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));
10362d61bbb3SSatish Balay         /* free up old matrix storage */
10372d61bbb3SSatish Balay         PetscFree(a->a);
10382d61bbb3SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
10392d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
10402d61bbb3SSatish Balay         a->singlemalloc = 1;
10412d61bbb3SSatish Balay 
10422d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10432d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
10443f1db9ecSBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
10452d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10462d61bbb3SSatish Balay         a->reallocs++;
10472d61bbb3SSatish Balay         a->nz++;
10482d61bbb3SSatish Balay       }
10492d61bbb3SSatish Balay       N = nrow++ - 1;
10502d61bbb3SSatish Balay       /* shift up all the later entries in this row */
10512d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
10522d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
10533f1db9ecSBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
10542d61bbb3SSatish Balay       }
10553f1db9ecSBarry Smith       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
10562d61bbb3SSatish Balay       rp[i]                      = bcol;
10572d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
10582d61bbb3SSatish Balay       noinsert1:;
10592d61bbb3SSatish Balay       low = i;
10602d61bbb3SSatish Balay     }
10612d61bbb3SSatish Balay     ailen[brow] = nrow;
10622d61bbb3SSatish Balay   }
10632d61bbb3SSatish Balay   PetscFunctionReturn(0);
10642d61bbb3SSatish Balay }
10652d61bbb3SSatish Balay 
10662d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
10672d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
10682d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
10697b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
10707b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
10712d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
10722d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
10732d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
10742d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
10752d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
10762d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
10772d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
10782d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
10792d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
10802d61bbb3SSatish Balay 
10812d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10822d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
10832d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
10842d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
10852d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
10862d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
1087*15091d37SBarry Smith extern int MatSolve_SeqBAIJ_6(Mat,Vec,Vec);
10882d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
10892d61bbb3SSatish Balay 
10902d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
10912d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
10922d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
10932d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
10942d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
10952d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
1096*15091d37SBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_6(Mat,Mat*);
10972d61bbb3SSatish Balay 
10982d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
10992d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
11002d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
11012d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
11022d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
1103*15091d37SBarry Smith extern int MatMult_SeqBAIJ_6(Mat,Vec,Vec);
11042d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
11052d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
11062d61bbb3SSatish Balay 
11072d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
11082d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
11092d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
11102d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
11112d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
1112*15091d37SBarry Smith extern int MatMultAdd_SeqBAIJ_6(Mat,Vec,Vec,Vec);
11132d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
11142d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
11152d61bbb3SSatish Balay 
11162d61bbb3SSatish Balay #undef __FUNC__
11172d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
11185ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
11192d61bbb3SSatish Balay {
11202d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
11212d61bbb3SSatish Balay   Mat         outA;
11222d61bbb3SSatish Balay   int         ierr;
1123667159a5SBarry Smith   PetscTruth  row_identity, col_identity;
11242d61bbb3SSatish Balay 
11252d61bbb3SSatish Balay   PetscFunctionBegin;
11265ef9f2a5SBarry Smith   if (info && info->fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1127667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity); CHKERRQ(ierr);
1128667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity); CHKERRQ(ierr);
1129667159a5SBarry Smith   if (!row_identity || !col_identity) {
1130667159a5SBarry Smith     SETERRQ(1,1,"Row and column permutations must be identity");
1131667159a5SBarry Smith   }
11322d61bbb3SSatish Balay 
11332d61bbb3SSatish Balay   outA          = inA;
11342d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11352d61bbb3SSatish Balay   a->row        = row;
11362d61bbb3SSatish Balay   a->col        = col;
11372d61bbb3SSatish Balay 
1138e51c0b9cSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1139e51c0b9cSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
11401e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1141e51c0b9cSSatish Balay 
11422d61bbb3SSatish Balay   if (!a->solve_work) {
11432d61bbb3SSatish Balay     a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
11442d61bbb3SSatish Balay     PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
11452d61bbb3SSatish Balay   }
11462d61bbb3SSatish Balay 
11472d61bbb3SSatish Balay   if (!a->diag) {
11482d61bbb3SSatish Balay     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
11492d61bbb3SSatish Balay   }
11502d61bbb3SSatish Balay   /*
1151*15091d37SBarry Smith       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1152*15091d37SBarry Smith       for ILU(0) factorization with natural ordering
11532d61bbb3SSatish Balay   */
1154*15091d37SBarry Smith   switch (a->bs) {
1155*15091d37SBarry Smith     case 2:
1156*15091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
1157*15091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering;
1158*15091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
1159*15091d37SBarry Smith     break;
1160*15091d37SBarry Smith   case 3:
1161*15091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
1162*15091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_3_NaturalOrdering;
1163*15091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
1164*15091d37SBarry Smith     break;
1165*15091d37SBarry Smith   case 4:
1166667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
1167f830108cSBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
1168*15091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
1169*15091d37SBarry Smith     break;
1170*15091d37SBarry Smith   case 5:
1171667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
1172667159a5SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
1173*15091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
1174*15091d37SBarry Smith     break;
1175*15091d37SBarry Smith   case 6:
1176*15091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
1177*15091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_6_NaturalOrdering;
1178*15091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
1179*15091d37SBarry Smith     break;
1180*15091d37SBarry Smith   case 7:
1181*15091d37SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
1182*15091d37SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_7_NaturalOrdering;
1183*15091d37SBarry Smith     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
1184*15091d37SBarry Smith     break;
11852d61bbb3SSatish Balay   }
11862d61bbb3SSatish Balay 
1187667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
1188667159a5SBarry Smith 
1189667159a5SBarry Smith 
11902d61bbb3SSatish Balay   PetscFunctionReturn(0);
11912d61bbb3SSatish Balay }
11922d61bbb3SSatish Balay #undef __FUNC__
1193d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1194ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1195ba4ca20aSSatish Balay {
1196ba4ca20aSSatish Balay   static int called = 0;
1197ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1198ba4ca20aSSatish Balay 
11993a40ed3dSBarry Smith   PetscFunctionBegin;
12003a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
120176be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
120276be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");
12033a40ed3dSBarry Smith   PetscFunctionReturn(0);
1204ba4ca20aSSatish Balay }
1205d9b7c43dSSatish Balay 
1206fb2e594dSBarry Smith EXTERN_C_BEGIN
120727a8da17SBarry Smith #undef __FUNC__
120827a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
120927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
121027a8da17SBarry Smith {
121127a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
121227a8da17SBarry Smith   int         i,nz,n;
121327a8da17SBarry Smith 
121427a8da17SBarry Smith   PetscFunctionBegin;
121527a8da17SBarry Smith   nz = baij->maxnz;
121627a8da17SBarry Smith   n  = baij->n;
121727a8da17SBarry Smith   for (i=0; i<nz; i++) {
121827a8da17SBarry Smith     baij->j[i] = indices[i];
121927a8da17SBarry Smith   }
122027a8da17SBarry Smith   baij->nz = nz;
122127a8da17SBarry Smith   for ( i=0; i<n; i++ ) {
122227a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
122327a8da17SBarry Smith   }
122427a8da17SBarry Smith 
122527a8da17SBarry Smith   PetscFunctionReturn(0);
122627a8da17SBarry Smith }
1227fb2e594dSBarry Smith EXTERN_C_END
122827a8da17SBarry Smith 
122927a8da17SBarry Smith #undef __FUNC__
123027a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices"
123127a8da17SBarry Smith /*@
123227a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
123327a8da17SBarry Smith        in the matrix.
123427a8da17SBarry Smith 
123527a8da17SBarry Smith   Input Parameters:
123627a8da17SBarry Smith +  mat - the SeqBAIJ matrix
123727a8da17SBarry Smith -  indices - the column indices
123827a8da17SBarry Smith 
1239*15091d37SBarry Smith   Level: advanced
1240*15091d37SBarry Smith 
124127a8da17SBarry Smith   Notes:
124227a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
124327a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
124427a8da17SBarry Smith   of the MatSetValues() operation.
124527a8da17SBarry Smith 
124627a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
124727a8da17SBarry Smith   MatCreateSeqBAIJ().
124827a8da17SBarry Smith 
124927a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
125027a8da17SBarry Smith 
125127a8da17SBarry Smith @*/
125227a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
125327a8da17SBarry Smith {
125427a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
125527a8da17SBarry Smith 
125627a8da17SBarry Smith   PetscFunctionBegin;
125727a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
125827a8da17SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);
125927a8da17SBarry Smith          CHKERRQ(ierr);
126027a8da17SBarry Smith   if (f) {
126127a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
126227a8da17SBarry Smith   } else {
126327a8da17SBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
126427a8da17SBarry Smith   }
126527a8da17SBarry Smith   PetscFunctionReturn(0);
126627a8da17SBarry Smith }
126727a8da17SBarry Smith 
12682593348eSBarry Smith /* -------------------------------------------------------------------*/
1269cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1270cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1271cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1272cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1273cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
1274cc2dc46cSBarry Smith        MatMultTrans_SeqBAIJ,
1275cc2dc46cSBarry Smith        MatMultTransAdd_SeqBAIJ,
1276cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1277cc2dc46cSBarry Smith        0,
1278cc2dc46cSBarry Smith        0,
1279cc2dc46cSBarry Smith        0,
1280cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1281cc2dc46cSBarry Smith        0,
1282b6490206SBarry Smith        0,
1283f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1284cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1285cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1286cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1287cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1288cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1289b6490206SBarry Smith        0,
1290cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1291cc2dc46cSBarry Smith        0,
1292cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1293cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1294cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1295cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1296cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1297cc2dc46cSBarry Smith        0,
1298cc2dc46cSBarry Smith        0,
1299cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1300cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1301cc2dc46cSBarry Smith        MatGetOwnershipRange_SeqBAIJ,
1302cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1303cc2dc46cSBarry Smith        0,
1304cc2dc46cSBarry Smith        0,
1305cc2dc46cSBarry Smith        0,
13062e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1307cc2dc46cSBarry Smith        0,
1308cc2dc46cSBarry Smith        0,
1309cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1310cc2dc46cSBarry Smith        0,
1311cc2dc46cSBarry Smith        0,
1312cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1313cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1314cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1315cc2dc46cSBarry Smith        0,
1316cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1317cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1318cc2dc46cSBarry Smith        0,
1319cc2dc46cSBarry Smith        0,
1320cc2dc46cSBarry Smith        0,
1321cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
13223b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
132392c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1324cc2dc46cSBarry Smith        0,
1325cc2dc46cSBarry Smith        0,
1326cc2dc46cSBarry Smith        0,
1327cc2dc46cSBarry Smith        0,
1328cc2dc46cSBarry Smith        0,
1329cc2dc46cSBarry Smith        0,
1330d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1331cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1332cc2dc46cSBarry Smith        0,
1333cc2dc46cSBarry Smith        0,
1334cc2dc46cSBarry Smith        MatGetMaps_Petsc};
13352593348eSBarry Smith 
13365615d1e5SSatish Balay #undef __FUNC__
13375615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
13382593348eSBarry Smith /*@C
133944cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
134044cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
134144cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
13422bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
13432bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
13442593348eSBarry Smith 
1345db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1346db81eaa0SLois Curfman McInnes 
13472593348eSBarry Smith    Input Parameters:
1348db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
1349b6490206SBarry Smith .  bs - size of block
13502593348eSBarry Smith .  m - number of rows
13512593348eSBarry Smith .  n - number of columns
1352b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1353db81eaa0SLois Curfman McInnes -  nzz - array containing the number of block nonzeros in the various block rows
13542bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
13552593348eSBarry Smith 
13562593348eSBarry Smith    Output Parameter:
13572593348eSBarry Smith .  A - the matrix
13582593348eSBarry Smith 
13590513a670SBarry Smith    Options Database Keys:
1360db81eaa0SLois Curfman McInnes .   -mat_no_unroll - uses code that does not unroll the loops in the
1361db81eaa0SLois Curfman McInnes                      block calculations (much slower)
1362db81eaa0SLois Curfman McInnes .    -mat_block_size - size of the blocks to use
13630513a670SBarry Smith 
1364*15091d37SBarry Smith    Level: intermediate
1365*15091d37SBarry Smith 
13662593348eSBarry Smith    Notes:
136744cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
13682593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
136944cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
13702593348eSBarry Smith 
13712593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
13722593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
13732593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
13746da5968aSLois Curfman McInnes    matrices.
13752593348eSBarry Smith 
1376db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
13772593348eSBarry Smith @*/
1378b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
13792593348eSBarry Smith {
13802593348eSBarry Smith   Mat         B;
1381b6490206SBarry Smith   Mat_SeqBAIJ *b;
13823b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
13833b2fbd54SBarry Smith 
13843a40ed3dSBarry Smith   PetscFunctionBegin;
13853b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1386a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1387b6490206SBarry Smith 
13880513a670SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
13890513a670SBarry Smith 
13903a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1391a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
13923a40ed3dSBarry Smith   }
13932593348eSBarry Smith 
13942593348eSBarry Smith   *A = 0;
13953f1db9ecSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView);
13962593348eSBarry Smith   PLogObjectCreate(B);
1397b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
139844cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
1399cc2dc46cSBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
14001a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
14011a6a6d98SBarry Smith   if (!flg) {
14027fc0212eSBarry Smith     switch (bs) {
14037fc0212eSBarry Smith     case 1:
1404f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1405f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1406f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1407f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
14087fc0212eSBarry Smith       break;
14094eeb42bcSBarry Smith     case 2:
1410f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1411f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1412f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1413f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
14146d84be18SBarry Smith       break;
1415f327dd97SBarry Smith     case 3:
1416f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1417f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1418f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1419f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
14204eeb42bcSBarry Smith       break;
14216d84be18SBarry Smith     case 4:
1422f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1423f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1424f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1425f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
14266d84be18SBarry Smith       break;
14276d84be18SBarry Smith     case 5:
1428f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1429f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1430f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1431f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
14326d84be18SBarry Smith       break;
1433*15091d37SBarry Smith     case 6:
1434*15091d37SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
1435*15091d37SBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_6;
1436*15091d37SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
1437*15091d37SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
1438*15091d37SBarry Smith       break;
143948e9ddb2SSatish Balay     case 7:
1440f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1441f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1442f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
144348e9ddb2SSatish Balay       break;
14447fc0212eSBarry Smith     }
14451a6a6d98SBarry Smith   }
1446e1311b90SBarry Smith   B->ops->destroy     = MatDestroy_SeqBAIJ;
1447e1311b90SBarry Smith   B->ops->view        = MatView_SeqBAIJ;
14482593348eSBarry Smith   B->factor           = 0;
14492593348eSBarry Smith   B->lupivotthreshold = 1.0;
145090f02eecSBarry Smith   B->mapping          = 0;
14512593348eSBarry Smith   b->row              = 0;
14522593348eSBarry Smith   b->col              = 0;
1453e51c0b9cSSatish Balay   b->icol             = 0;
14542593348eSBarry Smith   b->reallocs         = 0;
14552593348eSBarry Smith 
145644cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
145744cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1458a5ae1ecdSBarry Smith 
14597413bad6SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
14607413bad6SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1461a5ae1ecdSBarry Smith 
1462b6490206SBarry Smith   b->mbs     = mbs;
1463f2501298SSatish Balay   b->nbs     = nbs;
1464b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
14652593348eSBarry Smith   if (nnz == PETSC_NULL) {
1466119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
14672593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1468b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1469b6490206SBarry Smith     nz = nz*mbs;
14703a40ed3dSBarry Smith   } else {
14712593348eSBarry Smith     nz = 0;
1472b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
14732593348eSBarry Smith   }
14742593348eSBarry Smith 
14752593348eSBarry Smith   /* allocate the matrix space */
14763f1db9ecSBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
14773f1db9ecSBarry Smith   b->a  = (MatScalar *) PetscMalloc( len ); CHKPTRQ(b->a);
14783f1db9ecSBarry Smith   PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
14797e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
14802593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
14812593348eSBarry Smith   b->i  = b->j + nz;
14822593348eSBarry Smith   b->singlemalloc = 1;
14832593348eSBarry Smith 
1484de6a44a3SBarry Smith   b->i[0] = 0;
1485b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
14862593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
14872593348eSBarry Smith   }
14882593348eSBarry Smith 
1489b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1490b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1491f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1492b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
14932593348eSBarry Smith 
1494b6490206SBarry Smith   b->bs               = bs;
14959df24120SSatish Balay   b->bs2              = bs2;
1496b6490206SBarry Smith   b->mbs              = mbs;
14972593348eSBarry Smith   b->nz               = 0;
149818e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
14992593348eSBarry Smith   b->sorted           = 0;
15002593348eSBarry Smith   b->roworiented      = 1;
15012593348eSBarry Smith   b->nonew            = 0;
15022593348eSBarry Smith   b->diag             = 0;
15032593348eSBarry Smith   b->solve_work       = 0;
1504de6a44a3SBarry Smith   b->mult_work        = 0;
15052593348eSBarry Smith   b->spptr            = 0;
15064e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
15074e220ebcSLois Curfman McInnes 
15082593348eSBarry Smith   *A = B;
15092593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
15102593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
151127a8da17SBarry Smith 
151227a8da17SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
151327a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
151427a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
15153a40ed3dSBarry Smith   PetscFunctionReturn(0);
15162593348eSBarry Smith }
15172593348eSBarry Smith 
15185615d1e5SSatish Balay #undef __FUNC__
15192e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ"
15202e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
15212593348eSBarry Smith {
15222593348eSBarry Smith   Mat         C;
1523b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
152427a8da17SBarry Smith   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1525de6a44a3SBarry Smith 
15263a40ed3dSBarry Smith   PetscFunctionBegin;
1527a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
15282593348eSBarry Smith 
15292593348eSBarry Smith   *B = 0;
15303f1db9ecSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView);
15312593348eSBarry Smith   PLogObjectCreate(C);
1532b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
1533c3c66ee5SSatish Balay   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1534e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqBAIJ;
1535e1311b90SBarry Smith   C->ops->view       = MatView_SeqBAIJ;
15362593348eSBarry Smith   C->factor          = A->factor;
15372593348eSBarry Smith   c->row             = 0;
15382593348eSBarry Smith   c->col             = 0;
1539cac97260SSatish Balay   c->icol            = 0;
15402593348eSBarry Smith   C->assembled       = PETSC_TRUE;
15412593348eSBarry Smith 
154244cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
154344cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
154444cd7ae7SLois Curfman McInnes   C->M          = a->m;
154544cd7ae7SLois Curfman McInnes   C->N          = a->n;
154644cd7ae7SLois Curfman McInnes 
1547b6490206SBarry Smith   c->bs         = a->bs;
15489df24120SSatish Balay   c->bs2        = a->bs2;
1549b6490206SBarry Smith   c->mbs        = a->mbs;
1550b6490206SBarry Smith   c->nbs        = a->nbs;
15512593348eSBarry Smith 
1552b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1553b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1554b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
15552593348eSBarry Smith     c->imax[i] = a->imax[i];
15562593348eSBarry Smith     c->ilen[i] = a->ilen[i];
15572593348eSBarry Smith   }
15582593348eSBarry Smith 
15592593348eSBarry Smith   /* allocate the matrix space */
15602593348eSBarry Smith   c->singlemalloc = 1;
15613f1db9ecSBarry Smith   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
15623f1db9ecSBarry Smith   c->a  = (MatScalar *) PetscMalloc( len ); CHKPTRQ(c->a);
15637e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1564de6a44a3SBarry Smith   c->i  = c->j + nz;
1565b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1566b6490206SBarry Smith   if (mbs > 0) {
1567de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
15682e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
15693f1db9ecSBarry Smith       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
15702e8a6d31SBarry Smith     } else {
15713f1db9ecSBarry Smith       PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
15722593348eSBarry Smith     }
15732593348eSBarry Smith   }
15742593348eSBarry Smith 
1575f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
15762593348eSBarry Smith   c->sorted      = a->sorted;
15772593348eSBarry Smith   c->roworiented = a->roworiented;
15782593348eSBarry Smith   c->nonew       = a->nonew;
15792593348eSBarry Smith 
15802593348eSBarry Smith   if (a->diag) {
1581b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1582b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1583b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
15842593348eSBarry Smith       c->diag[i] = a->diag[i];
15852593348eSBarry Smith     }
158698305bb5SBarry Smith   } else c->diag        = 0;
15872593348eSBarry Smith   c->nz                 = a->nz;
15882593348eSBarry Smith   c->maxnz              = a->maxnz;
15892593348eSBarry Smith   c->solve_work         = 0;
15902593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
15917fc0212eSBarry Smith   c->mult_work          = 0;
15922593348eSBarry Smith   *B = C;
15937b2a1423SBarry Smith   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
15943a40ed3dSBarry Smith   PetscFunctionReturn(0);
15952593348eSBarry Smith }
15962593348eSBarry Smith 
15975615d1e5SSatish Balay #undef __FUNC__
15985615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
159919bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
16002593348eSBarry Smith {
1601b6490206SBarry Smith   Mat_SeqBAIJ  *a;
16022593348eSBarry Smith   Mat          B;
1603de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1604b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
160535aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1606a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1607b6490206SBarry Smith   Scalar       *aa;
160819bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
16092593348eSBarry Smith 
16103a40ed3dSBarry Smith   PetscFunctionBegin;
16111a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1612de6a44a3SBarry Smith   bs2  = bs*bs;
1613b6490206SBarry Smith 
16142593348eSBarry Smith   MPI_Comm_size(comm,&size);
1615cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
161690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
16170752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1618a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
16192593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
16202593348eSBarry Smith 
1621d64ed03dSBarry Smith   if (header[3] < 0) {
1622a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1623d64ed03dSBarry Smith   }
1624d64ed03dSBarry Smith 
1625a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
162635aab85fSBarry Smith 
162735aab85fSBarry Smith   /*
162835aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
162935aab85fSBarry Smith     divisible by the blocksize
163035aab85fSBarry Smith   */
1631b6490206SBarry Smith   mbs        = M/bs;
163235aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
163335aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
163435aab85fSBarry Smith   else                  mbs++;
163535aab85fSBarry Smith   if (extra_rows) {
1636537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
163735aab85fSBarry Smith   }
1638b6490206SBarry Smith 
16392593348eSBarry Smith   /* read in row lengths */
164035aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
16410752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
164235aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
16432593348eSBarry Smith 
1644b6490206SBarry Smith   /* read in column indices */
164535aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
16460752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr);
164735aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1648b6490206SBarry Smith 
1649b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1650b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1651b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
165235aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
165335aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
165435aab85fSBarry Smith   masked = mask + mbs;
1655b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1656b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
165735aab85fSBarry Smith     nmask = 0;
1658b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1659b6490206SBarry Smith       kmax = rowlengths[rowcount];
1660b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
166135aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
166235aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1663b6490206SBarry Smith       }
1664b6490206SBarry Smith       rowcount++;
1665b6490206SBarry Smith     }
166635aab85fSBarry Smith     browlengths[i] += nmask;
166735aab85fSBarry Smith     /* zero out the mask elements we set */
166835aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1669b6490206SBarry Smith   }
1670b6490206SBarry Smith 
16712593348eSBarry Smith   /* create our matrix */
16723eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
16732593348eSBarry Smith   B = *A;
1674b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
16752593348eSBarry Smith 
16762593348eSBarry Smith   /* set matrix "i" values */
1677de6a44a3SBarry Smith   a->i[0] = 0;
1678b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1679b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1680b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
16812593348eSBarry Smith   }
1682b6490206SBarry Smith   a->nz         = 0;
1683b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
16842593348eSBarry Smith 
1685b6490206SBarry Smith   /* read in nonzero values */
168635aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
16870752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr);
168835aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1689b6490206SBarry Smith 
1690b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1691b6490206SBarry Smith   nzcount = 0; jcount = 0;
1692b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1693b6490206SBarry Smith     nzcountb = nzcount;
169435aab85fSBarry Smith     nmask    = 0;
1695b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1696b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1697b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
169835aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
169935aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1700b6490206SBarry Smith       }
1701b6490206SBarry Smith       rowcount++;
1702b6490206SBarry Smith     }
1703de6a44a3SBarry Smith     /* sort the masked values */
170477c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1705de6a44a3SBarry Smith 
1706b6490206SBarry Smith     /* set "j" values into matrix */
1707b6490206SBarry Smith     maskcount = 1;
170835aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
170935aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1710de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1711b6490206SBarry Smith     }
1712b6490206SBarry Smith     /* set "a" values into matrix */
1713de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1714b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1715b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1716b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1717de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1718de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1719de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1720de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1721b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1722b6490206SBarry Smith       }
1723b6490206SBarry Smith     }
172435aab85fSBarry Smith     /* zero out the mask elements we set */
172535aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1726b6490206SBarry Smith   }
1727a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1728b6490206SBarry Smith 
1729b6490206SBarry Smith   PetscFree(rowlengths);
1730b6490206SBarry Smith   PetscFree(browlengths);
1731b6490206SBarry Smith   PetscFree(aa);
1732b6490206SBarry Smith   PetscFree(jj);
1733b6490206SBarry Smith   PetscFree(mask);
1734b6490206SBarry Smith 
1735b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1736de6a44a3SBarry Smith 
17379c01be13SBarry Smith   ierr = MatView_Private(B); CHKERRQ(ierr);
17383a40ed3dSBarry Smith   PetscFunctionReturn(0);
17392593348eSBarry Smith }
17402593348eSBarry Smith 
17412593348eSBarry Smith 
17422593348eSBarry Smith 
1743