xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 0cd7f59a393d874074aadc4ddaa1c8ce9a6fd403)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay /*
3a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
449b5e25fSSatish Balay   matrix storage format.
549b5e25fSSatish Balay */
6c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>         /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
8c6db04a5SJed Brown #include <petscblaslapack.h>
949b5e25fSSatish Balay 
10c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/relax.h>
1170dcbbb9SBarry Smith #define USESHORT
12c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/relax.h>
1370dcbbb9SBarry Smith 
146214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
15cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
166214f412SHong Zhang #endif
17b5b17502SBarry Smith 
1849b5e25fSSatish Balay /*
1949b5e25fSSatish Balay      Checks for missing diagonals
2049b5e25fSSatish Balay */
21ace3abfcSBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscBool  *missing,PetscInt *dd)
2249b5e25fSSatish Balay {
23045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
246849ba73SBarry Smith   PetscErrorCode ierr;
257734d3b5SMatthew G. Knepley   PetscInt       *diag,*ii = a->i,i;
2649b5e25fSSatish Balay 
2749b5e25fSSatish Balay   PetscFunctionBegin;
28045c9aa0SHong Zhang   ierr     = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
292af78befSBarry Smith   *missing = PETSC_FALSE;
307734d3b5SMatthew G. Knepley   if (A->rmap->n > 0 && !ii) {
31358d2f5dSShri Abhyankar     *missing = PETSC_TRUE;
32358d2f5dSShri Abhyankar     if (dd) *dd = 0;
33955c1f14SBarry Smith     ierr = PetscInfo(A,"Matrix has no entries therefore is missing diagonal\n");CHKERRQ(ierr);
34358d2f5dSShri Abhyankar   } else {
35358d2f5dSShri Abhyankar     diag = a->diag;
3649b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
377734d3b5SMatthew G. Knepley       if (diag[i] >= ii[i+1]) {
382af78befSBarry Smith         *missing = PETSC_TRUE;
392af78befSBarry Smith         if (dd) *dd = i;
402af78befSBarry Smith         break;
412af78befSBarry Smith       }
4249b5e25fSSatish Balay     }
43358d2f5dSShri Abhyankar   }
4449b5e25fSSatish Balay   PetscFunctionReturn(0);
4549b5e25fSSatish Balay }
4649b5e25fSSatish Balay 
47dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
4849b5e25fSSatish Balay {
49045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
506849ba73SBarry Smith   PetscErrorCode ierr;
5148dd3d27SHong Zhang   PetscInt       i,j;
5249b5e25fSSatish Balay 
5349b5e25fSSatish Balay   PetscFunctionBegin;
5409f38230SBarry Smith   if (!a->diag) {
55785e854fSJed Brown     ierr         = PetscMalloc1(a->mbs,&a->diag);CHKERRQ(ierr);
563bb1ff40SBarry Smith     ierr         = PetscLogObjectMemory((PetscObject)A,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
57c760cd28SBarry Smith     a->free_diag = PETSC_TRUE;
5809f38230SBarry Smith   }
5948dd3d27SHong Zhang   for (i=0; i<a->mbs; i++) {
6048dd3d27SHong Zhang     a->diag[i] = a->i[i+1];
6148dd3d27SHong Zhang     for (j=a->i[i]; j<a->i[i+1]; j++) {
6248dd3d27SHong Zhang       if (a->j[j] == i) {
6348dd3d27SHong Zhang         a->diag[i] = j;
6448dd3d27SHong Zhang         break;
6548dd3d27SHong Zhang       }
6648dd3d27SHong Zhang     }
6748dd3d27SHong Zhang   }
6849b5e25fSSatish Balay   PetscFunctionReturn(0);
6949b5e25fSSatish Balay }
7049b5e25fSSatish Balay 
711a83f524SJed Brown static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
7249b5e25fSSatish Balay {
73a6ece127SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
748f7157efSSatish Balay   PetscErrorCode ierr;
752462f5fdSStefano Zampini   PetscInt       i,j,n = a->mbs,nz = a->i[n],*tia,*tja,bs = A->rmap->bs,k,l,cnt;
762462f5fdSStefano Zampini   PetscInt       **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
7749b5e25fSSatish Balay 
7849b5e25fSSatish Balay   PetscFunctionBegin;
79d3e5a4abSHong Zhang   *nn = n;
80a1373b80SHong Zhang   if (!ia) PetscFunctionReturn(0);
812462f5fdSStefano Zampini   if (symmetric) {
822462f5fdSStefano Zampini     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,PETSC_FALSE,0,0,&tia,&tja);CHKERRQ(ierr);
832462f5fdSStefano Zampini     nz   = tia[n];
842462f5fdSStefano Zampini   } else {
852462f5fdSStefano Zampini     tia = a->i; tja = a->j;
862462f5fdSStefano Zampini   }
872462f5fdSStefano Zampini 
882462f5fdSStefano Zampini   if (!blockcompressed && bs > 1) {
892462f5fdSStefano Zampini     (*nn) *= bs;
908f7157efSSatish Balay     /* malloc & create the natural set of indices */
912462f5fdSStefano Zampini     ierr = PetscMalloc1((n+1)*bs,ia);CHKERRQ(ierr);
922462f5fdSStefano Zampini     if (n) {
932462f5fdSStefano Zampini       (*ia)[0] = oshift;
942462f5fdSStefano Zampini       for (j=1; j<bs; j++) {
952462f5fdSStefano Zampini         (*ia)[j] = (tia[1]-tia[0])*bs+(*ia)[j-1];
962462f5fdSStefano Zampini       }
972462f5fdSStefano Zampini     }
982462f5fdSStefano Zampini 
992462f5fdSStefano Zampini     for (i=1; i<n; i++) {
1002462f5fdSStefano Zampini       (*ia)[i*bs] = (tia[i]-tia[i-1])*bs + (*ia)[i*bs-1];
1012462f5fdSStefano Zampini       for (j=1; j<bs; j++) {
1022462f5fdSStefano Zampini         (*ia)[i*bs+j] = (tia[i+1]-tia[i])*bs + (*ia)[i*bs+j-1];
1032462f5fdSStefano Zampini       }
1042462f5fdSStefano Zampini     }
1052462f5fdSStefano Zampini     if (n) {
1062462f5fdSStefano Zampini       (*ia)[n*bs] = (tia[n]-tia[n-1])*bs + (*ia)[n*bs-1];
1072462f5fdSStefano Zampini     }
1082462f5fdSStefano Zampini 
1092462f5fdSStefano Zampini     if (inja) {
1102462f5fdSStefano Zampini       ierr = PetscMalloc1(nz*bs*bs,ja);CHKERRQ(ierr);
1112462f5fdSStefano Zampini       cnt = 0;
1122462f5fdSStefano Zampini       for (i=0; i<n; i++) {
1138f7157efSSatish Balay         for (j=0; j<bs; j++) {
1142462f5fdSStefano Zampini           for (k=tia[i]; k<tia[i+1]; k++) {
1152462f5fdSStefano Zampini             for (l=0; l<bs; l++) {
1162462f5fdSStefano Zampini               (*ja)[cnt++] = bs*tja[k] + l;
1178f7157efSSatish Balay             }
1188f7157efSSatish Balay           }
1198f7157efSSatish Balay         }
1208f7157efSSatish Balay       }
1218f7157efSSatish Balay     }
1222462f5fdSStefano Zampini 
1232462f5fdSStefano Zampini     if (symmetric) { /* deallocate memory allocated in MatToSymmetricIJ_SeqAIJ() */
1242462f5fdSStefano Zampini       ierr = PetscFree(tia);CHKERRQ(ierr);
1252462f5fdSStefano Zampini       ierr = PetscFree(tja);CHKERRQ(ierr);
1262462f5fdSStefano Zampini     }
1272462f5fdSStefano Zampini   } else if (oshift == 1) {
1282462f5fdSStefano Zampini     if (symmetric) {
1292462f5fdSStefano Zampini       nz = tia[A->rmap->n/bs];
1302462f5fdSStefano Zampini       /*  add 1 to i and j indices */
1312462f5fdSStefano Zampini       for (i=0; i<A->rmap->n/bs+1; i++) tia[i] = tia[i] + 1;
1322462f5fdSStefano Zampini       *ia = tia;
1332462f5fdSStefano Zampini       if (ja) {
1342462f5fdSStefano Zampini         for (i=0; i<nz; i++) tja[i] = tja[i] + 1;
1352462f5fdSStefano Zampini         *ja = tja;
1362462f5fdSStefano Zampini       }
1372462f5fdSStefano Zampini     } else {
1382462f5fdSStefano Zampini       nz = a->i[A->rmap->n/bs];
1392462f5fdSStefano Zampini       /* malloc space and  add 1 to i and j indices */
1402462f5fdSStefano Zampini       ierr = PetscMalloc1(A->rmap->n/bs+1,ia);CHKERRQ(ierr);
1412462f5fdSStefano Zampini       for (i=0; i<A->rmap->n/bs+1; i++) (*ia)[i] = a->i[i] + 1;
1422462f5fdSStefano Zampini       if (ja) {
1432462f5fdSStefano Zampini         ierr = PetscMalloc1(nz,ja);CHKERRQ(ierr);
1442462f5fdSStefano Zampini         for (i=0; i<nz; i++) (*ja)[i] = a->j[i] + 1;
1452462f5fdSStefano Zampini       }
1462462f5fdSStefano Zampini     }
1472462f5fdSStefano Zampini   } else {
1482462f5fdSStefano Zampini     *ia = tia;
1492462f5fdSStefano Zampini     if (ja) *ja = tja;
150a6ece127SHong Zhang   }
15149b5e25fSSatish Balay   PetscFunctionReturn(0);
15249b5e25fSSatish Balay }
15349b5e25fSSatish Balay 
1541a83f524SJed Brown static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
15549b5e25fSSatish Balay {
1568f7157efSSatish Balay   PetscErrorCode ierr;
157a6ece127SHong Zhang 
15849b5e25fSSatish Balay   PetscFunctionBegin;
15949b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
1602462f5fdSStefano Zampini   if ((!blockcompressed && A->rmap->bs > 1) || (symmetric || oshift == 1)) {
1612462f5fdSStefano Zampini     ierr = PetscFree(*ia);CHKERRQ(ierr);
1622462f5fdSStefano Zampini     if (ja) {ierr = PetscFree(*ja);CHKERRQ(ierr);}
163a6ece127SHong Zhang   }
164a6ece127SHong Zhang   PetscFunctionReturn(0);
16549b5e25fSSatish Balay }
16649b5e25fSSatish Balay 
167dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
16849b5e25fSSatish Balay {
16949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
170dfbe8321SBarry Smith   PetscErrorCode ierr;
17149b5e25fSSatish Balay 
17249b5e25fSSatish Balay   PetscFunctionBegin;
173a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
174d0f46423SBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
175a9f03627SSatish Balay #endif
176e6b907acSBarry Smith   ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr);
1777f53bb6cSHong Zhang   if (a->free_diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
1786bf464f9SBarry Smith   ierr = ISDestroy(&a->row);CHKERRQ(ierr);
1796bf464f9SBarry Smith   ierr = ISDestroy(&a->col);CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = ISDestroy(&a->icol);CHKERRQ(ierr);
181c31cb41cSBarry Smith   ierr = PetscFree(a->idiag);CHKERRQ(ierr);
182c31cb41cSBarry Smith   ierr = PetscFree(a->inode.size);CHKERRQ(ierr);
183c760cd28SBarry Smith   if (a->free_imax_ilen) {ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);}
18405b42c5fSBarry Smith   ierr = PetscFree(a->solve_work);CHKERRQ(ierr);
18541f059aeSBarry Smith   ierr = PetscFree(a->sor_work);CHKERRQ(ierr);
18605b42c5fSBarry Smith   ierr = PetscFree(a->solves_work);CHKERRQ(ierr);
18705b42c5fSBarry Smith   ierr = PetscFree(a->mult_work);CHKERRQ(ierr);
18805b42c5fSBarry Smith   ierr = PetscFree(a->saved_values);CHKERRQ(ierr);
1894da8f245SBarry Smith   if (a->free_jshort) {ierr = PetscFree(a->jshort);CHKERRQ(ierr);}
1901a3463dfSHong Zhang   ierr = PetscFree(a->inew);CHKERRQ(ierr);
1916bf464f9SBarry Smith   ierr = MatDestroy(&a->parent);CHKERRQ(ierr);
192bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
193901853e0SKris Buschelman 
194dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
195bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C",NULL);CHKERRQ(ierr);
196bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C",NULL);CHKERRQ(ierr);
197bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C",NULL);CHKERRQ(ierr);
198bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C",NULL);CHKERRQ(ierr);
199bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C",NULL);CHKERRQ(ierr);
200bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",NULL);CHKERRQ(ierr);
20138f409ebSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocationCSR_C",NULL);CHKERRQ(ierr);
2026214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
2036214f412SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_elemental_C",NULL);CHKERRQ(ierr);
2046214f412SHong Zhang #endif
20549b5e25fSSatish Balay   PetscFunctionReturn(0);
20649b5e25fSSatish Balay }
20749b5e25fSSatish Balay 
208ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
20949b5e25fSSatish Balay {
210045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
21163ba0a88SBarry Smith   PetscErrorCode ierr;
21249b5e25fSSatish Balay 
21349b5e25fSSatish Balay   PetscFunctionBegin;
2144d9d31abSKris Buschelman   switch (op) {
2154d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
2164e0d8c25SBarry Smith     a->roworiented = flg;
2174d9d31abSKris Buschelman     break;
218a9817697SBarry Smith   case MAT_KEEP_NONZERO_PATTERN:
219a9817697SBarry Smith     a->keepnonzeropattern = flg;
2204d9d31abSKris Buschelman     break;
221512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
222512a5fc5SBarry Smith     a->nonew = (flg ? 0 : 1);
2234d9d31abSKris Buschelman     break;
2244d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
2254e0d8c25SBarry Smith     a->nonew = (flg ? -1 : 0);
2264d9d31abSKris Buschelman     break;
2274d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
2284e0d8c25SBarry Smith     a->nonew = (flg ? -2 : 0);
2294d9d31abSKris Buschelman     break;
23028b2fa4aSMatthew Knepley   case MAT_UNUSED_NONZERO_LOCATION_ERR:
23128b2fa4aSMatthew Knepley     a->nounused = (flg ? -1 : 0);
23228b2fa4aSMatthew Knepley     break;
2334e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
2344d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
2354d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
236071fcb05SBarry Smith   case MAT_SORTED_FULL:
237290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
2384d9d31abSKris Buschelman     break;
2399a4540c5SBarry Smith   case MAT_HERMITIAN:
2400f2140c7SStefano Zampini #if defined(PETSC_USE_COMPLEX) /* MAT_HERMITIAN is a synonym for MAT_SYMMETRIC with reals */
241e32f2f54SBarry Smith     if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
2420bc54ff2SBarry Smith     if (A->cmap->n < 65536 && A->cmap->bs == 1) {
243eeffb40dSHong Zhang       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort;
2440bc54ff2SBarry Smith     } else if (A->cmap->bs == 1) {
245eeffb40dSHong Zhang       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian;
246e32f2f54SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
2470f2140c7SStefano Zampini #endif
248eeffb40dSHong Zhang     break;
2493d472b54SHong Zhang   case MAT_SPD:
2505021d80fSJed Brown     /* These options are handled directly by MatSetOption() */
2513d472b54SHong Zhang     break;
25277e54ba9SKris Buschelman   case MAT_SYMMETRIC:
25377e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
2549a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
255672ba085SHong Zhang   case MAT_STRUCTURE_ONLY:
2564dcd73b1SHong Zhang     /* These options are handled directly by MatSetOption() */
257290bbb0aSBarry Smith     break;
258941593c8SHong Zhang   case MAT_IGNORE_LOWER_TRIANGULAR:
2594e0d8c25SBarry Smith     a->ignore_ltriangular = flg;
260941593c8SHong Zhang     break;
261941593c8SHong Zhang   case MAT_ERROR_LOWER_TRIANGULAR:
2624e0d8c25SBarry Smith     a->ignore_ltriangular = flg;
26377e54ba9SKris Buschelman     break;
264f5edf698SHong Zhang   case MAT_GETROW_UPPERTRIANGULAR:
2654e0d8c25SBarry Smith     a->getrow_utriangular = flg;
266f5edf698SHong Zhang     break;
267c10200c1SHong Zhang   case MAT_SUBMAT_SINGLEIS:
268c10200c1SHong Zhang     break;
2694d9d31abSKris Buschelman   default:
270e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
27149b5e25fSSatish Balay   }
27249b5e25fSSatish Balay   PetscFunctionReturn(0);
27349b5e25fSSatish Balay }
27449b5e25fSSatish Balay 
27552768537SHong Zhang PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
27649b5e25fSSatish Balay {
27749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
2786849ba73SBarry Smith   PetscErrorCode ierr;
27949b5e25fSSatish Balay 
28049b5e25fSSatish Balay   PetscFunctionBegin;
281e32f2f54SBarry Smith   if (A && !a->getrow_utriangular) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");
28252768537SHong Zhang 
283f5edf698SHong Zhang   /* Get the upper triangular part of the row */
28452768537SHong Zhang   ierr = MatGetRow_SeqBAIJ_private(A,row,nz,idx,v,a->i,a->j,a->a);CHKERRQ(ierr);
28549b5e25fSSatish Balay   PetscFunctionReturn(0);
28649b5e25fSSatish Balay }
28749b5e25fSSatish Balay 
28813f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
28949b5e25fSSatish Balay {
290dfbe8321SBarry Smith   PetscErrorCode ierr;
29149b5e25fSSatish Balay 
29249b5e25fSSatish Balay   PetscFunctionBegin;
29305b42c5fSBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
29405b42c5fSBarry Smith   if (v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}
29549b5e25fSSatish Balay   PetscFunctionReturn(0);
29649b5e25fSSatish Balay }
29749b5e25fSSatish Balay 
298f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
299f5edf698SHong Zhang {
300f5edf698SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
301f5edf698SHong Zhang 
302f5edf698SHong Zhang   PetscFunctionBegin;
303f5edf698SHong Zhang   a->getrow_utriangular = PETSC_TRUE;
304f5edf698SHong Zhang   PetscFunctionReturn(0);
305f5edf698SHong Zhang }
306a323099bSStefano Zampini 
307f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
308f5edf698SHong Zhang {
309f5edf698SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
310f5edf698SHong Zhang 
311f5edf698SHong Zhang   PetscFunctionBegin;
312f5edf698SHong Zhang   a->getrow_utriangular = PETSC_FALSE;
313f5edf698SHong Zhang   PetscFunctionReturn(0);
314f5edf698SHong Zhang }
315f5edf698SHong Zhang 
316fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
31749b5e25fSSatish Balay {
318dfbe8321SBarry Smith   PetscErrorCode ierr;
3195fd66863SKarl Rupp 
32049b5e25fSSatish Balay   PetscFunctionBegin;
321cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
322999d9058SBarry Smith     ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
323cf37664fSBarry Smith   } else if (reuse == MAT_REUSE_MATRIX) {
324cf37664fSBarry Smith     ierr = MatCopy(A,*B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
325fc4dec0aSBarry Smith   }
3268115998fSBarry Smith   PetscFunctionReturn(0);
32749b5e25fSSatish Balay }
32849b5e25fSSatish Balay 
3297da1fb6eSBarry Smith PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
33049b5e25fSSatish Balay {
33149b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
332dfbe8321SBarry Smith   PetscErrorCode    ierr;
333d0f46423SBarry Smith   PetscInt          i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
334f3ef73ceSBarry Smith   PetscViewerFormat format;
335121deb67SSatish Balay   PetscInt          *diag;
33649b5e25fSSatish Balay 
33749b5e25fSSatish Balay   PetscFunctionBegin;
338b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
339456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
34077431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);CHKERRQ(ierr);
341fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
342d2507d54SMatthew Knepley     Mat        aij;
343ade3a672SBarry Smith     const char *matname;
344ade3a672SBarry Smith 
345d5f3da31SBarry Smith     if (A->factortype && bs>1) {
34670d5e725SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr);
34770d5e725SHong Zhang       PetscFunctionReturn(0);
34870d5e725SHong Zhang     }
349c9f458caSMatthew Knepley     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr);
350ade3a672SBarry Smith     ierr = PetscObjectGetName((PetscObject)A,&matname);CHKERRQ(ierr);
351ade3a672SBarry Smith     ierr = PetscObjectSetName((PetscObject)aij,matname);CHKERRQ(ierr);
352c9f458caSMatthew Knepley     ierr = MatView(aij,viewer);CHKERRQ(ierr);
3536bf464f9SBarry Smith     ierr = MatDestroy(&aij);CHKERRQ(ierr);
354fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
355d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
35649b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
35749b5e25fSSatish Balay       for (j=0; j<bs; j++) {
35877431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
35949b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
36049b5e25fSSatish Balay           for (l=0; l<bs; l++) {
36149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
36249b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
36357622a8eSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
36457622a8eSBarry Smith                                             (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
36549b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
36657622a8eSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
36757622a8eSBarry Smith                                             (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
36849b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
36957622a8eSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37049b5e25fSSatish Balay             }
37149b5e25fSSatish Balay #else
37249b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
37357622a8eSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
37449b5e25fSSatish Balay             }
37549b5e25fSSatish Balay #endif
37649b5e25fSSatish Balay           }
37749b5e25fSSatish Balay         }
378b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
37949b5e25fSSatish Balay       }
38049b5e25fSSatish Balay     }
381d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
382c1490034SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
383c1490034SHong Zhang     PetscFunctionReturn(0);
38449b5e25fSSatish Balay   } else {
385d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
3862c990fa1SHong Zhang     if (A->factortype) { /* for factored matrix */
3872c990fa1SHong Zhang       if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");
3882c990fa1SHong Zhang 
389121deb67SSatish Balay       diag=a->diag;
390121deb67SSatish Balay       for (i=0; i<a->mbs; i++) { /* for row block i */
3912c990fa1SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
3922c990fa1SHong Zhang         /* diagonal entry */
3932c990fa1SHong Zhang #if defined(PETSC_USE_COMPLEX)
3942c990fa1SHong Zhang         if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
39557622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),(double)PetscImaginaryPart(1.0/a->a[diag[i]]));CHKERRQ(ierr);
3962c990fa1SHong Zhang         } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
39757622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]),-(double)PetscImaginaryPart(1.0/a->a[diag[i]]));CHKERRQ(ierr);
3982c990fa1SHong Zhang         } else {
39957622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)PetscRealPart(1.0/a->a[diag[i]]));CHKERRQ(ierr);
4002c990fa1SHong Zhang         }
4012c990fa1SHong Zhang #else
4026712e2f1SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[diag[i]],(double)(1.0/a->a[diag[i]]));CHKERRQ(ierr);
4032c990fa1SHong Zhang #endif
4042c990fa1SHong Zhang         /* off-diagonal entries */
4052c990fa1SHong Zhang         for (k=a->i[i]; k<a->i[i+1]-1; k++) {
4062c990fa1SHong Zhang #if defined(PETSC_USE_COMPLEX)
407ca0704adSBarry Smith           if (PetscImaginaryPart(a->a[k]) > 0.0) {
40857622a8eSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),(double)PetscImaginaryPart(a->a[k]));CHKERRQ(ierr);
409ca0704adSBarry Smith           } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
41057622a8eSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k],(double)PetscRealPart(a->a[k]),-(double)PetscImaginaryPart(a->a[k]));CHKERRQ(ierr);
4112c990fa1SHong Zhang           } else {
41257622a8eSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k],(double)PetscRealPart(a->a[k]));CHKERRQ(ierr);
4132c990fa1SHong Zhang           }
4142c990fa1SHong Zhang #else
41557622a8eSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",a->j[k],(double)a->a[k]);CHKERRQ(ierr);
4162c990fa1SHong Zhang #endif
4172c990fa1SHong Zhang         }
4182c990fa1SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
4192c990fa1SHong Zhang       }
4202c990fa1SHong Zhang 
4212c990fa1SHong Zhang     } else { /* for non-factored matrix */
4220c74a584SJed Brown       for (i=0; i<a->mbs; i++) { /* for row block i */
4230c74a584SJed Brown         for (j=0; j<bs; j++) {   /* for row bs*i + j */
42477431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
4250c74a584SJed Brown           for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
4260c74a584SJed Brown             for (l=0; l<bs; l++) {            /* for column */
42749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
42849b5e25fSSatish Balay               if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
42957622a8eSBarry Smith                 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
43057622a8eSBarry Smith                                               (double)PetscRealPart(a->a[bs2*k + l*bs + j]),(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
43149b5e25fSSatish Balay               } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
43257622a8eSBarry Smith                 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
43357622a8eSBarry Smith                                               (double)PetscRealPart(a->a[bs2*k + l*bs + j]),-(double)PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
43449b5e25fSSatish Balay               } else {
43557622a8eSBarry Smith                 ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
43649b5e25fSSatish Balay               }
43749b5e25fSSatish Balay #else
43857622a8eSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,(double)a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
43949b5e25fSSatish Balay #endif
44049b5e25fSSatish Balay             }
44149b5e25fSSatish Balay           }
442b0a32e0cSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
44349b5e25fSSatish Balay         }
44449b5e25fSSatish Balay       }
4452c990fa1SHong Zhang     }
446d00279f6SBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
44749b5e25fSSatish Balay   }
448b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
44949b5e25fSSatish Balay   PetscFunctionReturn(0);
45049b5e25fSSatish Balay }
45149b5e25fSSatish Balay 
4529804daf3SBarry Smith #include <petscdraw.h>
4536849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
45449b5e25fSSatish Balay {
45549b5e25fSSatish Balay   Mat            A = (Mat) Aa;
45649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
4576849ba73SBarry Smith   PetscErrorCode ierr;
458d0f46423SBarry Smith   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
45949b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
46049b5e25fSSatish Balay   MatScalar      *aa;
461b0a32e0cSBarry Smith   PetscViewer    viewer;
46249b5e25fSSatish Balay 
46349b5e25fSSatish Balay   PetscFunctionBegin;
46449b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
465b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
46649b5e25fSSatish Balay 
46749b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
468383922c3SLisandro Dalcin 
469383922c3SLisandro Dalcin   ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
470383922c3SLisandro Dalcin   ierr = PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");CHKERRQ(ierr);
471383922c3SLisandro Dalcin   /* Blue for negative, Cyan for zero and  Red for positive */
472b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
47349b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
47449b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
475d0f46423SBarry Smith       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
47649b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
47749b5e25fSSatish Balay       aa  = a->a + j*bs2;
47849b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
48049b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
481b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
48249b5e25fSSatish Balay         }
48349b5e25fSSatish Balay       }
48449b5e25fSSatish Balay     }
48549b5e25fSSatish Balay   }
486b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
48749b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
48849b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
489d0f46423SBarry Smith       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
49049b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
49149b5e25fSSatish Balay       aa = a->a + j*bs2;
49249b5e25fSSatish Balay       for (k=0; k<bs; k++) {
49349b5e25fSSatish Balay         for (l=0; l<bs; l++) {
49449b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
495b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
49649b5e25fSSatish Balay         }
49749b5e25fSSatish Balay       }
49849b5e25fSSatish Balay     }
49949b5e25fSSatish Balay   }
500b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
50149b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
50249b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
503d0f46423SBarry Smith       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
50449b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
50549b5e25fSSatish Balay       aa = a->a + j*bs2;
50649b5e25fSSatish Balay       for (k=0; k<bs; k++) {
50749b5e25fSSatish Balay         for (l=0; l<bs; l++) {
50849b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
509b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
51049b5e25fSSatish Balay         }
51149b5e25fSSatish Balay       }
51249b5e25fSSatish Balay     }
51349b5e25fSSatish Balay   }
514383922c3SLisandro Dalcin   ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
51549b5e25fSSatish Balay   PetscFunctionReturn(0);
51649b5e25fSSatish Balay }
51749b5e25fSSatish Balay 
5186849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
51949b5e25fSSatish Balay {
520dfbe8321SBarry Smith   PetscErrorCode ierr;
52149b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,w,h;
522b0a32e0cSBarry Smith   PetscDraw      draw;
523ace3abfcSBarry Smith   PetscBool      isnull;
52449b5e25fSSatish Balay 
52549b5e25fSSatish Balay   PetscFunctionBegin;
526b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
527383922c3SLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
528383922c3SLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
52949b5e25fSSatish Balay 
530d0f46423SBarry Smith   xr   = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
53149b5e25fSSatish Balay   xr  += w;          yr += h;        xl = -w;     yl = -h;
532b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
533832b7cebSLisandro Dalcin   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
534b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
5350298fd71SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);CHKERRQ(ierr);
536832b7cebSLisandro Dalcin   ierr = PetscDrawSave(draw);CHKERRQ(ierr);
53749b5e25fSSatish Balay   PetscFunctionReturn(0);
53849b5e25fSSatish Balay }
53949b5e25fSSatish Balay 
540dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
54149b5e25fSSatish Balay {
542dfbe8321SBarry Smith   PetscErrorCode ierr;
543ace3abfcSBarry Smith   PetscBool      iascii,isdraw;
54408917f38SBarry Smith   FILE           *file = 0;
54549b5e25fSSatish Balay 
54649b5e25fSSatish Balay   PetscFunctionBegin;
547251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
548251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
54932077d6dSBarry Smith   if (iascii) {
55049b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
55149b5e25fSSatish Balay   } else if (isdraw) {
55249b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
55349b5e25fSSatish Balay   } else {
554a5e6ed63SBarry Smith     Mat        B;
555ade3a672SBarry Smith     const char *matname;
556ceb03754SKris Buschelman     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
557ade3a672SBarry Smith     ierr = PetscObjectGetName((PetscObject)A,&matname);CHKERRQ(ierr);
558ade3a672SBarry Smith     ierr = PetscObjectSetName((PetscObject)B,matname);CHKERRQ(ierr);
559a5e6ed63SBarry Smith     ierr = MatView(B,viewer);CHKERRQ(ierr);
5606bf464f9SBarry Smith     ierr = MatDestroy(&B);CHKERRQ(ierr);
56108917f38SBarry Smith     ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
56208917f38SBarry Smith     if (file) {
56308917f38SBarry Smith       fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
56408917f38SBarry Smith     }
56549b5e25fSSatish Balay   }
56649b5e25fSSatish Balay   PetscFunctionReturn(0);
56749b5e25fSSatish Balay }
56849b5e25fSSatish Balay 
56949b5e25fSSatish Balay 
57013f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
57149b5e25fSSatish Balay {
572045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
57313f74950SBarry Smith   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
57413f74950SBarry Smith   PetscInt     *ai = a->i,*ailen = a->ilen;
575d0f46423SBarry Smith   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
57697e567efSBarry Smith   MatScalar    *ap,*aa = a->a;
57749b5e25fSSatish Balay 
57849b5e25fSSatish Balay   PetscFunctionBegin;
57949b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
58049b5e25fSSatish Balay     row = im[k]; brow = row/bs;
581e32f2f54SBarry Smith     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
582e32f2f54SBarry Smith     if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
58349b5e25fSSatish Balay     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
58449b5e25fSSatish Balay     nrow = ailen[brow];
58549b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
586e32f2f54SBarry Smith       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
587e32f2f54SBarry Smith       if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
58849b5e25fSSatish Balay       col  = in[l];
58949b5e25fSSatish Balay       bcol = col/bs;
59049b5e25fSSatish Balay       cidx = col%bs;
59149b5e25fSSatish Balay       ridx = row%bs;
59249b5e25fSSatish Balay       high = nrow;
59349b5e25fSSatish Balay       low  = 0; /* assume unsorted */
59449b5e25fSSatish Balay       while (high-low > 5) {
59549b5e25fSSatish Balay         t = (low+high)/2;
59649b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
59749b5e25fSSatish Balay         else              low  = t;
59849b5e25fSSatish Balay       }
59949b5e25fSSatish Balay       for (i=low; i<high; i++) {
60049b5e25fSSatish Balay         if (rp[i] > bcol) break;
60149b5e25fSSatish Balay         if (rp[i] == bcol) {
60249b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
60349b5e25fSSatish Balay           goto finished;
60449b5e25fSSatish Balay         }
60549b5e25fSSatish Balay       }
60697e567efSBarry Smith       *v++ = 0.0;
60749b5e25fSSatish Balay finished:;
60849b5e25fSSatish Balay     }
60949b5e25fSSatish Balay   }
61049b5e25fSSatish Balay   PetscFunctionReturn(0);
61149b5e25fSSatish Balay }
61249b5e25fSSatish Balay 
61349b5e25fSSatish Balay 
61413f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
61549b5e25fSSatish Balay {
6160880e062SHong Zhang   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
6176849ba73SBarry Smith   PetscErrorCode    ierr;
618e2ee6c50SBarry Smith   PetscInt          *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
61913f74950SBarry Smith   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
620d0f46423SBarry Smith   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
621ace3abfcSBarry Smith   PetscBool         roworiented=a->roworiented;
622dd6ea824SBarry Smith   const PetscScalar *value     = v;
623f15d580aSBarry Smith   MatScalar         *ap,*aa = a->a,*bap;
6240880e062SHong Zhang 
62549b5e25fSSatish Balay   PetscFunctionBegin;
62626fbe8dcSKarl Rupp   if (roworiented) stepval = (n-1)*bs;
62726fbe8dcSKarl Rupp   else stepval = (m-1)*bs;
62826fbe8dcSKarl Rupp 
6290880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
6300880e062SHong Zhang     row = im[k];
6310880e062SHong Zhang     if (row < 0) continue;
6322515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
6332f7d4af7SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index row too large %D max %D",row,a->mbs-1);
6340880e062SHong Zhang #endif
6350880e062SHong Zhang     rp   = aj + ai[row];
6360880e062SHong Zhang     ap   = aa + bs2*ai[row];
6370880e062SHong Zhang     rmax = imax[row];
6380880e062SHong Zhang     nrow = ailen[row];
6390880e062SHong Zhang     low  = 0;
640818f2c47SBarry Smith     high = nrow;
6410880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
6420880e062SHong Zhang       if (in[l] < 0) continue;
6430880e062SHong Zhang       col = in[l];
6442515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
6452f7d4af7SBarry Smith       if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block index column too large %D max %D",col,a->nbs-1);
646b1823623SSatish Balay #endif
647b98bf0e1SJed Brown       if (col < row) {
64826fbe8dcSKarl Rupp         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
64926fbe8dcSKarl Rupp         else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
650b98bf0e1SJed Brown       }
65126fbe8dcSKarl Rupp       if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
65226fbe8dcSKarl Rupp       else value = v + l*(stepval+bs)*bs + k*bs;
65326fbe8dcSKarl Rupp 
65426fbe8dcSKarl Rupp       if (col <= lastcol) low = 0;
65526fbe8dcSKarl Rupp       else high = nrow;
65626fbe8dcSKarl Rupp 
657e2ee6c50SBarry Smith       lastcol = col;
6580880e062SHong Zhang       while (high-low > 7) {
6590880e062SHong Zhang         t = (low+high)/2;
6600880e062SHong Zhang         if (rp[t] > col) high = t;
6610880e062SHong Zhang         else             low  = t;
6620880e062SHong Zhang       }
6630880e062SHong Zhang       for (i=low; i<high; i++) {
6640880e062SHong Zhang         if (rp[i] > col) break;
6650880e062SHong Zhang         if (rp[i] == col) {
6660880e062SHong Zhang           bap = ap +  bs2*i;
6670880e062SHong Zhang           if (roworiented) {
6680880e062SHong Zhang             if (is == ADD_VALUES) {
6690880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6700880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6710880e062SHong Zhang                   bap[jj] += *value++;
6720880e062SHong Zhang                 }
6730880e062SHong Zhang               }
6740880e062SHong Zhang             } else {
6750880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6760880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6770880e062SHong Zhang                   bap[jj] = *value++;
6780880e062SHong Zhang                 }
6790880e062SHong Zhang                }
6800880e062SHong Zhang             }
6810880e062SHong Zhang           } else {
6820880e062SHong Zhang             if (is == ADD_VALUES) {
6830880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6840880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6850880e062SHong Zhang                   *bap++ += *value++;
6860880e062SHong Zhang                 }
6870880e062SHong Zhang               }
6880880e062SHong Zhang             } else {
6890880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6900880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6910880e062SHong Zhang                   *bap++  = *value++;
6920880e062SHong Zhang                 }
6930880e062SHong Zhang               }
6940880e062SHong Zhang             }
6950880e062SHong Zhang           }
6960880e062SHong Zhang           goto noinsert2;
6970880e062SHong Zhang         }
6980880e062SHong Zhang       }
6990880e062SHong Zhang       if (nonew == 1) goto noinsert2;
7002f7d4af7SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new block index nonzero block (%D, %D) in the matrix", row, col);
701fef13f97SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
702c03d1d03SSatish Balay       N = nrow++ - 1; high++;
7030880e062SHong Zhang       /* shift up all the later entries in this row */
704580bdb30SBarry Smith       ierr  = PetscArraymove(rp+i+1,rp+i,N-i+1);CHKERRQ(ierr);
705580bdb30SBarry Smith       ierr  = PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));CHKERRQ(ierr);
706580bdb30SBarry Smith       ierr  = PetscArrayzero(ap+bs2*i,bs2);CHKERRQ(ierr);
7070880e062SHong Zhang       rp[i] = col;
7080880e062SHong Zhang       bap   = ap +  bs2*i;
7090880e062SHong Zhang       if (roworiented) {
7100880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7110880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
7120880e062SHong Zhang             bap[jj] = *value++;
7130880e062SHong Zhang           }
7140880e062SHong Zhang         }
7150880e062SHong Zhang       } else {
7160880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7170880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
7180880e062SHong Zhang             *bap++ = *value++;
7190880e062SHong Zhang           }
7200880e062SHong Zhang         }
7210880e062SHong Zhang        }
7220880e062SHong Zhang     noinsert2:;
7230880e062SHong Zhang       low = i;
7240880e062SHong Zhang     }
7250880e062SHong Zhang     ailen[row] = nrow;
7260880e062SHong Zhang   }
7270880e062SHong Zhang   PetscFunctionReturn(0);
72849b5e25fSSatish Balay }
72949b5e25fSSatish Balay 
73064831d72SBarry Smith /*
73164831d72SBarry Smith     This is not yet used
73264831d72SBarry Smith */
7334108e4d5SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
7340def2e27SBarry Smith {
7350def2e27SBarry Smith   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
7360def2e27SBarry Smith   PetscErrorCode ierr;
7370def2e27SBarry Smith   const PetscInt *ai = a->i, *aj = a->j,*cols;
7380def2e27SBarry Smith   PetscInt       i   = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
739ace3abfcSBarry Smith   PetscBool      flag;
7400def2e27SBarry Smith 
7410def2e27SBarry Smith   PetscFunctionBegin;
742785e854fSJed Brown   ierr = PetscMalloc1(m,&ns);CHKERRQ(ierr);
7430def2e27SBarry Smith   while (i < m) {
7440def2e27SBarry Smith     nzx = ai[i+1] - ai[i];       /* Number of nonzeros */
7450def2e27SBarry Smith     /* Limits the number of elements in a node to 'a->inode.limit' */
7460def2e27SBarry Smith     for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
7470def2e27SBarry Smith       nzy = ai[j+1] - ai[j];
7480def2e27SBarry Smith       if (nzy != (nzx - j + i)) break;
749580bdb30SBarry Smith       ierr = PetscArraycmp(aj + ai[i] + j - i,aj + ai[j],nzy,&flag);CHKERRQ(ierr);
7500def2e27SBarry Smith       if (!flag) break;
7510def2e27SBarry Smith     }
7520def2e27SBarry Smith     ns[node_count++] = blk_size;
75326fbe8dcSKarl Rupp 
7540def2e27SBarry Smith     i = j;
7550def2e27SBarry Smith   }
7560def2e27SBarry Smith   if (!a->inode.size && m && node_count > .9*m) {
7570def2e27SBarry Smith     ierr = PetscFree(ns);CHKERRQ(ierr);
7580def2e27SBarry Smith     ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr);
7590def2e27SBarry Smith   } else {
7600def2e27SBarry Smith     a->inode.node_count = node_count;
76126fbe8dcSKarl Rupp 
762785e854fSJed Brown     ierr = PetscMalloc1(node_count,&a->inode.size);CHKERRQ(ierr);
7633bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,node_count*sizeof(PetscInt));CHKERRQ(ierr);
764580bdb30SBarry Smith     ierr = PetscArraycpy(a->inode.size,ns,node_count);CHKERRQ(ierr);
7650def2e27SBarry Smith     ierr = PetscFree(ns);CHKERRQ(ierr);
7660def2e27SBarry Smith     ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr);
7670def2e27SBarry Smith 
7680def2e27SBarry Smith     /* count collections of adjacent columns in each inode */
7690def2e27SBarry Smith     row = 0;
7700def2e27SBarry Smith     cnt = 0;
7710def2e27SBarry Smith     for (i=0; i<node_count; i++) {
7720def2e27SBarry Smith       cols = aj + ai[row] + a->inode.size[i];
7730def2e27SBarry Smith       nz   = ai[row+1] - ai[row] - a->inode.size[i];
7740def2e27SBarry Smith       for (j=1; j<nz; j++) {
77526fbe8dcSKarl Rupp         if (cols[j] != cols[j-1]+1) cnt++;
7760def2e27SBarry Smith       }
7770def2e27SBarry Smith       cnt++;
7780def2e27SBarry Smith       row += a->inode.size[i];
7790def2e27SBarry Smith     }
780785e854fSJed Brown     ierr = PetscMalloc1(2*cnt,&counts);CHKERRQ(ierr);
7810def2e27SBarry Smith     cnt  = 0;
7820def2e27SBarry Smith     row  = 0;
7830def2e27SBarry Smith     for (i=0; i<node_count; i++) {
7840def2e27SBarry Smith       cols = aj + ai[row] + a->inode.size[i];
7850def2e27SBarry Smith       counts[2*cnt] = cols[0];
7860def2e27SBarry Smith       nz   = ai[row+1] - ai[row] - a->inode.size[i];
7870def2e27SBarry Smith       cnt2 = 1;
7880def2e27SBarry Smith       for (j=1; j<nz; j++) {
7890def2e27SBarry Smith         if (cols[j] != cols[j-1]+1) {
7900def2e27SBarry Smith           counts[2*(cnt++)+1] = cnt2;
7910def2e27SBarry Smith           counts[2*cnt]       = cols[j];
7920def2e27SBarry Smith           cnt2 = 1;
7930def2e27SBarry Smith         } else cnt2++;
7940def2e27SBarry Smith       }
7950def2e27SBarry Smith       counts[2*(cnt++)+1] = cnt2;
7960def2e27SBarry Smith       row += a->inode.size[i];
7970def2e27SBarry Smith     }
79822d28d08SBarry Smith     ierr = PetscIntView(2*cnt,counts,0);CHKERRQ(ierr);
7990def2e27SBarry Smith   }
80038702af4SBarry Smith   PetscFunctionReturn(0);
80138702af4SBarry Smith }
80238702af4SBarry Smith 
803dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
80449b5e25fSSatish Balay {
80549b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
8066849ba73SBarry Smith   PetscErrorCode ierr;
8078f8f2f0dSBarry Smith   PetscInt       fshift = 0,i,*ai = a->i,*aj = a->j,*imax = a->imax;
808d0f46423SBarry Smith   PetscInt       m      = A->rmap->N,*ip,N,*ailen = a->ilen;
80913f74950SBarry Smith   PetscInt       mbs    = a->mbs,bs2 = a->bs2,rmax = 0;
81049b5e25fSSatish Balay   MatScalar      *aa    = a->a,*ap;
81149b5e25fSSatish Balay 
81249b5e25fSSatish Balay   PetscFunctionBegin;
81349b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
81449b5e25fSSatish Balay 
81549b5e25fSSatish Balay   if (m) rmax = ailen[0];
81649b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
81749b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
81849b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
81949b5e25fSSatish Balay     rmax    = PetscMax(rmax,ailen[i]);
82049b5e25fSSatish Balay     if (fshift) {
821580bdb30SBarry Smith       ip = aj + ai[i];
822580bdb30SBarry Smith       ap = aa + bs2*ai[i];
82349b5e25fSSatish Balay       N  = ailen[i];
824580bdb30SBarry Smith       ierr  = PetscArraymove(ip-fshift,ip,N);CHKERRQ(ierr);
825580bdb30SBarry Smith       ierr  = PetscArraymove(ap-bs2*fshift,ap,bs2*N);CHKERRQ(ierr);
82649b5e25fSSatish Balay     }
82749b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
82849b5e25fSSatish Balay   }
82949b5e25fSSatish Balay   if (mbs) {
83049b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
83149b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
83249b5e25fSSatish Balay   }
83349b5e25fSSatish Balay   /* reset ilen and imax for each row */
83449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
83549b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
83649b5e25fSSatish Balay   }
8376c6c5352SBarry Smith   a->nz = ai[mbs];
83849b5e25fSSatish Balay 
839b424e231SHong Zhang   /* diagonals may have moved, reset it */
840b424e231SHong Zhang   if (a->diag) {
841580bdb30SBarry Smith     ierr = PetscArraycpy(a->diag,ai,mbs);CHKERRQ(ierr);
84249b5e25fSSatish Balay   }
84326fbe8dcSKarl Rupp   if (fshift && a->nounused == -1) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D block size %D, %D unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2);
84426fbe8dcSKarl Rupp 
845d0f46423SBarry Smith   ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap->N,A->rmap->bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr);
846ae15b995SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr);
847ae15b995SBarry Smith   ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr);
84826fbe8dcSKarl Rupp 
8498e58a170SBarry Smith   A->info.mallocs    += a->reallocs;
85049b5e25fSSatish Balay   a->reallocs         = 0;
85149b5e25fSSatish Balay   A->info.nz_unneeded = (PetscReal)fshift*bs2;
852061b2667SBarry Smith   a->idiagvalid       = PETSC_FALSE;
8534dcd73b1SHong Zhang   a->rmax             = rmax;
85438702af4SBarry Smith 
85538702af4SBarry Smith   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
85644e1c64aSLisandro Dalcin     if (a->jshort && a->free_jshort) {
85717803ae8SHong Zhang       /* when matrix data structure is changed, previous jshort must be replaced */
85817803ae8SHong Zhang       ierr = PetscFree(a->jshort);CHKERRQ(ierr);
85917803ae8SHong Zhang     }
860785e854fSJed Brown     ierr = PetscMalloc1(a->i[A->rmap->n],&a->jshort);CHKERRQ(ierr);
8613bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,a->i[A->rmap->n]*sizeof(unsigned short));CHKERRQ(ierr);
86238702af4SBarry Smith     for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
86338702af4SBarry Smith     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
86441f059aeSBarry Smith     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
8654da8f245SBarry Smith     a->free_jshort = PETSC_TRUE;
86638702af4SBarry Smith   }
86749b5e25fSSatish Balay   PetscFunctionReturn(0);
86849b5e25fSSatish Balay }
86949b5e25fSSatish Balay 
87049b5e25fSSatish Balay /*
87149b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
87249b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
87349b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
87449b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
87549b5e25fSSatish Balay */
87613f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
87749b5e25fSSatish Balay {
87813f74950SBarry Smith   PetscInt  i,j,k,row;
879ace3abfcSBarry Smith   PetscBool flg;
88049b5e25fSSatish Balay 
88149b5e25fSSatish Balay   PetscFunctionBegin;
88249b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
88349b5e25fSSatish Balay     row = idx[i];
88449b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
88549b5e25fSSatish Balay       sizes[j] = 1;
88649b5e25fSSatish Balay       i++;
88749b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
88849b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
88949b5e25fSSatish Balay       i++;
89049b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
89149b5e25fSSatish Balay       flg = PETSC_TRUE;
89249b5e25fSSatish Balay       for (k=1; k<bs; k++) {
89349b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
89449b5e25fSSatish Balay           flg = PETSC_FALSE;
89549b5e25fSSatish Balay           break;
89649b5e25fSSatish Balay         }
89749b5e25fSSatish Balay       }
898abc0a331SBarry Smith       if (flg) { /* No break in the bs */
89949b5e25fSSatish Balay         sizes[j] = bs;
90049b5e25fSSatish Balay         i       += bs;
90149b5e25fSSatish Balay       } else {
90249b5e25fSSatish Balay         sizes[j] = 1;
90349b5e25fSSatish Balay         i++;
90449b5e25fSSatish Balay       }
90549b5e25fSSatish Balay     }
90649b5e25fSSatish Balay   }
90749b5e25fSSatish Balay   *bs_max = j;
90849b5e25fSSatish Balay   PetscFunctionReturn(0);
90949b5e25fSSatish Balay }
91049b5e25fSSatish Balay 
91149b5e25fSSatish Balay 
91249b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
91349b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
91449b5e25fSSatish Balay */
91549b5e25fSSatish Balay 
91613f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
91749b5e25fSSatish Balay {
91849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
9196849ba73SBarry Smith   PetscErrorCode ierr;
920e2ee6c50SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
92113f74950SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
922d0f46423SBarry Smith   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
92313f74950SBarry Smith   PetscInt       ridx,cidx,bs2=a->bs2;
92449b5e25fSSatish Balay   MatScalar      *ap,value,*aa=a->a,*bap;
92549b5e25fSSatish Balay 
92649b5e25fSSatish Balay   PetscFunctionBegin;
92749b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
92849b5e25fSSatish Balay     row  = im[k];       /* row number */
92949b5e25fSSatish Balay     brow = row/bs;      /* block row number */
93049b5e25fSSatish Balay     if (row < 0) continue;
9312515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
932e32f2f54SBarry Smith     if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
93349b5e25fSSatish Balay #endif
93449b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
93549b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
93649b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
93749b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
93849b5e25fSSatish Balay     low  = 0;
9398509e838SStefano Zampini     high = nrow;
94049b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
94149b5e25fSSatish Balay       if (in[l] < 0) continue;
9422515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
9435ca5a24fSStefano Zampini       if (in[l] >= A->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->N-1);
94449b5e25fSSatish Balay #endif
94549b5e25fSSatish Balay       col  = in[l];
94649b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
94749b5e25fSSatish Balay 
948941593c8SHong Zhang       if (brow > bcol) {
94926fbe8dcSKarl Rupp         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
95026fbe8dcSKarl Rupp         else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
951941593c8SHong Zhang       }
952f4989cb3SHong Zhang 
95349b5e25fSSatish Balay       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
9548549e402SHong Zhang       if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
95549b5e25fSSatish Balay         /* element value a(k,l) */
95626fbe8dcSKarl Rupp         if (roworiented) value = v[l + k*n];
95726fbe8dcSKarl Rupp         else value = v[k + l*m];
95849b5e25fSSatish Balay 
95949b5e25fSSatish Balay         /* move pointer bap to a(k,l) quickly and add/insert value */
96026fbe8dcSKarl Rupp         if (col <= lastcol) low = 0;
9618509e838SStefano Zampini         else high = nrow;
9628509e838SStefano Zampini 
963e2ee6c50SBarry Smith         lastcol = col;
96449b5e25fSSatish Balay         while (high-low > 7) {
96549b5e25fSSatish Balay           t = (low+high)/2;
96649b5e25fSSatish Balay           if (rp[t] > bcol) high = t;
96749b5e25fSSatish Balay           else              low  = t;
96849b5e25fSSatish Balay         }
96949b5e25fSSatish Balay         for (i=low; i<high; i++) {
97049b5e25fSSatish Balay           if (rp[i] > bcol) break;
97149b5e25fSSatish Balay           if (rp[i] == bcol) {
97249b5e25fSSatish Balay             bap = ap +  bs2*i + bs*cidx + ridx;
97349b5e25fSSatish Balay             if (is == ADD_VALUES) *bap += value;
97449b5e25fSSatish Balay             else                  *bap  = value;
9758549e402SHong Zhang             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9768549e402SHong Zhang             if (brow == bcol && ridx < cidx) {
9778549e402SHong Zhang               bap = ap +  bs2*i + bs*ridx + cidx;
9788549e402SHong Zhang               if (is == ADD_VALUES) *bap += value;
9798549e402SHong Zhang               else                  *bap  = value;
9808549e402SHong Zhang             }
98149b5e25fSSatish Balay             goto noinsert1;
98249b5e25fSSatish Balay           }
98349b5e25fSSatish Balay         }
98449b5e25fSSatish Balay 
98549b5e25fSSatish Balay         if (nonew == 1) goto noinsert1;
986e32f2f54SBarry Smith         if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
987fef13f97SBarry Smith         MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
98849b5e25fSSatish Balay 
989c03d1d03SSatish Balay         N = nrow++ - 1; high++;
99049b5e25fSSatish Balay         /* shift up all the later entries in this row */
991580bdb30SBarry Smith         ierr = PetscArraymove(rp+i+1,rp+i,N-i+1);CHKERRQ(ierr);
992580bdb30SBarry Smith         ierr = PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));CHKERRQ(ierr);
993580bdb30SBarry Smith         ierr = PetscArrayzero(ap+bs2*i,bs2);CHKERRQ(ierr);
99449b5e25fSSatish Balay         rp[i]                      = bcol;
99549b5e25fSSatish Balay         ap[bs2*i + bs*cidx + ridx] = value;
9968509e838SStefano Zampini         /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9978509e838SStefano Zampini         if (brow == bcol && ridx < cidx) {
9988509e838SStefano Zampini           ap[bs2*i + bs*ridx + cidx] = value;
9998509e838SStefano Zampini         }
1000e56f5c9eSBarry Smith         A->nonzerostate++;
100149b5e25fSSatish Balay noinsert1:;
100249b5e25fSSatish Balay         low = i;
10038549e402SHong Zhang       }
100449b5e25fSSatish Balay     }   /* end of loop over added columns */
100549b5e25fSSatish Balay     ailen[brow] = nrow;
100649b5e25fSSatish Balay   }   /* end of loop over added rows */
100749b5e25fSSatish Balay   PetscFunctionReturn(0);
100849b5e25fSSatish Balay }
100949b5e25fSSatish Balay 
10100481f469SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
101149b5e25fSSatish Balay {
10124ccecd49SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
101349b5e25fSSatish Balay   Mat            outA;
1014dfbe8321SBarry Smith   PetscErrorCode ierr;
1015ace3abfcSBarry Smith   PetscBool      row_identity;
101649b5e25fSSatish Balay 
101749b5e25fSSatish Balay   PetscFunctionBegin;
1018e32f2f54SBarry Smith   if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1019c84f5b01SHong Zhang   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1020e32f2f54SBarry Smith   if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1021e32f2f54SBarry Smith   if (inA->rmap->bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */
1022c84f5b01SHong Zhang 
102349b5e25fSSatish Balay   outA            = inA;
1024d5f3da31SBarry Smith   inA->factortype = MAT_FACTOR_ICC;
1025f6224b95SHong Zhang   ierr = PetscFree(inA->solvertype);CHKERRQ(ierr);
1026f6224b95SHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&inA->solvertype);CHKERRQ(ierr);
102749b5e25fSSatish Balay 
10281a3463dfSHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
1029d595f711SHong Zhang   ierr = MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);CHKERRQ(ierr);
103049b5e25fSSatish Balay 
1031c3122656SLisandro Dalcin   ierr   = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
10326bf464f9SBarry Smith   ierr   = ISDestroy(&a->row);CHKERRQ(ierr);
1033c84f5b01SHong Zhang   a->row = row;
1034c3122656SLisandro Dalcin   ierr   = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
10356bf464f9SBarry Smith   ierr   = ISDestroy(&a->col);CHKERRQ(ierr);
1036c84f5b01SHong Zhang   a->col = row;
1037c84f5b01SHong Zhang 
1038c84f5b01SHong Zhang   /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
1039c84f5b01SHong Zhang   if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);}
10403bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)inA,(PetscObject)a->icol);CHKERRQ(ierr);
104149b5e25fSSatish Balay 
104249b5e25fSSatish Balay   if (!a->solve_work) {
1043854ce69bSBarry Smith     ierr = PetscMalloc1(inA->rmap->N+inA->rmap->bs,&a->solve_work);CHKERRQ(ierr);
10443bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));CHKERRQ(ierr);
104549b5e25fSSatish Balay   }
104649b5e25fSSatish Balay 
1047719d5645SBarry Smith   ierr = MatCholeskyFactorNumeric(outA,inA,info);CHKERRQ(ierr);
104849b5e25fSSatish Balay   PetscFunctionReturn(0);
104949b5e25fSSatish Balay }
1050950f1e5bSHong Zhang 
10517087cfbeSBarry Smith PetscErrorCode  MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
105249b5e25fSSatish Balay {
1053045c9aa0SHong Zhang   Mat_SeqSBAIJ   *baij = (Mat_SeqSBAIJ*)mat->data;
105413f74950SBarry Smith   PetscInt       i,nz,n;
10557827cd58SJed Brown   PetscErrorCode ierr;
105649b5e25fSSatish Balay 
105749b5e25fSSatish Balay   PetscFunctionBegin;
10586c6c5352SBarry Smith   nz = baij->maxnz;
1059d0f46423SBarry Smith   n  = mat->cmap->n;
106026fbe8dcSKarl Rupp   for (i=0; i<nz; i++) baij->j[i] = indices[i];
106126fbe8dcSKarl Rupp 
10626c6c5352SBarry Smith   baij->nz = nz;
106326fbe8dcSKarl Rupp   for (i=0; i<n; i++) baij->ilen[i] = baij->imax[i];
106426fbe8dcSKarl Rupp 
10657827cd58SJed Brown   ierr = MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
106649b5e25fSSatish Balay   PetscFunctionReturn(0);
106749b5e25fSSatish Balay }
106849b5e25fSSatish Balay 
106949b5e25fSSatish Balay /*@
107019585528SSatish Balay   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
107149b5e25fSSatish Balay   in the matrix.
107249b5e25fSSatish Balay 
107349b5e25fSSatish Balay   Input Parameters:
107419585528SSatish Balay   +  mat     - the SeqSBAIJ matrix
107549b5e25fSSatish Balay   -  indices - the column indices
107649b5e25fSSatish Balay 
107749b5e25fSSatish Balay   Level: advanced
107849b5e25fSSatish Balay 
107949b5e25fSSatish Balay   Notes:
108049b5e25fSSatish Balay   This can be called if you have precomputed the nonzero structure of the
108149b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
108249b5e25fSSatish Balay   of the MatSetValues() operation.
108349b5e25fSSatish Balay 
108449b5e25fSSatish Balay   You MUST have set the correct numbers of nonzeros per row in the call to
1085d1be2dadSMatthew Knepley   MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
108649b5e25fSSatish Balay 
1087ab9f2c04SSatish Balay   MUST be called before any calls to MatSetValues()
108849b5e25fSSatish Balay 
1089ab9f2c04SSatish Balay   .seealso: MatCreateSeqSBAIJ
109049b5e25fSSatish Balay @*/
10917087cfbeSBarry Smith PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
109249b5e25fSSatish Balay {
10934ac538c5SBarry Smith   PetscErrorCode ierr;
109449b5e25fSSatish Balay 
109549b5e25fSSatish Balay   PetscFunctionBegin;
10960700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
10974482741eSBarry Smith   PetscValidPointer(indices,2);
10984ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));CHKERRQ(ierr);
109949b5e25fSSatish Balay   PetscFunctionReturn(0);
110049b5e25fSSatish Balay }
110149b5e25fSSatish Balay 
11023c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
11033c896bc6SHong Zhang {
11043c896bc6SHong Zhang   PetscErrorCode ierr;
11054c7a3774SStefano Zampini   PetscBool      isbaij;
11063c896bc6SHong Zhang 
11073c896bc6SHong Zhang   PetscFunctionBegin;
11084c7a3774SStefano Zampini   ierr = PetscObjectTypeCompareAny((PetscObject)B,&isbaij,MATSEQSBAIJ,MATMPISBAIJ,"");CHKERRQ(ierr);
11094c7a3774SStefano Zampini   if (!isbaij) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)B)->type_name);
11104c7a3774SStefano Zampini   /* If the two matrices have the same copy implementation and nonzero pattern, use fast copy. */
11114c7a3774SStefano Zampini   if (str == SAME_NONZERO_PATTERN && A->ops->copy == B->ops->copy) {
11123c896bc6SHong Zhang     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
11133c896bc6SHong Zhang     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
11143c896bc6SHong Zhang 
11154c7a3774SStefano Zampini     if (a->i[a->mbs] != b->i[b->mbs]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
11164c7a3774SStefano Zampini     if (a->mbs != b->mbs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of rows in two matrices are different");
11174c7a3774SStefano Zampini     if (a->bs2 != b->bs2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Different block size");
1118580bdb30SBarry Smith     ierr = PetscArraycpy(b->a,a->a,a->bs2*a->i[a->mbs]);CHKERRQ(ierr);
1119cdc753b6SBarry Smith     ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr);
11203c896bc6SHong Zhang   } else {
1121f5edf698SHong Zhang     ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
11223c896bc6SHong Zhang     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1123f5edf698SHong Zhang     ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
11243c896bc6SHong Zhang   }
11253c896bc6SHong Zhang   PetscFunctionReturn(0);
11263c896bc6SHong Zhang }
11273c896bc6SHong Zhang 
11284994cf47SJed Brown PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1129273d9f13SBarry Smith {
1130dfbe8321SBarry Smith   PetscErrorCode ierr;
1131273d9f13SBarry Smith 
1132273d9f13SBarry Smith   PetscFunctionBegin;
1133367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0);CHKERRQ(ierr);
1134273d9f13SBarry Smith   PetscFunctionReturn(0);
1135273d9f13SBarry Smith }
1136273d9f13SBarry Smith 
11378c778c55SBarry Smith PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1138a6ece127SHong Zhang {
1139a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
11405fd66863SKarl Rupp 
1141a6ece127SHong Zhang   PetscFunctionBegin;
1142a6ece127SHong Zhang   *array = a->a;
1143a6ece127SHong Zhang   PetscFunctionReturn(0);
1144a6ece127SHong Zhang }
1145a6ece127SHong Zhang 
11468c778c55SBarry Smith PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1147a6ece127SHong Zhang {
1148a6ece127SHong Zhang   PetscFunctionBegin;
1149a6ece127SHong Zhang   PetscFunctionReturn(0);
1150a6ece127SHong Zhang }
1151a6ece127SHong Zhang 
115252768537SHong Zhang PetscErrorCode MatAXPYGetPreallocation_SeqSBAIJ(Mat Y,Mat X,PetscInt *nnz)
115352768537SHong Zhang {
1154b264fe52SHong Zhang   PetscInt       bs = Y->rmap->bs,mbs = Y->rmap->N/bs;
115552768537SHong Zhang   Mat_SeqSBAIJ   *x = (Mat_SeqSBAIJ*)X->data;
115652768537SHong Zhang   Mat_SeqSBAIJ   *y = (Mat_SeqSBAIJ*)Y->data;
1157b264fe52SHong Zhang   PetscErrorCode ierr;
115852768537SHong Zhang 
115952768537SHong Zhang   PetscFunctionBegin;
116052768537SHong Zhang   /* Set the number of nonzeros in the new matrix */
1161b264fe52SHong Zhang   ierr = MatAXPYGetPreallocation_SeqX_private(mbs,x->i,x->j,y->i,y->j,nnz);CHKERRQ(ierr);
116252768537SHong Zhang   PetscFunctionReturn(0);
116352768537SHong Zhang }
116452768537SHong Zhang 
1165f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
116642ee4b1aSHong Zhang {
116742ee4b1aSHong Zhang   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1168dfbe8321SBarry Smith   PetscErrorCode ierr;
116931ce2d13SHong Zhang   PetscInt       bs=Y->rmap->bs,bs2=bs*bs;
1170e838b9e7SJed Brown   PetscBLASInt   one = 1;
117142ee4b1aSHong Zhang 
117242ee4b1aSHong Zhang   PetscFunctionBegin;
117342ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1174f4df32b1SMatthew Knepley     PetscScalar  alpha = a;
1175c5df96a5SBarry Smith     PetscBLASInt bnz;
1176c5df96a5SBarry Smith     ierr = PetscBLASIntCast(x->nz*bs2,&bnz);CHKERRQ(ierr);
11778b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1178a3fa217bSJose E. Roman     ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);
1179ab784542SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1180ab784542SHong Zhang     ierr = MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
1181ab784542SHong Zhang     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
1182ab784542SHong Zhang     ierr = MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);CHKERRQ(ierr);
118342ee4b1aSHong Zhang   } else {
118452768537SHong Zhang     Mat      B;
118552768537SHong Zhang     PetscInt *nnz;
118652768537SHong Zhang     if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1187f5edf698SHong Zhang     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
118852768537SHong Zhang     ierr = MatGetRowUpperTriangular(Y);CHKERRQ(ierr);
118952768537SHong Zhang     ierr = PetscMalloc1(Y->rmap->N,&nnz);CHKERRQ(ierr);
119052768537SHong Zhang     ierr = MatCreate(PetscObjectComm((PetscObject)Y),&B);CHKERRQ(ierr);
119152768537SHong Zhang     ierr = PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);CHKERRQ(ierr);
119252768537SHong Zhang     ierr = MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);CHKERRQ(ierr);
119352768537SHong Zhang     ierr = MatSetBlockSizesFromMats(B,Y,Y);CHKERRQ(ierr);
11944c7a3774SStefano Zampini     ierr = MatSetType(B,((PetscObject)Y)->type_name);CHKERRQ(ierr);
119552768537SHong Zhang     ierr = MatAXPYGetPreallocation_SeqSBAIJ(Y,X,nnz);CHKERRQ(ierr);
119652768537SHong Zhang     ierr = MatSeqSBAIJSetPreallocation(B,bs,0,nnz);CHKERRQ(ierr);
119752768537SHong Zhang 
119852768537SHong Zhang     ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr);
119952768537SHong Zhang 
120028be2f97SBarry Smith     ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr);
120152768537SHong Zhang     ierr = PetscFree(nnz);CHKERRQ(ierr);
1202f5edf698SHong Zhang     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
120352768537SHong Zhang     ierr = MatRestoreRowUpperTriangular(Y);CHKERRQ(ierr);
120442ee4b1aSHong Zhang   }
120542ee4b1aSHong Zhang   PetscFunctionReturn(0);
120642ee4b1aSHong Zhang }
120742ee4b1aSHong Zhang 
1208ace3abfcSBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1209efcf0fc3SBarry Smith {
1210efcf0fc3SBarry Smith   PetscFunctionBegin;
1211efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1212efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1213efcf0fc3SBarry Smith }
1214efcf0fc3SBarry Smith 
1215ace3abfcSBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool  *flg)
1216efcf0fc3SBarry Smith {
1217efcf0fc3SBarry Smith   PetscFunctionBegin;
1218efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1219efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1220efcf0fc3SBarry Smith }
1221efcf0fc3SBarry Smith 
1222ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1223efcf0fc3SBarry Smith {
1224efcf0fc3SBarry Smith   PetscFunctionBegin;
1225efcf0fc3SBarry Smith   *flg = PETSC_FALSE;
1226efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1227efcf0fc3SBarry Smith }
1228efcf0fc3SBarry Smith 
122999cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
123099cafbc1SBarry Smith {
123199cafbc1SBarry Smith   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
123299cafbc1SBarry Smith   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1233dd6ea824SBarry Smith   MatScalar    *aa = a->a;
123499cafbc1SBarry Smith 
123599cafbc1SBarry Smith   PetscFunctionBegin;
123699cafbc1SBarry Smith   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
123799cafbc1SBarry Smith   PetscFunctionReturn(0);
123899cafbc1SBarry Smith }
123999cafbc1SBarry Smith 
124099cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
124199cafbc1SBarry Smith {
124299cafbc1SBarry Smith   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
124399cafbc1SBarry Smith   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1244dd6ea824SBarry Smith   MatScalar    *aa = a->a;
124599cafbc1SBarry Smith 
124699cafbc1SBarry Smith   PetscFunctionBegin;
124799cafbc1SBarry Smith   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
124899cafbc1SBarry Smith   PetscFunctionReturn(0);
124999cafbc1SBarry Smith }
125099cafbc1SBarry Smith 
12513bededecSBarry Smith PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
12523bededecSBarry Smith {
12533bededecSBarry Smith   Mat_SeqSBAIJ      *baij=(Mat_SeqSBAIJ*)A->data;
12543bededecSBarry Smith   PetscErrorCode    ierr;
12553bededecSBarry Smith   PetscInt          i,j,k,count;
12563bededecSBarry Smith   PetscInt          bs   =A->rmap->bs,bs2=baij->bs2,row,col;
12573bededecSBarry Smith   PetscScalar       zero = 0.0;
12583bededecSBarry Smith   MatScalar         *aa;
12593bededecSBarry Smith   const PetscScalar *xx;
12603bededecSBarry Smith   PetscScalar       *bb;
126156777dd2SBarry Smith   PetscBool         *zeroed,vecs = PETSC_FALSE;
12623bededecSBarry Smith 
12633bededecSBarry Smith   PetscFunctionBegin;
12643bededecSBarry Smith   /* fix right hand side if needed */
12653bededecSBarry Smith   if (x && b) {
12663bededecSBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
12673bededecSBarry Smith     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
126856777dd2SBarry Smith     vecs = PETSC_TRUE;
12693bededecSBarry Smith   }
12703bededecSBarry Smith 
12713bededecSBarry Smith   /* zero the columns */
12721795a4d1SJed Brown   ierr = PetscCalloc1(A->rmap->n,&zeroed);CHKERRQ(ierr);
12733bededecSBarry Smith   for (i=0; i<is_n; i++) {
12743bededecSBarry Smith     if (is_idx[i] < 0 || is_idx[i] >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range",is_idx[i]);
12753bededecSBarry Smith     zeroed[is_idx[i]] = PETSC_TRUE;
12763bededecSBarry Smith   }
127756777dd2SBarry Smith   if (vecs) {
127856777dd2SBarry Smith     for (i=0; i<A->rmap->N; i++) {
127956777dd2SBarry Smith       row = i/bs;
128056777dd2SBarry Smith       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
128156777dd2SBarry Smith         for (k=0; k<bs; k++) {
128256777dd2SBarry Smith           col = bs*baij->j[j] + k;
128356777dd2SBarry Smith           if (col <= i) continue;
128456777dd2SBarry Smith           aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
128526fbe8dcSKarl Rupp           if (!zeroed[i] && zeroed[col]) bb[i]   -= aa[0]*xx[col];
128626fbe8dcSKarl Rupp           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
128756777dd2SBarry Smith         }
128856777dd2SBarry Smith       }
128956777dd2SBarry Smith     }
129026fbe8dcSKarl Rupp     for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
129156777dd2SBarry Smith   }
129256777dd2SBarry Smith 
12933bededecSBarry Smith   for (i=0; i<A->rmap->N; i++) {
12943bededecSBarry Smith     if (!zeroed[i]) {
12953bededecSBarry Smith       row = i/bs;
12963bededecSBarry Smith       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
12973bededecSBarry Smith         for (k=0; k<bs; k++) {
12983bededecSBarry Smith           col = bs*baij->j[j] + k;
12993bededecSBarry Smith           if (zeroed[col]) {
13003bededecSBarry Smith             aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
13013bededecSBarry Smith             aa[0] = 0.0;
13023bededecSBarry Smith           }
13033bededecSBarry Smith         }
13043bededecSBarry Smith       }
13053bededecSBarry Smith     }
13063bededecSBarry Smith   }
13073bededecSBarry Smith   ierr = PetscFree(zeroed);CHKERRQ(ierr);
130856777dd2SBarry Smith   if (vecs) {
130956777dd2SBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
131056777dd2SBarry Smith     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
131156777dd2SBarry Smith   }
13123bededecSBarry Smith 
13133bededecSBarry Smith   /* zero the rows */
13143bededecSBarry Smith   for (i=0; i<is_n; i++) {
13153bededecSBarry Smith     row   = is_idx[i];
13163bededecSBarry Smith     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
13173bededecSBarry Smith     aa    = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
13183bededecSBarry Smith     for (k=0; k<count; k++) {
13193bededecSBarry Smith       aa[0] =  zero;
13203bededecSBarry Smith       aa   += bs;
13213bededecSBarry Smith     }
13223bededecSBarry Smith     if (diag != 0.0) {
13233bededecSBarry Smith       ierr = (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);CHKERRQ(ierr);
13243bededecSBarry Smith     }
13253bededecSBarry Smith   }
13263bededecSBarry Smith   ierr = MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13273bededecSBarry Smith   PetscFunctionReturn(0);
13283bededecSBarry Smith }
13293bededecSBarry Smith 
13307d68702bSBarry Smith PetscErrorCode MatShift_SeqSBAIJ(Mat Y,PetscScalar a)
13317d68702bSBarry Smith {
13327d68702bSBarry Smith   PetscErrorCode ierr;
13337d68702bSBarry Smith   Mat_SeqSBAIJ    *aij = (Mat_SeqSBAIJ*)Y->data;
13347d68702bSBarry Smith 
13357d68702bSBarry Smith   PetscFunctionBegin;
13366f33a894SBarry Smith   if (!Y->preallocated || !aij->nz) {
13377d68702bSBarry Smith     ierr = MatSeqSBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL);CHKERRQ(ierr);
13387d68702bSBarry Smith   }
13397d68702bSBarry Smith   ierr = MatShift_Basic(Y,a);CHKERRQ(ierr);
13407d68702bSBarry Smith   PetscFunctionReturn(0);
13417d68702bSBarry Smith }
13427d68702bSBarry Smith 
134349b5e25fSSatish Balay /* -------------------------------------------------------------------*/
13443964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
134549b5e25fSSatish Balay                                        MatGetRow_SeqSBAIJ,
134649b5e25fSSatish Balay                                        MatRestoreRow_SeqSBAIJ,
134749b5e25fSSatish Balay                                        MatMult_SeqSBAIJ_N,
134897304618SKris Buschelman                                /*  4*/ MatMultAdd_SeqSBAIJ_N,
1349431c96f7SBarry Smith                                        MatMult_SeqSBAIJ_N,       /* transpose versions are same as non-transpose versions */
1350e005ede5SBarry Smith                                        MatMultAdd_SeqSBAIJ_N,
1351db4efbfdSBarry Smith                                        0,
135249b5e25fSSatish Balay                                        0,
135349b5e25fSSatish Balay                                        0,
135497304618SKris Buschelman                                /* 10*/ 0,
135549b5e25fSSatish Balay                                        0,
1356c078aec8SLisandro Dalcin                                        MatCholeskyFactor_SeqSBAIJ,
135741f059aeSBarry Smith                                        MatSOR_SeqSBAIJ,
135849b5e25fSSatish Balay                                        MatTranspose_SeqSBAIJ,
135997304618SKris Buschelman                                /* 15*/ MatGetInfo_SeqSBAIJ,
136049b5e25fSSatish Balay                                        MatEqual_SeqSBAIJ,
136149b5e25fSSatish Balay                                        MatGetDiagonal_SeqSBAIJ,
136249b5e25fSSatish Balay                                        MatDiagonalScale_SeqSBAIJ,
136349b5e25fSSatish Balay                                        MatNorm_SeqSBAIJ,
136497304618SKris Buschelman                                /* 20*/ 0,
136549b5e25fSSatish Balay                                        MatAssemblyEnd_SeqSBAIJ,
136649b5e25fSSatish Balay                                        MatSetOption_SeqSBAIJ,
136749b5e25fSSatish Balay                                        MatZeroEntries_SeqSBAIJ,
1368d519adbfSMatthew Knepley                                /* 24*/ 0,
136949b5e25fSSatish Balay                                        0,
137049b5e25fSSatish Balay                                        0,
1371db4efbfdSBarry Smith                                        0,
1372db4efbfdSBarry Smith                                        0,
13734994cf47SJed Brown                                /* 29*/ MatSetUp_SeqSBAIJ,
1374c464158bSHong Zhang                                        0,
1375db4efbfdSBarry Smith                                        0,
13768c778c55SBarry Smith                                        0,
13778c778c55SBarry Smith                                        0,
1378d519adbfSMatthew Knepley                                /* 34*/ MatDuplicate_SeqSBAIJ,
1379719d5645SBarry Smith                                        0,
1380719d5645SBarry Smith                                        0,
138149b5e25fSSatish Balay                                        0,
1382c84f5b01SHong Zhang                                        MatICCFactor_SeqSBAIJ,
1383d519adbfSMatthew Knepley                                /* 39*/ MatAXPY_SeqSBAIJ,
13847dae84e0SHong Zhang                                        MatCreateSubMatrices_SeqSBAIJ,
138549b5e25fSSatish Balay                                        MatIncreaseOverlap_SeqSBAIJ,
138649b5e25fSSatish Balay                                        MatGetValues_SeqSBAIJ,
13873c896bc6SHong Zhang                                        MatCopy_SeqSBAIJ,
1388d519adbfSMatthew Knepley                                /* 44*/ 0,
138949b5e25fSSatish Balay                                        MatScale_SeqSBAIJ,
13907d68702bSBarry Smith                                        MatShift_SeqSBAIJ,
139149b5e25fSSatish Balay                                        0,
13923bededecSBarry Smith                                        MatZeroRowsColumns_SeqSBAIJ,
1393f73d5cc4SBarry Smith                                /* 49*/ 0,
139449b5e25fSSatish Balay                                        MatGetRowIJ_SeqSBAIJ,
139549b5e25fSSatish Balay                                        MatRestoreRowIJ_SeqSBAIJ,
139649b5e25fSSatish Balay                                        0,
139749b5e25fSSatish Balay                                        0,
1398d519adbfSMatthew Knepley                                /* 54*/ 0,
139949b5e25fSSatish Balay                                        0,
140049b5e25fSSatish Balay                                        0,
140149b5e25fSSatish Balay                                        0,
140249b5e25fSSatish Balay                                        MatSetValuesBlocked_SeqSBAIJ,
14037dae84e0SHong Zhang                                /* 59*/ MatCreateSubMatrix_SeqSBAIJ,
140449b5e25fSSatish Balay                                        0,
140549b5e25fSSatish Balay                                        0,
1406357abbc8SBarry Smith                                        0,
1407d959ec07SHong Zhang                                        0,
1408d519adbfSMatthew Knepley                                /* 64*/ 0,
1409d959ec07SHong Zhang                                        0,
1410d959ec07SHong Zhang                                        0,
1411d959ec07SHong Zhang                                        0,
1412d959ec07SHong Zhang                                        0,
1413d519adbfSMatthew Knepley                                /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
14143e0d88b5SBarry Smith                                        0,
14153e0d88b5SBarry Smith                                        0,
14163e0d88b5SBarry Smith                                        0,
14173e0d88b5SBarry Smith                                        0,
1418d519adbfSMatthew Knepley                                /* 74*/ 0,
14193e0d88b5SBarry Smith                                        0,
14203e0d88b5SBarry Smith                                        0,
14213e0d88b5SBarry Smith                                        0,
14223e0d88b5SBarry Smith                                        0,
1423d519adbfSMatthew Knepley                                /* 79*/ 0,
14243e0d88b5SBarry Smith                                        0,
14253e0d88b5SBarry Smith                                        0,
142697304618SKris Buschelman                                        MatGetInertia_SeqSBAIJ,
14275bba2384SShri Abhyankar                                        MatLoad_SeqSBAIJ,
1428d519adbfSMatthew Knepley                                /* 84*/ MatIsSymmetric_SeqSBAIJ,
1429865e5f61SKris Buschelman                                        MatIsHermitian_SeqSBAIJ,
1430efcf0fc3SBarry Smith                                        MatIsStructurallySymmetric_SeqSBAIJ,
1431865e5f61SKris Buschelman                                        0,
1432865e5f61SKris Buschelman                                        0,
1433d519adbfSMatthew Knepley                                /* 89*/ 0,
1434865e5f61SKris Buschelman                                        0,
1435865e5f61SKris Buschelman                                        0,
1436865e5f61SKris Buschelman                                        0,
1437865e5f61SKris Buschelman                                        0,
1438d519adbfSMatthew Knepley                                /* 94*/ 0,
1439865e5f61SKris Buschelman                                        0,
1440865e5f61SKris Buschelman                                        0,
144199cafbc1SBarry Smith                                        0,
144299cafbc1SBarry Smith                                        0,
1443d519adbfSMatthew Knepley                                /* 99*/ 0,
144499cafbc1SBarry Smith                                        0,
144599cafbc1SBarry Smith                                        0,
144699cafbc1SBarry Smith                                        0,
144799cafbc1SBarry Smith                                        0,
1448d519adbfSMatthew Knepley                                /*104*/ 0,
144999cafbc1SBarry Smith                                        MatRealPart_SeqSBAIJ,
1450f5edf698SHong Zhang                                        MatImaginaryPart_SeqSBAIJ,
1451f5edf698SHong Zhang                                        MatGetRowUpperTriangular_SeqSBAIJ,
14522af78befSBarry Smith                                        MatRestoreRowUpperTriangular_SeqSBAIJ,
1453d519adbfSMatthew Knepley                                /*109*/ 0,
14542af78befSBarry Smith                                        0,
14552af78befSBarry Smith                                        0,
14562af78befSBarry Smith                                        0,
1457547795f9SHong Zhang                                        MatMissingDiagonal_SeqSBAIJ,
1458547795f9SHong Zhang                                /*114*/ 0,
1459547795f9SHong Zhang                                        0,
1460547795f9SHong Zhang                                        0,
1461547795f9SHong Zhang                                        0,
1462547795f9SHong Zhang                                        0,
1463547795f9SHong Zhang                                /*119*/ 0,
1464547795f9SHong Zhang                                        0,
14652f480046SShri Abhyankar                                        0,
14663964eb88SJed Brown                                        0,
14673964eb88SJed Brown                                        0,
14683964eb88SJed Brown                                /*124*/ 0,
14693964eb88SJed Brown                                        0,
14703964eb88SJed Brown                                        0,
14713964eb88SJed Brown                                        0,
14723964eb88SJed Brown                                        0,
14733964eb88SJed Brown                                /*129*/ 0,
14743964eb88SJed Brown                                        0,
14753964eb88SJed Brown                                        0,
14763964eb88SJed Brown                                        0,
14773964eb88SJed Brown                                        0,
14783964eb88SJed Brown                                /*134*/ 0,
14793964eb88SJed Brown                                        0,
14803964eb88SJed Brown                                        0,
14813964eb88SJed Brown                                        0,
14823964eb88SJed Brown                                        0,
148346533700Sstefano_zampini                                /*139*/ MatSetBlockSizes_Default,
1484f9426fe0SMark Adams                                        0,
148559f5e6ceSHong Zhang                                        0,
148659f5e6ceSHong Zhang                                        0,
148759f5e6ceSHong Zhang                                        0,
148859f5e6ceSHong Zhang                                 /*144*/MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ
148999cafbc1SBarry Smith };
1490be1d678aSKris Buschelman 
14917087cfbeSBarry Smith PetscErrorCode  MatStoreValues_SeqSBAIJ(Mat mat)
149249b5e25fSSatish Balay {
14934afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1494d0f46423SBarry Smith   PetscInt       nz   = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1495dfbe8321SBarry Smith   PetscErrorCode ierr;
149649b5e25fSSatish Balay 
149749b5e25fSSatish Balay   PetscFunctionBegin;
1498e7e72b3dSBarry Smith   if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
149949b5e25fSSatish Balay 
150049b5e25fSSatish Balay   /* allocate space for values if not already there */
150149b5e25fSSatish Balay   if (!aij->saved_values) {
1502854ce69bSBarry Smith     ierr = PetscMalloc1(nz+1,&aij->saved_values);CHKERRQ(ierr);
150349b5e25fSSatish Balay   }
150449b5e25fSSatish Balay 
150549b5e25fSSatish Balay   /* copy values over */
1506580bdb30SBarry Smith   ierr = PetscArraycpy(aij->saved_values,aij->a,nz);CHKERRQ(ierr);
150749b5e25fSSatish Balay   PetscFunctionReturn(0);
150849b5e25fSSatish Balay }
150949b5e25fSSatish Balay 
15107087cfbeSBarry Smith PetscErrorCode  MatRetrieveValues_SeqSBAIJ(Mat mat)
151149b5e25fSSatish Balay {
15124afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
15136849ba73SBarry Smith   PetscErrorCode ierr;
1514d0f46423SBarry Smith   PetscInt       nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
151549b5e25fSSatish Balay 
151649b5e25fSSatish Balay   PetscFunctionBegin;
1517e7e72b3dSBarry Smith   if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1518e7e72b3dSBarry Smith   if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
151949b5e25fSSatish Balay 
152049b5e25fSSatish Balay   /* copy values over */
1521580bdb30SBarry Smith   ierr = PetscArraycpy(aij->a,aij->saved_values,nz);CHKERRQ(ierr);
152249b5e25fSSatish Balay   PetscFunctionReturn(0);
152349b5e25fSSatish Balay }
152449b5e25fSSatish Balay 
1525367daffbSBarry Smith static PetscErrorCode  MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
152649b5e25fSSatish Balay {
1527c464158bSHong Zhang   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
15286849ba73SBarry Smith   PetscErrorCode ierr;
15294dcd73b1SHong Zhang   PetscInt       i,mbs,nbs,bs2;
15302576faa2SJed Brown   PetscBool      skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;
153149b5e25fSSatish Balay 
153249b5e25fSSatish Balay   PetscFunctionBegin;
15332576faa2SJed Brown   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1534db4efbfdSBarry Smith 
153533d57670SJed Brown   ierr = MatSetBlockSize(B,PetscAbs(bs));CHKERRQ(ierr);
153626283091SBarry Smith   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
153726283091SBarry Smith   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
1538e02043d6SBarry Smith   ierr = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr);
1539899cda47SBarry Smith 
154021940c7eSstefano_zampini   B->preallocated = PETSC_TRUE;
154121940c7eSstefano_zampini 
1542d0f46423SBarry Smith   mbs = B->rmap->N/bs;
15434dcd73b1SHong Zhang   nbs = B->cmap->n/bs;
154449b5e25fSSatish Balay   bs2 = bs*bs;
154549b5e25fSSatish Balay 
15464dcd73b1SHong Zhang   if (mbs*bs != B->rmap->N || nbs*bs!=B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
154749b5e25fSSatish Balay 
1548ab93d7beSBarry Smith   if (nz == MAT_SKIP_ALLOCATION) {
1549ab93d7beSBarry Smith     skipallocation = PETSC_TRUE;
1550ab93d7beSBarry Smith     nz             = 0;
1551ab93d7beSBarry Smith   }
1552ab93d7beSBarry Smith 
1553435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1554e32f2f54SBarry Smith   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
155549b5e25fSSatish Balay   if (nnz) {
155649b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
1557e32f2f54SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
1558de64b629SHong Zhang       if (nnz[i] > nbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D block rowlength %D",i,nnz[i],nbs);
155949b5e25fSSatish Balay     }
156049b5e25fSSatish Balay   }
156149b5e25fSSatish Balay 
1562db4efbfdSBarry Smith   B->ops->mult             = MatMult_SeqSBAIJ_N;
1563db4efbfdSBarry Smith   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1564db4efbfdSBarry Smith   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1565db4efbfdSBarry Smith   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;
156626fbe8dcSKarl Rupp 
1567c5929fdfSBarry Smith   ierr  = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);CHKERRQ(ierr);
156849b5e25fSSatish Balay   if (!flg) {
156949b5e25fSSatish Balay     switch (bs) {
157049b5e25fSSatish Balay     case 1:
157149b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_1;
157249b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1573431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1574431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
157549b5e25fSSatish Balay       break;
157649b5e25fSSatish Balay     case 2:
157749b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_2;
157849b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1579431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1580431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
158149b5e25fSSatish Balay       break;
158249b5e25fSSatish Balay     case 3:
158349b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_3;
158449b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1585431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1586431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
158749b5e25fSSatish Balay       break;
158849b5e25fSSatish Balay     case 4:
158949b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_4;
159049b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1591431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1592431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
159349b5e25fSSatish Balay       break;
159449b5e25fSSatish Balay     case 5:
159549b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_5;
159649b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1597431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1598431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
159949b5e25fSSatish Balay       break;
160049b5e25fSSatish Balay     case 6:
160149b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_6;
160249b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1603431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1604431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
160549b5e25fSSatish Balay       break;
160649b5e25fSSatish Balay     case 7:
1607de53e5efSHong Zhang       B->ops->mult             = MatMult_SeqSBAIJ_7;
160849b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1609431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1610431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
161149b5e25fSSatish Balay       break;
161249b5e25fSSatish Balay     }
161349b5e25fSSatish Balay   }
161449b5e25fSSatish Balay 
161549b5e25fSSatish Balay   b->mbs = mbs;
16164dcd73b1SHong Zhang   b->nbs = nbs;
1617ab93d7beSBarry Smith   if (!skipallocation) {
16182ee49352SLisandro Dalcin     if (!b->imax) {
1619dcca6d9dSJed Brown       ierr = PetscMalloc2(mbs,&b->imax,mbs,&b->ilen);CHKERRQ(ierr);
162026fbe8dcSKarl Rupp 
1621c760cd28SBarry Smith       b->free_imax_ilen = PETSC_TRUE;
162226fbe8dcSKarl Rupp 
16233bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)B,2*mbs*sizeof(PetscInt));CHKERRQ(ierr);
16242ee49352SLisandro Dalcin     }
162549b5e25fSSatish Balay     if (!nnz) {
1626435da068SBarry Smith       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
162749b5e25fSSatish Balay       else if (nz <= 0) nz = 1;
16285d2a9ed1SStefano Zampini       nz = PetscMin(nbs,nz);
162926fbe8dcSKarl Rupp       for (i=0; i<mbs; i++) b->imax[i] = nz;
1630153ea458SHong Zhang       nz = nz*mbs; /* total nz */
163149b5e25fSSatish Balay     } else {
1632c73702f5SBarry Smith       PetscInt64 nz64 = 0;
1633c73702f5SBarry Smith       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz64 += nnz[i];}
1634c73702f5SBarry Smith       ierr = PetscIntCast(nz64,&nz);CHKERRQ(ierr);
163549b5e25fSSatish Balay     }
16362ee49352SLisandro Dalcin     /* b->ilen will count nonzeros in each block row so far. */
163726fbe8dcSKarl Rupp     for (i=0; i<mbs; i++) b->ilen[i] = 0;
16386c6c5352SBarry Smith     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
163949b5e25fSSatish Balay 
164049b5e25fSSatish Balay     /* allocate the matrix space */
16412ee49352SLisandro Dalcin     ierr = MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);CHKERRQ(ierr);
1642dcca6d9dSJed Brown     ierr = PetscMalloc3(bs2*nz,&b->a,nz,&b->j,B->rmap->N+1,&b->i);CHKERRQ(ierr);
16433bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
1644580bdb30SBarry Smith     ierr = PetscArrayzero(b->a,nz*bs2);CHKERRQ(ierr);
1645580bdb30SBarry Smith     ierr = PetscArrayzero(b->j,nz);CHKERRQ(ierr);
164626fbe8dcSKarl Rupp 
164749b5e25fSSatish Balay     b->singlemalloc = PETSC_TRUE;
164849b5e25fSSatish Balay 
164949b5e25fSSatish Balay     /* pointer to beginning of each row */
1650e60cf9a0SBarry Smith     b->i[0] = 0;
165126fbe8dcSKarl Rupp     for (i=1; i<mbs+1; i++) b->i[i] = b->i[i-1] + b->imax[i-1];
165226fbe8dcSKarl Rupp 
1653e6b907acSBarry Smith     b->free_a  = PETSC_TRUE;
1654e6b907acSBarry Smith     b->free_ij = PETSC_TRUE;
1655e811da20SHong Zhang   } else {
1656e6b907acSBarry Smith     b->free_a  = PETSC_FALSE;
1657e6b907acSBarry Smith     b->free_ij = PETSC_FALSE;
1658ab93d7beSBarry Smith   }
165949b5e25fSSatish Balay 
166049b5e25fSSatish Balay   b->bs2     = bs2;
16616c6c5352SBarry Smith   b->nz      = 0;
1662b32cb4a7SJed Brown   b->maxnz   = nz;
166316cdd363SHong Zhang   b->inew    = 0;
166416cdd363SHong Zhang   b->jnew    = 0;
166516cdd363SHong Zhang   b->anew    = 0;
166616cdd363SHong Zhang   b->a2anew  = 0;
16671a3463dfSHong Zhang   b->permute = PETSC_FALSE;
1668cb7b82ddSBarry Smith 
1669cb7b82ddSBarry Smith   B->was_assembled = PETSC_FALSE;
1670cb7b82ddSBarry Smith   B->assembled     = PETSC_FALSE;
16712576faa2SJed Brown   if (realalloc) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);}
1672c464158bSHong Zhang   PetscFunctionReturn(0);
1673c464158bSHong Zhang }
1674153ea458SHong Zhang 
167538f409ebSLisandro Dalcin PetscErrorCode MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[], const PetscScalar V[])
167638f409ebSLisandro Dalcin {
1677*0cd7f59aSBarry Smith   PetscInt       i,j,m,nz,anz, nz_max=0,*nnz;
167838f409ebSLisandro Dalcin   PetscScalar    *values=0;
167938f409ebSLisandro Dalcin   PetscBool      roworiented = ((Mat_SeqSBAIJ*)B->data)->roworiented;
168038f409ebSLisandro Dalcin   PetscErrorCode ierr;
1681*0cd7f59aSBarry Smith 
168238f409ebSLisandro Dalcin   PetscFunctionBegin;
168338f409ebSLisandro Dalcin   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
168438f409ebSLisandro Dalcin   ierr   = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr);
168538f409ebSLisandro Dalcin   ierr   = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr);
168638f409ebSLisandro Dalcin   ierr   = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
168738f409ebSLisandro Dalcin   ierr   = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
168838f409ebSLisandro Dalcin   ierr   = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr);
168938f409ebSLisandro Dalcin   m      = B->rmap->n/bs;
169038f409ebSLisandro Dalcin 
169138f409ebSLisandro Dalcin   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1692854ce69bSBarry Smith   ierr = PetscMalloc1(m+1,&nnz);CHKERRQ(ierr);
169338f409ebSLisandro Dalcin   for (i=0; i<m; i++) {
169438f409ebSLisandro Dalcin     nz = ii[i+1] - ii[i];
169538f409ebSLisandro Dalcin     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D has a negative number of columns %D",i,nz);
1696*0cd7f59aSBarry Smith     anz = 0;
1697*0cd7f59aSBarry Smith     for (j=0; j<nz; j++) {
1698*0cd7f59aSBarry Smith       /* count only values on the diagonal or above */
1699*0cd7f59aSBarry Smith       if (jj[ii[i] + j] >= i) {
1700*0cd7f59aSBarry Smith         anz = nz - j;
1701*0cd7f59aSBarry Smith         break;
1702*0cd7f59aSBarry Smith       }
1703*0cd7f59aSBarry Smith     }
1704*0cd7f59aSBarry Smith     nz_max = PetscMax(nz_max,anz);
1705*0cd7f59aSBarry Smith     nnz[i] = anz;
170638f409ebSLisandro Dalcin   }
170738f409ebSLisandro Dalcin   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,nnz);CHKERRQ(ierr);
170838f409ebSLisandro Dalcin   ierr = PetscFree(nnz);CHKERRQ(ierr);
170938f409ebSLisandro Dalcin 
171038f409ebSLisandro Dalcin   values = (PetscScalar*)V;
171138f409ebSLisandro Dalcin   if (!values) {
17121795a4d1SJed Brown     ierr = PetscCalloc1(bs*bs*nz_max,&values);CHKERRQ(ierr);
171338f409ebSLisandro Dalcin   }
171438f409ebSLisandro Dalcin   for (i=0; i<m; i++) {
171538f409ebSLisandro Dalcin     PetscInt          ncols  = ii[i+1] - ii[i];
171638f409ebSLisandro Dalcin     const PetscInt    *icols = jj + ii[i];
171738f409ebSLisandro Dalcin     if (!roworiented || bs == 1) {
171838f409ebSLisandro Dalcin       const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
171938f409ebSLisandro Dalcin       ierr = MatSetValuesBlocked_SeqSBAIJ(B,1,&i,ncols,icols,svals,INSERT_VALUES);CHKERRQ(ierr);
172038f409ebSLisandro Dalcin     } else {
172138f409ebSLisandro Dalcin       for (j=0; j<ncols; j++) {
172238f409ebSLisandro Dalcin         const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
172338f409ebSLisandro Dalcin         ierr = MatSetValuesBlocked_SeqSBAIJ(B,1,&i,1,&icols[j],svals,INSERT_VALUES);CHKERRQ(ierr);
172438f409ebSLisandro Dalcin       }
172538f409ebSLisandro Dalcin     }
172638f409ebSLisandro Dalcin   }
172738f409ebSLisandro Dalcin   if (!V) { ierr = PetscFree(values);CHKERRQ(ierr); }
172838f409ebSLisandro Dalcin   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
172938f409ebSLisandro Dalcin   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
173038f409ebSLisandro Dalcin   ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
173138f409ebSLisandro Dalcin   PetscFunctionReturn(0);
173238f409ebSLisandro Dalcin }
173338f409ebSLisandro Dalcin 
1734db4efbfdSBarry Smith /*
1735db4efbfdSBarry Smith    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1736db4efbfdSBarry Smith */
1737ace3abfcSBarry Smith PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1738db4efbfdSBarry Smith {
1739db4efbfdSBarry Smith   PetscErrorCode ierr;
1740ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1741db4efbfdSBarry Smith   PetscInt       bs  = B->rmap->bs;
1742db4efbfdSBarry Smith 
1743db4efbfdSBarry Smith   PetscFunctionBegin;
1744c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);CHKERRQ(ierr);
1745db4efbfdSBarry Smith   if (flg) bs = 8;
1746db4efbfdSBarry Smith 
1747db4efbfdSBarry Smith   if (!natural) {
1748db4efbfdSBarry Smith     switch (bs) {
1749db4efbfdSBarry Smith     case 1:
1750d595f711SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1751db4efbfdSBarry Smith       break;
1752db4efbfdSBarry Smith     case 2:
1753db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1754db4efbfdSBarry Smith       break;
1755db4efbfdSBarry Smith     case 3:
1756db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1757db4efbfdSBarry Smith       break;
1758db4efbfdSBarry Smith     case 4:
1759db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1760db4efbfdSBarry Smith       break;
1761db4efbfdSBarry Smith     case 5:
1762db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1763db4efbfdSBarry Smith       break;
1764db4efbfdSBarry Smith     case 6:
1765db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1766db4efbfdSBarry Smith       break;
1767db4efbfdSBarry Smith     case 7:
1768db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1769db4efbfdSBarry Smith       break;
1770db4efbfdSBarry Smith     default:
1771db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1772db4efbfdSBarry Smith       break;
1773db4efbfdSBarry Smith     }
1774db4efbfdSBarry Smith   } else {
1775db4efbfdSBarry Smith     switch (bs) {
1776db4efbfdSBarry Smith     case 1:
1777d595f711SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1778db4efbfdSBarry Smith       break;
1779db4efbfdSBarry Smith     case 2:
1780db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1781db4efbfdSBarry Smith       break;
1782db4efbfdSBarry Smith     case 3:
1783db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1784db4efbfdSBarry Smith       break;
1785db4efbfdSBarry Smith     case 4:
1786db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1787db4efbfdSBarry Smith       break;
1788db4efbfdSBarry Smith     case 5:
1789db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1790db4efbfdSBarry Smith       break;
1791db4efbfdSBarry Smith     case 6:
1792db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1793db4efbfdSBarry Smith       break;
1794db4efbfdSBarry Smith     case 7:
1795db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1796db4efbfdSBarry Smith       break;
1797db4efbfdSBarry Smith     default:
1798db4efbfdSBarry Smith       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1799db4efbfdSBarry Smith       break;
1800db4efbfdSBarry Smith     }
1801db4efbfdSBarry Smith   }
1802db4efbfdSBarry Smith   PetscFunctionReturn(0);
1803db4efbfdSBarry Smith }
1804db4efbfdSBarry Smith 
1805cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1806cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1807d769727bSBarry Smith 
1808cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
18095c9eb25fSBarry Smith {
1810d0f46423SBarry Smith   PetscInt       n = A->rmap->n;
18115c9eb25fSBarry Smith   PetscErrorCode ierr;
18125c9eb25fSBarry Smith 
18135c9eb25fSBarry Smith   PetscFunctionBegin;
18140e92d65fSHong Zhang #if defined(PETSC_USE_COMPLEX)
18150e92d65fSHong Zhang   if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
18160e92d65fSHong Zhang #endif
1817ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
18185c9eb25fSBarry Smith   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
18195c9eb25fSBarry Smith   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
18205c9eb25fSBarry Smith     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
18210298fd71SBarry Smith     ierr = MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
182226fbe8dcSKarl Rupp 
18237b056e98SHong Zhang     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1824c6d0d4f0SHong Zhang     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
1825e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
182600c67f3bSHong Zhang 
1827d5f3da31SBarry Smith   (*B)->factortype = ftype;
182800c67f3bSHong Zhang   ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr);
182900c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);CHKERRQ(ierr);
18305c9eb25fSBarry Smith   PetscFunctionReturn(0);
18315c9eb25fSBarry Smith }
18325c9eb25fSBarry Smith 
18338397e458SBarry Smith /*@C
18348397e458SBarry Smith    MatSeqSBAIJGetArray - gives access to the array where the data for a MATSEQSBAIJ matrix is stored
18358397e458SBarry Smith 
18368397e458SBarry Smith    Not Collective
18378397e458SBarry Smith 
18388397e458SBarry Smith    Input Parameter:
18398397e458SBarry Smith .  mat - a MATSEQSBAIJ matrix
18408397e458SBarry Smith 
18418397e458SBarry Smith    Output Parameter:
18428397e458SBarry Smith .   array - pointer to the data
18438397e458SBarry Smith 
18448397e458SBarry Smith    Level: intermediate
18458397e458SBarry Smith 
18468397e458SBarry Smith .seealso: MatSeqSBAIJRestoreArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray()
18478397e458SBarry Smith @*/
18488397e458SBarry Smith PetscErrorCode  MatSeqSBAIJGetArray(Mat A,PetscScalar **array)
18498397e458SBarry Smith {
18508397e458SBarry Smith   PetscErrorCode ierr;
18518397e458SBarry Smith 
18528397e458SBarry Smith   PetscFunctionBegin;
18538397e458SBarry Smith   ierr = PetscUseMethod(A,"MatSeqSBAIJGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr);
18548397e458SBarry Smith   PetscFunctionReturn(0);
18558397e458SBarry Smith }
18568397e458SBarry Smith 
18578397e458SBarry Smith /*@C
18588397e458SBarry Smith    MatSeqSBAIJRestoreArray - returns access to the array where the data for a MATSEQSBAIJ matrix is stored obtained by MatSeqSBAIJGetArray()
18598397e458SBarry Smith 
18608397e458SBarry Smith    Not Collective
18618397e458SBarry Smith 
18628397e458SBarry Smith    Input Parameters:
1863a2b725a8SWilliam Gropp +  mat - a MATSEQSBAIJ matrix
1864a2b725a8SWilliam Gropp -  array - pointer to the data
18658397e458SBarry Smith 
18668397e458SBarry Smith    Level: intermediate
18678397e458SBarry Smith 
18688397e458SBarry Smith .seealso: MatSeqSBAIJGetArray(), MatSeqAIJGetArray(), MatSeqAIJRestoreArray()
18698397e458SBarry Smith @*/
18708397e458SBarry Smith PetscErrorCode  MatSeqSBAIJRestoreArray(Mat A,PetscScalar **array)
18718397e458SBarry Smith {
18728397e458SBarry Smith   PetscErrorCode ierr;
18738397e458SBarry Smith 
18748397e458SBarry Smith   PetscFunctionBegin;
18758397e458SBarry Smith   ierr = PetscUseMethod(A,"MatSeqSBAIJRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr);
18768397e458SBarry Smith   PetscFunctionReturn(0);
18778397e458SBarry Smith }
18788397e458SBarry Smith 
18790bad9183SKris Buschelman /*MC
1880fafad747SKris Buschelman   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
18810bad9183SKris Buschelman   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
18820bad9183SKris Buschelman 
1883828413b8SBarry Smith   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
188471dad5bbSBarry Smith   can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd()
1885828413b8SBarry Smith 
18860bad9183SKris Buschelman   Options Database Keys:
18870bad9183SKris Buschelman   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
18880bad9183SKris Buschelman 
188995452b02SPatrick Sanan   Notes:
189095452b02SPatrick Sanan     By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
189171dad5bbSBarry Smith      stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
189271dad5bbSBarry Smith      the options database -mat_ignore_lower_triangular false it will generate an error if you try to set a value in the lower triangular portion.
189371dad5bbSBarry Smith 
189471dad5bbSBarry Smith 
18950bad9183SKris Buschelman   Level: beginner
18960bad9183SKris Buschelman 
18970bad9183SKris Buschelman   .seealso: MatCreateSeqSBAIJ
18980bad9183SKris Buschelman M*/
18998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1900a23d5eceSKris Buschelman {
1901a23d5eceSKris Buschelman   Mat_SeqSBAIJ   *b;
1902dfbe8321SBarry Smith   PetscErrorCode ierr;
190313f74950SBarry Smith   PetscMPIInt    size;
1904ace3abfcSBarry Smith   PetscBool      no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;
1905a23d5eceSKris Buschelman 
1906a23d5eceSKris Buschelman   PetscFunctionBegin;
1907ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr);
1908e32f2f54SBarry Smith   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1909a23d5eceSKris Buschelman 
1910b00a9115SJed Brown   ierr    = PetscNewLog(B,&b);CHKERRQ(ierr);
1911a23d5eceSKris Buschelman   B->data = (void*)b;
1912a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
191326fbe8dcSKarl Rupp 
1914a23d5eceSKris Buschelman   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1915a23d5eceSKris Buschelman   B->ops->view       = MatView_SeqSBAIJ;
1916a23d5eceSKris Buschelman   b->row             = 0;
1917a23d5eceSKris Buschelman   b->icol            = 0;
1918a23d5eceSKris Buschelman   b->reallocs        = 0;
1919a23d5eceSKris Buschelman   b->saved_values    = 0;
19200def2e27SBarry Smith   b->inode.limit     = 5;
19210def2e27SBarry Smith   b->inode.max_limit = 5;
1922a23d5eceSKris Buschelman 
1923a23d5eceSKris Buschelman   b->roworiented        = PETSC_TRUE;
1924a23d5eceSKris Buschelman   b->nonew              = 0;
1925a23d5eceSKris Buschelman   b->diag               = 0;
1926a23d5eceSKris Buschelman   b->solve_work         = 0;
1927a23d5eceSKris Buschelman   b->mult_work          = 0;
1928a23d5eceSKris Buschelman   B->spptr              = 0;
1929f2cbd3d5SJed Brown   B->info.nz_unneeded   = (PetscReal)b->maxnz*b->bs2;
1930a9817697SBarry Smith   b->keepnonzeropattern = PETSC_FALSE;
1931a23d5eceSKris Buschelman 
1932a23d5eceSKris Buschelman   b->inew    = 0;
1933a23d5eceSKris Buschelman   b->jnew    = 0;
1934a23d5eceSKris Buschelman   b->anew    = 0;
1935a23d5eceSKris Buschelman   b->a2anew  = 0;
1936a23d5eceSKris Buschelman   b->permute = PETSC_FALSE;
1937a23d5eceSKris Buschelman 
193871dad5bbSBarry Smith   b->ignore_ltriangular = PETSC_TRUE;
193926fbe8dcSKarl Rupp 
1940c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,NULL);CHKERRQ(ierr);
1941941593c8SHong Zhang 
1942f5edf698SHong Zhang   b->getrow_utriangular = PETSC_FALSE;
194326fbe8dcSKarl Rupp 
1944c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)B)->options,((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,NULL);CHKERRQ(ierr);
1945f5edf698SHong Zhang 
19468397e458SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJGetArray_C",MatSeqSBAIJGetArray_SeqSBAIJ);CHKERRQ(ierr);
19478397e458SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJRestoreArray_C",MatSeqSBAIJRestoreArray_SeqSBAIJ);CHKERRQ(ierr);
1948bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1949bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1950bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1951bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1952bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1953bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
195438f409ebSLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocationCSR_C",MatSeqSBAIJSetPreallocationCSR_SeqSBAIJ);CHKERRQ(ierr);
19556214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
19566214f412SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_elemental_C",MatConvert_SeqSBAIJ_Elemental);CHKERRQ(ierr);
19576214f412SHong Zhang #endif
195823ce1328SBarry Smith 
195923ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
196023ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
196123ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
196223ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
19639899f194SHong Zhang   B->symmetric_eternal          = PETSC_TRUE;
196426fbe8dcSKarl Rupp 
196513647f61SHong Zhang   B->hermitian                  = PETSC_FALSE;
196613647f61SHong Zhang   B->hermitian_set              = PETSC_FALSE;
196713647f61SHong Zhang 
196817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr);
19690def2e27SBarry Smith 
1970ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");CHKERRQ(ierr);
19710298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);CHKERRQ(ierr);
197226fbe8dcSKarl Rupp   if (no_unroll) {
197326fbe8dcSKarl Rupp     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);
197426fbe8dcSKarl Rupp   }
19750298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);CHKERRQ(ierr);
197626fbe8dcSKarl Rupp   if (no_inode) {
197726fbe8dcSKarl Rupp     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);
197826fbe8dcSKarl Rupp   }
19790298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);CHKERRQ(ierr);
19800def2e27SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1981ace3abfcSBarry Smith   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
19820def2e27SBarry Smith   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1983a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1984a23d5eceSKris Buschelman }
1985a23d5eceSKris Buschelman 
1986a23d5eceSKris Buschelman /*@C
1987a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1988a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1989a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1990a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1991a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1992a23d5eceSKris Buschelman 
1993a23d5eceSKris Buschelman    Collective on Mat
1994a23d5eceSKris Buschelman 
1995a23d5eceSKris Buschelman    Input Parameters:
19961c4f3114SJed Brown +  B - the symmetric matrix
1997bb7ae925SBarry Smith .  bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
1998bb7ae925SBarry Smith           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
1999a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
2000a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
20010298fd71SBarry Smith          diagonal portion of each block (possibly different for each block row) or NULL
2002a23d5eceSKris Buschelman 
2003a23d5eceSKris Buschelman    Options Database Keys:
2004a2b725a8SWilliam Gropp +   -mat_no_unroll - uses code that does not unroll the loops in the
2005a23d5eceSKris Buschelman                      block calculations (much slower)
2006a2b725a8SWilliam Gropp -   -mat_block_size - size of the blocks to use (only works if a negative bs is passed in
2007a23d5eceSKris Buschelman 
2008a23d5eceSKris Buschelman    Level: intermediate
2009a23d5eceSKris Buschelman 
2010a23d5eceSKris Buschelman    Notes:
2011a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
20120298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2013a7f22e61SSatish Balay    allocation.  See Users-Manual: ch_mat for details.
2014a23d5eceSKris Buschelman 
2015aa95bbe8SBarry Smith    You can call MatGetInfo() to get information on how effective the preallocation was;
2016aa95bbe8SBarry Smith    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2017aa95bbe8SBarry Smith    You can also run with the option -info and look for messages with the string
2018aa95bbe8SBarry Smith    malloc in them to see if additional memory allocation was needed.
2019aa95bbe8SBarry Smith 
202049a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
202149a6f317SBarry Smith 
202249a6f317SBarry Smith 
202369b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2024a23d5eceSKris Buschelman @*/
20257087cfbeSBarry Smith PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
202613f74950SBarry Smith {
20274ac538c5SBarry Smith   PetscErrorCode ierr;
2028a23d5eceSKris Buschelman 
2029a23d5eceSKris Buschelman   PetscFunctionBegin;
20306ba663aaSJed Brown   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
20316ba663aaSJed Brown   PetscValidType(B,1);
20326ba663aaSJed Brown   PetscValidLogicalCollectiveInt(B,bs,2);
20334ac538c5SBarry Smith   ierr = PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr);
2034a23d5eceSKris Buschelman   PetscFunctionReturn(0);
2035a23d5eceSKris Buschelman }
203649b5e25fSSatish Balay 
203738f409ebSLisandro Dalcin /*@C
2038664954b6SBarry Smith    MatSeqSBAIJSetPreallocationCSR - Creates a sparse parallel matrix in SBAIJ format using the given nonzero structure and (optional) numerical values
203938f409ebSLisandro Dalcin 
204038f409ebSLisandro Dalcin    Input Parameters:
20411c4f3114SJed Brown +  B - the matrix
2042eab78319SHong Zhang .  bs - size of block, the blocks are ALWAYS square.
204338f409ebSLisandro Dalcin .  i - the indices into j for the start of each local row (starts with zero)
204438f409ebSLisandro Dalcin .  j - the column indices for each local row (starts with zero) these must be sorted for each row
204538f409ebSLisandro Dalcin -  v - optional values in the matrix
204638f409ebSLisandro Dalcin 
2047664954b6SBarry Smith    Level: advanced
204838f409ebSLisandro Dalcin 
204938f409ebSLisandro Dalcin    Notes:
205038f409ebSLisandro Dalcin    The order of the entries in values is specified by the MatOption MAT_ROW_ORIENTED.  For example, C programs
205138f409ebSLisandro Dalcin    may want to use the default MAT_ROW_ORIENTED=PETSC_TRUE and use an array v[nnz][bs][bs] where the second index is
205238f409ebSLisandro Dalcin    over rows within a block and the last index is over columns within a block row.  Fortran programs will likely set
205338f409ebSLisandro Dalcin    MAT_ROW_ORIENTED=PETSC_FALSE and use a Fortran array v(bs,bs,nnz) in which the first index is over rows within a
205438f409ebSLisandro Dalcin    block column and the second index is over columns within a block.
205538f409ebSLisandro Dalcin 
2056*0cd7f59aSBarry Smith    Any entries below the diagaonl in are ignored
2057*0cd7f59aSBarry Smith 
2058*0cd7f59aSBarry Smith    Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
2059*0cd7f59aSBarry Smith    and usually the numerical values as well
2060664954b6SBarry Smith 
206138f409ebSLisandro Dalcin .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValuesBlocked(), MatSeqSBAIJSetPreallocation(), MATSEQSBAIJ
206238f409ebSLisandro Dalcin @*/
206338f409ebSLisandro Dalcin PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
206438f409ebSLisandro Dalcin {
206538f409ebSLisandro Dalcin   PetscErrorCode ierr;
206638f409ebSLisandro Dalcin 
206738f409ebSLisandro Dalcin   PetscFunctionBegin;
206838f409ebSLisandro Dalcin   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
206938f409ebSLisandro Dalcin   PetscValidType(B,1);
207038f409ebSLisandro Dalcin   PetscValidLogicalCollectiveInt(B,bs,2);
207138f409ebSLisandro Dalcin   ierr = PetscTryMethod(B,"MatSeqSBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));CHKERRQ(ierr);
207238f409ebSLisandro Dalcin   PetscFunctionReturn(0);
207338f409ebSLisandro Dalcin }
207438f409ebSLisandro Dalcin 
2075c464158bSHong Zhang /*@C
2076c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
2077c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
2078c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
2079c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
2080c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
208149b5e25fSSatish Balay 
2082d083f849SBarry Smith    Collective
2083c464158bSHong Zhang 
2084c464158bSHong Zhang    Input Parameters:
2085c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
2086bb7ae925SBarry Smith .  bs - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2087bb7ae925SBarry Smith           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2088c464158bSHong Zhang .  m - number of rows, or number of columns
2089c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
2090744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
20910298fd71SBarry Smith          diagonal portion of each block (possibly different for each block row) or NULL
2092c464158bSHong Zhang 
2093c464158bSHong Zhang    Output Parameter:
2094c464158bSHong Zhang .  A - the symmetric matrix
2095c464158bSHong Zhang 
2096c464158bSHong Zhang    Options Database Keys:
2097a2b725a8SWilliam Gropp +   -mat_no_unroll - uses code that does not unroll the loops in the
2098c464158bSHong Zhang                      block calculations (much slower)
2099a2b725a8SWilliam Gropp -   -mat_block_size - size of the blocks to use
2100c464158bSHong Zhang 
2101c464158bSHong Zhang    Level: intermediate
2102c464158bSHong Zhang 
2103175b88e8SBarry Smith    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2104f6f02116SRichard Tran Mills    MatXXXXSetPreallocation() paradigm instead of this routine directly.
2105175b88e8SBarry Smith    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2106175b88e8SBarry Smith 
2107c464158bSHong Zhang    Notes:
21086d6d819aSHong Zhang    The number of rows and columns must be divisible by blocksize.
21096d6d819aSHong Zhang    This matrix type does not support complex Hermitian operation.
2110c464158bSHong Zhang 
2111c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
21120298fd71SBarry Smith    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2113a7f22e61SSatish Balay    allocation.  See Users-Manual: ch_mat for details.
2114c464158bSHong Zhang 
211549a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
211649a6f317SBarry Smith 
211769b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2118c464158bSHong Zhang @*/
21197087cfbeSBarry Smith PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2120c464158bSHong Zhang {
2121dfbe8321SBarry Smith   PetscErrorCode ierr;
2122c464158bSHong Zhang 
2123c464158bSHong Zhang   PetscFunctionBegin;
2124f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2125f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
2126c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
2127367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr);
212849b5e25fSSatish Balay   PetscFunctionReturn(0);
212949b5e25fSSatish Balay }
213049b5e25fSSatish Balay 
2131dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
213249b5e25fSSatish Balay {
213349b5e25fSSatish Balay   Mat            C;
213449b5e25fSSatish Balay   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
21356849ba73SBarry Smith   PetscErrorCode ierr;
2136b40805acSSatish Balay   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
213749b5e25fSSatish Balay 
213849b5e25fSSatish Balay   PetscFunctionBegin;
2139e32f2f54SBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix");
214049b5e25fSSatish Balay 
214149b5e25fSSatish Balay   *B   = 0;
2142ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
2143d0f46423SBarry Smith   ierr = MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);CHKERRQ(ierr);
21444c7a3774SStefano Zampini   ierr = MatSetBlockSizesFromMats(C,A,A);CHKERRQ(ierr);
21458e9a0fb8SHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
2146692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
2147692f9cbeSHong Zhang 
2148273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
2149d5f3da31SBarry Smith   C->factortype         = A->factortype;
215049b5e25fSSatish Balay   c->row                = 0;
215149b5e25fSSatish Balay   c->icol               = 0;
215249b5e25fSSatish Balay   c->saved_values       = 0;
2153a9817697SBarry Smith   c->keepnonzeropattern = a->keepnonzeropattern;
215449b5e25fSSatish Balay   C->assembled          = PETSC_TRUE;
215549b5e25fSSatish Balay 
21561e1e43feSBarry Smith   ierr   = PetscLayoutReference(A->rmap,&C->rmap);CHKERRQ(ierr);
21571e1e43feSBarry Smith   ierr   = PetscLayoutReference(A->cmap,&C->cmap);CHKERRQ(ierr);
215849b5e25fSSatish Balay   c->bs2 = a->bs2;
215949b5e25fSSatish Balay   c->mbs = a->mbs;
216049b5e25fSSatish Balay   c->nbs = a->nbs;
216149b5e25fSSatish Balay 
2162c760cd28SBarry Smith   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2163c760cd28SBarry Smith     c->imax           = a->imax;
2164c760cd28SBarry Smith     c->ilen           = a->ilen;
2165c760cd28SBarry Smith     c->free_imax_ilen = PETSC_FALSE;
2166c760cd28SBarry Smith   } else {
2167dcca6d9dSJed Brown     ierr = PetscMalloc2((mbs+1),&c->imax,(mbs+1),&c->ilen);CHKERRQ(ierr);
21683bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)C,2*(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
216949b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
217049b5e25fSSatish Balay       c->imax[i] = a->imax[i];
217149b5e25fSSatish Balay       c->ilen[i] = a->ilen[i];
217249b5e25fSSatish Balay     }
2173c760cd28SBarry Smith     c->free_imax_ilen = PETSC_TRUE;
2174c760cd28SBarry Smith   }
217549b5e25fSSatish Balay 
217649b5e25fSSatish Balay   /* allocate the matrix space */
21774da8f245SBarry Smith   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2178785e854fSJed Brown     ierr            = PetscMalloc1(bs2*nz,&c->a);CHKERRQ(ierr);
21793bb1ff40SBarry Smith     ierr            = PetscLogObjectMemory((PetscObject)C,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
218044e1c64aSLisandro Dalcin     c->i            = a->i;
218144e1c64aSLisandro Dalcin     c->j            = a->j;
21824da8f245SBarry Smith     c->singlemalloc = PETSC_FALSE;
218344e1c64aSLisandro Dalcin     c->free_a       = PETSC_TRUE;
21844da8f245SBarry Smith     c->free_ij      = PETSC_FALSE;
21854da8f245SBarry Smith     c->parent       = A;
21864da8f245SBarry Smith     ierr            = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
21874da8f245SBarry Smith     ierr            = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
21884da8f245SBarry Smith     ierr            = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
21894da8f245SBarry Smith   } else {
2190dcca6d9dSJed Brown     ierr            = PetscMalloc3(bs2*nz,&c->a,nz,&c->j,mbs+1,&c->i);CHKERRQ(ierr);
2191580bdb30SBarry Smith     ierr            = PetscArraycpy(c->i,a->i,mbs+1);CHKERRQ(ierr);
21923bb1ff40SBarry Smith     ierr            = PetscLogObjectMemory((PetscObject)C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr);
21934da8f245SBarry Smith     c->singlemalloc = PETSC_TRUE;
219444e1c64aSLisandro Dalcin     c->free_a       = PETSC_TRUE;
21954da8f245SBarry Smith     c->free_ij      = PETSC_TRUE;
21964da8f245SBarry Smith   }
219749b5e25fSSatish Balay   if (mbs > 0) {
21984da8f245SBarry Smith     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2199580bdb30SBarry Smith       ierr = PetscArraycpy(c->j,a->j,nz);CHKERRQ(ierr);
22004da8f245SBarry Smith     }
220149b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
2202580bdb30SBarry Smith       ierr = PetscArraycpy(c->a,a->a,bs2*nz);CHKERRQ(ierr);
220349b5e25fSSatish Balay     } else {
2204580bdb30SBarry Smith       ierr = PetscArrayzero(c->a,bs2*nz);CHKERRQ(ierr);
220549b5e25fSSatish Balay     }
2206a1c3900fSBarry Smith     if (a->jshort) {
220744e1c64aSLisandro Dalcin       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
220844e1c64aSLisandro Dalcin       /* if the parent matrix is reassembled, this child matrix will never notice */
2209785e854fSJed Brown       ierr = PetscMalloc1(nz,&c->jshort);CHKERRQ(ierr);
22103bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)C,nz*sizeof(unsigned short));CHKERRQ(ierr);
2211580bdb30SBarry Smith       ierr = PetscArraycpy(c->jshort,a->jshort,nz);CHKERRQ(ierr);
221226fbe8dcSKarl Rupp 
22134da8f245SBarry Smith       c->free_jshort = PETSC_TRUE;
22144da8f245SBarry Smith     }
2215a1c3900fSBarry Smith   }
221649b5e25fSSatish Balay 
221749b5e25fSSatish Balay   c->roworiented = a->roworiented;
221849b5e25fSSatish Balay   c->nonew       = a->nonew;
221949b5e25fSSatish Balay 
222049b5e25fSSatish Balay   if (a->diag) {
2221c760cd28SBarry Smith     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2222c760cd28SBarry Smith       c->diag      = a->diag;
2223c760cd28SBarry Smith       c->free_diag = PETSC_FALSE;
2224c760cd28SBarry Smith     } else {
2225785e854fSJed Brown       ierr = PetscMalloc1(mbs,&c->diag);CHKERRQ(ierr);
22263bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)C,mbs*sizeof(PetscInt));CHKERRQ(ierr);
222726fbe8dcSKarl Rupp       for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2228c760cd28SBarry Smith       c->free_diag = PETSC_TRUE;
2229c760cd28SBarry Smith     }
223044e1c64aSLisandro Dalcin   }
22316c6c5352SBarry Smith   c->nz         = a->nz;
2232f2cbd3d5SJed Brown   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
223349b5e25fSSatish Balay   c->solve_work = 0;
223449b5e25fSSatish Balay   c->mult_work  = 0;
223526fbe8dcSKarl Rupp 
223649b5e25fSSatish Balay   *B   = C;
2237140e18c1SBarry Smith   ierr = PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);CHKERRQ(ierr);
223849b5e25fSSatish Balay   PetscFunctionReturn(0);
223949b5e25fSSatish Balay }
224049b5e25fSSatish Balay 
2241112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
22422f480046SShri Abhyankar {
22432f480046SShri Abhyankar   Mat_SeqSBAIJ   *a;
22442f480046SShri Abhyankar   PetscErrorCode ierr;
22452f480046SShri Abhyankar   int            fd;
22462f480046SShri Abhyankar   PetscMPIInt    size;
22473059b6faSBarry Smith   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs = newmat->rmap->bs;
22482f480046SShri Abhyankar   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
22492f480046SShri Abhyankar   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
22502f480046SShri Abhyankar   PetscInt       *masked,nmask,tmp,bs2,ishift;
22512f480046SShri Abhyankar   PetscScalar    *aa;
2252ce94432eSBarry Smith   MPI_Comm       comm;
22537f489da9SVaclav Hapla   PetscBool      isbinary;
22542f480046SShri Abhyankar 
22552f480046SShri Abhyankar   PetscFunctionBegin;
22567f489da9SVaclav Hapla   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
22577f489da9SVaclav Hapla   if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)newmat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newmat)->type_name);
22587f489da9SVaclav Hapla 
2259c98fd787SBarry Smith   /* force binary viewer to load .info file if it has not yet done so */
2260c98fd787SBarry Smith   ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr);
2261ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
2262c5929fdfSBarry Smith   ierr = PetscOptionsGetInt(((PetscObject)newmat)->options,((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);CHKERRQ(ierr);
22633059b6faSBarry Smith   if (bs < 0) bs = 1;
22642f480046SShri Abhyankar   bs2  = bs*bs;
22652f480046SShri Abhyankar 
22662f480046SShri Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
22672f480046SShri Abhyankar   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
22682f480046SShri Abhyankar   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
22699860990eSLisandro Dalcin   ierr = PetscBinaryRead(fd,header,4,NULL,PETSC_INT);CHKERRQ(ierr);
22702f480046SShri Abhyankar   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
22712f480046SShri Abhyankar   M = header[1]; N = header[2]; nz = header[3];
22722f480046SShri Abhyankar 
22732f480046SShri Abhyankar   if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
22742f480046SShri Abhyankar 
22752f480046SShri Abhyankar   if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
22762f480046SShri Abhyankar 
22772f480046SShri Abhyankar   /*
22782f480046SShri Abhyankar      This code adds extra rows to make sure the number of rows is
22792f480046SShri Abhyankar     divisible by the blocksize
22802f480046SShri Abhyankar   */
22812f480046SShri Abhyankar   mbs        = M/bs;
22822f480046SShri Abhyankar   extra_rows = bs - M + bs*(mbs);
22832f480046SShri Abhyankar   if (extra_rows == bs) extra_rows = 0;
22842f480046SShri Abhyankar   else                  mbs++;
22852f480046SShri Abhyankar   if (extra_rows) {
22862f480046SShri Abhyankar     ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr);
22872f480046SShri Abhyankar   }
22882f480046SShri Abhyankar 
22892f480046SShri Abhyankar   /* Set global sizes if not already set */
22902f480046SShri Abhyankar   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
22912f480046SShri Abhyankar     ierr = MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);CHKERRQ(ierr);
22922f480046SShri Abhyankar   } else { /* Check if the matrix global sizes are correct */
22932f480046SShri Abhyankar     ierr = MatGetSize(newmat,&rows,&cols);CHKERRQ(ierr);
22942f480046SShri Abhyankar     if (M != rows ||  N != cols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols);
22952f480046SShri Abhyankar   }
22962f480046SShri Abhyankar 
22972f480046SShri Abhyankar   /* read in row lengths */
2298854ce69bSBarry Smith   ierr = PetscMalloc1(M+extra_rows,&rowlengths);CHKERRQ(ierr);
22999860990eSLisandro Dalcin   ierr = PetscBinaryRead(fd,rowlengths,M,NULL,PETSC_INT);CHKERRQ(ierr);
23002f480046SShri Abhyankar   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
23012f480046SShri Abhyankar 
23022f480046SShri Abhyankar   /* read in column indices */
2303854ce69bSBarry Smith   ierr = PetscMalloc1(nz+extra_rows,&jj);CHKERRQ(ierr);
23049860990eSLisandro Dalcin   ierr = PetscBinaryRead(fd,jj,nz,NULL,PETSC_INT);CHKERRQ(ierr);
23052f480046SShri Abhyankar   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
23062f480046SShri Abhyankar 
23072f480046SShri Abhyankar   /* loop over row lengths determining block row lengths */
2308071fcb05SBarry Smith   ierr     = PetscCalloc2(mbs,&s_browlengths,mbs,&mask);CHKERRQ(ierr);
2309071fcb05SBarry Smith   ierr     = PetscMalloc1(mbs,&masked);CHKERRQ(ierr);
23102f480046SShri Abhyankar   rowcount = 0;
23112f480046SShri Abhyankar   nzcount  = 0;
23122f480046SShri Abhyankar   for (i=0; i<mbs; i++) {
23132f480046SShri Abhyankar     nmask = 0;
23142f480046SShri Abhyankar     for (j=0; j<bs; j++) {
23152f480046SShri Abhyankar       kmax = rowlengths[rowcount];
23162f480046SShri Abhyankar       for (k=0; k<kmax; k++) {
23172f480046SShri Abhyankar         tmp = jj[nzcount++]/bs;   /* block col. index */
23182f480046SShri Abhyankar         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
23192f480046SShri Abhyankar       }
23202f480046SShri Abhyankar       rowcount++;
23212f480046SShri Abhyankar     }
23222f480046SShri Abhyankar     s_browlengths[i] += nmask;
23232f480046SShri Abhyankar 
23242f480046SShri Abhyankar     /* zero out the mask elements we set */
23252f480046SShri Abhyankar     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
23262f480046SShri Abhyankar   }
23272f480046SShri Abhyankar 
23282f480046SShri Abhyankar   /* Do preallocation */
2329367daffbSBarry Smith   ierr = MatSeqSBAIJSetPreallocation(newmat,bs,0,s_browlengths);CHKERRQ(ierr);
23302f480046SShri Abhyankar   a    = (Mat_SeqSBAIJ*)newmat->data;
23312f480046SShri Abhyankar 
23322f480046SShri Abhyankar   /* set matrix "i" values */
23332f480046SShri Abhyankar   a->i[0] = 0;
23342f480046SShri Abhyankar   for (i=1; i<= mbs; i++) {
23352f480046SShri Abhyankar     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
23362f480046SShri Abhyankar     a->ilen[i-1] = s_browlengths[i-1];
23372f480046SShri Abhyankar   }
23382f480046SShri Abhyankar   a->nz = a->i[mbs];
23392f480046SShri Abhyankar 
23402f480046SShri Abhyankar   /* read in nonzero values */
2341854ce69bSBarry Smith   ierr = PetscMalloc1(nz+extra_rows,&aa);CHKERRQ(ierr);
23429860990eSLisandro Dalcin   ierr = PetscBinaryRead(fd,aa,nz,NULL,PETSC_SCALAR);CHKERRQ(ierr);
23432f480046SShri Abhyankar   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
23442f480046SShri Abhyankar 
23452f480046SShri Abhyankar   /* set "a" and "j" values into matrix */
23462f480046SShri Abhyankar   nzcount = 0; jcount = 0;
23472f480046SShri Abhyankar   for (i=0; i<mbs; i++) {
23482f480046SShri Abhyankar     nzcountb = nzcount;
23492f480046SShri Abhyankar     nmask    = 0;
23502f480046SShri Abhyankar     for (j=0; j<bs; j++) {
23512f480046SShri Abhyankar       kmax = rowlengths[i*bs+j];
23522f480046SShri Abhyankar       for (k=0; k<kmax; k++) {
23532f480046SShri Abhyankar         tmp = jj[nzcount++]/bs; /* block col. index */
23542f480046SShri Abhyankar         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
23552f480046SShri Abhyankar       }
23562f480046SShri Abhyankar     }
23572f480046SShri Abhyankar     /* sort the masked values */
23582f480046SShri Abhyankar     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
23592f480046SShri Abhyankar 
23602f480046SShri Abhyankar     /* set "j" values into matrix */
23612f480046SShri Abhyankar     maskcount = 1;
23622f480046SShri Abhyankar     for (j=0; j<nmask; j++) {
23632f480046SShri Abhyankar       a->j[jcount++]  = masked[j];
23642f480046SShri Abhyankar       mask[masked[j]] = maskcount++;
23652f480046SShri Abhyankar     }
23662f480046SShri Abhyankar 
23672f480046SShri Abhyankar     /* set "a" values into matrix */
23682f480046SShri Abhyankar     ishift = bs2*a->i[i];
23692f480046SShri Abhyankar     for (j=0; j<bs; j++) {
23702f480046SShri Abhyankar       kmax = rowlengths[i*bs+j];
23712f480046SShri Abhyankar       for (k=0; k<kmax; k++) {
23722f480046SShri Abhyankar         tmp = jj[nzcountb]/bs;        /* block col. index */
23732f480046SShri Abhyankar         if (tmp >= i) {
23742f480046SShri Abhyankar           block     = mask[tmp] - 1;
23752f480046SShri Abhyankar           point     = jj[nzcountb] - bs*tmp;
23762f480046SShri Abhyankar           idx       = ishift + bs2*block + j + bs*point;
23772f480046SShri Abhyankar           a->a[idx] = aa[nzcountb];
23782f480046SShri Abhyankar         }
23792f480046SShri Abhyankar         nzcountb++;
23802f480046SShri Abhyankar       }
23812f480046SShri Abhyankar     }
23822f480046SShri Abhyankar     /* zero out the mask elements we set */
23832f480046SShri Abhyankar     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
23842f480046SShri Abhyankar   }
23852f480046SShri Abhyankar   if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
23862f480046SShri Abhyankar 
23872f480046SShri Abhyankar   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
2388071fcb05SBarry Smith   ierr = PetscFree2(s_browlengths,mask);CHKERRQ(ierr);
23892f480046SShri Abhyankar   ierr = PetscFree(aa);CHKERRQ(ierr);
23902f480046SShri Abhyankar   ierr = PetscFree(jj);CHKERRQ(ierr);
2391071fcb05SBarry Smith   ierr = PetscFree(masked);CHKERRQ(ierr);
23922f480046SShri Abhyankar 
23932f480046SShri Abhyankar   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23942f480046SShri Abhyankar   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
23952f480046SShri Abhyankar   PetscFunctionReturn(0);
23962f480046SShri Abhyankar }
23972f480046SShri Abhyankar 
2398c75a6043SHong Zhang /*@
2399c75a6043SHong Zhang      MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2400c75a6043SHong Zhang               (upper triangular entries in CSR format) provided by the user.
2401c75a6043SHong Zhang 
2402d083f849SBarry Smith      Collective
2403c75a6043SHong Zhang 
2404c75a6043SHong Zhang    Input Parameters:
2405c75a6043SHong Zhang +  comm - must be an MPI communicator of size 1
2406c75a6043SHong Zhang .  bs - size of block
2407c75a6043SHong Zhang .  m - number of rows
2408c75a6043SHong Zhang .  n - number of columns
2409483a2f95SBarry Smith .  i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix
2410c75a6043SHong Zhang .  j - column indices
2411c75a6043SHong Zhang -  a - matrix values
2412c75a6043SHong Zhang 
2413c75a6043SHong Zhang    Output Parameter:
2414c75a6043SHong Zhang .  mat - the matrix
2415c75a6043SHong Zhang 
2416dfb205c3SBarry Smith    Level: advanced
2417c75a6043SHong Zhang 
2418c75a6043SHong Zhang    Notes:
2419c75a6043SHong Zhang        The i, j, and a arrays are not copied by this routine, the user must free these arrays
2420c75a6043SHong Zhang     once the matrix is destroyed
2421c75a6043SHong Zhang 
2422c75a6043SHong Zhang        You cannot set new nonzero locations into this matrix, that will generate an error.
2423c75a6043SHong Zhang 
2424c75a6043SHong Zhang        The i and j indices are 0 based
2425c75a6043SHong Zhang 
2426dfb205c3SBarry Smith        When block size is greater than 1 the matrix values must be stored using the SBAIJ storage format (see the SBAIJ code to determine this). For block size of 1
2427dfb205c3SBarry Smith        it is the regular CSR format excluding the lower triangular elements.
2428dfb205c3SBarry Smith 
242969b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ()
2430c75a6043SHong Zhang 
2431c75a6043SHong Zhang @*/
2432c3c607ccSBarry Smith PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],Mat *mat)
2433c75a6043SHong Zhang {
2434c75a6043SHong Zhang   PetscErrorCode ierr;
2435c75a6043SHong Zhang   PetscInt       ii;
2436c75a6043SHong Zhang   Mat_SeqSBAIJ   *sbaij;
2437c75a6043SHong Zhang 
2438c75a6043SHong Zhang   PetscFunctionBegin;
2439e32f2f54SBarry Smith   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
244041096f02SStefano Zampini   if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2441c75a6043SHong Zhang 
2442c75a6043SHong Zhang   ierr  = MatCreate(comm,mat);CHKERRQ(ierr);
2443c75a6043SHong Zhang   ierr  = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr);
2444c75a6043SHong Zhang   ierr  = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr);
2445367daffbSBarry Smith   ierr  = MatSeqSBAIJSetPreallocation(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr);
2446c75a6043SHong Zhang   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2447dcca6d9dSJed Brown   ierr  = PetscMalloc2(m,&sbaij->imax,m,&sbaij->ilen);CHKERRQ(ierr);
24483bb1ff40SBarry Smith   ierr  = PetscLogObjectMemory((PetscObject)*mat,2*m*sizeof(PetscInt));CHKERRQ(ierr);
2449c75a6043SHong Zhang 
2450c75a6043SHong Zhang   sbaij->i = i;
2451c75a6043SHong Zhang   sbaij->j = j;
2452c75a6043SHong Zhang   sbaij->a = a;
245326fbe8dcSKarl Rupp 
2454c75a6043SHong Zhang   sbaij->singlemalloc   = PETSC_FALSE;
2455c75a6043SHong Zhang   sbaij->nonew          = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2456e6b907acSBarry Smith   sbaij->free_a         = PETSC_FALSE;
2457e6b907acSBarry Smith   sbaij->free_ij        = PETSC_FALSE;
2458ddf7884eSMatthew Knepley   sbaij->free_imax_ilen = PETSC_TRUE;
2459c75a6043SHong Zhang 
2460c75a6043SHong Zhang   for (ii=0; ii<m; ii++) {
2461c75a6043SHong Zhang     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2462c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG)
2463e32f2f54SBarry Smith     if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
2464c75a6043SHong Zhang #endif
2465c75a6043SHong Zhang   }
2466c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG)
2467c75a6043SHong Zhang   for (ii=0; ii<sbaij->i[m]; ii++) {
2468e32f2f54SBarry Smith     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2469e32f2f54SBarry Smith     if (j[ii] > n - 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
2470c75a6043SHong Zhang   }
2471c75a6043SHong Zhang #endif
2472c75a6043SHong Zhang 
2473c75a6043SHong Zhang   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2474c75a6043SHong Zhang   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2475c75a6043SHong Zhang   PetscFunctionReturn(0);
2476c75a6043SHong Zhang }
2477d06b337dSHong Zhang 
247859f5e6ceSHong Zhang PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqSBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
247959f5e6ceSHong Zhang {
248059f5e6ceSHong Zhang   PetscErrorCode ierr;
24818761c3d6SHong Zhang   PetscMPIInt    size;
248259f5e6ceSHong Zhang 
248359f5e6ceSHong Zhang   PetscFunctionBegin;
24848761c3d6SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
24858761c3d6SHong Zhang   if (size == 1 && scall == MAT_REUSE_MATRIX) {
24868761c3d6SHong Zhang     ierr = MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
24878761c3d6SHong Zhang   } else {
248859f5e6ceSHong Zhang     ierr = MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(comm,inmat,n,scall,outmat);CHKERRQ(ierr);
24898761c3d6SHong Zhang   }
249059f5e6ceSHong Zhang   PetscFunctionReturn(0);
249159f5e6ceSHong Zhang }
2492