xref: /petsc/src/mat/impls/baij/seq/baij.c (revision e4fda26c8cf8e28b761bdf9c05f26aafb81d2f8c)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e4fda26cSSatish Balay static char vcid[] = "$Id: baij.c,v 1.161 1999/02/15 17:07:34 bsmith Exp balay $";
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;
1702d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERROR)   a->nonew       = -1;
1712d61bbb3SSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) 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 
3442d61bbb3SSatish Balay 
3452d61bbb3SSatish Balay 
3463b2fbd54SBarry Smith 
3475615d1e5SSatish Balay #undef __FUNC__
348d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Binary"
349b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
3502593348eSBarry Smith {
351b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3529df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
353b6490206SBarry Smith   Scalar      *aa;
3542593348eSBarry Smith 
3553a40ed3dSBarry Smith   PetscFunctionBegin;
35690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
3572593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
3582593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
3593b2fbd54SBarry Smith 
3602593348eSBarry Smith   col_lens[1] = a->m;
3612593348eSBarry Smith   col_lens[2] = a->n;
3627e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3632593348eSBarry Smith 
3642593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
365b6490206SBarry Smith   count = 0;
366b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
367b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
368b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
369b6490206SBarry Smith     }
3702593348eSBarry Smith   }
3710752156aSBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr);
3722593348eSBarry Smith   PetscFree(col_lens);
3732593348eSBarry Smith 
3742593348eSBarry Smith   /* store column indices (zero start index) */
37566963ce1SSatish Balay   jj = (int *) PetscMalloc( (a->nz+1)*bs2*sizeof(int) ); CHKPTRQ(jj);
376b6490206SBarry Smith   count = 0;
377b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
378b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
379b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
380b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
381b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3822593348eSBarry Smith         }
3832593348eSBarry Smith       }
384b6490206SBarry Smith     }
385b6490206SBarry Smith   }
3860752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0); CHKERRQ(ierr);
387b6490206SBarry Smith   PetscFree(jj);
3882593348eSBarry Smith 
3892593348eSBarry Smith   /* store nonzero values */
39066963ce1SSatish Balay   aa = (Scalar *) PetscMalloc((a->nz+1)*bs2*sizeof(Scalar)); CHKPTRQ(aa);
391b6490206SBarry Smith   count = 0;
392b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
393b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
394b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
395b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
3967e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
397b6490206SBarry Smith         }
398b6490206SBarry Smith       }
399b6490206SBarry Smith     }
400b6490206SBarry Smith   }
4010752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0); CHKERRQ(ierr);
402b6490206SBarry Smith   PetscFree(aa);
4033a40ed3dSBarry Smith   PetscFunctionReturn(0);
4042593348eSBarry Smith }
4052593348eSBarry Smith 
4065615d1e5SSatish Balay #undef __FUNC__
407d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ_ASCII"
408b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
4092593348eSBarry Smith {
410b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4119df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
4122593348eSBarry Smith   FILE        *fd;
4132593348eSBarry Smith   char        *outputname;
4142593348eSBarry Smith 
4153a40ed3dSBarry Smith   PetscFunctionBegin;
41690ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
41777ed5343SBarry Smith   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
41890ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
419639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
4200ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);
4210ef38995SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
4227b2a1423SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Socket format not supported");
4230ef38995SBarry Smith   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
42444cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
42544cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
42644cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
42744cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
42844cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
4293a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
4300ef38995SBarry Smith             if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
43144cd7ae7SLois Curfman McInnes               fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
432e20fef11SSatish Balay                       PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
4330ef38995SBarry Smith             } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0 && PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
434766eeae4SLois 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 (PetscReal(a->a[bs2*k + l*bs + j]) != 0.0) {
437e20fef11SSatish Balay               fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
4380ef38995SBarry Smith             }
43944cd7ae7SLois Curfman McInnes #else
4400ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
44144cd7ae7SLois Curfman McInnes               fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
4420ef38995SBarry Smith             }
44344cd7ae7SLois Curfman McInnes #endif
44444cd7ae7SLois Curfman McInnes           }
44544cd7ae7SLois Curfman McInnes         }
44644cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
44744cd7ae7SLois Curfman McInnes       }
44844cd7ae7SLois Curfman McInnes     }
4490ef38995SBarry Smith   } else {
450b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
451b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
452b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
453b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
454b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
4553a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
456e20fef11SSatish Balay             if (PetscImaginary(a->a[bs2*k + l*bs + j]) > 0.0) {
45788685aaeSLois Curfman McInnes               fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
458e20fef11SSatish Balay                 PetscReal(a->a[bs2*k + l*bs + j]),PetscImaginary(a->a[bs2*k + l*bs + j]));
4590ef38995SBarry Smith             } else if (PetscImaginary(a->a[bs2*k + l*bs + j]) < 0.0) {
460766eeae4SLois 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 {
463e20fef11SSatish Balay               fprintf(fd," %d %g ",bs*a->j[k]+l,PetscReal(a->a[bs2*k + l*bs + j]));
46488685aaeSLois Curfman McInnes             }
46588685aaeSLois Curfman McInnes #else
4667e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
46788685aaeSLois Curfman McInnes #endif
4682593348eSBarry Smith           }
4692593348eSBarry Smith         }
4702593348eSBarry Smith         fprintf(fd,"\n");
4712593348eSBarry Smith       }
4722593348eSBarry Smith     }
473b6490206SBarry Smith   }
4742593348eSBarry Smith   fflush(fd);
4753a40ed3dSBarry Smith   PetscFunctionReturn(0);
4762593348eSBarry Smith }
4772593348eSBarry Smith 
4785615d1e5SSatish Balay #undef __FUNC__
47977ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw_Zoom"
48077ed5343SBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(Draw draw,void *Aa)
4813270192aSSatish Balay {
48277ed5343SBarry Smith   Mat          A = (Mat) Aa;
4833270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
4847dce120fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
485fae8c767SSatish Balay   double       xl,yl,xr,yr,x_l,x_r,y_l,y_r;
4863f1db9ecSBarry Smith   MatScalar    *aa;
48777ed5343SBarry Smith   MPI_Comm     comm;
48877ed5343SBarry Smith   Viewer       viewer;
4893270192aSSatish Balay 
4903a40ed3dSBarry Smith   PetscFunctionBegin;
49177ed5343SBarry Smith   /*
49277ed5343SBarry Smith       This is nasty. If this is called from an originally parallel matrix
49377ed5343SBarry Smith    then all processes call this, but only the first has the matrix so the
49477ed5343SBarry Smith    rest should return immediately.
49577ed5343SBarry Smith   */
49677ed5343SBarry Smith   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
49777ed5343SBarry Smith   MPI_Comm_rank(comm,&rank);
49877ed5343SBarry Smith   if (rank) PetscFunctionReturn(0);
4993270192aSSatish Balay 
50077ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr);
50177ed5343SBarry Smith 
50277ed5343SBarry Smith   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr);
50377ed5343SBarry Smith 
5043270192aSSatish Balay   /* loop over matrix elements drawing boxes */
5053270192aSSatish Balay   color = DRAW_BLUE;
5063270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5073270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5083270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5093270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5103270192aSSatish Balay       aa = a->a + j*bs2;
5113270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5123270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5133270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
514b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5153270192aSSatish Balay         }
5163270192aSSatish Balay       }
5173270192aSSatish Balay     }
5183270192aSSatish Balay   }
5193270192aSSatish Balay   color = DRAW_CYAN;
5203270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5213270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5223270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5233270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5243270192aSSatish Balay       aa = a->a + j*bs2;
5253270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5263270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5273270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
528b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5293270192aSSatish Balay         }
5303270192aSSatish Balay       }
5313270192aSSatish Balay     }
5323270192aSSatish Balay   }
5333270192aSSatish Balay 
5343270192aSSatish Balay   color = DRAW_RED;
5353270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
5363270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
5373270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
5383270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5393270192aSSatish Balay       aa = a->a + j*bs2;
5403270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
5413270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
5423270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
543b34f160eSSatish Balay           DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);
5443270192aSSatish Balay         }
5453270192aSSatish Balay       }
5463270192aSSatish Balay     }
5473270192aSSatish Balay   }
54877ed5343SBarry Smith   PetscFunctionReturn(0);
54977ed5343SBarry Smith }
5503270192aSSatish Balay 
55177ed5343SBarry Smith #undef __FUNC__
55277ed5343SBarry Smith #define __FUNC__ "MatView_SeqBAIJ_Draw"
55377ed5343SBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
55477ed5343SBarry Smith {
55577ed5343SBarry Smith   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
5567dce120fSSatish Balay   int          ierr;
5577dce120fSSatish Balay   double       xl,yl,xr,yr,w,h;
55877ed5343SBarry Smith   Draw         draw;
55977ed5343SBarry Smith   PetscTruth   isnull;
5603270192aSSatish Balay 
56177ed5343SBarry Smith   PetscFunctionBegin;
56277ed5343SBarry Smith 
56377ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
56477ed5343SBarry Smith   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
56577ed5343SBarry Smith 
56677ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
56777ed5343SBarry Smith   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
56877ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
5693270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
57077ed5343SBarry Smith   ierr = DrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A); CHKERRQ(ierr);
57177ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5723a40ed3dSBarry Smith   PetscFunctionReturn(0);
5733270192aSSatish Balay }
5743270192aSSatish Balay 
5755615d1e5SSatish Balay #undef __FUNC__
576d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqBAIJ"
577e1311b90SBarry Smith int MatView_SeqBAIJ(Mat A,Viewer viewer)
5782593348eSBarry Smith {
57919bcc07fSBarry Smith   ViewerType  vtype;
58019bcc07fSBarry Smith   int         ierr;
5812593348eSBarry Smith 
5823a40ed3dSBarry Smith   PetscFunctionBegin;
58319bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5847b2a1423SBarry Smith   if (PetscTypeCompare(vtype,SOCKET_VIEWER)) {
5857b2a1423SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported");
5863f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){
5873a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
5883f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) {
5893a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
5903f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) {
5913a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
5925cd90555SBarry Smith   } else {
5935cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported by PETSc object");
5942593348eSBarry Smith   }
5953a40ed3dSBarry Smith   PetscFunctionReturn(0);
5962593348eSBarry Smith }
597b6490206SBarry Smith 
598cd0e1443SSatish Balay 
5995615d1e5SSatish Balay #undef __FUNC__
6002d61bbb3SSatish Balay #define __FUNC__ "MatGetValues_SeqBAIJ"
6012d61bbb3SSatish Balay int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
602cd0e1443SSatish Balay {
603cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6042d61bbb3SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
6052d61bbb3SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
6062d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6073f1db9ecSBarry Smith   MatScalar  *ap, *aa = a->a, zero = 0.0;
608cd0e1443SSatish Balay 
6093a40ed3dSBarry Smith   PetscFunctionBegin;
6102d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
611cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
612a8c6a408SBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
613a8c6a408SBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
6142d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6152c3acbe9SBarry Smith     nrow = ailen[brow];
6162d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
617a8c6a408SBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
618a8c6a408SBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
6192d61bbb3SSatish Balay       col  = in[l] ;
6202d61bbb3SSatish Balay       bcol = col/bs;
6212d61bbb3SSatish Balay       cidx = col%bs;
6222d61bbb3SSatish Balay       ridx = row%bs;
6232d61bbb3SSatish Balay       high = nrow;
6242d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6252d61bbb3SSatish Balay       while (high-low > 5) {
626cd0e1443SSatish Balay         t = (low+high)/2;
627cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
628cd0e1443SSatish Balay         else             low  = t;
629cd0e1443SSatish Balay       }
630cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
631cd0e1443SSatish Balay         if (rp[i] > bcol) break;
632cd0e1443SSatish Balay         if (rp[i] == bcol) {
6332d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6342d61bbb3SSatish Balay           goto finished;
635cd0e1443SSatish Balay         }
636cd0e1443SSatish Balay       }
6372d61bbb3SSatish Balay       *v++ = zero;
6382d61bbb3SSatish Balay       finished:;
639cd0e1443SSatish Balay     }
640cd0e1443SSatish Balay   }
6413a40ed3dSBarry Smith   PetscFunctionReturn(0);
642cd0e1443SSatish Balay }
643cd0e1443SSatish Balay 
6442d61bbb3SSatish Balay 
6455615d1e5SSatish Balay #undef __FUNC__
64605a5f48eSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqBAIJ"
64792c4ed94SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
64892c4ed94SBarry Smith {
64992c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6508a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
651dd9472c6SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen, roworiented=a->roworiented;
652dd9472c6SBarry Smith   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
6533f1db9ecSBarry Smith   Scalar      *value = v;
6543f1db9ecSBarry Smith   MatScalar   *ap,*aa=a->a,*bap;
65592c4ed94SBarry Smith 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
6570e324ae4SSatish Balay   if (roworiented) {
6580e324ae4SSatish Balay     stepval = (n-1)*bs;
6590e324ae4SSatish Balay   } else {
6600e324ae4SSatish Balay     stepval = (m-1)*bs;
6610e324ae4SSatish Balay   }
66292c4ed94SBarry Smith   for ( k=0; k<m; k++ ) { /* loop over added rows */
66392c4ed94SBarry Smith     row  = im[k];
6645ef9f2a5SBarry Smith     if (row < 0) continue;
6653a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
666a8c6a408SBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
66792c4ed94SBarry Smith #endif
66892c4ed94SBarry Smith     rp   = aj + ai[row];
66992c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
67092c4ed94SBarry Smith     rmax = imax[row];
67192c4ed94SBarry Smith     nrow = ailen[row];
67292c4ed94SBarry Smith     low  = 0;
67392c4ed94SBarry Smith     for ( l=0; l<n; l++ ) { /* loop over added columns */
6745ef9f2a5SBarry Smith       if (in[l] < 0) continue;
6753a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g)
676a8c6a408SBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
67792c4ed94SBarry Smith #endif
67892c4ed94SBarry Smith       col = in[l];
67992c4ed94SBarry Smith       if (roworiented) {
6800e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
6810e324ae4SSatish Balay       } else {
6820e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
68392c4ed94SBarry Smith       }
68492c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
68592c4ed94SBarry Smith       while (high-low > 7) {
68692c4ed94SBarry Smith         t = (low+high)/2;
68792c4ed94SBarry Smith         if (rp[t] > col) high = t;
68892c4ed94SBarry Smith         else             low  = t;
68992c4ed94SBarry Smith       }
69092c4ed94SBarry Smith       for ( i=low; i<high; i++ ) {
69192c4ed94SBarry Smith         if (rp[i] > col) break;
69292c4ed94SBarry Smith         if (rp[i] == col) {
6938a84c255SSatish Balay           bap  = ap +  bs2*i;
6940e324ae4SSatish Balay           if (roworiented) {
6958a84c255SSatish Balay             if (is == ADD_VALUES) {
696dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
697dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
6988a84c255SSatish Balay                   bap[jj] += *value++;
699dd9472c6SBarry Smith                 }
700dd9472c6SBarry Smith               }
7010e324ae4SSatish Balay             } else {
702dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
703dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs ) {
7040e324ae4SSatish Balay                   bap[jj] = *value++;
7058a84c255SSatish Balay                 }
706dd9472c6SBarry Smith               }
707dd9472c6SBarry Smith             }
7080e324ae4SSatish Balay           } else {
7090e324ae4SSatish Balay             if (is == ADD_VALUES) {
710dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
711dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7120e324ae4SSatish Balay                   *bap++ += *value++;
713dd9472c6SBarry Smith                 }
714dd9472c6SBarry Smith               }
7150e324ae4SSatish Balay             } else {
716dd9472c6SBarry Smith               for ( ii=0; ii<bs; ii++,value+=stepval ) {
717dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++ ) {
7180e324ae4SSatish Balay                   *bap++  = *value++;
7190e324ae4SSatish Balay                 }
7208a84c255SSatish Balay               }
721dd9472c6SBarry Smith             }
722dd9472c6SBarry Smith           }
723f1241b54SBarry Smith           goto noinsert2;
72492c4ed94SBarry Smith         }
72592c4ed94SBarry Smith       }
72689280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
727a8c6a408SBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
72892c4ed94SBarry Smith       if (nrow >= rmax) {
72992c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
73092c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
7313f1db9ecSBarry Smith         MatScalar *new_a;
73292c4ed94SBarry Smith 
733a8c6a408SBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
73489280ab3SLois Curfman McInnes 
73592c4ed94SBarry Smith         /* malloc new storage space */
7363f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
7373f1db9ecSBarry Smith         new_a   = (MatScalar *) PetscMalloc( len ); CHKPTRQ(new_a);
73892c4ed94SBarry Smith         new_j   = (int *) (new_a + bs2*new_nz);
73992c4ed94SBarry Smith         new_i   = new_j + new_nz;
74092c4ed94SBarry Smith 
74192c4ed94SBarry Smith         /* copy over old data into new slots */
74292c4ed94SBarry Smith         for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];}
74392c4ed94SBarry Smith         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
74492c4ed94SBarry Smith         PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));
74592c4ed94SBarry Smith         len = (new_nz - CHUNKSIZE - ai[row] - nrow);
746dd9472c6SBarry Smith         PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));
7473f1db9ecSBarry Smith         PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));
7483f1db9ecSBarry Smith         PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));
7493f1db9ecSBarry Smith         PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));
75092c4ed94SBarry Smith         /* free up old matrix storage */
75192c4ed94SBarry Smith         PetscFree(a->a);
75292c4ed94SBarry Smith         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
75392c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
75492c4ed94SBarry Smith         a->singlemalloc = 1;
75592c4ed94SBarry Smith 
75692c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
75792c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
7583f1db9ecSBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
75992c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
76092c4ed94SBarry Smith         a->reallocs++;
76192c4ed94SBarry Smith         a->nz++;
76292c4ed94SBarry Smith       }
76392c4ed94SBarry Smith       N = nrow++ - 1;
76492c4ed94SBarry Smith       /* shift up all the later entries in this row */
76592c4ed94SBarry Smith       for ( ii=N; ii>=i; ii-- ) {
76692c4ed94SBarry Smith         rp[ii+1] = rp[ii];
7673f1db9ecSBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
76892c4ed94SBarry Smith       }
7693f1db9ecSBarry Smith       if (N >= i) PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
77092c4ed94SBarry Smith       rp[i] = col;
7718a84c255SSatish Balay       bap   = ap +  bs2*i;
7720e324ae4SSatish Balay       if (roworiented) {
773dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval) {
774dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs ) {
7750e324ae4SSatish Balay             bap[jj] = *value++;
776dd9472c6SBarry Smith           }
777dd9472c6SBarry Smith         }
7780e324ae4SSatish Balay       } else {
779dd9472c6SBarry Smith         for ( ii=0; ii<bs; ii++,value+=stepval ) {
780dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++ ) {
7810e324ae4SSatish Balay             *bap++  = *value++;
7820e324ae4SSatish Balay           }
783dd9472c6SBarry Smith         }
784dd9472c6SBarry Smith       }
785f1241b54SBarry Smith       noinsert2:;
78692c4ed94SBarry Smith       low = i;
78792c4ed94SBarry Smith     }
78892c4ed94SBarry Smith     ailen[row] = nrow;
78992c4ed94SBarry Smith   }
7903a40ed3dSBarry Smith   PetscFunctionReturn(0);
79192c4ed94SBarry Smith }
79292c4ed94SBarry Smith 
793f2501298SSatish Balay 
7945615d1e5SSatish Balay #undef __FUNC__
7955615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqBAIJ"
7968f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
797584200bdSSatish Balay {
798584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
799584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
800a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
80143ee02c3SBarry Smith   int        mbs = a->mbs, bs2 = a->bs2,rmax = 0;
8023f1db9ecSBarry Smith   MatScalar  *aa = a->a, *ap;
803584200bdSSatish Balay 
8043a40ed3dSBarry Smith   PetscFunctionBegin;
8053a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
806584200bdSSatish Balay 
80743ee02c3SBarry Smith   if (m) rmax = ailen[0];
808584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
809584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
810584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
811d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
812584200bdSSatish Balay     if (fshift) {
813a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
814584200bdSSatish Balay       N = ailen[i];
815584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
816584200bdSSatish Balay         ip[j-fshift] = ip[j];
8173f1db9ecSBarry Smith         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));
818584200bdSSatish Balay       }
819584200bdSSatish Balay     }
820584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
821584200bdSSatish Balay   }
822584200bdSSatish Balay   if (mbs) {
823584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
824584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
825584200bdSSatish Balay   }
826584200bdSSatish Balay   /* reset ilen and imax for each row */
827584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
828584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
829584200bdSSatish Balay   }
830a7c10996SSatish Balay   a->nz = ai[mbs];
831584200bdSSatish Balay 
832584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
833584200bdSSatish Balay   if (fshift && a->diag) {
834584200bdSSatish Balay     PetscFree(a->diag);
835584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
836584200bdSSatish Balay     a->diag = 0;
837584200bdSSatish Balay   }
8383d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
839219d9a1aSLois Curfman McInnes            m,a->n,a->bs,fshift*bs2,a->nz*bs2);
8403d2013a6SLois Curfman McInnes   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",
841584200bdSSatish Balay            a->reallocs);
842d402145bSBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
843e2f3b5e9SSatish Balay   a->reallocs          = 0;
8444e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
8454e220ebcSLois Curfman McInnes 
8463a40ed3dSBarry Smith   PetscFunctionReturn(0);
847584200bdSSatish Balay }
848584200bdSSatish Balay 
8495a838052SSatish Balay 
850bea157c4SSatish Balay 
851bea157c4SSatish Balay /*
852bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
853bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
854bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
855bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
856bea157c4SSatish Balay */
8575615d1e5SSatish Balay #undef __FUNC__
858bea157c4SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ_Check_Blocks"
859bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
860d9b7c43dSSatish Balay {
861bea157c4SSatish Balay   int i,j,k,row;
862bea157c4SSatish Balay   PetscTruth flg;
8633a40ed3dSBarry Smith 
864bea157c4SSatish Balay   /*   PetscFunctionBegin;*/
865bea157c4SSatish Balay   for ( i=0,j=0; i<n; j++ ) {
866bea157c4SSatish Balay     row = idx[i];
867bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
868bea157c4SSatish Balay       sizes[j] = 1;
869bea157c4SSatish Balay       i++;
870*e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
871bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
872bea157c4SSatish Balay       i++;
873bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
874bea157c4SSatish Balay       flg = PETSC_TRUE;
875bea157c4SSatish Balay       for ( k=1; k<bs; k++ ) {
876bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
877bea157c4SSatish Balay           flg = PETSC_FALSE;
878bea157c4SSatish Balay           break;
879d9b7c43dSSatish Balay         }
880bea157c4SSatish Balay       }
881bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
882bea157c4SSatish Balay         sizes[j] = bs;
883bea157c4SSatish Balay         i+= bs;
884bea157c4SSatish Balay       } else {
885bea157c4SSatish Balay         sizes[j] = 1;
886bea157c4SSatish Balay         i++;
887bea157c4SSatish Balay       }
888bea157c4SSatish Balay     }
889bea157c4SSatish Balay   }
890bea157c4SSatish Balay   *bs_max = j;
8913a40ed3dSBarry Smith   PetscFunctionReturn(0);
892d9b7c43dSSatish Balay }
893d9b7c43dSSatish Balay 
8945615d1e5SSatish Balay #undef __FUNC__
8955615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqBAIJ"
8968f6be9afSLois Curfman McInnes int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
897d9b7c43dSSatish Balay {
898d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
899b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
900bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
9013f1db9ecSBarry Smith   Scalar      zero = 0.0;
9023f1db9ecSBarry Smith   MatScalar   *aa;
903d9b7c43dSSatish Balay 
9043a40ed3dSBarry Smith   PetscFunctionBegin;
905d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
906d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
907d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
908d9b7c43dSSatish Balay 
909bea157c4SSatish Balay   /* allocate memory for rows,sizes */
910bea157c4SSatish Balay   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int)); CHKPTRQ(rows);
911bea157c4SSatish Balay   sizes = rows + is_n;
912bea157c4SSatish Balay 
913bea157c4SSatish Balay   /* initialize copy IS valurs to rows, and sort them */
914bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
915bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows); CHKERRQ(ierr);
916bea157c4SSatish Balay   ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max); CHKERRQ(ierr);
917b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
918bea157c4SSatish Balay 
919bea157c4SSatish Balay   for ( i=0,j=0; i<bs_max; j+=sizes[i],i++ ) {
920bea157c4SSatish Balay     row   = rows[j];
921b6815fffSBarry Smith     if (row < 0 || row > baij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row);
922bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
923bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
924bea157c4SSatish Balay     if (sizes[i] == bs) {
925bea157c4SSatish Balay       if (diag) {
926bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
927bea157c4SSatish Balay           baij->ilen[row/bs] = 1;
928bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
929bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar)); CHKERRQ(ierr);
930a07cd24cSSatish Balay         }
931bea157c4SSatish Balay         /* Now insert all the diagoanl values for this bs */
932bea157c4SSatish Balay         for ( k=0; k<bs; k++ ) {
933bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
934bea157c4SSatish Balay         }
935bea157c4SSatish Balay       } else { /* (!diag) */
936bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
937bea157c4SSatish Balay       } /* end (!diag) */
938bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
939bea157c4SSatish Balay #if defined (USE_PETSC_DEBUG)
940bea157c4SSatish Balay       if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1");
941bea157c4SSatish Balay #endif
942bea157c4SSatish Balay       for ( k=0; k<count; k++ ) {
943d9b7c43dSSatish Balay         aa[0] = zero;
944d9b7c43dSSatish Balay         aa+=bs;
945d9b7c43dSSatish Balay       }
946d9b7c43dSSatish Balay       if (diag) {
947f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
948d9b7c43dSSatish Balay       }
949d9b7c43dSSatish Balay     }
950bea157c4SSatish Balay   }
951bea157c4SSatish Balay 
952bea157c4SSatish Balay   PetscFree(rows);
9539a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
9543a40ed3dSBarry Smith   PetscFunctionReturn(0);
955d9b7c43dSSatish Balay }
9561c351548SSatish Balay 
9575615d1e5SSatish Balay #undef __FUNC__
9582d61bbb3SSatish Balay #define __FUNC__ "MatSetValues_SeqBAIJ"
9592d61bbb3SSatish Balay int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
9602d61bbb3SSatish Balay {
9612d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9622d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
9632d61bbb3SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
9642d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
9652d61bbb3SSatish Balay   int         ridx,cidx,bs2=a->bs2;
9663f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
9672d61bbb3SSatish Balay 
9682d61bbb3SSatish Balay   PetscFunctionBegin;
9692d61bbb3SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
9702d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
9715ef9f2a5SBarry Smith     if (row < 0) continue;
9722d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9732d61bbb3SSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
9742d61bbb3SSatish Balay #endif
9752d61bbb3SSatish Balay     rp   = aj + ai[brow];
9762d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
9772d61bbb3SSatish Balay     rmax = imax[brow];
9782d61bbb3SSatish Balay     nrow = ailen[brow];
9792d61bbb3SSatish Balay     low  = 0;
9802d61bbb3SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
9815ef9f2a5SBarry Smith       if (in[l] < 0) continue;
9822d61bbb3SSatish Balay #if defined(USE_PETSC_BOPT_g)
9832d61bbb3SSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
9842d61bbb3SSatish Balay #endif
9852d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
9862d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
9872d61bbb3SSatish Balay       if (roworiented) {
9885ef9f2a5SBarry Smith         value = v[l + k*n];
9892d61bbb3SSatish Balay       } else {
9902d61bbb3SSatish Balay         value = v[k + l*m];
9912d61bbb3SSatish Balay       }
9922d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
9932d61bbb3SSatish Balay       while (high-low > 7) {
9942d61bbb3SSatish Balay         t = (low+high)/2;
9952d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
9962d61bbb3SSatish Balay         else              low  = t;
9972d61bbb3SSatish Balay       }
9982d61bbb3SSatish Balay       for ( i=low; i<high; i++ ) {
9992d61bbb3SSatish Balay         if (rp[i] > bcol) break;
10002d61bbb3SSatish Balay         if (rp[i] == bcol) {
10012d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
10022d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
10032d61bbb3SSatish Balay           else                  *bap  = value;
10042d61bbb3SSatish Balay           goto noinsert1;
10052d61bbb3SSatish Balay         }
10062d61bbb3SSatish Balay       }
10072d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
10082d61bbb3SSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
10092d61bbb3SSatish Balay       if (nrow >= rmax) {
10102d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
10112d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
10123f1db9ecSBarry Smith         MatScalar *new_a;
10132d61bbb3SSatish Balay 
10142d61bbb3SSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
10152d61bbb3SSatish Balay 
10162d61bbb3SSatish Balay         /* Malloc new storage space */
10173f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
10183f1db9ecSBarry Smith         new_a   = (MatScalar *) PetscMalloc( len ); CHKPTRQ(new_a);
10192d61bbb3SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
10202d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10212d61bbb3SSatish Balay 
10222d61bbb3SSatish Balay         /* copy over old data into new slots */
10232d61bbb3SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
10242d61bbb3SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
10252d61bbb3SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
10262d61bbb3SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
10272d61bbb3SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
10282d61bbb3SSatish Balay                                                            len*sizeof(int));
10293f1db9ecSBarry Smith         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));
10303f1db9ecSBarry Smith         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));
10312d61bbb3SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
10323f1db9ecSBarry Smith                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));
10332d61bbb3SSatish Balay         /* free up old matrix storage */
10342d61bbb3SSatish Balay         PetscFree(a->a);
10352d61bbb3SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
10362d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
10372d61bbb3SSatish Balay         a->singlemalloc = 1;
10382d61bbb3SSatish Balay 
10392d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
10402d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
10413f1db9ecSBarry Smith         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
10422d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
10432d61bbb3SSatish Balay         a->reallocs++;
10442d61bbb3SSatish Balay         a->nz++;
10452d61bbb3SSatish Balay       }
10462d61bbb3SSatish Balay       N = nrow++ - 1;
10472d61bbb3SSatish Balay       /* shift up all the later entries in this row */
10482d61bbb3SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
10492d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
10503f1db9ecSBarry Smith         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));
10512d61bbb3SSatish Balay       }
10523f1db9ecSBarry Smith       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));
10532d61bbb3SSatish Balay       rp[i]                      = bcol;
10542d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
10552d61bbb3SSatish Balay       noinsert1:;
10562d61bbb3SSatish Balay       low = i;
10572d61bbb3SSatish Balay     }
10582d61bbb3SSatish Balay     ailen[brow] = nrow;
10592d61bbb3SSatish Balay   }
10602d61bbb3SSatish Balay   PetscFunctionReturn(0);
10612d61bbb3SSatish Balay }
10622d61bbb3SSatish Balay 
10632d61bbb3SSatish Balay extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
10642d61bbb3SSatish Balay extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
10652d61bbb3SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
10667b2a1423SBarry Smith extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
10677b2a1423SBarry Smith extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
10682d61bbb3SSatish Balay extern int MatMultTrans_SeqBAIJ(Mat,Vec,Vec);
10692d61bbb3SSatish Balay extern int MatMultTransAdd_SeqBAIJ(Mat,Vec,Vec,Vec);
10702d61bbb3SSatish Balay extern int MatScale_SeqBAIJ(Scalar*,Mat);
10712d61bbb3SSatish Balay extern int MatNorm_SeqBAIJ(Mat,NormType,double *);
10722d61bbb3SSatish Balay extern int MatEqual_SeqBAIJ(Mat,Mat, PetscTruth*);
10732d61bbb3SSatish Balay extern int MatGetDiagonal_SeqBAIJ(Mat,Vec);
10742d61bbb3SSatish Balay extern int MatDiagonalScale_SeqBAIJ(Mat,Vec,Vec);
10752d61bbb3SSatish Balay extern int MatGetInfo_SeqBAIJ(Mat,MatInfoType,MatInfo *);
10762d61bbb3SSatish Balay extern int MatZeroEntries_SeqBAIJ(Mat);
10772d61bbb3SSatish Balay 
10782d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
10792d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
10802d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
10812d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
10822d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
10832d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
10842d61bbb3SSatish Balay extern int MatSolve_SeqBAIJ_7(Mat,Vec,Vec);
10852d61bbb3SSatish Balay 
10862d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
10872d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
10882d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
10892d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
10902d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
10912d61bbb3SSatish Balay extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
10922d61bbb3SSatish Balay 
10932d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_1(Mat,Vec,Vec);
10942d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_2(Mat,Vec,Vec);
10952d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_3(Mat,Vec,Vec);
10962d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_4(Mat,Vec,Vec);
10972d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_5(Mat,Vec,Vec);
10982d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_7(Mat,Vec,Vec);
10992d61bbb3SSatish Balay extern int MatMult_SeqBAIJ_N(Mat,Vec,Vec);
11002d61bbb3SSatish Balay 
11012d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_1(Mat,Vec,Vec,Vec);
11022d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_2(Mat,Vec,Vec,Vec);
11032d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_3(Mat,Vec,Vec,Vec);
11042d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_4(Mat,Vec,Vec,Vec);
11052d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_5(Mat,Vec,Vec,Vec);
11062d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_7(Mat,Vec,Vec,Vec);
11072d61bbb3SSatish Balay extern int MatMultAdd_SeqBAIJ_N(Mat,Vec,Vec,Vec);
11082d61bbb3SSatish Balay 
11092d61bbb3SSatish Balay #undef __FUNC__
11102d61bbb3SSatish Balay #define __FUNC__ "MatILUFactor_SeqBAIJ"
11115ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
11122d61bbb3SSatish Balay {
11132d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
11142d61bbb3SSatish Balay   Mat         outA;
11152d61bbb3SSatish Balay   int         ierr;
1116667159a5SBarry Smith   PetscTruth  row_identity, col_identity;
11172d61bbb3SSatish Balay 
11182d61bbb3SSatish Balay   PetscFunctionBegin;
11195ef9f2a5SBarry Smith   if (info && info->fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported");
1120667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity); CHKERRQ(ierr);
1121667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity); CHKERRQ(ierr);
1122667159a5SBarry Smith   if (!row_identity || !col_identity) {
1123667159a5SBarry Smith     SETERRQ(1,1,"Row and column permutations must be identity");
1124667159a5SBarry Smith   }
11252d61bbb3SSatish Balay 
11262d61bbb3SSatish Balay   outA          = inA;
11272d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11282d61bbb3SSatish Balay   a->row        = row;
11292d61bbb3SSatish Balay   a->col        = col;
11302d61bbb3SSatish Balay 
1131e51c0b9cSSatish Balay   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
1132e51c0b9cSSatish Balay   ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr);
11331e4dd532SSatish Balay   PLogObjectParent(inA,a->icol);
1134e51c0b9cSSatish Balay 
11352d61bbb3SSatish Balay   if (!a->solve_work) {
11362d61bbb3SSatish Balay     a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
11372d61bbb3SSatish Balay     PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
11382d61bbb3SSatish Balay   }
11392d61bbb3SSatish Balay 
11402d61bbb3SSatish Balay   if (!a->diag) {
11412d61bbb3SSatish Balay     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
11422d61bbb3SSatish Balay   }
11432d61bbb3SSatish Balay   /*
1144667159a5SBarry Smith       Blocksize 4 and 5 have a special faster factorization/solver for ILU(0) factorization
11452d61bbb3SSatish Balay     with natural ordering
11462d61bbb3SSatish Balay   */
11472d61bbb3SSatish Balay   if (a->bs == 4) {
1148667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
1149f830108cSBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
1150667159a5SBarry Smith     PLogInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");
1151667159a5SBarry Smith   } else if (a->bs == 5) {
1152667159a5SBarry Smith     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
1153667159a5SBarry Smith     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
1154667159a5SBarry Smith     PLogInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");
11552d61bbb3SSatish Balay   }
11562d61bbb3SSatish Balay 
1157667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
1158667159a5SBarry Smith 
1159667159a5SBarry Smith 
11602d61bbb3SSatish Balay   PetscFunctionReturn(0);
11612d61bbb3SSatish Balay }
11622d61bbb3SSatish Balay #undef __FUNC__
1163d4bb536fSBarry Smith #define __FUNC__ "MatPrintHelp_SeqBAIJ"
1164ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1165ba4ca20aSSatish Balay {
1166ba4ca20aSSatish Balay   static int called = 0;
1167ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1168ba4ca20aSSatish Balay 
11693a40ed3dSBarry Smith   PetscFunctionBegin;
11703a40ed3dSBarry Smith   if (called) {PetscFunctionReturn(0);} else called = 1;
117176be9ce4SBarry Smith   (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
117276be9ce4SBarry Smith   (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");
11733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1174ba4ca20aSSatish Balay }
1175d9b7c43dSSatish Balay 
1176fb2e594dSBarry Smith EXTERN_C_BEGIN
117727a8da17SBarry Smith #undef __FUNC__
117827a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
117927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
118027a8da17SBarry Smith {
118127a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
118227a8da17SBarry Smith   int         i,nz,n;
118327a8da17SBarry Smith 
118427a8da17SBarry Smith   PetscFunctionBegin;
118527a8da17SBarry Smith   nz = baij->maxnz;
118627a8da17SBarry Smith   n  = baij->n;
118727a8da17SBarry Smith   for (i=0; i<nz; i++) {
118827a8da17SBarry Smith     baij->j[i] = indices[i];
118927a8da17SBarry Smith   }
119027a8da17SBarry Smith   baij->nz = nz;
119127a8da17SBarry Smith   for ( i=0; i<n; i++ ) {
119227a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
119327a8da17SBarry Smith   }
119427a8da17SBarry Smith 
119527a8da17SBarry Smith   PetscFunctionReturn(0);
119627a8da17SBarry Smith }
1197fb2e594dSBarry Smith EXTERN_C_END
119827a8da17SBarry Smith 
119927a8da17SBarry Smith #undef __FUNC__
120027a8da17SBarry Smith #define __FUNC__ "MatSeqBAIJSetColumnIndices"
120127a8da17SBarry Smith /*@
120227a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
120327a8da17SBarry Smith        in the matrix.
120427a8da17SBarry Smith 
120527a8da17SBarry Smith   Input Parameters:
120627a8da17SBarry Smith +  mat - the SeqBAIJ matrix
120727a8da17SBarry Smith -  indices - the column indices
120827a8da17SBarry Smith 
120927a8da17SBarry Smith   Notes:
121027a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
121127a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
121227a8da17SBarry Smith   of the MatSetValues() operation.
121327a8da17SBarry Smith 
121427a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
121527a8da17SBarry Smith   MatCreateSeqBAIJ().
121627a8da17SBarry Smith 
121727a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
121827a8da17SBarry Smith 
121927a8da17SBarry Smith @*/
122027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
122127a8da17SBarry Smith {
122227a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
122327a8da17SBarry Smith 
122427a8da17SBarry Smith   PetscFunctionBegin;
122527a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
122627a8da17SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void **)&f);
122727a8da17SBarry Smith          CHKERRQ(ierr);
122827a8da17SBarry Smith   if (f) {
122927a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
123027a8da17SBarry Smith   } else {
123127a8da17SBarry Smith     SETERRQ(1,1,"Wrong type of matrix to set column indices");
123227a8da17SBarry Smith   }
123327a8da17SBarry Smith   PetscFunctionReturn(0);
123427a8da17SBarry Smith }
123527a8da17SBarry Smith 
12362593348eSBarry Smith /* -------------------------------------------------------------------*/
1237cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1238cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1239cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1240cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1241cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
1242cc2dc46cSBarry Smith        MatMultTrans_SeqBAIJ,
1243cc2dc46cSBarry Smith        MatMultTransAdd_SeqBAIJ,
1244cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1245cc2dc46cSBarry Smith        0,
1246cc2dc46cSBarry Smith        0,
1247cc2dc46cSBarry Smith        0,
1248cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1249cc2dc46cSBarry Smith        0,
1250b6490206SBarry Smith        0,
1251f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1252cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1253cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1254cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1255cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1256cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1257b6490206SBarry Smith        0,
1258cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1259cc2dc46cSBarry Smith        0,
1260cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1261cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1262cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1263cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1264cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1265cc2dc46cSBarry Smith        0,
1266cc2dc46cSBarry Smith        0,
1267cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1268cc2dc46cSBarry Smith        MatGetSize_SeqBAIJ,
1269cc2dc46cSBarry Smith        MatGetOwnershipRange_SeqBAIJ,
1270cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1271cc2dc46cSBarry Smith        0,
1272cc2dc46cSBarry Smith        0,
1273cc2dc46cSBarry Smith        0,
12742e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1275cc2dc46cSBarry Smith        0,
1276cc2dc46cSBarry Smith        0,
1277cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1278cc2dc46cSBarry Smith        0,
1279cc2dc46cSBarry Smith        0,
1280cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1281cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1282cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1283cc2dc46cSBarry Smith        0,
1284cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1285cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1286cc2dc46cSBarry Smith        0,
1287cc2dc46cSBarry Smith        0,
1288cc2dc46cSBarry Smith        0,
1289cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
12903b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
129192c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1292cc2dc46cSBarry Smith        0,
1293cc2dc46cSBarry Smith        0,
1294cc2dc46cSBarry Smith        0,
1295cc2dc46cSBarry Smith        0,
1296cc2dc46cSBarry Smith        0,
1297cc2dc46cSBarry Smith        0,
1298d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1299cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1300cc2dc46cSBarry Smith        0,
1301cc2dc46cSBarry Smith        0,
1302cc2dc46cSBarry Smith        MatGetMaps_Petsc};
13032593348eSBarry Smith 
13045615d1e5SSatish Balay #undef __FUNC__
13055615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqBAIJ"
13062593348eSBarry Smith /*@C
130744cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
130844cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
130944cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
13102bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
13112bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
13122593348eSBarry Smith 
1313db81eaa0SLois Curfman McInnes    Collective on MPI_Comm
1314db81eaa0SLois Curfman McInnes 
13152593348eSBarry Smith    Input Parameters:
1316db81eaa0SLois Curfman McInnes +  comm - MPI communicator, set to PETSC_COMM_SELF
1317b6490206SBarry Smith .  bs - size of block
13182593348eSBarry Smith .  m - number of rows
13192593348eSBarry Smith .  n - number of columns
1320b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1321db81eaa0SLois Curfman McInnes -  nzz - array containing the number of block nonzeros in the various block rows
13222bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
13232593348eSBarry Smith 
13242593348eSBarry Smith    Output Parameter:
13252593348eSBarry Smith .  A - the matrix
13262593348eSBarry Smith 
13270513a670SBarry Smith    Options Database Keys:
1328db81eaa0SLois Curfman McInnes .   -mat_no_unroll - uses code that does not unroll the loops in the
1329db81eaa0SLois Curfman McInnes                      block calculations (much slower)
1330db81eaa0SLois Curfman McInnes .    -mat_block_size - size of the blocks to use
13310513a670SBarry Smith 
13322593348eSBarry Smith    Notes:
133344cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
13342593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
133544cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
13362593348eSBarry Smith 
13372593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
13382593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
13392593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
13406da5968aSLois Curfman McInnes    matrices.
13412593348eSBarry Smith 
1342027ccd11SLois Curfman McInnes    Level: intermediate
1343027ccd11SLois Curfman McInnes 
1344db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
13452593348eSBarry Smith @*/
1346b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
13472593348eSBarry Smith {
13482593348eSBarry Smith   Mat         B;
1349b6490206SBarry Smith   Mat_SeqBAIJ *b;
13503b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
13513b2fbd54SBarry Smith 
13523a40ed3dSBarry Smith   PetscFunctionBegin;
13533b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
1354a8c6a408SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
1355b6490206SBarry Smith 
13560513a670SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
13570513a670SBarry Smith 
13583a40ed3dSBarry Smith   if (mbs*bs!=m || nbs*bs!=n) {
1359a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
13603a40ed3dSBarry Smith   }
13612593348eSBarry Smith 
13622593348eSBarry Smith   *A = 0;
13633f1db9ecSBarry Smith   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",comm,MatDestroy,MatView);
13642593348eSBarry Smith   PLogObjectCreate(B);
1365b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
136644cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
1367cc2dc46cSBarry Smith   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));
13681a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
13691a6a6d98SBarry Smith   if (!flg) {
13707fc0212eSBarry Smith     switch (bs) {
13717fc0212eSBarry Smith     case 1:
1372f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1373f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_1;
1374f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1375f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
13767fc0212eSBarry Smith       break;
13774eeb42bcSBarry Smith     case 2:
1378f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1379f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_2;
1380f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1381f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
13826d84be18SBarry Smith       break;
1383f327dd97SBarry Smith     case 3:
1384f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1385f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_3;
1386f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1387f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
13884eeb42bcSBarry Smith       break;
13896d84be18SBarry Smith     case 4:
1390f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1391f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_4;
1392f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1393f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
13946d84be18SBarry Smith       break;
13956d84be18SBarry Smith     case 5:
1396f830108cSBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1397f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_5;
1398f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1399f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
14006d84be18SBarry Smith       break;
140148e9ddb2SSatish Balay     case 7:
1402f830108cSBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1403f830108cSBarry Smith       B->ops->solve           = MatSolve_SeqBAIJ_7;
1404f830108cSBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
140548e9ddb2SSatish Balay       break;
14067fc0212eSBarry Smith     }
14071a6a6d98SBarry Smith   }
1408e1311b90SBarry Smith   B->ops->destroy     = MatDestroy_SeqBAIJ;
1409e1311b90SBarry Smith   B->ops->view        = MatView_SeqBAIJ;
14102593348eSBarry Smith   B->factor           = 0;
14112593348eSBarry Smith   B->lupivotthreshold = 1.0;
141290f02eecSBarry Smith   B->mapping          = 0;
14132593348eSBarry Smith   b->row              = 0;
14142593348eSBarry Smith   b->col              = 0;
1415e51c0b9cSSatish Balay   b->icol             = 0;
14162593348eSBarry Smith   b->reallocs         = 0;
14172593348eSBarry Smith 
141844cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
141944cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1420a5ae1ecdSBarry Smith 
14217413bad6SBarry Smith   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
14227413bad6SBarry Smith   ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr);
1423a5ae1ecdSBarry Smith 
1424b6490206SBarry Smith   b->mbs     = mbs;
1425f2501298SSatish Balay   b->nbs     = nbs;
1426b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
14272593348eSBarry Smith   if (nnz == PETSC_NULL) {
1428119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
14292593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1430b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1431b6490206SBarry Smith     nz = nz*mbs;
14323a40ed3dSBarry Smith   } else {
14332593348eSBarry Smith     nz = 0;
1434b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
14352593348eSBarry Smith   }
14362593348eSBarry Smith 
14372593348eSBarry Smith   /* allocate the matrix space */
14383f1db9ecSBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
14393f1db9ecSBarry Smith   b->a  = (MatScalar *) PetscMalloc( len ); CHKPTRQ(b->a);
14403f1db9ecSBarry Smith   PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));
14417e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
14422593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
14432593348eSBarry Smith   b->i  = b->j + nz;
14442593348eSBarry Smith   b->singlemalloc = 1;
14452593348eSBarry Smith 
1446de6a44a3SBarry Smith   b->i[0] = 0;
1447b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
14482593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
14492593348eSBarry Smith   }
14502593348eSBarry Smith 
1451b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1452b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1453f09e8eb9SSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1454b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
14552593348eSBarry Smith 
1456b6490206SBarry Smith   b->bs               = bs;
14579df24120SSatish Balay   b->bs2              = bs2;
1458b6490206SBarry Smith   b->mbs              = mbs;
14592593348eSBarry Smith   b->nz               = 0;
146018e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
14612593348eSBarry Smith   b->sorted           = 0;
14622593348eSBarry Smith   b->roworiented      = 1;
14632593348eSBarry Smith   b->nonew            = 0;
14642593348eSBarry Smith   b->diag             = 0;
14652593348eSBarry Smith   b->solve_work       = 0;
1466de6a44a3SBarry Smith   b->mult_work        = 0;
14672593348eSBarry Smith   b->spptr            = 0;
14684e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
14694e220ebcSLois Curfman McInnes 
14702593348eSBarry Smith   *A = B;
14712593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
14722593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
147327a8da17SBarry Smith 
147427a8da17SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
147527a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
147627a8da17SBarry Smith                                      (void*)MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
14773a40ed3dSBarry Smith   PetscFunctionReturn(0);
14782593348eSBarry Smith }
14792593348eSBarry Smith 
14805615d1e5SSatish Balay #undef __FUNC__
14812e8a6d31SBarry Smith #define __FUNC__ "MatDuplicate_SeqBAIJ"
14822e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
14832593348eSBarry Smith {
14842593348eSBarry Smith   Mat         C;
1485b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
148627a8da17SBarry Smith   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1487de6a44a3SBarry Smith 
14883a40ed3dSBarry Smith   PetscFunctionBegin;
1489a8c6a408SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
14902593348eSBarry Smith 
14912593348eSBarry Smith   *B = 0;
14923f1db9ecSBarry Smith   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQBAIJ,"Mat",A->comm,MatDestroy,MatView);
14932593348eSBarry Smith   PLogObjectCreate(C);
1494b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
1495c3c66ee5SSatish Balay   PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));
1496e1311b90SBarry Smith   C->ops->destroy    = MatDestroy_SeqBAIJ;
1497e1311b90SBarry Smith   C->ops->view       = MatView_SeqBAIJ;
14982593348eSBarry Smith   C->factor          = A->factor;
14992593348eSBarry Smith   c->row             = 0;
15002593348eSBarry Smith   c->col             = 0;
15012593348eSBarry Smith   C->assembled       = PETSC_TRUE;
15022593348eSBarry Smith 
150344cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
150444cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
150544cd7ae7SLois Curfman McInnes   C->M          = a->m;
150644cd7ae7SLois Curfman McInnes   C->N          = a->n;
150744cd7ae7SLois Curfman McInnes 
1508b6490206SBarry Smith   c->bs         = a->bs;
15099df24120SSatish Balay   c->bs2        = a->bs2;
1510b6490206SBarry Smith   c->mbs        = a->mbs;
1511b6490206SBarry Smith   c->nbs        = a->nbs;
15122593348eSBarry Smith 
1513b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1514b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1515b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
15162593348eSBarry Smith     c->imax[i] = a->imax[i];
15172593348eSBarry Smith     c->ilen[i] = a->ilen[i];
15182593348eSBarry Smith   }
15192593348eSBarry Smith 
15202593348eSBarry Smith   /* allocate the matrix space */
15212593348eSBarry Smith   c->singlemalloc = 1;
15223f1db9ecSBarry Smith   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
15233f1db9ecSBarry Smith   c->a  = (MatScalar *) PetscMalloc( len ); CHKPTRQ(c->a);
15247e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1525de6a44a3SBarry Smith   c->i  = c->j + nz;
1526b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1527b6490206SBarry Smith   if (mbs > 0) {
1528de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
15292e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
15303f1db9ecSBarry Smith       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));
15312e8a6d31SBarry Smith     } else {
15323f1db9ecSBarry Smith       PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));
15332593348eSBarry Smith     }
15342593348eSBarry Smith   }
15352593348eSBarry Smith 
1536f09e8eb9SSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
15372593348eSBarry Smith   c->sorted      = a->sorted;
15382593348eSBarry Smith   c->roworiented = a->roworiented;
15392593348eSBarry Smith   c->nonew       = a->nonew;
15402593348eSBarry Smith 
15412593348eSBarry Smith   if (a->diag) {
1542b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1543b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1544b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
15452593348eSBarry Smith       c->diag[i] = a->diag[i];
15462593348eSBarry Smith     }
15472593348eSBarry Smith   }
15482593348eSBarry Smith   else c->diag          = 0;
15492593348eSBarry Smith   c->nz                 = a->nz;
15502593348eSBarry Smith   c->maxnz              = a->maxnz;
15512593348eSBarry Smith   c->solve_work         = 0;
15522593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
15537fc0212eSBarry Smith   c->mult_work          = 0;
15542593348eSBarry Smith   *B = C;
15557b2a1423SBarry Smith   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
15563a40ed3dSBarry Smith   PetscFunctionReturn(0);
15572593348eSBarry Smith }
15582593348eSBarry Smith 
15595615d1e5SSatish Balay #undef __FUNC__
15605615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqBAIJ"
156119bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
15622593348eSBarry Smith {
1563b6490206SBarry Smith   Mat_SeqBAIJ  *a;
15642593348eSBarry Smith   Mat          B;
1565de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1566b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
156735aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1568a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1569b6490206SBarry Smith   Scalar       *aa;
157019bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
15712593348eSBarry Smith 
15723a40ed3dSBarry Smith   PetscFunctionBegin;
15731a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1574de6a44a3SBarry Smith   bs2  = bs*bs;
1575b6490206SBarry Smith 
15762593348eSBarry Smith   MPI_Comm_size(comm,&size);
1577cc0fae11SBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
157890ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
15790752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr);
1580a8c6a408SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
15812593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
15822593348eSBarry Smith 
1583d64ed03dSBarry Smith   if (header[3] < 0) {
1584a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqBAIJ");
1585d64ed03dSBarry Smith   }
1586d64ed03dSBarry Smith 
1587a8c6a408SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
158835aab85fSBarry Smith 
158935aab85fSBarry Smith   /*
159035aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
159135aab85fSBarry Smith     divisible by the blocksize
159235aab85fSBarry Smith   */
1593b6490206SBarry Smith   mbs        = M/bs;
159435aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
159535aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
159635aab85fSBarry Smith   else                  mbs++;
159735aab85fSBarry Smith   if (extra_rows) {
1598537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
159935aab85fSBarry Smith   }
1600b6490206SBarry Smith 
16012593348eSBarry Smith   /* read in row lengths */
160235aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
16030752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr);
160435aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
16052593348eSBarry Smith 
1606b6490206SBarry Smith   /* read in column indices */
160735aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
16080752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT); CHKERRQ(ierr);
160935aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1610b6490206SBarry Smith 
1611b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1612b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1613b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
161435aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
161535aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
161635aab85fSBarry Smith   masked = mask + mbs;
1617b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1618b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
161935aab85fSBarry Smith     nmask = 0;
1620b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1621b6490206SBarry Smith       kmax = rowlengths[rowcount];
1622b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
162335aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
162435aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1625b6490206SBarry Smith       }
1626b6490206SBarry Smith       rowcount++;
1627b6490206SBarry Smith     }
162835aab85fSBarry Smith     browlengths[i] += nmask;
162935aab85fSBarry Smith     /* zero out the mask elements we set */
163035aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1631b6490206SBarry Smith   }
1632b6490206SBarry Smith 
16332593348eSBarry Smith   /* create our matrix */
16343eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
16352593348eSBarry Smith   B = *A;
1636b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
16372593348eSBarry Smith 
16382593348eSBarry Smith   /* set matrix "i" values */
1639de6a44a3SBarry Smith   a->i[0] = 0;
1640b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1641b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1642b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
16432593348eSBarry Smith   }
1644b6490206SBarry Smith   a->nz         = 0;
1645b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
16462593348eSBarry Smith 
1647b6490206SBarry Smith   /* read in nonzero values */
164835aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
16490752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR); CHKERRQ(ierr);
165035aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1651b6490206SBarry Smith 
1652b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1653b6490206SBarry Smith   nzcount = 0; jcount = 0;
1654b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1655b6490206SBarry Smith     nzcountb = nzcount;
165635aab85fSBarry Smith     nmask    = 0;
1657b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1658b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1659b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
166035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
166135aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1662b6490206SBarry Smith       }
1663b6490206SBarry Smith       rowcount++;
1664b6490206SBarry Smith     }
1665de6a44a3SBarry Smith     /* sort the masked values */
166677c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1667de6a44a3SBarry Smith 
1668b6490206SBarry Smith     /* set "j" values into matrix */
1669b6490206SBarry Smith     maskcount = 1;
167035aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
167135aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1672de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1673b6490206SBarry Smith     }
1674b6490206SBarry Smith     /* set "a" values into matrix */
1675de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1676b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1677b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1678b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1679de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1680de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1681de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1682de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1683b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1684b6490206SBarry Smith       }
1685b6490206SBarry Smith     }
168635aab85fSBarry Smith     /* zero out the mask elements we set */
168735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1688b6490206SBarry Smith   }
1689a8c6a408SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
1690b6490206SBarry Smith 
1691b6490206SBarry Smith   PetscFree(rowlengths);
1692b6490206SBarry Smith   PetscFree(browlengths);
1693b6490206SBarry Smith   PetscFree(aa);
1694b6490206SBarry Smith   PetscFree(jj);
1695b6490206SBarry Smith   PetscFree(mask);
1696b6490206SBarry Smith 
1697b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1698de6a44a3SBarry Smith 
16999c01be13SBarry Smith   ierr = MatView_Private(B); CHKERRQ(ierr);
17003a40ed3dSBarry Smith   PetscFunctionReturn(0);
17012593348eSBarry Smith }
17022593348eSBarry Smith 
17032593348eSBarry Smith 
17042593348eSBarry Smith 
1705