xref: /petsc/src/mat/impls/sbaij/mpi/mpisbaij.c (revision 99a7f59ec9f9ccee3cdf9fd88ce16d7f72105d85)
1a30f8f8cSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h>    /*I "petscmat.h" I*/
3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
4c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
5c6db04a5SJed Brown #include <petscblaslapack.h>
6a30f8f8cSSatish Balay 
76214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
8cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
96214f412SHong Zhang #endif
10d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
11d24d4204SJose E. Roman PETSC_INTERN PetscErrorCode MatConvert_SBAIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
12d24d4204SJose E. Roman #endif
13b147fbf3SStefano Zampini 
14b147fbf3SStefano Zampini /* This could be moved to matimpl.h */
15b147fbf3SStefano Zampini static PetscErrorCode MatPreallocateWithMats_Private(Mat B, PetscInt nm, Mat X[], PetscBool symm[], PetscBool fill)
16b147fbf3SStefano Zampini {
17b147fbf3SStefano Zampini   Mat            preallocator;
18b147fbf3SStefano Zampini   PetscInt       r,rstart,rend;
19b147fbf3SStefano Zampini   PetscInt       bs,i,m,n,M,N;
20b147fbf3SStefano Zampini   PetscBool      cong = PETSC_TRUE;
21b147fbf3SStefano Zampini 
22b147fbf3SStefano Zampini   PetscFunctionBegin;
23b147fbf3SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
24b147fbf3SStefano Zampini   PetscValidLogicalCollectiveInt(B,nm,2);
25b147fbf3SStefano Zampini   for (i = 0; i < nm; i++) {
26b147fbf3SStefano Zampini     PetscValidHeaderSpecific(X[i],MAT_CLASSID,3);
279566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCompare(B->rmap,X[i]->rmap,&cong));
285f80ce2aSJacob Faibussowitsch     PetscCheck(cong,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for different layouts");
29b147fbf3SStefano Zampini   }
30b147fbf3SStefano Zampini   PetscValidLogicalCollectiveBool(B,fill,5);
319566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(B,&bs));
329566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B,&M,&N));
339566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(B,&m,&n));
349566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)B),&preallocator));
359566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator,MATPREALLOCATOR));
369566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(preallocator,bs));
379566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator,m,n,M,N));
389566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
399566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(preallocator,&rstart,&rend));
40b147fbf3SStefano Zampini   for (r = rstart; r < rend; ++r) {
41b147fbf3SStefano Zampini     PetscInt          ncols;
42b147fbf3SStefano Zampini     const PetscInt    *row;
43b147fbf3SStefano Zampini     const PetscScalar *vals;
44b147fbf3SStefano Zampini 
45b147fbf3SStefano Zampini     for (i = 0; i < nm; i++) {
469566063dSJacob Faibussowitsch       PetscCall(MatGetRow(X[i],r,&ncols,&row,&vals));
479566063dSJacob Faibussowitsch       PetscCall(MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES));
48b147fbf3SStefano Zampini       if (symm && symm[i]) {
499566063dSJacob Faibussowitsch         PetscCall(MatSetValues(preallocator,ncols,row,1,&r,vals,INSERT_VALUES));
50b147fbf3SStefano Zampini       }
519566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(X[i],r,&ncols,&row,&vals));
52b147fbf3SStefano Zampini     }
53b147fbf3SStefano Zampini   }
549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
569566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator,fill,B));
579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
58b147fbf3SStefano Zampini   PetscFunctionReturn(0);
59b147fbf3SStefano Zampini }
60b147fbf3SStefano Zampini 
6128d58a37SPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
62b147fbf3SStefano Zampini {
63b147fbf3SStefano Zampini   Mat            B;
64b147fbf3SStefano Zampini   PetscInt       r;
65b147fbf3SStefano Zampini 
66b147fbf3SStefano Zampini   PetscFunctionBegin;
67b147fbf3SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
6828d58a37SPierre Jolivet     PetscBool symm = PETSC_TRUE,isdense;
69b147fbf3SStefano Zampini     PetscInt  bs;
70b147fbf3SStefano Zampini 
719566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B));
729566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N));
739566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,newtype));
749566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A,&bs));
759566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSize(B,bs));
769566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(B->rmap));
779566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(B->cmap));
789566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&isdense,MATSEQDENSE,MATMPIDENSE,MATSEQDENSECUDA,""));
7928d58a37SPierre Jolivet     if (!isdense) {
809566063dSJacob Faibussowitsch       PetscCall(MatGetRowUpperTriangular(A));
819566063dSJacob Faibussowitsch       PetscCall(MatPreallocateWithMats_Private(B,1,&A,&symm,PETSC_TRUE));
829566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowUpperTriangular(A));
8328d58a37SPierre Jolivet     } else {
849566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B));
8528d58a37SPierre Jolivet     }
8628d58a37SPierre Jolivet   } else {
8728d58a37SPierre Jolivet     B    = *newmat;
889566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(B));
8928d58a37SPierre Jolivet   }
90b147fbf3SStefano Zampini 
919566063dSJacob Faibussowitsch   PetscCall(MatGetRowUpperTriangular(A));
92b147fbf3SStefano Zampini   for (r = A->rmap->rstart; r < A->rmap->rend; r++) {
93b147fbf3SStefano Zampini     PetscInt          ncols;
94b147fbf3SStefano Zampini     const PetscInt    *row;
95b147fbf3SStefano Zampini     const PetscScalar *vals;
96b147fbf3SStefano Zampini 
979566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,r,&ncols,&row,&vals));
989566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,1,&r,ncols,row,vals,INSERT_VALUES));
99eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
100eb1ec7c1SStefano Zampini     if (A->hermitian) {
101eb1ec7c1SStefano Zampini       PetscInt i;
102eb1ec7c1SStefano Zampini       for (i = 0; i < ncols; i++) {
1039566063dSJacob Faibussowitsch         PetscCall(MatSetValue(B,row[i],r,PetscConj(vals[i]),INSERT_VALUES));
104eb1ec7c1SStefano Zampini       }
105eb1ec7c1SStefano Zampini     } else {
1069566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B,ncols,row,1,&r,vals,INSERT_VALUES));
107eb1ec7c1SStefano Zampini     }
108eb1ec7c1SStefano Zampini #else
1099566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B,ncols,row,1,&r,vals,INSERT_VALUES));
110eb1ec7c1SStefano Zampini #endif
1119566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,r,&ncols,&row,&vals));
112b147fbf3SStefano Zampini   }
1139566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowUpperTriangular(A));
1149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
116b147fbf3SStefano Zampini 
117b147fbf3SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
1189566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A,&B));
119b147fbf3SStefano Zampini   } else {
120b147fbf3SStefano Zampini     *newmat = B;
121b147fbf3SStefano Zampini   }
122b147fbf3SStefano Zampini   PetscFunctionReturn(0);
123b147fbf3SStefano Zampini }
124b147fbf3SStefano Zampini 
1257087cfbeSBarry Smith PetscErrorCode  MatStoreValues_MPISBAIJ(Mat mat)
126a30f8f8cSSatish Balay {
127f3566a2aSHong Zhang   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;
128a30f8f8cSSatish Balay 
129a30f8f8cSSatish Balay   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(aij->A));
1319566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(aij->B));
132a30f8f8cSSatish Balay   PetscFunctionReturn(0);
133a30f8f8cSSatish Balay }
134a30f8f8cSSatish Balay 
1357087cfbeSBarry Smith PetscErrorCode  MatRetrieveValues_MPISBAIJ(Mat mat)
136a30f8f8cSSatish Balay {
137f3566a2aSHong Zhang   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;
138a30f8f8cSSatish Balay 
139a30f8f8cSSatish Balay   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(aij->A));
1419566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(aij->B));
142a30f8f8cSSatish Balay   PetscFunctionReturn(0);
143a30f8f8cSSatish Balay }
144a30f8f8cSSatish Balay 
145d40312a9SBarry Smith #define  MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,orow,ocol)      \
146a30f8f8cSSatish Balay   { \
147a30f8f8cSSatish Balay     brow = row/bs;  \
148a30f8f8cSSatish Balay     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
149a30f8f8cSSatish Balay     rmax = aimax[brow]; nrow = ailen[brow]; \
150a30f8f8cSSatish Balay     bcol = col/bs; \
151a30f8f8cSSatish Balay     ridx = row % bs; cidx = col % bs; \
152a30f8f8cSSatish Balay     low  = 0; high = nrow; \
153a30f8f8cSSatish Balay     while (high-low > 3) { \
154a30f8f8cSSatish Balay       t = (low+high)/2; \
155a30f8f8cSSatish Balay       if (rp[t] > bcol) high = t; \
156a30f8f8cSSatish Balay       else              low  = t; \
157a30f8f8cSSatish Balay     } \
158a30f8f8cSSatish Balay     for (_i=low; _i<high; _i++) { \
159a30f8f8cSSatish Balay       if (rp[_i] > bcol) break; \
160a30f8f8cSSatish Balay       if (rp[_i] == bcol) { \
161a30f8f8cSSatish Balay         bap = ap + bs2*_i + bs*cidx + ridx; \
162a30f8f8cSSatish Balay         if (addv == ADD_VALUES) *bap += value;  \
163a30f8f8cSSatish Balay         else                    *bap  = value;  \
164a30f8f8cSSatish Balay         goto a_noinsert; \
165a30f8f8cSSatish Balay       } \
166a30f8f8cSSatish Balay     } \
167a30f8f8cSSatish Balay     if (a->nonew == 1) goto a_noinsert; \
1685f80ce2aSJacob Faibussowitsch     PetscCheck(a->nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
169fef13f97SBarry Smith     MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
170a30f8f8cSSatish Balay     N = nrow++ - 1;  \
171a30f8f8cSSatish Balay     /* shift up all the later entries in this row */ \
1729566063dSJacob Faibussowitsch     PetscCall(PetscArraymove(rp+_i+1,rp+_i,N-_i+1)); \
1739566063dSJacob Faibussowitsch     PetscCall(PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1))); \
1749566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ap+bs2*_i,bs2));  \
175a30f8f8cSSatish Balay     rp[_i]                      = bcol;  \
176a30f8f8cSSatish Balay     ap[bs2*_i + bs*cidx + ridx] = value;  \
177e56f5c9eSBarry Smith     A->nonzerostate++;\
178a30f8f8cSSatish Balay a_noinsert:; \
179a30f8f8cSSatish Balay     ailen[brow] = nrow; \
180a30f8f8cSSatish Balay   }
181e5e170daSBarry Smith 
182d40312a9SBarry Smith #define  MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,orow,ocol) \
183a30f8f8cSSatish Balay   { \
184a30f8f8cSSatish Balay     brow = row/bs;  \
185a30f8f8cSSatish Balay     rp   = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
186a30f8f8cSSatish Balay     rmax = bimax[brow]; nrow = bilen[brow]; \
187a30f8f8cSSatish Balay     bcol = col/bs; \
188a30f8f8cSSatish Balay     ridx = row % bs; cidx = col % bs; \
189a30f8f8cSSatish Balay     low  = 0; high = nrow; \
190a30f8f8cSSatish Balay     while (high-low > 3) { \
191a30f8f8cSSatish Balay       t = (low+high)/2; \
192a30f8f8cSSatish Balay       if (rp[t] > bcol) high = t; \
193a30f8f8cSSatish Balay       else              low  = t; \
194a30f8f8cSSatish Balay     } \
195a30f8f8cSSatish Balay     for (_i=low; _i<high; _i++) { \
196a30f8f8cSSatish Balay       if (rp[_i] > bcol) break; \
197a30f8f8cSSatish Balay       if (rp[_i] == bcol) { \
198a30f8f8cSSatish Balay         bap = ap + bs2*_i + bs*cidx + ridx; \
199a30f8f8cSSatish Balay         if (addv == ADD_VALUES) *bap += value;  \
200a30f8f8cSSatish Balay         else                    *bap  = value;  \
201a30f8f8cSSatish Balay         goto b_noinsert; \
202a30f8f8cSSatish Balay       } \
203a30f8f8cSSatish Balay     } \
204a30f8f8cSSatish Balay     if (b->nonew == 1) goto b_noinsert; \
2055f80ce2aSJacob Faibussowitsch     PetscCheck(b->nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
206fef13f97SBarry Smith     MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
207a30f8f8cSSatish Balay     N = nrow++ - 1;  \
208a30f8f8cSSatish Balay     /* shift up all the later entries in this row */ \
2099566063dSJacob Faibussowitsch     PetscCall(PetscArraymove(rp+_i+1,rp+_i,N-_i+1)); \
2109566063dSJacob Faibussowitsch     PetscCall(PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1))); \
2119566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(ap+bs2*_i,bs2)); \
212a30f8f8cSSatish Balay     rp[_i]                      = bcol;  \
213a30f8f8cSSatish Balay     ap[bs2*_i + bs*cidx + ridx] = value;  \
214e56f5c9eSBarry Smith     B->nonzerostate++;\
215a30f8f8cSSatish Balay b_noinsert:; \
216a30f8f8cSSatish Balay     bilen[brow] = nrow; \
217a30f8f8cSSatish Balay   }
218a30f8f8cSSatish Balay 
219a30f8f8cSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
220476417e5SBarry Smith    Any a(i,j) with i>j input by user is ingored or generates an error
221a30f8f8cSSatish Balay */
222dd6ea824SBarry Smith PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
223a30f8f8cSSatish Balay {
224a30f8f8cSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
225a30f8f8cSSatish Balay   MatScalar      value;
226ace3abfcSBarry Smith   PetscBool      roworiented = baij->roworiented;
2271302d50aSBarry Smith   PetscInt       i,j,row,col;
228d0f46423SBarry Smith   PetscInt       rstart_orig=mat->rmap->rstart;
229d0f46423SBarry Smith   PetscInt       rend_orig  =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
230d0f46423SBarry Smith   PetscInt       cend_orig  =mat->cmap->rend,bs=mat->rmap->bs;
231a30f8f8cSSatish Balay 
232a30f8f8cSSatish Balay   /* Some Variables required in the macro */
233a30f8f8cSSatish Balay   Mat          A     = baij->A;
234a30f8f8cSSatish Balay   Mat_SeqSBAIJ *a    = (Mat_SeqSBAIJ*)(A)->data;
2351302d50aSBarry Smith   PetscInt     *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
236a30f8f8cSSatish Balay   MatScalar    *aa   =a->a;
237a30f8f8cSSatish Balay 
238a30f8f8cSSatish Balay   Mat         B     = baij->B;
239a30f8f8cSSatish Balay   Mat_SeqBAIJ *b    = (Mat_SeqBAIJ*)(B)->data;
2401302d50aSBarry Smith   PetscInt    *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
241a30f8f8cSSatish Balay   MatScalar   *ba   =b->a;
242a30f8f8cSSatish Balay 
2431302d50aSBarry Smith   PetscInt  *rp,ii,nrow,_i,rmax,N,brow,bcol;
2441302d50aSBarry Smith   PetscInt  low,high,t,ridx,cidx,bs2=a->bs2;
245a30f8f8cSSatish Balay   MatScalar *ap,*bap;
246a30f8f8cSSatish Balay 
247a30f8f8cSSatish Balay   /* for stash */
2480298fd71SBarry Smith   PetscInt  n_loc, *in_loc = NULL;
2490298fd71SBarry Smith   MatScalar *v_loc = NULL;
250a30f8f8cSSatish Balay 
251a30f8f8cSSatish Balay   PetscFunctionBegin;
252a30f8f8cSSatish Balay   if (!baij->donotstash) {
25359ffdab8SBarry Smith     if (n > baij->n_loc) {
2549566063dSJacob Faibussowitsch       PetscCall(PetscFree(baij->in_loc));
2559566063dSJacob Faibussowitsch       PetscCall(PetscFree(baij->v_loc));
2569566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n,&baij->in_loc));
2579566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n,&baij->v_loc));
25826fbe8dcSKarl Rupp 
25959ffdab8SBarry Smith       baij->n_loc = n;
26059ffdab8SBarry Smith     }
26159ffdab8SBarry Smith     in_loc = baij->in_loc;
26259ffdab8SBarry Smith     v_loc  = baij->v_loc;
263a30f8f8cSSatish Balay   }
264a30f8f8cSSatish Balay 
265a30f8f8cSSatish Balay   for (i=0; i<m; i++) {
266a30f8f8cSSatish Balay     if (im[i] < 0) continue;
2675f80ce2aSJacob Faibussowitsch     PetscCheck(im[i] < mat->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,im[i],mat->rmap->N-1);
268a30f8f8cSSatish Balay     if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
269a30f8f8cSSatish Balay       row = im[i] - rstart_orig;              /* local row index */
270a30f8f8cSSatish Balay       for (j=0; j<n; j++) {
27101b2bd88SHong Zhang         if (im[i]/bs > in[j]/bs) {
27201b2bd88SHong Zhang           if (a->ignore_ltriangular) {
27301b2bd88SHong Zhang             continue;    /* ignore lower triangular blocks */
27426fbe8dcSKarl 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)");
27501b2bd88SHong Zhang         }
276a30f8f8cSSatish Balay         if (in[j] >= cstart_orig && in[j] < cend_orig) {  /* diag entry (A) */
277a30f8f8cSSatish Balay           col  = in[j] - cstart_orig;         /* local col index */
278a30f8f8cSSatish Balay           brow = row/bs; bcol = col/bs;
279a30f8f8cSSatish Balay           if (brow > bcol) continue;  /* ignore lower triangular blocks of A */
280db4deed7SKarl Rupp           if (roworiented) value = v[i*n+j];
281db4deed7SKarl Rupp           else             value = v[i+j*m];
282d40312a9SBarry Smith           MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,im[i],in[j]);
2839566063dSJacob Faibussowitsch           /* PetscCall(MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv)); */
284f7d195e4SLawrence Mitchell         } else if (in[j] < 0) {
285f7d195e4SLawrence Mitchell           continue;
286f7d195e4SLawrence Mitchell         } else {
287f7d195e4SLawrence Mitchell           PetscCheck(in[j] < mat->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,in[j],mat->cmap->N-1);
288f7d195e4SLawrence Mitchell           /* off-diag entry (B) */
289a30f8f8cSSatish Balay           if (mat->was_assembled) {
290a30f8f8cSSatish Balay             if (!baij->colmap) {
2919566063dSJacob Faibussowitsch               PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
292a30f8f8cSSatish Balay             }
293a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
2949566063dSJacob Faibussowitsch             PetscCall(PetscTableFind(baij->colmap,in[j]/bs + 1,&col));
29571730473SSatish Balay             col  = col - 1;
296a30f8f8cSSatish Balay #else
29771730473SSatish Balay             col = baij->colmap[in[j]/bs] - 1;
298a30f8f8cSSatish Balay #endif
299a30f8f8cSSatish Balay             if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
3009566063dSJacob Faibussowitsch               PetscCall(MatDisAssemble_MPISBAIJ(mat));
301a30f8f8cSSatish Balay               col  =  in[j];
302a30f8f8cSSatish Balay               /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
303a30f8f8cSSatish Balay               B    = baij->B;
304a30f8f8cSSatish Balay               b    = (Mat_SeqBAIJ*)(B)->data;
305a30f8f8cSSatish Balay               bimax= b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
306a30f8f8cSSatish Balay               ba   = b->a;
30771730473SSatish Balay             } else col += in[j]%bs;
308a30f8f8cSSatish Balay           } else col = in[j];
309db4deed7SKarl Rupp           if (roworiented) value = v[i*n+j];
310db4deed7SKarl Rupp           else             value = v[i+j*m];
311d40312a9SBarry Smith           MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,im[i],in[j]);
3129566063dSJacob Faibussowitsch           /* PetscCall(MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv)); */
313a30f8f8cSSatish Balay         }
314a30f8f8cSSatish Balay       }
315a30f8f8cSSatish Balay     } else {  /* off processor entry */
3165f80ce2aSJacob Faibussowitsch       PetscCheck(!mat->nooffprocentries,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
317a30f8f8cSSatish Balay       if (!baij->donotstash) {
3185080c13bSMatthew G Knepley         mat->assembled = PETSC_FALSE;
319a30f8f8cSSatish Balay         n_loc          = 0;
320a30f8f8cSSatish Balay         for (j=0; j<n; j++) {
321f65c83cfSHong Zhang           if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
322a30f8f8cSSatish Balay           in_loc[n_loc] = in[j];
323a30f8f8cSSatish Balay           if (roworiented) {
324a30f8f8cSSatish Balay             v_loc[n_loc] = v[i*n+j];
325a30f8f8cSSatish Balay           } else {
326a30f8f8cSSatish Balay             v_loc[n_loc] = v[j*m+i];
327a30f8f8cSSatish Balay           }
328a30f8f8cSSatish Balay           n_loc++;
329a30f8f8cSSatish Balay         }
3309566063dSJacob Faibussowitsch         PetscCall(MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE));
331a30f8f8cSSatish Balay       }
332a30f8f8cSSatish Balay     }
333a30f8f8cSSatish Balay   }
334a30f8f8cSSatish Balay   PetscFunctionReturn(0);
335a30f8f8cSSatish Balay }
336a30f8f8cSSatish Balay 
3379fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetValuesBlocked_SeqSBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
33836bd2089SBarry Smith {
33936bd2089SBarry Smith   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
34036bd2089SBarry Smith   PetscInt          *rp,low,high,t,ii,jj,nrow,i,rmax,N;
34136bd2089SBarry Smith   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
34236bd2089SBarry Smith   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
34336bd2089SBarry Smith   PetscBool         roworiented=a->roworiented;
34436bd2089SBarry Smith   const PetscScalar *value     = v;
34536bd2089SBarry Smith   MatScalar         *ap,*aa = a->a,*bap;
34636bd2089SBarry Smith 
34736bd2089SBarry Smith   PetscFunctionBegin;
34836bd2089SBarry Smith   if (col < row) {
34936bd2089SBarry Smith     if (a->ignore_ltriangular) PetscFunctionReturn(0); /* ignore lower triangular block */
35036bd2089SBarry Smith     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)");
35136bd2089SBarry Smith   }
35236bd2089SBarry Smith   rp   = aj + ai[row];
35336bd2089SBarry Smith   ap   = aa + bs2*ai[row];
35436bd2089SBarry Smith   rmax = imax[row];
35536bd2089SBarry Smith   nrow = ailen[row];
35636bd2089SBarry Smith   value = v;
35736bd2089SBarry Smith   low   = 0;
35836bd2089SBarry Smith   high  = nrow;
35936bd2089SBarry Smith 
36036bd2089SBarry Smith   while (high-low > 7) {
36136bd2089SBarry Smith     t = (low+high)/2;
36236bd2089SBarry Smith     if (rp[t] > col) high = t;
36336bd2089SBarry Smith     else             low  = t;
36436bd2089SBarry Smith   }
36536bd2089SBarry Smith   for (i=low; i<high; i++) {
36636bd2089SBarry Smith     if (rp[i] > col) break;
36736bd2089SBarry Smith     if (rp[i] == col) {
36836bd2089SBarry Smith       bap = ap +  bs2*i;
36936bd2089SBarry Smith       if (roworiented) {
37036bd2089SBarry Smith         if (is == ADD_VALUES) {
37136bd2089SBarry Smith           for (ii=0; ii<bs; ii++) {
37236bd2089SBarry Smith             for (jj=ii; jj<bs2; jj+=bs) {
37336bd2089SBarry Smith               bap[jj] += *value++;
37436bd2089SBarry Smith             }
37536bd2089SBarry Smith           }
37636bd2089SBarry Smith         } else {
37736bd2089SBarry Smith           for (ii=0; ii<bs; ii++) {
37836bd2089SBarry Smith             for (jj=ii; jj<bs2; jj+=bs) {
37936bd2089SBarry Smith               bap[jj] = *value++;
38036bd2089SBarry Smith             }
38136bd2089SBarry Smith           }
38236bd2089SBarry Smith         }
38336bd2089SBarry Smith       } else {
38436bd2089SBarry Smith         if (is == ADD_VALUES) {
38536bd2089SBarry Smith           for (ii=0; ii<bs; ii++) {
38636bd2089SBarry Smith             for (jj=0; jj<bs; jj++) {
38736bd2089SBarry Smith               *bap++ += *value++;
38836bd2089SBarry Smith             }
38936bd2089SBarry Smith           }
39036bd2089SBarry Smith         } else {
39136bd2089SBarry Smith           for (ii=0; ii<bs; ii++) {
39236bd2089SBarry Smith             for (jj=0; jj<bs; jj++) {
39336bd2089SBarry Smith               *bap++  = *value++;
39436bd2089SBarry Smith             }
39536bd2089SBarry Smith           }
39636bd2089SBarry Smith         }
39736bd2089SBarry Smith       }
39836bd2089SBarry Smith       goto noinsert2;
39936bd2089SBarry Smith     }
40036bd2089SBarry Smith   }
40136bd2089SBarry Smith   if (nonew == 1) goto noinsert2;
4025f80ce2aSJacob Faibussowitsch   PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new block index nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", orow, ocol);
40336bd2089SBarry Smith   MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
40436bd2089SBarry Smith   N = nrow++ - 1; high++;
40536bd2089SBarry Smith   /* shift up all the later entries in this row */
4069566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(rp+i+1,rp+i,N-i+1));
4079566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1)));
40836bd2089SBarry Smith   rp[i] = col;
40936bd2089SBarry Smith   bap   = ap +  bs2*i;
41036bd2089SBarry Smith   if (roworiented) {
41136bd2089SBarry Smith     for (ii=0; ii<bs; ii++) {
41236bd2089SBarry Smith       for (jj=ii; jj<bs2; jj+=bs) {
41336bd2089SBarry Smith         bap[jj] = *value++;
41436bd2089SBarry Smith       }
41536bd2089SBarry Smith     }
41636bd2089SBarry Smith   } else {
41736bd2089SBarry Smith     for (ii=0; ii<bs; ii++) {
41836bd2089SBarry Smith       for (jj=0; jj<bs; jj++) {
41936bd2089SBarry Smith         *bap++ = *value++;
42036bd2089SBarry Smith       }
42136bd2089SBarry Smith     }
42236bd2089SBarry Smith   }
42336bd2089SBarry Smith   noinsert2:;
42436bd2089SBarry Smith   ailen[row] = nrow;
42536bd2089SBarry Smith   PetscFunctionReturn(0);
42636bd2089SBarry Smith }
42736bd2089SBarry Smith 
42836bd2089SBarry Smith /*
42936bd2089SBarry Smith    This routine is exactly duplicated in mpibaij.c
43036bd2089SBarry Smith */
4319fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetValuesBlocked_SeqBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
43236bd2089SBarry Smith {
43336bd2089SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
43436bd2089SBarry Smith   PetscInt          *rp,low,high,t,ii,jj,nrow,i,rmax,N;
43536bd2089SBarry Smith   PetscInt          *imax=a->imax,*ai=a->i,*ailen=a->ilen;
43636bd2089SBarry Smith   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
43736bd2089SBarry Smith   PetscBool         roworiented=a->roworiented;
43836bd2089SBarry Smith   const PetscScalar *value     = v;
43936bd2089SBarry Smith   MatScalar         *ap,*aa = a->a,*bap;
44036bd2089SBarry Smith 
44136bd2089SBarry Smith   PetscFunctionBegin;
44236bd2089SBarry Smith   rp   = aj + ai[row];
44336bd2089SBarry Smith   ap   = aa + bs2*ai[row];
44436bd2089SBarry Smith   rmax = imax[row];
44536bd2089SBarry Smith   nrow = ailen[row];
44636bd2089SBarry Smith   low  = 0;
44736bd2089SBarry Smith   high = nrow;
44836bd2089SBarry Smith   value = v;
44936bd2089SBarry Smith   while (high-low > 7) {
45036bd2089SBarry Smith     t = (low+high)/2;
45136bd2089SBarry Smith     if (rp[t] > col) high = t;
45236bd2089SBarry Smith     else             low  = t;
45336bd2089SBarry Smith   }
45436bd2089SBarry Smith   for (i=low; i<high; i++) {
45536bd2089SBarry Smith     if (rp[i] > col) break;
45636bd2089SBarry Smith     if (rp[i] == col) {
45736bd2089SBarry Smith       bap = ap +  bs2*i;
45836bd2089SBarry Smith       if (roworiented) {
45936bd2089SBarry Smith         if (is == ADD_VALUES) {
46036bd2089SBarry Smith           for (ii=0; ii<bs; ii++) {
46136bd2089SBarry Smith             for (jj=ii; jj<bs2; jj+=bs) {
46236bd2089SBarry Smith               bap[jj] += *value++;
46336bd2089SBarry Smith             }
46436bd2089SBarry Smith           }
46536bd2089SBarry Smith         } else {
46636bd2089SBarry Smith           for (ii=0; ii<bs; ii++) {
46736bd2089SBarry Smith             for (jj=ii; jj<bs2; jj+=bs) {
46836bd2089SBarry Smith               bap[jj] = *value++;
46936bd2089SBarry Smith             }
47036bd2089SBarry Smith           }
47136bd2089SBarry Smith         }
47236bd2089SBarry Smith       } else {
47336bd2089SBarry Smith         if (is == ADD_VALUES) {
47436bd2089SBarry Smith           for (ii=0; ii<bs; ii++,value+=bs) {
47536bd2089SBarry Smith             for (jj=0; jj<bs; jj++) {
47636bd2089SBarry Smith               bap[jj] += value[jj];
47736bd2089SBarry Smith             }
47836bd2089SBarry Smith             bap += bs;
47936bd2089SBarry Smith           }
48036bd2089SBarry Smith         } else {
48136bd2089SBarry Smith           for (ii=0; ii<bs; ii++,value+=bs) {
48236bd2089SBarry Smith             for (jj=0; jj<bs; jj++) {
48336bd2089SBarry Smith               bap[jj]  = value[jj];
48436bd2089SBarry Smith             }
48536bd2089SBarry Smith             bap += bs;
48636bd2089SBarry Smith           }
48736bd2089SBarry Smith         }
48836bd2089SBarry Smith       }
48936bd2089SBarry Smith       goto noinsert2;
49036bd2089SBarry Smith     }
49136bd2089SBarry Smith   }
49236bd2089SBarry Smith   if (nonew == 1) goto noinsert2;
4935f80ce2aSJacob Faibussowitsch   PetscCheck(nonew != -1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new global block indexed nonzero block (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", orow, ocol);
49436bd2089SBarry Smith   MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
49536bd2089SBarry Smith   N = nrow++ - 1; high++;
49636bd2089SBarry Smith   /* shift up all the later entries in this row */
4979566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(rp+i+1,rp+i,N-i+1));
4989566063dSJacob Faibussowitsch   PetscCall(PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1)));
49936bd2089SBarry Smith   rp[i] = col;
50036bd2089SBarry Smith   bap   = ap +  bs2*i;
50136bd2089SBarry Smith   if (roworiented) {
50236bd2089SBarry Smith     for (ii=0; ii<bs; ii++) {
50336bd2089SBarry Smith       for (jj=ii; jj<bs2; jj+=bs) {
50436bd2089SBarry Smith         bap[jj] = *value++;
50536bd2089SBarry Smith       }
50636bd2089SBarry Smith     }
50736bd2089SBarry Smith   } else {
50836bd2089SBarry Smith     for (ii=0; ii<bs; ii++) {
50936bd2089SBarry Smith       for (jj=0; jj<bs; jj++) {
51036bd2089SBarry Smith         *bap++ = *value++;
51136bd2089SBarry Smith       }
51236bd2089SBarry Smith     }
51336bd2089SBarry Smith   }
51436bd2089SBarry Smith   noinsert2:;
51536bd2089SBarry Smith   ailen[row] = nrow;
51636bd2089SBarry Smith   PetscFunctionReturn(0);
51736bd2089SBarry Smith }
51836bd2089SBarry Smith 
51936bd2089SBarry Smith /*
52036bd2089SBarry Smith     This routine could be optimized by removing the need for the block copy below and passing stride information
52136bd2089SBarry Smith   to the above inline routines; similarly in MatSetValuesBlocked_MPIBAIJ()
52236bd2089SBarry Smith */
523dd6ea824SBarry Smith PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
524a30f8f8cSSatish Balay {
5250880e062SHong Zhang   Mat_MPISBAIJ    *baij = (Mat_MPISBAIJ*)mat->data;
526f15d580aSBarry Smith   const MatScalar *value;
527f15d580aSBarry Smith   MatScalar       *barray     =baij->barray;
528ace3abfcSBarry Smith   PetscBool       roworiented = baij->roworiented,ignore_ltriangular = ((Mat_SeqSBAIJ*)baij->A->data)->ignore_ltriangular;
529899cda47SBarry Smith   PetscInt        i,j,ii,jj,row,col,rstart=baij->rstartbs;
530476417e5SBarry Smith   PetscInt        rend=baij->rendbs,cstart=baij->cstartbs,stepval;
531476417e5SBarry Smith   PetscInt        cend=baij->cendbs,bs=mat->rmap->bs,bs2=baij->bs2;
5320880e062SHong Zhang 
533a30f8f8cSSatish Balay   PetscFunctionBegin;
5340880e062SHong Zhang   if (!barray) {
5359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs2,&barray));
5360880e062SHong Zhang     baij->barray = barray;
5370880e062SHong Zhang   }
5380880e062SHong Zhang 
5390880e062SHong Zhang   if (roworiented) {
5400880e062SHong Zhang     stepval = (n-1)*bs;
5410880e062SHong Zhang   } else {
5420880e062SHong Zhang     stepval = (m-1)*bs;
5430880e062SHong Zhang   }
5440880e062SHong Zhang   for (i=0; i<m; i++) {
5450880e062SHong Zhang     if (im[i] < 0) continue;
5466bdcaf15SBarry Smith     PetscCheck(im[i] < baij->Mbs,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block indexed row too large %" PetscInt_FMT " max %" PetscInt_FMT,im[i],baij->Mbs-1);
5470880e062SHong Zhang     if (im[i] >= rstart && im[i] < rend) {
5480880e062SHong Zhang       row = im[i] - rstart;
5490880e062SHong Zhang       for (j=0; j<n; j++) {
550f3f98c53SJed Brown         if (im[i] > in[j]) {
551f3f98c53SJed Brown           if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
552e32f2f54SBarry Smith           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)");
553f3f98c53SJed Brown         }
5540880e062SHong Zhang         /* If NumCol = 1 then a copy is not required */
5550880e062SHong Zhang         if ((roworiented) && (n == 1)) {
556f15d580aSBarry Smith           barray = (MatScalar*) v + i*bs2;
5570880e062SHong Zhang         } else if ((!roworiented) && (m == 1)) {
558f15d580aSBarry Smith           barray = (MatScalar*) v + j*bs2;
5590880e062SHong Zhang         } else { /* Here a copy is required */
5600880e062SHong Zhang           if (roworiented) {
5610880e062SHong Zhang             value = v + i*(stepval+bs)*bs + j*bs;
5620880e062SHong Zhang           } else {
5630880e062SHong Zhang             value = v + j*(stepval+bs)*bs + i*bs;
5640880e062SHong Zhang           }
5650880e062SHong Zhang           for (ii=0; ii<bs; ii++,value+=stepval) {
5660880e062SHong Zhang             for (jj=0; jj<bs; jj++) {
5670880e062SHong Zhang               *barray++ = *value++;
5680880e062SHong Zhang             }
5690880e062SHong Zhang           }
5700880e062SHong Zhang           barray -=bs2;
5710880e062SHong Zhang         }
5720880e062SHong Zhang 
5730880e062SHong Zhang         if (in[j] >= cstart && in[j] < cend) {
5740880e062SHong Zhang           col  = in[j] - cstart;
5759566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlocked_SeqSBAIJ_Inlined(baij->A,row,col,barray,addv,im[i],in[j]));
576f7d195e4SLawrence Mitchell         } else if (in[j] < 0) {
577f7d195e4SLawrence Mitchell           continue;
578f7d195e4SLawrence Mitchell         } else {
579f7d195e4SLawrence Mitchell           PetscCheck(in[j] < baij->Nbs,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block indexed column too large %" PetscInt_FMT " max %" PetscInt_FMT,in[j],baij->Nbs-1);
5800880e062SHong Zhang           if (mat->was_assembled) {
5810880e062SHong Zhang             if (!baij->colmap) {
5829566063dSJacob Faibussowitsch               PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
5830880e062SHong Zhang             }
5840880e062SHong Zhang 
5852515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
5860880e062SHong Zhang #if defined(PETSC_USE_CTABLE)
5871302d50aSBarry Smith             { PetscInt data;
5889566063dSJacob Faibussowitsch               PetscCall(PetscTableFind(baij->colmap,in[j]+1,&data));
58908401ef6SPierre Jolivet               PetscCheck((data - 1) % bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
5900880e062SHong Zhang             }
5910880e062SHong Zhang #else
59208401ef6SPierre Jolivet             PetscCheck((baij->colmap[in[j]] - 1) % bs == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
5930880e062SHong Zhang #endif
5940880e062SHong Zhang #endif
5950880e062SHong Zhang #if defined(PETSC_USE_CTABLE)
5969566063dSJacob Faibussowitsch             PetscCall(PetscTableFind(baij->colmap,in[j]+1,&col));
5970880e062SHong Zhang             col  = (col - 1)/bs;
5980880e062SHong Zhang #else
5990880e062SHong Zhang             col = (baij->colmap[in[j]] - 1)/bs;
6000880e062SHong Zhang #endif
6010880e062SHong Zhang             if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
6029566063dSJacob Faibussowitsch               PetscCall(MatDisAssemble_MPISBAIJ(mat));
6030880e062SHong Zhang               col  = in[j];
6040880e062SHong Zhang             }
60526fbe8dcSKarl Rupp           } else col = in[j];
6069566063dSJacob Faibussowitsch           PetscCall(MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B,row,col,barray,addv,im[i],in[j]));
6070880e062SHong Zhang         }
6080880e062SHong Zhang       }
6090880e062SHong Zhang     } else {
6105f80ce2aSJacob Faibussowitsch       PetscCheck(!mat->nooffprocentries,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process block indexed row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
6110880e062SHong Zhang       if (!baij->donotstash) {
6120880e062SHong Zhang         if (roworiented) {
6139566063dSJacob Faibussowitsch           PetscCall(MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i));
6140880e062SHong Zhang         } else {
6159566063dSJacob Faibussowitsch           PetscCall(MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i));
6160880e062SHong Zhang         }
6170880e062SHong Zhang       }
6180880e062SHong Zhang     }
6190880e062SHong Zhang   }
6200880e062SHong Zhang   PetscFunctionReturn(0);
621a30f8f8cSSatish Balay }
622a30f8f8cSSatish Balay 
6231302d50aSBarry Smith PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
624a30f8f8cSSatish Balay {
625f3566a2aSHong Zhang   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
626d0f46423SBarry Smith   PetscInt       bs       = mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
627d0f46423SBarry Smith   PetscInt       bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;
628a30f8f8cSSatish Balay 
629a30f8f8cSSatish Balay   PetscFunctionBegin;
630a30f8f8cSSatish Balay   for (i=0; i<m; i++) {
63154c59aa7SJacob Faibussowitsch     if (idxm[i] < 0) continue; /* negative row */
63254c59aa7SJacob Faibussowitsch     PetscCheck(idxm[i] < mat->rmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,idxm[i],mat->rmap->N-1);
633a30f8f8cSSatish Balay     if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
634a30f8f8cSSatish Balay       row = idxm[i] - bsrstart;
635a30f8f8cSSatish Balay       for (j=0; j<n; j++) {
63654c59aa7SJacob Faibussowitsch         if (idxn[j] < 0) continue; /* negative column */
63754c59aa7SJacob Faibussowitsch         PetscCheck(idxn[j] < mat->cmap->N,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,idxn[j],mat->cmap->N-1);
638a30f8f8cSSatish Balay         if (idxn[j] >= bscstart && idxn[j] < bscend) {
639a30f8f8cSSatish Balay           col  = idxn[j] - bscstart;
6409566063dSJacob Faibussowitsch           PetscCall(MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j));
641a30f8f8cSSatish Balay         } else {
642a30f8f8cSSatish Balay           if (!baij->colmap) {
6439566063dSJacob Faibussowitsch             PetscCall(MatCreateColmap_MPIBAIJ_Private(mat));
644a30f8f8cSSatish Balay           }
645a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
6469566063dSJacob Faibussowitsch           PetscCall(PetscTableFind(baij->colmap,idxn[j]/bs+1,&data));
647a30f8f8cSSatish Balay           data--;
648a30f8f8cSSatish Balay #else
649a30f8f8cSSatish Balay           data = baij->colmap[idxn[j]/bs]-1;
650a30f8f8cSSatish Balay #endif
651a30f8f8cSSatish Balay           if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
652a30f8f8cSSatish Balay           else {
653a30f8f8cSSatish Balay             col  = data + idxn[j]%bs;
6549566063dSJacob Faibussowitsch             PetscCall(MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j));
655a30f8f8cSSatish Balay           }
656a30f8f8cSSatish Balay         }
657a30f8f8cSSatish Balay       }
658f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
659a30f8f8cSSatish Balay   }
660a30f8f8cSSatish Balay   PetscFunctionReturn(0);
661a30f8f8cSSatish Balay }
662a30f8f8cSSatish Balay 
663dfbe8321SBarry Smith PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
664a30f8f8cSSatish Balay {
665a30f8f8cSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
666a30f8f8cSSatish Balay   PetscReal      sum[2],*lnorm2;
667a30f8f8cSSatish Balay 
668a30f8f8cSSatish Balay   PetscFunctionBegin;
669a30f8f8cSSatish Balay   if (baij->size == 1) {
6709566063dSJacob Faibussowitsch     PetscCall(MatNorm(baij->A,type,norm));
671a30f8f8cSSatish Balay   } else {
672a30f8f8cSSatish Balay     if (type == NORM_FROBENIUS) {
6739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2,&lnorm2));
6749566063dSJacob Faibussowitsch       PetscCall(MatNorm(baij->A,type,lnorm2));
675a30f8f8cSSatish Balay       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++;            /* squar power of norm(A) */
6769566063dSJacob Faibussowitsch       PetscCall(MatNorm(baij->B,type,lnorm2));
677a30f8f8cSSatish Balay       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--;             /* squar power of norm(B) */
6781c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(lnorm2,sum,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat)));
6798f1a2a5eSBarry Smith       *norm   = PetscSqrtReal(sum[0] + 2*sum[1]);
6809566063dSJacob Faibussowitsch       PetscCall(PetscFree(lnorm2));
6810b8dc8d2SHong Zhang     } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
6820b8dc8d2SHong Zhang       Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
6830b8dc8d2SHong Zhang       Mat_SeqBAIJ  *bmat=(Mat_SeqBAIJ*)baij->B->data;
6840b8dc8d2SHong Zhang       PetscReal    *rsum,*rsum2,vabs;
685899cda47SBarry Smith       PetscInt     *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
686d0f46423SBarry Smith       PetscInt     brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
6870b8dc8d2SHong Zhang       MatScalar    *v;
6880b8dc8d2SHong Zhang 
6899566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(mat->cmap->N,&rsum,mat->cmap->N,&rsum2));
6909566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(rsum,mat->cmap->N));
6910b8dc8d2SHong Zhang       /* Amat */
6920b8dc8d2SHong Zhang       v = amat->a; jj = amat->j;
6930b8dc8d2SHong Zhang       for (brow=0; brow<mbs; brow++) {
6940b8dc8d2SHong Zhang         grow = bs*(rstart + brow);
6950b8dc8d2SHong Zhang         nz   = amat->i[brow+1] - amat->i[brow];
6960b8dc8d2SHong Zhang         for (bcol=0; bcol<nz; bcol++) {
6970b8dc8d2SHong Zhang           gcol = bs*(rstart + *jj); jj++;
6980b8dc8d2SHong Zhang           for (col=0; col<bs; col++) {
6990b8dc8d2SHong Zhang             for (row=0; row<bs; row++) {
7000b8dc8d2SHong Zhang               vabs            = PetscAbsScalar(*v); v++;
7010b8dc8d2SHong Zhang               rsum[gcol+col] += vabs;
7020b8dc8d2SHong Zhang               /* non-diagonal block */
7030b8dc8d2SHong Zhang               if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
7040b8dc8d2SHong Zhang             }
7050b8dc8d2SHong Zhang           }
7060b8dc8d2SHong Zhang         }
7079566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(nz*bs*bs));
7080b8dc8d2SHong Zhang       }
7090b8dc8d2SHong Zhang       /* Bmat */
7100b8dc8d2SHong Zhang       v = bmat->a; jj = bmat->j;
7110b8dc8d2SHong Zhang       for (brow=0; brow<mbs; brow++) {
7120b8dc8d2SHong Zhang         grow = bs*(rstart + brow);
7130b8dc8d2SHong Zhang         nz = bmat->i[brow+1] - bmat->i[brow];
7140b8dc8d2SHong Zhang         for (bcol=0; bcol<nz; bcol++) {
7150b8dc8d2SHong Zhang           gcol = bs*garray[*jj]; jj++;
7160b8dc8d2SHong Zhang           for (col=0; col<bs; col++) {
7170b8dc8d2SHong Zhang             for (row=0; row<bs; row++) {
7180b8dc8d2SHong Zhang               vabs            = PetscAbsScalar(*v); v++;
7190b8dc8d2SHong Zhang               rsum[gcol+col] += vabs;
7200b8dc8d2SHong Zhang               rsum[grow+row] += vabs;
7210b8dc8d2SHong Zhang             }
7220b8dc8d2SHong Zhang           }
7230b8dc8d2SHong Zhang         }
7249566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(nz*bs*bs));
7250b8dc8d2SHong Zhang       }
7261c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat)));
7270b8dc8d2SHong Zhang       *norm = 0.0;
728d0f46423SBarry Smith       for (col=0; col<mat->cmap->N; col++) {
7290b8dc8d2SHong Zhang         if (rsum2[col] > *norm) *norm = rsum2[col];
7300b8dc8d2SHong Zhang       }
7319566063dSJacob Faibussowitsch       PetscCall(PetscFree2(rsum,rsum2));
732f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
733a30f8f8cSSatish Balay   }
734a30f8f8cSSatish Balay   PetscFunctionReturn(0);
735a30f8f8cSSatish Balay }
736a30f8f8cSSatish Balay 
737dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
738a30f8f8cSSatish Balay {
739a30f8f8cSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
7401302d50aSBarry Smith   PetscInt       nstash,reallocs;
741a30f8f8cSSatish Balay 
742a30f8f8cSSatish Balay   PetscFunctionBegin;
74326fbe8dcSKarl Rupp   if (baij->donotstash || mat->nooffprocentries) PetscFunctionReturn(0);
744a30f8f8cSSatish Balay 
7459566063dSJacob Faibussowitsch   PetscCall(MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range));
7469566063dSJacob Faibussowitsch   PetscCall(MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs));
7479566063dSJacob Faibussowitsch   PetscCall(MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs));
7489566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat,"Stash has %" PetscInt_FMT " entries,uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs));
7499566063dSJacob Faibussowitsch   PetscCall(MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs));
7509566063dSJacob Faibussowitsch   PetscCall(PetscInfo(mat,"Block-Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs));
751a30f8f8cSSatish Balay   PetscFunctionReturn(0);
752a30f8f8cSSatish Balay }
753a30f8f8cSSatish Balay 
754dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
755a30f8f8cSSatish Balay {
756a30f8f8cSSatish Balay   Mat_MPISBAIJ   *baij=(Mat_MPISBAIJ*)mat->data;
757a30f8f8cSSatish Balay   Mat_SeqSBAIJ   *a   =(Mat_SeqSBAIJ*)baij->A->data;
75813f74950SBarry Smith   PetscInt       i,j,rstart,ncols,flg,bs2=baij->bs2;
759e44c0bd4SBarry Smith   PetscInt       *row,*col;
760ace3abfcSBarry Smith   PetscBool      other_disassembled;
76113f74950SBarry Smith   PetscMPIInt    n;
762ace3abfcSBarry Smith   PetscBool      r1,r2,r3;
763a30f8f8cSSatish Balay   MatScalar      *val;
764a30f8f8cSSatish Balay 
76591c97fd4SSatish Balay   /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
766a30f8f8cSSatish Balay   PetscFunctionBegin;
7674cb17eb5SBarry Smith   if (!baij->donotstash &&  !mat->nooffprocentries) {
768a30f8f8cSSatish Balay     while (1) {
7699566063dSJacob Faibussowitsch       PetscCall(MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg));
770a30f8f8cSSatish Balay       if (!flg) break;
771a30f8f8cSSatish Balay 
772a30f8f8cSSatish Balay       for (i=0; i<n;) {
773a30f8f8cSSatish Balay         /* Now identify the consecutive vals belonging to the same row */
77426fbe8dcSKarl Rupp         for (j=i,rstart=row[j]; j<n; j++) {
77526fbe8dcSKarl Rupp           if (row[j] != rstart) break;
77626fbe8dcSKarl Rupp         }
777a30f8f8cSSatish Balay         if (j < n) ncols = j-i;
778a30f8f8cSSatish Balay         else       ncols = n-i;
779a30f8f8cSSatish Balay         /* Now assemble all these values with a single function call */
7809566063dSJacob Faibussowitsch         PetscCall(MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode));
781a30f8f8cSSatish Balay         i    = j;
782a30f8f8cSSatish Balay       }
783a30f8f8cSSatish Balay     }
7849566063dSJacob Faibussowitsch     PetscCall(MatStashScatterEnd_Private(&mat->stash));
785a30f8f8cSSatish Balay     /* Now process the block-stash. Since the values are stashed column-oriented,
786a30f8f8cSSatish Balay        set the roworiented flag to column oriented, and after MatSetValues()
787a30f8f8cSSatish Balay        restore the original flags */
788a30f8f8cSSatish Balay     r1 = baij->roworiented;
789a30f8f8cSSatish Balay     r2 = a->roworiented;
79091c97fd4SSatish Balay     r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
79126fbe8dcSKarl Rupp 
792a30f8f8cSSatish Balay     baij->roworiented = PETSC_FALSE;
793a30f8f8cSSatish Balay     a->roworiented    = PETSC_FALSE;
79426fbe8dcSKarl Rupp 
79591c97fd4SSatish Balay     ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
796a30f8f8cSSatish Balay     while (1) {
7979566063dSJacob Faibussowitsch       PetscCall(MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg));
798a30f8f8cSSatish Balay       if (!flg) break;
799a30f8f8cSSatish Balay 
800a30f8f8cSSatish Balay       for (i=0; i<n;) {
801a30f8f8cSSatish Balay         /* Now identify the consecutive vals belonging to the same row */
80226fbe8dcSKarl Rupp         for (j=i,rstart=row[j]; j<n; j++) {
80326fbe8dcSKarl Rupp           if (row[j] != rstart) break;
80426fbe8dcSKarl Rupp         }
805a30f8f8cSSatish Balay         if (j < n) ncols = j-i;
806a30f8f8cSSatish Balay         else       ncols = n-i;
8079566063dSJacob Faibussowitsch         PetscCall(MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode));
808a30f8f8cSSatish Balay         i    = j;
809a30f8f8cSSatish Balay       }
810a30f8f8cSSatish Balay     }
8119566063dSJacob Faibussowitsch     PetscCall(MatStashScatterEnd_Private(&mat->bstash));
81226fbe8dcSKarl Rupp 
813a30f8f8cSSatish Balay     baij->roworiented = r1;
814a30f8f8cSSatish Balay     a->roworiented    = r2;
81526fbe8dcSKarl Rupp 
81691c97fd4SSatish Balay     ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
817a30f8f8cSSatish Balay   }
818a30f8f8cSSatish Balay 
8199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(baij->A,mode));
8209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(baij->A,mode));
821a30f8f8cSSatish Balay 
822a30f8f8cSSatish Balay   /* determine if any processor has disassembled, if so we must
8236aad120cSJose E. Roman      also disassemble ourselves, in order that we may reassemble. */
824a30f8f8cSSatish Balay   /*
825a30f8f8cSSatish Balay      if nonzero structure of submatrix B cannot change then we know that
826a30f8f8cSSatish Balay      no processor disassembled thus we can skip this stuff
827a30f8f8cSSatish Balay   */
828a30f8f8cSSatish Balay   if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
8291c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat)));
830a30f8f8cSSatish Balay     if (mat->was_assembled && !other_disassembled) {
8319566063dSJacob Faibussowitsch       PetscCall(MatDisAssemble_MPISBAIJ(mat));
832a30f8f8cSSatish Balay     }
833a30f8f8cSSatish Balay   }
834a30f8f8cSSatish Balay 
835a30f8f8cSSatish Balay   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
8369566063dSJacob Faibussowitsch     PetscCall(MatSetUpMultiply_MPISBAIJ(mat)); /* setup Mvctx and sMvctx */
837a30f8f8cSSatish Balay   }
8389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(baij->B,mode));
8399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(baij->B,mode));
840a30f8f8cSSatish Balay 
8419566063dSJacob Faibussowitsch   PetscCall(PetscFree2(baij->rowvalues,baij->rowindices));
84226fbe8dcSKarl Rupp 
843f4259b30SLisandro Dalcin   baij->rowvalues = NULL;
8444f9cfa9eSBarry Smith 
8454f9cfa9eSBarry Smith   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
8464f9cfa9eSBarry Smith   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
847e56f5c9eSBarry Smith     PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
8481c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat)));
849e56f5c9eSBarry Smith   }
850a30f8f8cSSatish Balay   PetscFunctionReturn(0);
851a30f8f8cSSatish Balay }
852a30f8f8cSSatish Balay 
853dd6ea824SBarry Smith extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
8549804daf3SBarry Smith #include <petscdraw.h>
8556849ba73SBarry Smith static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
856a30f8f8cSSatish Balay {
857a30f8f8cSSatish Balay   Mat_MPISBAIJ      *baij = (Mat_MPISBAIJ*)mat->data;
858d0f46423SBarry Smith   PetscInt          bs   = mat->rmap->bs;
8597da1fb6eSBarry Smith   PetscMPIInt       rank = baij->rank;
860ace3abfcSBarry Smith   PetscBool         iascii,isdraw;
861b0a32e0cSBarry Smith   PetscViewer       sviewer;
862f3ef73ceSBarry Smith   PetscViewerFormat format;
863a30f8f8cSSatish Balay 
864a30f8f8cSSatish Balay   PetscFunctionBegin;
8659566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
8669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
86732077d6dSBarry Smith   if (iascii) {
8689566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
869456192e2SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
870a30f8f8cSSatish Balay       MatInfo info;
8719566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank));
8729566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(mat,MAT_LOCAL,&info));
8739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
8749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " bs %" PetscInt_FMT " mem %g\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(double)info.memory));
8759566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(baij->A,MAT_LOCAL,&info));
8769566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %" PetscInt_FMT " \n",rank,(PetscInt)info.nz_used));
8779566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(baij->B,MAT_LOCAL,&info));
8789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %" PetscInt_FMT " \n",rank,(PetscInt)info.nz_used));
8799566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
8809566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
8819566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n"));
8829566063dSJacob Faibussowitsch       PetscCall(VecScatterView(baij->Mvctx,viewer));
883a30f8f8cSSatish Balay       PetscFunctionReturn(0);
884fb9695e5SSatish Balay     } else if (format == PETSC_VIEWER_ASCII_INFO) {
8859566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  block size is %" PetscInt_FMT "\n",bs));
886a30f8f8cSSatish Balay       PetscFunctionReturn(0);
887c1490034SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
888c1490034SHong Zhang       PetscFunctionReturn(0);
889a30f8f8cSSatish Balay     }
890a30f8f8cSSatish Balay   }
891a30f8f8cSSatish Balay 
892a30f8f8cSSatish Balay   if (isdraw) {
893b0a32e0cSBarry Smith     PetscDraw draw;
894ace3abfcSBarry Smith     PetscBool isnull;
8959566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
8969566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw,&isnull));
89745f3bb6eSLisandro Dalcin     if (isnull) PetscFunctionReturn(0);
898a30f8f8cSSatish Balay   }
899a30f8f8cSSatish Balay 
9007da1fb6eSBarry Smith   {
901a30f8f8cSSatish Balay     /* assemble the entire matrix onto first processor. */
902a30f8f8cSSatish Balay     Mat          A;
90365d70643SHong Zhang     Mat_SeqSBAIJ *Aloc;
90465d70643SHong Zhang     Mat_SeqBAIJ  *Bloc;
905d0f46423SBarry Smith     PetscInt     M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
906a30f8f8cSSatish Balay     MatScalar    *a;
9073e219373SBarry Smith     const char   *matname;
908a30f8f8cSSatish Balay 
909f204ca49SKris Buschelman     /* Should this be the same type as mat? */
9109566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat),&A));
911dd400576SPatrick Sanan     if (rank == 0) {
9129566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A,M,N,M,N));
913a30f8f8cSSatish Balay     } else {
9149566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A,0,0,M,N));
915a30f8f8cSSatish Balay     }
9169566063dSJacob Faibussowitsch     PetscCall(MatSetType(A,MATMPISBAIJ));
9179566063dSJacob Faibussowitsch     PetscCall(MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL));
9189566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE));
9199566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)A));
920a30f8f8cSSatish Balay 
921a30f8f8cSSatish Balay     /* copy over the A part */
92265d70643SHong Zhang     Aloc = (Mat_SeqSBAIJ*)baij->A->data;
923a30f8f8cSSatish Balay     ai   = Aloc->i; aj = Aloc->j; a = Aloc->a;
9249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs,&rvals));
925a30f8f8cSSatish Balay 
926a30f8f8cSSatish Balay     for (i=0; i<mbs; i++) {
927e9f7bc9eSHong Zhang       rvals[0] = bs*(baij->rstartbs + i);
92826fbe8dcSKarl Rupp       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
929a30f8f8cSSatish Balay       for (j=ai[i]; j<ai[i+1]; j++) {
930e9f7bc9eSHong Zhang         col = (baij->cstartbs+aj[j])*bs;
931a30f8f8cSSatish Balay         for (k=0; k<bs; k++) {
9329566063dSJacob Faibussowitsch           PetscCall(MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES));
93326fbe8dcSKarl Rupp           col++;
93426fbe8dcSKarl Rupp           a += bs;
935a30f8f8cSSatish Balay         }
936a30f8f8cSSatish Balay       }
937a30f8f8cSSatish Balay     }
938a30f8f8cSSatish Balay     /* copy over the B part */
93965d70643SHong Zhang     Bloc = (Mat_SeqBAIJ*)baij->B->data;
94065d70643SHong Zhang     ai   = Bloc->i; aj = Bloc->j; a = Bloc->a;
941a30f8f8cSSatish Balay     for (i=0; i<mbs; i++) {
942e9f7bc9eSHong Zhang 
943e9f7bc9eSHong Zhang       rvals[0] = bs*(baij->rstartbs + i);
94426fbe8dcSKarl Rupp       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
945a30f8f8cSSatish Balay       for (j=ai[i]; j<ai[i+1]; j++) {
946a30f8f8cSSatish Balay         col = baij->garray[aj[j]]*bs;
947a30f8f8cSSatish Balay         for (k=0; k<bs; k++) {
9489566063dSJacob Faibussowitsch           PetscCall(MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES));
94926fbe8dcSKarl Rupp           col++;
95026fbe8dcSKarl Rupp           a += bs;
951a30f8f8cSSatish Balay         }
952a30f8f8cSSatish Balay       }
953a30f8f8cSSatish Balay     }
9549566063dSJacob Faibussowitsch     PetscCall(PetscFree(rvals));
9559566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
9569566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
957a30f8f8cSSatish Balay     /*
958a30f8f8cSSatish Balay        Everyone has to call to draw the matrix since the graphics waits are
959b0a32e0cSBarry Smith        synchronized across all processors that share the PetscDraw object
960a30f8f8cSSatish Balay     */
9619566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
9629566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)mat,&matname));
963dd400576SPatrick Sanan     if (rank == 0) {
9649566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,matname));
9659566063dSJacob Faibussowitsch       PetscCall(MatView_SeqSBAIJ(((Mat_MPISBAIJ*)(A->data))->A,sviewer));
966a30f8f8cSSatish Balay     }
9679566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
9689566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
9699566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
970a30f8f8cSSatish Balay   }
971a30f8f8cSSatish Balay   PetscFunctionReturn(0);
972a30f8f8cSSatish Balay }
973a30f8f8cSSatish Balay 
974618cc2edSLisandro Dalcin /* Used for both MPIBAIJ and MPISBAIJ matrices */
975618cc2edSLisandro Dalcin #define MatView_MPISBAIJ_Binary MatView_MPIBAIJ_Binary
976d1654148SHong Zhang 
977dfbe8321SBarry Smith PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
978a30f8f8cSSatish Balay {
979ace3abfcSBarry Smith   PetscBool      iascii,isdraw,issocket,isbinary;
980a30f8f8cSSatish Balay 
981a30f8f8cSSatish Balay   PetscFunctionBegin;
9829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
9839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
9849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket));
9859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
986d1654148SHong Zhang   if (iascii || isdraw || issocket) {
9879566063dSJacob Faibussowitsch     PetscCall(MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer));
9881baa6e33SBarry Smith   } else if (isbinary) PetscCall(MatView_MPISBAIJ_Binary(mat,viewer));
989a30f8f8cSSatish Balay   PetscFunctionReturn(0);
990a30f8f8cSSatish Balay }
991a30f8f8cSSatish Balay 
992dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
993a30f8f8cSSatish Balay {
994a30f8f8cSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
995a30f8f8cSSatish Balay 
996a30f8f8cSSatish Balay   PetscFunctionBegin;
997a30f8f8cSSatish Balay #if defined(PETSC_USE_LOG)
998c0aa6a63SJacob Faibussowitsch   PetscLogObjectState((PetscObject)mat,"Rows=%" PetscInt_FMT ",Cols=%" PetscInt_FMT,mat->rmap->N,mat->cmap->N);
999a30f8f8cSSatish Balay #endif
10009566063dSJacob Faibussowitsch   PetscCall(MatStashDestroy_Private(&mat->stash));
10019566063dSJacob Faibussowitsch   PetscCall(MatStashDestroy_Private(&mat->bstash));
10029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&baij->A));
10039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&baij->B));
1004a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
10059566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&baij->colmap));
1006a30f8f8cSSatish Balay #else
10079566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->colmap));
1008a30f8f8cSSatish Balay #endif
10099566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->garray));
10109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->lvec));
10119566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&baij->Mvctx));
10129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->slvec0));
10139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->slvec0b));
10149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->slvec1));
10159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->slvec1a));
10169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->slvec1b));
10179566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&baij->sMvctx));
10189566063dSJacob Faibussowitsch   PetscCall(PetscFree2(baij->rowvalues,baij->rowindices));
10199566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->barray));
10209566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->hd));
10219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->diag));
10229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->bb1));
10239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->xx1));
1024ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
10259566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->setvaluescopy));
1026a30f8f8cSSatish Balay #endif
10279566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->in_loc));
10289566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->v_loc));
10299566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->rangebs));
10309566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
1031901853e0SKris Buschelman 
10329566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat,NULL));
10339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL));
10349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL));
10359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL));
10369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocationCSR_C",NULL));
10376214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
10389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_elemental_C",NULL));
10396214f412SHong Zhang #endif
1040d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
10419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_scalapack_C",NULL));
1042d24d4204SJose E. Roman #endif
10439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpiaij_C",NULL));
10449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpibaij_C",NULL));
1045a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1046a30f8f8cSSatish Balay }
1047a30f8f8cSSatish Balay 
1048547795f9SHong Zhang PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
1049547795f9SHong Zhang {
1050547795f9SHong Zhang   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1051eb1ec7c1SStefano Zampini   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
10526de40e93SBarry Smith   PetscScalar       *from;
10536de40e93SBarry Smith   const PetscScalar *x;
1054547795f9SHong Zhang 
1055547795f9SHong Zhang   PetscFunctionBegin;
1056547795f9SHong Zhang   /* diagonal part */
10579566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->mult)(a->A,xx,a->slvec1a));
10589566063dSJacob Faibussowitsch   PetscCall(VecSet(a->slvec1b,0.0));
1059547795f9SHong Zhang 
1060547795f9SHong Zhang   /* subdiagonal part */
10615f80ce2aSJacob Faibussowitsch   PetscCheck(a->B->ops->multhermitiantranspose,PetscObjectComm((PetscObject)a->B),PETSC_ERR_SUP,"Not for type %s",((PetscObject)a->B)->type_name);
10629566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b));
1063547795f9SHong Zhang 
1064547795f9SHong Zhang   /* copy x into the vec slvec0 */
10659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a->slvec0,&from));
10669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx,&x));
1067547795f9SHong Zhang 
10689566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(from,x,bs*mbs));
10699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a->slvec0,&from));
10709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx,&x));
1071547795f9SHong Zhang 
10729566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
10739566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
1074547795f9SHong Zhang   /* supperdiagonal part */
10759566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy));
1076547795f9SHong Zhang   PetscFunctionReturn(0);
1077547795f9SHong Zhang }
1078547795f9SHong Zhang 
1079dfbe8321SBarry Smith PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
1080a9d4b620SHong Zhang {
1081a9d4b620SHong Zhang   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1082eb1ec7c1SStefano Zampini   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1083d9ca1df4SBarry Smith   PetscScalar       *from;
1084d9ca1df4SBarry Smith   const PetscScalar *x;
1085a9d4b620SHong Zhang 
1086a9d4b620SHong Zhang   PetscFunctionBegin;
1087a9d4b620SHong Zhang   /* diagonal part */
10889566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->mult)(a->A,xx,a->slvec1a));
10899566063dSJacob Faibussowitsch   PetscCall(VecSet(a->slvec1b,0.0));
1090a9d4b620SHong Zhang 
1091a9d4b620SHong Zhang   /* subdiagonal part */
10929566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B,xx,a->slvec0b));
1093fc165ae2SBarry Smith 
1094a9d4b620SHong Zhang   /* copy x into the vec slvec0 */
10959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a->slvec0,&from));
10969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx,&x));
1097a9d4b620SHong Zhang 
10989566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(from,x,bs*mbs));
10999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a->slvec0,&from));
11009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx,&x));
1101fc165ae2SBarry Smith 
11029566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
11039566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
1104a9d4b620SHong Zhang   /* supperdiagonal part */
11059566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy));
1106a9d4b620SHong Zhang   PetscFunctionReturn(0);
1107a9d4b620SHong Zhang }
1108a9d4b620SHong Zhang 
1109eb1ec7c1SStefano Zampini PetscErrorCode MatMultAdd_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy,Vec zz)
1110eb1ec7c1SStefano Zampini {
1111eb1ec7c1SStefano Zampini   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1112eb1ec7c1SStefano Zampini   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1113eb1ec7c1SStefano Zampini   PetscScalar       *from,zero=0.0;
1114eb1ec7c1SStefano Zampini   const PetscScalar *x;
1115eb1ec7c1SStefano Zampini 
1116eb1ec7c1SStefano Zampini   PetscFunctionBegin;
1117eb1ec7c1SStefano Zampini   /* diagonal part */
11189566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a));
11199566063dSJacob Faibussowitsch   PetscCall(VecSet(a->slvec1b,zero));
1120eb1ec7c1SStefano Zampini 
1121eb1ec7c1SStefano Zampini   /* subdiagonal part */
11225f80ce2aSJacob Faibussowitsch   PetscCheck(a->B->ops->multhermitiantranspose,PetscObjectComm((PetscObject)a->B),PETSC_ERR_SUP,"Not for type %s",((PetscObject)a->B)->type_name);
11239566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b));
1124eb1ec7c1SStefano Zampini 
1125eb1ec7c1SStefano Zampini   /* copy x into the vec slvec0 */
11269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a->slvec0,&from));
11279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx,&x));
11289566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(from,x,bs*mbs));
11299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a->slvec0,&from));
1130eb1ec7c1SStefano Zampini 
11319566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
11329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx,&x));
11339566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
1134eb1ec7c1SStefano Zampini 
1135eb1ec7c1SStefano Zampini   /* supperdiagonal part */
11369566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz));
1137eb1ec7c1SStefano Zampini   PetscFunctionReturn(0);
1138eb1ec7c1SStefano Zampini }
1139eb1ec7c1SStefano Zampini 
1140dfbe8321SBarry Smith PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1141a30f8f8cSSatish Balay {
1142de8b6608SHong Zhang   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1143d0f46423SBarry Smith   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1144d9ca1df4SBarry Smith   PetscScalar       *from,zero=0.0;
1145d9ca1df4SBarry Smith   const PetscScalar *x;
1146a9d4b620SHong Zhang 
1147a9d4b620SHong Zhang   PetscFunctionBegin;
1148a9d4b620SHong Zhang   /* diagonal part */
11499566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a));
11509566063dSJacob Faibussowitsch   PetscCall(VecSet(a->slvec1b,zero));
1151a9d4b620SHong Zhang 
1152a9d4b620SHong Zhang   /* subdiagonal part */
11539566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B,xx,a->slvec0b));
1154a9d4b620SHong Zhang 
1155a9d4b620SHong Zhang   /* copy x into the vec slvec0 */
11569566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a->slvec0,&from));
11579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx,&x));
11589566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(from,x,bs*mbs));
11599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a->slvec0,&from));
1160a9d4b620SHong Zhang 
11619566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
11629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx,&x));
11639566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD));
1164a9d4b620SHong Zhang 
1165a9d4b620SHong Zhang   /* supperdiagonal part */
11669566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz));
1167a9d4b620SHong Zhang   PetscFunctionReturn(0);
1168a9d4b620SHong Zhang }
1169a9d4b620SHong Zhang 
1170a30f8f8cSSatish Balay /*
1171a30f8f8cSSatish Balay   This only works correctly for square matrices where the subblock A->A is the
1172a30f8f8cSSatish Balay    diagonal block
1173a30f8f8cSSatish Balay */
1174dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1175a30f8f8cSSatish Balay {
1176a30f8f8cSSatish Balay   Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1177a30f8f8cSSatish Balay 
1178a30f8f8cSSatish Balay   PetscFunctionBegin;
117908401ef6SPierre Jolivet   /* PetscCheck(a->rmap->N == a->cmap->N,PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
11809566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(a->A,v));
1181a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1182a30f8f8cSSatish Balay }
1183a30f8f8cSSatish Balay 
1184f4df32b1SMatthew Knepley PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1185a30f8f8cSSatish Balay {
1186a30f8f8cSSatish Balay   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1187a30f8f8cSSatish Balay 
1188a30f8f8cSSatish Balay   PetscFunctionBegin;
11899566063dSJacob Faibussowitsch   PetscCall(MatScale(a->A,aa));
11909566063dSJacob Faibussowitsch   PetscCall(MatScale(a->B,aa));
1191a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1192a30f8f8cSSatish Balay }
1193a30f8f8cSSatish Balay 
11941302d50aSBarry Smith PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1195a30f8f8cSSatish Balay {
1196d0d4cfc2SHong Zhang   Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)matin->data;
1197d0d4cfc2SHong Zhang   PetscScalar  *vworkA,*vworkB,**pvA,**pvB,*v_p;
1198d0f46423SBarry Smith   PetscInt     bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1199d0f46423SBarry Smith   PetscInt     nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1200899cda47SBarry Smith   PetscInt    *cmap,*idx_p,cstart = mat->rstartbs;
1201d0d4cfc2SHong Zhang 
1202a30f8f8cSSatish Balay   PetscFunctionBegin;
12035f80ce2aSJacob Faibussowitsch   PetscCheck(!mat->getrowactive,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
1204d0d4cfc2SHong Zhang   mat->getrowactive = PETSC_TRUE;
1205d0d4cfc2SHong Zhang 
1206d0d4cfc2SHong Zhang   if (!mat->rowvalues && (idx || v)) {
1207d0d4cfc2SHong Zhang     /*
1208d0d4cfc2SHong Zhang         allocate enough space to hold information from the longest row.
1209d0d4cfc2SHong Zhang     */
1210d0d4cfc2SHong Zhang     Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1211d0d4cfc2SHong Zhang     Mat_SeqBAIJ  *Ba = (Mat_SeqBAIJ*)mat->B->data;
1212d0d4cfc2SHong Zhang     PetscInt     max = 1,mbs = mat->mbs,tmp;
1213d0d4cfc2SHong Zhang     for (i=0; i<mbs; i++) {
1214d0d4cfc2SHong Zhang       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
121526fbe8dcSKarl Rupp       if (max < tmp) max = tmp;
1216d0d4cfc2SHong Zhang     }
12179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices));
1218d0d4cfc2SHong Zhang   }
1219d0d4cfc2SHong Zhang 
12205f80ce2aSJacob Faibussowitsch   PetscCheck(row >= brstart && row < brend,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1221d0d4cfc2SHong Zhang   lrow = row - brstart;  /* local row index */
1222d0d4cfc2SHong Zhang 
1223d0d4cfc2SHong Zhang   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1224f4259b30SLisandro Dalcin   if (!v)   {pvA = NULL; pvB = NULL;}
1225f4259b30SLisandro Dalcin   if (!idx) {pcA = NULL; if (!v) pcB = NULL;}
12269566063dSJacob Faibussowitsch   PetscCall((*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA));
12279566063dSJacob Faibussowitsch   PetscCall((*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB));
1228d0d4cfc2SHong Zhang   nztot = nzA + nzB;
1229d0d4cfc2SHong Zhang 
1230d0d4cfc2SHong Zhang   cmap = mat->garray;
1231d0d4cfc2SHong Zhang   if (v  || idx) {
1232d0d4cfc2SHong Zhang     if (nztot) {
1233d0d4cfc2SHong Zhang       /* Sort by increasing column numbers, assuming A and B already sorted */
1234d0d4cfc2SHong Zhang       PetscInt imark = -1;
1235d0d4cfc2SHong Zhang       if (v) {
1236d0d4cfc2SHong Zhang         *v = v_p = mat->rowvalues;
1237d0d4cfc2SHong Zhang         for (i=0; i<nzB; i++) {
1238d0d4cfc2SHong Zhang           if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1239d0d4cfc2SHong Zhang           else break;
1240d0d4cfc2SHong Zhang         }
1241d0d4cfc2SHong Zhang         imark = i;
1242d0d4cfc2SHong Zhang         for (i=0; i<nzA; i++)     v_p[imark+i] = vworkA[i];
1243d0d4cfc2SHong Zhang         for (i=imark; i<nzB; i++) v_p[nzA+i]   = vworkB[i];
1244d0d4cfc2SHong Zhang       }
1245d0d4cfc2SHong Zhang       if (idx) {
1246d0d4cfc2SHong Zhang         *idx = idx_p = mat->rowindices;
1247d0d4cfc2SHong Zhang         if (imark > -1) {
1248d0d4cfc2SHong Zhang           for (i=0; i<imark; i++) {
1249d0d4cfc2SHong Zhang             idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1250d0d4cfc2SHong Zhang           }
1251d0d4cfc2SHong Zhang         } else {
1252d0d4cfc2SHong Zhang           for (i=0; i<nzB; i++) {
125326fbe8dcSKarl Rupp             if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1254d0d4cfc2SHong Zhang             else break;
1255d0d4cfc2SHong Zhang           }
1256d0d4cfc2SHong Zhang           imark = i;
1257d0d4cfc2SHong Zhang         }
1258d0d4cfc2SHong Zhang         for (i=0; i<nzA; i++)     idx_p[imark+i] = cstart*bs + cworkA[i];
1259d0d4cfc2SHong Zhang         for (i=imark; i<nzB; i++) idx_p[nzA+i]   = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1260d0d4cfc2SHong Zhang       }
1261d0d4cfc2SHong Zhang     } else {
1262f4259b30SLisandro Dalcin       if (idx) *idx = NULL;
1263f4259b30SLisandro Dalcin       if (v)   *v   = NULL;
1264d0d4cfc2SHong Zhang     }
1265d0d4cfc2SHong Zhang   }
1266d0d4cfc2SHong Zhang   *nz  = nztot;
12679566063dSJacob Faibussowitsch   PetscCall((*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA));
12689566063dSJacob Faibussowitsch   PetscCall((*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB));
1269a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1270a30f8f8cSSatish Balay }
1271a30f8f8cSSatish Balay 
12721302d50aSBarry Smith PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1273a30f8f8cSSatish Balay {
1274a30f8f8cSSatish Balay   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1275a30f8f8cSSatish Balay 
1276a30f8f8cSSatish Balay   PetscFunctionBegin;
12775f80ce2aSJacob Faibussowitsch   PetscCheck(baij->getrowactive,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1278a30f8f8cSSatish Balay   baij->getrowactive = PETSC_FALSE;
1279a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1280a30f8f8cSSatish Balay }
1281a30f8f8cSSatish Balay 
1282d0d4cfc2SHong Zhang PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1283d0d4cfc2SHong Zhang {
1284d0d4cfc2SHong Zhang   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1285d0d4cfc2SHong Zhang   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1286d0d4cfc2SHong Zhang 
1287d0d4cfc2SHong Zhang   PetscFunctionBegin;
1288d0d4cfc2SHong Zhang   aA->getrow_utriangular = PETSC_TRUE;
1289d0d4cfc2SHong Zhang   PetscFunctionReturn(0);
1290d0d4cfc2SHong Zhang }
1291d0d4cfc2SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1292d0d4cfc2SHong Zhang {
1293d0d4cfc2SHong Zhang   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1294d0d4cfc2SHong Zhang   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1295d0d4cfc2SHong Zhang 
1296d0d4cfc2SHong Zhang   PetscFunctionBegin;
1297d0d4cfc2SHong Zhang   aA->getrow_utriangular = PETSC_FALSE;
1298d0d4cfc2SHong Zhang   PetscFunctionReturn(0);
1299d0d4cfc2SHong Zhang }
1300d0d4cfc2SHong Zhang 
13012726fb6dSPierre Jolivet PetscErrorCode MatConjugate_MPISBAIJ(Mat mat)
13022726fb6dSPierre Jolivet {
13035f80ce2aSJacob Faibussowitsch   PetscFunctionBegin;
13045f80ce2aSJacob Faibussowitsch   if (PetscDefined(USE_COMPLEX)) {
13052726fb6dSPierre Jolivet     Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)mat->data;
13062726fb6dSPierre Jolivet 
13079566063dSJacob Faibussowitsch     PetscCall(MatConjugate(a->A));
13089566063dSJacob Faibussowitsch     PetscCall(MatConjugate(a->B));
13095f80ce2aSJacob Faibussowitsch   }
13102726fb6dSPierre Jolivet   PetscFunctionReturn(0);
13112726fb6dSPierre Jolivet }
13122726fb6dSPierre Jolivet 
131399cafbc1SBarry Smith PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
131499cafbc1SBarry Smith {
131599cafbc1SBarry Smith   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
131699cafbc1SBarry Smith 
131799cafbc1SBarry Smith   PetscFunctionBegin;
13189566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->A));
13199566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->B));
132099cafbc1SBarry Smith   PetscFunctionReturn(0);
132199cafbc1SBarry Smith }
132299cafbc1SBarry Smith 
132399cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
132499cafbc1SBarry Smith {
132599cafbc1SBarry Smith   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
132699cafbc1SBarry Smith 
132799cafbc1SBarry Smith   PetscFunctionBegin;
13289566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->A));
13299566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->B));
133099cafbc1SBarry Smith   PetscFunctionReturn(0);
133199cafbc1SBarry Smith }
133299cafbc1SBarry Smith 
13337dae84e0SHong Zhang /* Check if isrow is a subset of iscol_local, called by MatCreateSubMatrix_MPISBAIJ()
133436032a97SHong Zhang    Input: isrow       - distributed(parallel),
133536032a97SHong Zhang           iscol_local - locally owned (seq)
133636032a97SHong Zhang */
133736032a97SHong Zhang PetscErrorCode ISEqual_private(IS isrow,IS iscol_local,PetscBool  *flg)
13388f46ffcaSHong Zhang {
13398f46ffcaSHong Zhang   PetscInt       sz1,sz2,*a1,*a2,i,j,k,nmatch;
13408f46ffcaSHong Zhang   const PetscInt *ptr1,*ptr2;
134136032a97SHong Zhang 
134236032a97SHong Zhang   PetscFunctionBegin;
13439566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow,&sz1));
13449566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscol_local,&sz2));
13451098a8e8SHong Zhang   if (sz1 > sz2) {
13461098a8e8SHong Zhang     *flg = PETSC_FALSE;
13471098a8e8SHong Zhang     PetscFunctionReturn(0);
13481098a8e8SHong Zhang   }
13498f46ffcaSHong Zhang 
13509566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow,&ptr1));
13519566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscol_local,&ptr2));
13528f46ffcaSHong Zhang 
13539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(sz1,&a1));
13549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(sz2,&a2));
13559566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a1,ptr1,sz1));
13569566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a2,ptr2,sz2));
13579566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(sz1,a1));
13589566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(sz2,a2));
13598f46ffcaSHong Zhang 
13608f46ffcaSHong Zhang   nmatch=0;
13618f46ffcaSHong Zhang   k     = 0;
13628f46ffcaSHong Zhang   for (i=0; i<sz1; i++) {
13638f46ffcaSHong Zhang     for (j=k; j<sz2; j++) {
13648f46ffcaSHong Zhang       if (a1[i] == a2[j]) {
13658f46ffcaSHong Zhang         k = j; nmatch++;
13668f46ffcaSHong Zhang         break;
13678f46ffcaSHong Zhang       }
13688f46ffcaSHong Zhang     }
13698f46ffcaSHong Zhang   }
13709566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow,&ptr1));
13719566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscol_local,&ptr2));
13729566063dSJacob Faibussowitsch   PetscCall(PetscFree(a1));
13739566063dSJacob Faibussowitsch   PetscCall(PetscFree(a2));
13741098a8e8SHong Zhang   if (nmatch < sz1) {
13751098a8e8SHong Zhang     *flg = PETSC_FALSE;
13761098a8e8SHong Zhang   } else {
13771098a8e8SHong Zhang     *flg = PETSC_TRUE;
13781098a8e8SHong Zhang   }
137936032a97SHong Zhang   PetscFunctionReturn(0);
13808f46ffcaSHong Zhang }
138136032a97SHong Zhang 
13827dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrix_MPISBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
138336032a97SHong Zhang {
138436032a97SHong Zhang   IS             iscol_local;
138536032a97SHong Zhang   PetscInt       csize;
138636032a97SHong Zhang   PetscBool      isequal;
138736032a97SHong Zhang 
138836032a97SHong Zhang   PetscFunctionBegin;
13899566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscol,&csize));
139036032a97SHong Zhang   if (call == MAT_REUSE_MATRIX) {
13919566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local));
13925f80ce2aSJacob Faibussowitsch     PetscCheck(iscol_local,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Submatrix passed in was not used before, cannot reuse");
139336032a97SHong Zhang   } else {
1394068661f9SToby Isaac     PetscBool issorted;
1395068661f9SToby Isaac 
13969566063dSJacob Faibussowitsch     PetscCall(ISAllGather(iscol,&iscol_local));
13979566063dSJacob Faibussowitsch     PetscCall(ISEqual_private(isrow,iscol_local,&isequal));
13989566063dSJacob Faibussowitsch     PetscCall(ISSorted(iscol_local, &issorted));
13995f80ce2aSJacob Faibussowitsch     PetscCheck(isequal && issorted,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow and be sorted");
14008f46ffcaSHong Zhang   }
14018f46ffcaSHong Zhang 
14027dae84e0SHong Zhang   /* now call MatCreateSubMatrix_MPIBAIJ() */
14039566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat));
14048f46ffcaSHong Zhang   if (call == MAT_INITIAL_MATRIX) {
14059566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local));
14069566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&iscol_local));
14078f46ffcaSHong Zhang   }
14088f46ffcaSHong Zhang   PetscFunctionReturn(0);
14098f46ffcaSHong Zhang }
14108f46ffcaSHong Zhang 
1411dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1412a30f8f8cSSatish Balay {
1413a30f8f8cSSatish Balay   Mat_MPISBAIJ   *l = (Mat_MPISBAIJ*)A->data;
1414a30f8f8cSSatish Balay 
1415a30f8f8cSSatish Balay   PetscFunctionBegin;
14169566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->A));
14179566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->B));
1418a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1419a30f8f8cSSatish Balay }
1420a30f8f8cSSatish Balay 
1421dfbe8321SBarry Smith PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1422a30f8f8cSSatish Balay {
1423a30f8f8cSSatish Balay   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)matin->data;
1424a30f8f8cSSatish Balay   Mat            A  = a->A,B = a->B;
14253966268fSBarry Smith   PetscLogDouble isend[5],irecv[5];
1426a30f8f8cSSatish Balay 
1427a30f8f8cSSatish Balay   PetscFunctionBegin;
1428d0f46423SBarry Smith   info->block_size = (PetscReal)matin->rmap->bs;
142926fbe8dcSKarl Rupp 
14309566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A,MAT_LOCAL,info));
143126fbe8dcSKarl Rupp 
1432a30f8f8cSSatish Balay   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1433a30f8f8cSSatish Balay   isend[3] = info->memory;  isend[4] = info->mallocs;
143426fbe8dcSKarl Rupp 
14359566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(B,MAT_LOCAL,info));
143626fbe8dcSKarl Rupp 
1437a30f8f8cSSatish Balay   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1438a30f8f8cSSatish Balay   isend[3] += info->memory;  isend[4] += info->mallocs;
1439a30f8f8cSSatish Balay   if (flag == MAT_LOCAL) {
1440a30f8f8cSSatish Balay     info->nz_used      = isend[0];
1441a30f8f8cSSatish Balay     info->nz_allocated = isend[1];
1442a30f8f8cSSatish Balay     info->nz_unneeded  = isend[2];
1443a30f8f8cSSatish Balay     info->memory       = isend[3];
1444a30f8f8cSSatish Balay     info->mallocs      = isend[4];
1445a30f8f8cSSatish Balay   } else if (flag == MAT_GLOBAL_MAX) {
14461c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)matin)));
144726fbe8dcSKarl Rupp 
1448a30f8f8cSSatish Balay     info->nz_used      = irecv[0];
1449a30f8f8cSSatish Balay     info->nz_allocated = irecv[1];
1450a30f8f8cSSatish Balay     info->nz_unneeded  = irecv[2];
1451a30f8f8cSSatish Balay     info->memory       = irecv[3];
1452a30f8f8cSSatish Balay     info->mallocs      = irecv[4];
1453a30f8f8cSSatish Balay   } else if (flag == MAT_GLOBAL_SUM) {
14541c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)matin)));
145526fbe8dcSKarl Rupp 
1456a30f8f8cSSatish Balay     info->nz_used      = irecv[0];
1457a30f8f8cSSatish Balay     info->nz_allocated = irecv[1];
1458a30f8f8cSSatish Balay     info->nz_unneeded  = irecv[2];
1459a30f8f8cSSatish Balay     info->memory       = irecv[3];
1460a30f8f8cSSatish Balay     info->mallocs      = irecv[4];
146198921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1462a30f8f8cSSatish Balay   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
1463a30f8f8cSSatish Balay   info->fill_ratio_needed = 0;
1464a30f8f8cSSatish Balay   info->factor_mallocs    = 0;
1465a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1466a30f8f8cSSatish Balay }
1467a30f8f8cSSatish Balay 
1468ace3abfcSBarry Smith PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1469a30f8f8cSSatish Balay {
1470a30f8f8cSSatish Balay   Mat_MPISBAIJ   *a  = (Mat_MPISBAIJ*)A->data;
1471d0d4cfc2SHong Zhang   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;
1472a30f8f8cSSatish Balay 
1473a30f8f8cSSatish Balay   PetscFunctionBegin;
1474e98b92d7SKris Buschelman   switch (op) {
1475512a5fc5SBarry Smith   case MAT_NEW_NONZERO_LOCATIONS:
1476e98b92d7SKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
147728b2fa4aSMatthew Knepley   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1478a9817697SBarry Smith   case MAT_KEEP_NONZERO_PATTERN:
1479c10200c1SHong Zhang   case MAT_SUBMAT_SINGLEIS:
1480e98b92d7SKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
148143674050SBarry Smith     MatCheckPreallocated(A,1);
14829566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
14839566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B,op,flg));
1484e98b92d7SKris Buschelman     break;
1485e98b92d7SKris Buschelman   case MAT_ROW_ORIENTED:
148643674050SBarry Smith     MatCheckPreallocated(A,1);
14874e0d8c25SBarry Smith     a->roworiented = flg;
148826fbe8dcSKarl Rupp 
14899566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
14909566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B,op,flg));
1491e98b92d7SKris Buschelman     break;
14928c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
1493071fcb05SBarry Smith   case MAT_SORTED_FULL:
14949566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op]));
1495e98b92d7SKris Buschelman     break;
1496e98b92d7SKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
14974e0d8c25SBarry Smith     a->donotstash = flg;
1498e98b92d7SKris Buschelman     break;
1499e98b92d7SKris Buschelman   case MAT_USE_HASH_TABLE:
15004e0d8c25SBarry Smith     a->ht_flag = flg;
1501e98b92d7SKris Buschelman     break;
15029a4540c5SBarry Smith   case MAT_HERMITIAN:
150343674050SBarry Smith     MatCheckPreallocated(A,1);
15049566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
15050f2140c7SStefano Zampini #if defined(PETSC_USE_COMPLEX)
1506eb1ec7c1SStefano Zampini     if (flg) { /* need different mat-vec ops */
1507547795f9SHong Zhang       A->ops->mult             = MatMult_MPISBAIJ_Hermitian;
1508eb1ec7c1SStefano Zampini       A->ops->multadd          = MatMultAdd_MPISBAIJ_Hermitian;
1509eb1ec7c1SStefano Zampini       A->ops->multtranspose    = NULL;
1510eb1ec7c1SStefano Zampini       A->ops->multtransposeadd = NULL;
1511eb1ec7c1SStefano Zampini       A->symmetric = PETSC_FALSE;
1512eb1ec7c1SStefano Zampini     }
15130f2140c7SStefano Zampini #endif
1514eeffb40dSHong Zhang     break;
1515ffa07934SHong Zhang   case MAT_SPD:
151677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
151743674050SBarry Smith     MatCheckPreallocated(A,1);
15189566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
1519eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
1520eb1ec7c1SStefano Zampini     if (flg) { /* restore to use default mat-vec ops */
1521eb1ec7c1SStefano Zampini       A->ops->mult             = MatMult_MPISBAIJ;
1522eb1ec7c1SStefano Zampini       A->ops->multadd          = MatMultAdd_MPISBAIJ;
1523eb1ec7c1SStefano Zampini       A->ops->multtranspose    = MatMult_MPISBAIJ;
1524eb1ec7c1SStefano Zampini       A->ops->multtransposeadd = MatMultAdd_MPISBAIJ;
1525eb1ec7c1SStefano Zampini     }
1526eb1ec7c1SStefano Zampini #endif
1527eeffb40dSHong Zhang     break;
152877e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
152943674050SBarry Smith     MatCheckPreallocated(A,1);
15309566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A,op,flg));
1531eeffb40dSHong Zhang     break;
15329a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
15335f80ce2aSJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
15349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op]));
153577e54ba9SKris Buschelman     break;
1536d0d4cfc2SHong Zhang   case MAT_IGNORE_LOWER_TRIANGULAR:
15374e0d8c25SBarry Smith     aA->ignore_ltriangular = flg;
1538d0d4cfc2SHong Zhang     break;
1539d0d4cfc2SHong Zhang   case MAT_ERROR_LOWER_TRIANGULAR:
15404e0d8c25SBarry Smith     aA->ignore_ltriangular = flg;
1541d0d4cfc2SHong Zhang     break;
1542d0d4cfc2SHong Zhang   case MAT_GETROW_UPPERTRIANGULAR:
15434e0d8c25SBarry Smith     aA->getrow_utriangular = flg;
1544d0d4cfc2SHong Zhang     break;
1545e98b92d7SKris Buschelman   default:
154698921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1547a30f8f8cSSatish Balay   }
1548a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1549a30f8f8cSSatish Balay }
1550a30f8f8cSSatish Balay 
1551fc4dec0aSBarry Smith PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1552a30f8f8cSSatish Balay {
1553a30f8f8cSSatish Balay   PetscFunctionBegin;
1554cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX) {
15559566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(A,MAT_COPY_VALUES,B));
1556cf37664fSBarry Smith   }  else if (reuse == MAT_REUSE_MATRIX) {
15579566063dSJacob Faibussowitsch     PetscCall(MatCopy(A,*B,SAME_NONZERO_PATTERN));
1558fc4dec0aSBarry Smith   }
15598115998fSBarry Smith   PetscFunctionReturn(0);
1560a30f8f8cSSatish Balay }
1561a30f8f8cSSatish Balay 
1562dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1563a30f8f8cSSatish Balay {
1564a30f8f8cSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
1565a30f8f8cSSatish Balay   Mat            a     = baij->A, b=baij->B;
15665e90f9d9SHong Zhang   PetscInt       nv,m,n;
1567ace3abfcSBarry Smith   PetscBool      flg;
1568a30f8f8cSSatish Balay 
1569a30f8f8cSSatish Balay   PetscFunctionBegin;
1570a30f8f8cSSatish Balay   if (ll != rr) {
15719566063dSJacob Faibussowitsch     PetscCall(VecEqual(ll,rr,&flg));
15725f80ce2aSJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same");
1573a30f8f8cSSatish Balay   }
1574b3bf805bSHong Zhang   if (!ll) PetscFunctionReturn(0);
1575b3bf805bSHong Zhang 
15769566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat,&m,&n));
15775f80ce2aSJacob Faibussowitsch   PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %" PetscInt_FMT " %" PetscInt_FMT " must be same",m,n);
1578b3bf805bSHong Zhang 
15799566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(rr,&nv));
15805f80ce2aSJacob Faibussowitsch   PetscCheck(nv==n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
15815e90f9d9SHong Zhang 
15829566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD));
15835e90f9d9SHong Zhang 
15845e90f9d9SHong Zhang   /* left diagonalscale the off-diagonal part */
15859566063dSJacob Faibussowitsch   PetscCall((*b->ops->diagonalscale)(b,ll,NULL));
15865e90f9d9SHong Zhang 
15875e90f9d9SHong Zhang   /* scale the diagonal part */
15889566063dSJacob Faibussowitsch   PetscCall((*a->ops->diagonalscale)(a,ll,rr));
1589a30f8f8cSSatish Balay 
15905e90f9d9SHong Zhang   /* right diagonalscale the off-diagonal part */
15919566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD));
15929566063dSJacob Faibussowitsch   PetscCall((*b->ops->diagonalscale)(b,NULL,baij->lvec));
1593a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1594a30f8f8cSSatish Balay }
1595a30f8f8cSSatish Balay 
1596dfbe8321SBarry Smith PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1597a30f8f8cSSatish Balay {
1598f3566a2aSHong Zhang   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1599a30f8f8cSSatish Balay 
1600a30f8f8cSSatish Balay   PetscFunctionBegin;
16019566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(a->A));
1602a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1603a30f8f8cSSatish Balay }
1604a30f8f8cSSatish Balay 
16056849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);
1606a30f8f8cSSatish Balay 
1607ace3abfcSBarry Smith PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool  *flag)
1608a30f8f8cSSatish Balay {
1609a30f8f8cSSatish Balay   Mat_MPISBAIJ   *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1610a30f8f8cSSatish Balay   Mat            a,b,c,d;
1611ace3abfcSBarry Smith   PetscBool      flg;
1612a30f8f8cSSatish Balay 
1613a30f8f8cSSatish Balay   PetscFunctionBegin;
1614a30f8f8cSSatish Balay   a = matA->A; b = matA->B;
1615a30f8f8cSSatish Balay   c = matB->A; d = matB->B;
1616a30f8f8cSSatish Balay 
16179566063dSJacob Faibussowitsch   PetscCall(MatEqual(a,c,&flg));
1618abc0a331SBarry Smith   if (flg) {
16199566063dSJacob Faibussowitsch     PetscCall(MatEqual(b,d,&flg));
1620a30f8f8cSSatish Balay   }
16211c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A)));
1622a30f8f8cSSatish Balay   PetscFunctionReturn(0);
1623a30f8f8cSSatish Balay }
1624a30f8f8cSSatish Balay 
16253c896bc6SHong Zhang PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
16263c896bc6SHong Zhang {
16274c7a3774SStefano Zampini   PetscBool      isbaij;
16283c896bc6SHong Zhang 
16293c896bc6SHong Zhang   PetscFunctionBegin;
16309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)B,&isbaij,MATSEQSBAIJ,MATMPISBAIJ,""));
16315f80ce2aSJacob Faibussowitsch   PetscCheck(isbaij,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)B)->type_name);
16323c896bc6SHong Zhang   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
16333c896bc6SHong Zhang   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
16349566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(A));
16359566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A,B,str));
16369566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(A));
16373c896bc6SHong Zhang   } else {
16384c7a3774SStefano Zampini     Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
16394c7a3774SStefano Zampini     Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;
16404c7a3774SStefano Zampini 
16419566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->A,b->A,str));
16429566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->B,b->B,str));
16433c896bc6SHong Zhang   }
16449566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)B));
16453c896bc6SHong Zhang   PetscFunctionReturn(0);
16463c896bc6SHong Zhang }
16473c896bc6SHong Zhang 
16484994cf47SJed Brown PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1649273d9f13SBarry Smith {
1650273d9f13SBarry Smith   PetscFunctionBegin;
16519566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
1652273d9f13SBarry Smith   PetscFunctionReturn(0);
1653273d9f13SBarry Smith }
1654a5e6ed63SBarry Smith 
16554fe895cdSHong Zhang PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
16564fe895cdSHong Zhang {
16574fe895cdSHong Zhang   Mat_MPISBAIJ   *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
16584fe895cdSHong Zhang   PetscBLASInt   bnz,one=1;
16594fe895cdSHong Zhang   Mat_SeqSBAIJ   *xa,*ya;
16604fe895cdSHong Zhang   Mat_SeqBAIJ    *xb,*yb;
16614fe895cdSHong Zhang 
16624fe895cdSHong Zhang   PetscFunctionBegin;
16634fe895cdSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
16644fe895cdSHong Zhang     PetscScalar alpha = a;
16654fe895cdSHong Zhang     xa   = (Mat_SeqSBAIJ*)xx->A->data;
16664fe895cdSHong Zhang     ya   = (Mat_SeqSBAIJ*)yy->A->data;
16679566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(xa->nz,&bnz));
16688b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
16694fe895cdSHong Zhang     xb   = (Mat_SeqBAIJ*)xx->B->data;
16704fe895cdSHong Zhang     yb   = (Mat_SeqBAIJ*)yy->B->data;
16719566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(xb->nz,&bnz));
16728b83055fSJed Brown     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
16739566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateIncrease((PetscObject)Y));
1674ab784542SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
16759566063dSJacob Faibussowitsch     PetscCall(MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE));
16769566063dSJacob Faibussowitsch     PetscCall(MatAXPY_Basic(Y,a,X,str));
16779566063dSJacob Faibussowitsch     PetscCall(MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE));
16784fe895cdSHong Zhang   } else {
16794de5dceeSHong Zhang     Mat      B;
16804de5dceeSHong Zhang     PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
16815f80ce2aSJacob Faibussowitsch     PetscCheck(bs == X->rmap->bs,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
16829566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(X));
16839566063dSJacob Faibussowitsch     PetscCall(MatGetRowUpperTriangular(Y));
16849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(yy->A->rmap->N,&nnz_d));
16859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(yy->B->rmap->N,&nnz_o));
16869566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)Y),&B));
16879566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name));
16889566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N));
16899566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(B,Y,Y));
16909566063dSJacob Faibussowitsch     PetscCall(MatSetType(B,MATMPISBAIJ));
16919566063dSJacob Faibussowitsch     PetscCall(MatAXPYGetPreallocation_SeqSBAIJ(yy->A,xx->A,nnz_d));
16929566063dSJacob Faibussowitsch     PetscCall(MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o));
16939566063dSJacob Faibussowitsch     PetscCall(MatMPISBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o));
16949566063dSJacob Faibussowitsch     PetscCall(MatAXPY_BasicWithPreallocation(B,Y,a,X,str));
16959566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(Y,&B));
16969566063dSJacob Faibussowitsch     PetscCall(PetscFree(nnz_d));
16979566063dSJacob Faibussowitsch     PetscCall(PetscFree(nnz_o));
16989566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(X));
16999566063dSJacob Faibussowitsch     PetscCall(MatRestoreRowUpperTriangular(Y));
17004fe895cdSHong Zhang   }
17014fe895cdSHong Zhang   PetscFunctionReturn(0);
17024fe895cdSHong Zhang }
17034fe895cdSHong Zhang 
17047dae84e0SHong Zhang PetscErrorCode MatCreateSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1705a5e6ed63SBarry Smith {
17061302d50aSBarry Smith   PetscInt       i;
1707afebec48SHong Zhang   PetscBool      flg;
1708a5e6ed63SBarry Smith 
17096849ba73SBarry Smith   PetscFunctionBegin;
17109566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B)); /* B[] are sbaij matrices */
1711a5e6ed63SBarry Smith   for (i=0; i<n; i++) {
17129566063dSJacob Faibussowitsch     PetscCall(ISEqual(irow[i],icol[i],&flg));
1713afebec48SHong Zhang     if (!flg) {
17149566063dSJacob Faibussowitsch       PetscCall(MatSeqSBAIJZeroOps_Private(*B[i]));
1715a5e6ed63SBarry Smith     }
17164dcd73b1SHong Zhang   }
1717a5e6ed63SBarry Smith   PetscFunctionReturn(0);
1718a5e6ed63SBarry Smith }
1719a5e6ed63SBarry Smith 
17207d68702bSBarry Smith PetscErrorCode MatShift_MPISBAIJ(Mat Y,PetscScalar a)
17217d68702bSBarry Smith {
17227d68702bSBarry Smith   Mat_MPISBAIJ    *maij = (Mat_MPISBAIJ*)Y->data;
17236f33a894SBarry Smith   Mat_SeqSBAIJ    *aij = (Mat_SeqSBAIJ*)maij->A->data;
17247d68702bSBarry Smith 
17257d68702bSBarry Smith   PetscFunctionBegin;
17266f33a894SBarry Smith   if (!Y->preallocated) {
17279566063dSJacob Faibussowitsch     PetscCall(MatMPISBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL));
17286f33a894SBarry Smith   } else if (!aij->nz) {
1729b83222d8SBarry Smith     PetscInt nonew = aij->nonew;
17309566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL));
1731b83222d8SBarry Smith     aij->nonew = nonew;
17327d68702bSBarry Smith   }
17339566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y,a));
17347d68702bSBarry Smith   PetscFunctionReturn(0);
17357d68702bSBarry Smith }
17367d68702bSBarry Smith 
17373b49f96aSBarry Smith PetscErrorCode MatMissingDiagonal_MPISBAIJ(Mat A,PetscBool  *missing,PetscInt *d)
17383b49f96aSBarry Smith {
17393b49f96aSBarry Smith   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
17403b49f96aSBarry Smith 
17413b49f96aSBarry Smith   PetscFunctionBegin;
17425f80ce2aSJacob Faibussowitsch   PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
17439566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(a->A,missing,d));
17443b49f96aSBarry Smith   if (d) {
17453b49f96aSBarry Smith     PetscInt rstart;
17469566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A,&rstart,NULL));
17473b49f96aSBarry Smith     *d += rstart/A->rmap->bs;
17483b49f96aSBarry Smith 
17493b49f96aSBarry Smith   }
17503b49f96aSBarry Smith   PetscFunctionReturn(0);
17513b49f96aSBarry Smith }
17523b49f96aSBarry Smith 
1753a5b7ff6bSBarry Smith PetscErrorCode  MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1754a5b7ff6bSBarry Smith {
1755a5b7ff6bSBarry Smith   PetscFunctionBegin;
1756a5b7ff6bSBarry Smith   *a = ((Mat_MPISBAIJ*)A->data)->A;
1757a5b7ff6bSBarry Smith   PetscFunctionReturn(0);
1758a5b7ff6bSBarry Smith }
17593b49f96aSBarry Smith 
1760a30f8f8cSSatish Balay /* -------------------------------------------------------------------*/
17613964eb88SJed Brown static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1762a30f8f8cSSatish Balay                                        MatGetRow_MPISBAIJ,
1763a30f8f8cSSatish Balay                                        MatRestoreRow_MPISBAIJ,
1764a9d4b620SHong Zhang                                        MatMult_MPISBAIJ,
176597304618SKris Buschelman                                /*  4*/ MatMultAdd_MPISBAIJ,
1766431c96f7SBarry Smith                                        MatMult_MPISBAIJ,       /* transpose versions are same as non-transpose */
1767431c96f7SBarry Smith                                        MatMultAdd_MPISBAIJ,
1768f4259b30SLisandro Dalcin                                        NULL,
1769f4259b30SLisandro Dalcin                                        NULL,
1770f4259b30SLisandro Dalcin                                        NULL,
1771f4259b30SLisandro Dalcin                                /* 10*/ NULL,
1772f4259b30SLisandro Dalcin                                        NULL,
1773f4259b30SLisandro Dalcin                                        NULL,
177441f059aeSBarry Smith                                        MatSOR_MPISBAIJ,
1775a30f8f8cSSatish Balay                                        MatTranspose_MPISBAIJ,
177697304618SKris Buschelman                                /* 15*/ MatGetInfo_MPISBAIJ,
1777a30f8f8cSSatish Balay                                        MatEqual_MPISBAIJ,
1778a30f8f8cSSatish Balay                                        MatGetDiagonal_MPISBAIJ,
1779a30f8f8cSSatish Balay                                        MatDiagonalScale_MPISBAIJ,
1780a30f8f8cSSatish Balay                                        MatNorm_MPISBAIJ,
178197304618SKris Buschelman                                /* 20*/ MatAssemblyBegin_MPISBAIJ,
1782a30f8f8cSSatish Balay                                        MatAssemblyEnd_MPISBAIJ,
1783a30f8f8cSSatish Balay                                        MatSetOption_MPISBAIJ,
1784a30f8f8cSSatish Balay                                        MatZeroEntries_MPISBAIJ,
1785f4259b30SLisandro Dalcin                                /* 24*/ NULL,
1786f4259b30SLisandro Dalcin                                        NULL,
1787f4259b30SLisandro Dalcin                                        NULL,
1788f4259b30SLisandro Dalcin                                        NULL,
1789f4259b30SLisandro Dalcin                                        NULL,
17904994cf47SJed Brown                                /* 29*/ MatSetUp_MPISBAIJ,
1791f4259b30SLisandro Dalcin                                        NULL,
1792f4259b30SLisandro Dalcin                                        NULL,
1793a5b7ff6bSBarry Smith                                        MatGetDiagonalBlock_MPISBAIJ,
1794f4259b30SLisandro Dalcin                                        NULL,
1795d519adbfSMatthew Knepley                                /* 34*/ MatDuplicate_MPISBAIJ,
1796f4259b30SLisandro Dalcin                                        NULL,
1797f4259b30SLisandro Dalcin                                        NULL,
1798f4259b30SLisandro Dalcin                                        NULL,
1799f4259b30SLisandro Dalcin                                        NULL,
1800d519adbfSMatthew Knepley                                /* 39*/ MatAXPY_MPISBAIJ,
18017dae84e0SHong Zhang                                        MatCreateSubMatrices_MPISBAIJ,
1802d94109b8SHong Zhang                                        MatIncreaseOverlap_MPISBAIJ,
1803a30f8f8cSSatish Balay                                        MatGetValues_MPISBAIJ,
18043c896bc6SHong Zhang                                        MatCopy_MPISBAIJ,
1805f4259b30SLisandro Dalcin                                /* 44*/ NULL,
1806a30f8f8cSSatish Balay                                        MatScale_MPISBAIJ,
18077d68702bSBarry Smith                                        MatShift_MPISBAIJ,
1808f4259b30SLisandro Dalcin                                        NULL,
1809f4259b30SLisandro Dalcin                                        NULL,
1810f4259b30SLisandro Dalcin                                /* 49*/ NULL,
1811f4259b30SLisandro Dalcin                                        NULL,
1812f4259b30SLisandro Dalcin                                        NULL,
1813f4259b30SLisandro Dalcin                                        NULL,
1814f4259b30SLisandro Dalcin                                        NULL,
1815f4259b30SLisandro Dalcin                                /* 54*/ NULL,
1816f4259b30SLisandro Dalcin                                        NULL,
1817a30f8f8cSSatish Balay                                        MatSetUnfactored_MPISBAIJ,
1818f4259b30SLisandro Dalcin                                        NULL,
1819a30f8f8cSSatish Balay                                        MatSetValuesBlocked_MPISBAIJ,
18207dae84e0SHong Zhang                                /* 59*/ MatCreateSubMatrix_MPISBAIJ,
1821f4259b30SLisandro Dalcin                                        NULL,
1822f4259b30SLisandro Dalcin                                        NULL,
1823f4259b30SLisandro Dalcin                                        NULL,
1824f4259b30SLisandro Dalcin                                        NULL,
1825f4259b30SLisandro Dalcin                                /* 64*/ NULL,
1826f4259b30SLisandro Dalcin                                        NULL,
1827f4259b30SLisandro Dalcin                                        NULL,
1828f4259b30SLisandro Dalcin                                        NULL,
1829f4259b30SLisandro Dalcin                                        NULL,
1830d519adbfSMatthew Knepley                                /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1831f4259b30SLisandro Dalcin                                        NULL,
183228d58a37SPierre Jolivet                                        MatConvert_MPISBAIJ_Basic,
1833f4259b30SLisandro Dalcin                                        NULL,
1834f4259b30SLisandro Dalcin                                        NULL,
1835f4259b30SLisandro Dalcin                                /* 74*/ NULL,
1836f4259b30SLisandro Dalcin                                        NULL,
1837f4259b30SLisandro Dalcin                                        NULL,
1838f4259b30SLisandro Dalcin                                        NULL,
1839f4259b30SLisandro Dalcin                                        NULL,
1840f4259b30SLisandro Dalcin                                /* 79*/ NULL,
1841f4259b30SLisandro Dalcin                                        NULL,
1842f4259b30SLisandro Dalcin                                        NULL,
1843f4259b30SLisandro Dalcin                                        NULL,
18445bba2384SShri Abhyankar                                        MatLoad_MPISBAIJ,
1845f4259b30SLisandro Dalcin                                /* 84*/ NULL,
1846f4259b30SLisandro Dalcin                                        NULL,
1847f4259b30SLisandro Dalcin                                        NULL,
1848f4259b30SLisandro Dalcin                                        NULL,
1849f4259b30SLisandro Dalcin                                        NULL,
1850f4259b30SLisandro Dalcin                                /* 89*/ NULL,
1851f4259b30SLisandro Dalcin                                        NULL,
1852f4259b30SLisandro Dalcin                                        NULL,
1853f4259b30SLisandro Dalcin                                        NULL,
1854f4259b30SLisandro Dalcin                                        NULL,
1855f4259b30SLisandro Dalcin                                /* 94*/ NULL,
1856f4259b30SLisandro Dalcin                                        NULL,
1857f4259b30SLisandro Dalcin                                        NULL,
1858f4259b30SLisandro Dalcin                                        NULL,
1859f4259b30SLisandro Dalcin                                        NULL,
1860f4259b30SLisandro Dalcin                                /* 99*/ NULL,
1861f4259b30SLisandro Dalcin                                        NULL,
1862f4259b30SLisandro Dalcin                                        NULL,
18632726fb6dSPierre Jolivet                                        MatConjugate_MPISBAIJ,
1864f4259b30SLisandro Dalcin                                        NULL,
1865f4259b30SLisandro Dalcin                                /*104*/ NULL,
186699cafbc1SBarry Smith                                        MatRealPart_MPISBAIJ,
1867d0d4cfc2SHong Zhang                                        MatImaginaryPart_MPISBAIJ,
1868d0d4cfc2SHong Zhang                                        MatGetRowUpperTriangular_MPISBAIJ,
186995936485SShri Abhyankar                                        MatRestoreRowUpperTriangular_MPISBAIJ,
1870f4259b30SLisandro Dalcin                                /*109*/ NULL,
1871f4259b30SLisandro Dalcin                                        NULL,
1872f4259b30SLisandro Dalcin                                        NULL,
1873f4259b30SLisandro Dalcin                                        NULL,
18743b49f96aSBarry Smith                                        MatMissingDiagonal_MPISBAIJ,
1875f4259b30SLisandro Dalcin                                /*114*/ NULL,
1876f4259b30SLisandro Dalcin                                        NULL,
1877f4259b30SLisandro Dalcin                                        NULL,
1878f4259b30SLisandro Dalcin                                        NULL,
1879f4259b30SLisandro Dalcin                                        NULL,
1880f4259b30SLisandro Dalcin                                /*119*/ NULL,
1881f4259b30SLisandro Dalcin                                        NULL,
1882f4259b30SLisandro Dalcin                                        NULL,
1883f4259b30SLisandro Dalcin                                        NULL,
1884f4259b30SLisandro Dalcin                                        NULL,
1885f4259b30SLisandro Dalcin                                /*124*/ NULL,
1886f4259b30SLisandro Dalcin                                        NULL,
1887f4259b30SLisandro Dalcin                                        NULL,
1888f4259b30SLisandro Dalcin                                        NULL,
1889f4259b30SLisandro Dalcin                                        NULL,
1890f4259b30SLisandro Dalcin                                /*129*/ NULL,
1891f4259b30SLisandro Dalcin                                        NULL,
1892f4259b30SLisandro Dalcin                                        NULL,
1893f4259b30SLisandro Dalcin                                        NULL,
1894f4259b30SLisandro Dalcin                                        NULL,
1895f4259b30SLisandro Dalcin                                /*134*/ NULL,
1896f4259b30SLisandro Dalcin                                        NULL,
1897f4259b30SLisandro Dalcin                                        NULL,
1898f4259b30SLisandro Dalcin                                        NULL,
1899f4259b30SLisandro Dalcin                                        NULL,
190046533700Sstefano_zampini                                /*139*/ MatSetBlockSizes_Default,
1901f4259b30SLisandro Dalcin                                        NULL,
1902f4259b30SLisandro Dalcin                                        NULL,
1903f4259b30SLisandro Dalcin                                        NULL,
1904f4259b30SLisandro Dalcin                                        NULL,
1905d70f29a3SPierre Jolivet                                 /*144*/MatCreateMPIMatConcatenateSeqMat_MPISBAIJ,
1906d70f29a3SPierre Jolivet                                        NULL,
1907d70f29a3SPierre Jolivet                                        NULL,
1908*99a7f59eSMark Adams                                        NULL,
1909*99a7f59eSMark Adams                                        NULL,
1910d70f29a3SPierre Jolivet                                        NULL
191199cafbc1SBarry Smith };
1912a30f8f8cSSatish Balay 
1913b2573a8aSBarry Smith PetscErrorCode  MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
1914a23d5eceSKris Buschelman {
1915476417e5SBarry Smith   Mat_MPISBAIJ   *b = (Mat_MPISBAIJ*)B->data;
1916535b19f3SBarry Smith   PetscInt       i,mbs,Mbs;
19175d2a9ed1SStefano Zampini   PetscMPIInt    size;
1918a23d5eceSKris Buschelman 
1919a23d5eceSKris Buschelman   PetscFunctionBegin;
19209566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(B,PetscAbs(bs)));
19219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
19229566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
19239566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(B->rmap,&bs));
19245f80ce2aSJacob Faibussowitsch   PetscCheck(B->rmap->N <= B->cmap->N,PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MPISBAIJ matrix cannot have more rows %" PetscInt_FMT " than columns %" PetscInt_FMT,B->rmap->N,B->cmap->N);
19255f80ce2aSJacob Faibussowitsch   PetscCheck(B->rmap->n <= B->cmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"MPISBAIJ matrix cannot have more local rows %" PetscInt_FMT " than columns %" PetscInt_FMT,B->rmap->n,B->cmap->n);
1926899cda47SBarry Smith 
1927d0f46423SBarry Smith   mbs = B->rmap->n/bs;
1928d0f46423SBarry Smith   Mbs = B->rmap->N/bs;
19295f80ce2aSJacob Faibussowitsch   PetscCheck(mbs*bs == B->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"No of local rows %" PetscInt_FMT " must be divisible by blocksize %" PetscInt_FMT,B->rmap->N,bs);
1930a23d5eceSKris Buschelman 
1931d0f46423SBarry Smith   B->rmap->bs = bs;
1932a23d5eceSKris Buschelman   b->bs2      = bs*bs;
1933a23d5eceSKris Buschelman   b->mbs      = mbs;
1934a23d5eceSKris Buschelman   b->Mbs      = Mbs;
1935de64b629SHong Zhang   b->nbs      = B->cmap->n/bs;
1936de64b629SHong Zhang   b->Nbs      = B->cmap->N/bs;
1937a23d5eceSKris Buschelman 
1938a23d5eceSKris Buschelman   for (i=0; i<=b->size; i++) {
1939d0f46423SBarry Smith     b->rangebs[i] = B->rmap->range[i]/bs;
1940a23d5eceSKris Buschelman   }
1941d0f46423SBarry Smith   b->rstartbs = B->rmap->rstart/bs;
1942d0f46423SBarry Smith   b->rendbs   = B->rmap->rend/bs;
1943a23d5eceSKris Buschelman 
1944d0f46423SBarry Smith   b->cstartbs = B->cmap->rstart/bs;
1945d0f46423SBarry Smith   b->cendbs   = B->cmap->rend/bs;
1946a23d5eceSKris Buschelman 
1947cb7b82ddSBarry Smith #if defined(PETSC_USE_CTABLE)
19489566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&b->colmap));
1949cb7b82ddSBarry Smith #else
19509566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->colmap));
1951cb7b82ddSBarry Smith #endif
19529566063dSJacob Faibussowitsch   PetscCall(PetscFree(b->garray));
19539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->lvec));
19549566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&b->Mvctx));
19559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->slvec0));
19569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->slvec0b));
19579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->slvec1));
19589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->slvec1a));
19599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b->slvec1b));
19609566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&b->sMvctx));
1961cb7b82ddSBarry Smith 
1962cb7b82ddSBarry Smith   /* Because the B will have been resized we simply destroy it and create a new one each time */
19639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B),&size));
19649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&b->B));
19659566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&b->B));
19669566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(b->B,B->rmap->n,size > 1 ? B->cmap->N : 0,B->rmap->n,size > 1 ? B->cmap->N : 0));
19679566063dSJacob Faibussowitsch   PetscCall(MatSetType(b->B,MATSEQBAIJ));
19689566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)b->B));
1969cb7b82ddSBarry Smith 
1970526dfc15SBarry Smith   if (!B->preallocated) {
19719566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF,&b->A));
19729566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n));
19739566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->A,MATSEQSBAIJ));
19749566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)B,(PetscObject)b->A));
19759566063dSJacob Faibussowitsch     PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash));
1976526dfc15SBarry Smith   }
1977a23d5eceSKris Buschelman 
19789566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz));
19799566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz));
198026fbe8dcSKarl Rupp 
1981526dfc15SBarry Smith   B->preallocated  = PETSC_TRUE;
1982cb7b82ddSBarry Smith   B->was_assembled = PETSC_FALSE;
1983cb7b82ddSBarry Smith   B->assembled     = PETSC_FALSE;
1984a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1985a23d5eceSKris Buschelman }
1986a23d5eceSKris Buschelman 
1987dfb205c3SBarry Smith PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
1988dfb205c3SBarry Smith {
198902106b30SBarry Smith   PetscInt       m,rstart,cend;
1990f4259b30SLisandro Dalcin   PetscInt       i,j,d,nz,bd, nz_max=0,*d_nnz=NULL,*o_nnz=NULL;
1991f4259b30SLisandro Dalcin   const PetscInt *JJ    =NULL;
1992f4259b30SLisandro Dalcin   PetscScalar    *values=NULL;
1993bb80cfbbSStefano Zampini   PetscBool      roworiented = ((Mat_MPISBAIJ*)B->data)->roworiented;
19943bd0feecSPierre Jolivet   PetscBool      nooffprocentries;
1995dfb205c3SBarry Smith 
1996dfb205c3SBarry Smith   PetscFunctionBegin;
19975f80ce2aSJacob Faibussowitsch   PetscCheck(bs >= 1,PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %" PetscInt_FMT,bs);
19989566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(B->rmap,bs));
19999566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(B->cmap,bs));
20009566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
20019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
20029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetBlockSize(B->rmap,&bs));
2003dfb205c3SBarry Smith   m      = B->rmap->n/bs;
2004dfb205c3SBarry Smith   rstart = B->rmap->rstart/bs;
2005dfb205c3SBarry Smith   cend   = B->cmap->rend/bs;
2006dfb205c3SBarry Smith 
20075f80ce2aSJacob Faibussowitsch   PetscCheck(!ii[0],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %" PetscInt_FMT,ii[0]);
20089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m,&d_nnz,m,&o_nnz));
2009dfb205c3SBarry Smith   for (i=0; i<m; i++) {
2010dfb205c3SBarry Smith     nz = ii[i+1] - ii[i];
20115f80ce2aSJacob Faibussowitsch     PetscCheck(nz >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %" PetscInt_FMT " has a negative number of columns %" PetscInt_FMT,i,nz);
20120cd7f59aSBarry Smith     /* count the ones on the diagonal and above, split into diagonal and off diagonal portions. */
2013dfb205c3SBarry Smith     JJ     = jj + ii[i];
20140cd7f59aSBarry Smith     bd     = 0;
2015dfb205c3SBarry Smith     for (j=0; j<nz; j++) {
20160cd7f59aSBarry Smith       if (*JJ >= i + rstart) break;
2017dfb205c3SBarry Smith       JJ++;
20180cd7f59aSBarry Smith       bd++;
2019dfb205c3SBarry Smith     }
2020dfb205c3SBarry Smith     d  = 0;
2021dfb205c3SBarry Smith     for (; j<nz; j++) {
2022dfb205c3SBarry Smith       if (*JJ++ >= cend) break;
2023dfb205c3SBarry Smith       d++;
2024dfb205c3SBarry Smith     }
2025dfb205c3SBarry Smith     d_nnz[i] = d;
20260cd7f59aSBarry Smith     o_nnz[i] = nz - d - bd;
20270cd7f59aSBarry Smith     nz       = nz - bd;
20280cd7f59aSBarry Smith     nz_max = PetscMax(nz_max,nz);
2029dfb205c3SBarry Smith   }
20309566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz));
20319566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
20329566063dSJacob Faibussowitsch   PetscCall(PetscFree2(d_nnz,o_nnz));
2033dfb205c3SBarry Smith 
2034dfb205c3SBarry Smith   values = (PetscScalar*)V;
2035dfb205c3SBarry Smith   if (!values) {
20369566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(bs*bs*nz_max,&values));
2037dfb205c3SBarry Smith   }
2038dfb205c3SBarry Smith   for (i=0; i<m; i++) {
2039dfb205c3SBarry Smith     PetscInt          row    = i + rstart;
2040dfb205c3SBarry Smith     PetscInt          ncols  = ii[i+1] - ii[i];
2041dfb205c3SBarry Smith     const PetscInt    *icols = jj + ii[i];
2042bb80cfbbSStefano Zampini     if (bs == 1 || !roworiented) {         /* block ordering matches the non-nested layout of MatSetValues so we can insert entire rows */
2043dfb205c3SBarry Smith       const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
20449566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES));
2045bb80cfbbSStefano Zampini     } else {                    /* block ordering does not match so we can only insert one block at a time. */
2046bb80cfbbSStefano Zampini       PetscInt j;
20470cd7f59aSBarry Smith       for (j=0; j<ncols; j++) {
20480cd7f59aSBarry Smith         const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
20499566063dSJacob Faibussowitsch         PetscCall(MatSetValuesBlocked_MPISBAIJ(B,1,&row,1,&icols[j],svals,INSERT_VALUES));
20500cd7f59aSBarry Smith       }
20510cd7f59aSBarry Smith     }
2052dfb205c3SBarry Smith   }
2053dfb205c3SBarry Smith 
20549566063dSJacob Faibussowitsch   if (!V) PetscCall(PetscFree(values));
20553bd0feecSPierre Jolivet   nooffprocentries    = B->nooffprocentries;
20563bd0feecSPierre Jolivet   B->nooffprocentries = PETSC_TRUE;
20579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
20589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
20593bd0feecSPierre Jolivet   B->nooffprocentries = nooffprocentries;
20603bd0feecSPierre Jolivet 
20619566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
2062dfb205c3SBarry Smith   PetscFunctionReturn(0);
2063dfb205c3SBarry Smith }
2064dfb205c3SBarry Smith 
20650bad9183SKris Buschelman /*MC
2066fafad747SKris Buschelman    MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2067828413b8SBarry Smith    based on block compressed sparse row format.  Only the upper triangular portion of the "diagonal" portion of
2068828413b8SBarry Smith    the matrix is stored.
2069828413b8SBarry Smith 
2070828413b8SBarry Smith    For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
2071828413b8SBarry Smith    can call MatSetOption(Mat, MAT_HERMITIAN);
20720bad9183SKris Buschelman 
20730bad9183SKris Buschelman    Options Database Keys:
20740bad9183SKris Buschelman . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
20750bad9183SKris Buschelman 
2076476417e5SBarry Smith    Notes:
2077476417e5SBarry Smith      The number of rows in the matrix must be less than or equal to the number of columns. Similarly the number of rows in the
2078476417e5SBarry Smith      diagonal portion of the matrix of each process has to less than or equal the number of columns.
2079476417e5SBarry Smith 
20800bad9183SKris Buschelman    Level: beginner
20810bad9183SKris Buschelman 
2082db781477SPatrick Sanan .seealso: `MatCreateBAIJ()`, `MATSEQSBAIJ`, `MatType`
20830bad9183SKris Buschelman M*/
20840bad9183SKris Buschelman 
20858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2086b5df2d14SHong Zhang {
2087b5df2d14SHong Zhang   Mat_MPISBAIJ   *b;
208894ae4db5SBarry Smith   PetscBool      flg = PETSC_FALSE;
2089b5df2d14SHong Zhang 
2090b5df2d14SHong Zhang   PetscFunctionBegin;
20919566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&b));
2092b0a32e0cSBarry Smith   B->data = (void*)b;
20939566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)));
2094b5df2d14SHong Zhang 
2095b5df2d14SHong Zhang   B->ops->destroy = MatDestroy_MPISBAIJ;
2096b5df2d14SHong Zhang   B->ops->view    = MatView_MPISBAIJ;
2097b5df2d14SHong Zhang   B->assembled    = PETSC_FALSE;
2098b5df2d14SHong Zhang   B->insertmode   = NOT_SET_VALUES;
209926fbe8dcSKarl Rupp 
21009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank));
21019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size));
2102b5df2d14SHong Zhang 
2103b5df2d14SHong Zhang   /* build local table of row and column ownerships */
21049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(b->size+2,&b->rangebs));
2105b5df2d14SHong Zhang 
2106b5df2d14SHong Zhang   /* build cache for off array entries formed */
21079566063dSJacob Faibussowitsch   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash));
210826fbe8dcSKarl Rupp 
2109b5df2d14SHong Zhang   b->donotstash  = PETSC_FALSE;
21100298fd71SBarry Smith   b->colmap      = NULL;
21110298fd71SBarry Smith   b->garray      = NULL;
2112b5df2d14SHong Zhang   b->roworiented = PETSC_TRUE;
2113b5df2d14SHong Zhang 
2114b5df2d14SHong Zhang   /* stuff used in block assembly */
2115f4259b30SLisandro Dalcin   b->barray = NULL;
2116b5df2d14SHong Zhang 
2117b5df2d14SHong Zhang   /* stuff used for matrix vector multiply */
2118f4259b30SLisandro Dalcin   b->lvec    = NULL;
2119f4259b30SLisandro Dalcin   b->Mvctx   = NULL;
2120f4259b30SLisandro Dalcin   b->slvec0  = NULL;
2121f4259b30SLisandro Dalcin   b->slvec0b = NULL;
2122f4259b30SLisandro Dalcin   b->slvec1  = NULL;
2123f4259b30SLisandro Dalcin   b->slvec1a = NULL;
2124f4259b30SLisandro Dalcin   b->slvec1b = NULL;
2125f4259b30SLisandro Dalcin   b->sMvctx  = NULL;
2126b5df2d14SHong Zhang 
2127b5df2d14SHong Zhang   /* stuff for MatGetRow() */
2128f4259b30SLisandro Dalcin   b->rowindices   = NULL;
2129f4259b30SLisandro Dalcin   b->rowvalues    = NULL;
2130b5df2d14SHong Zhang   b->getrowactive = PETSC_FALSE;
2131b5df2d14SHong Zhang 
2132b5df2d14SHong Zhang   /* hash table stuff */
2133f4259b30SLisandro Dalcin   b->ht           = NULL;
2134f4259b30SLisandro Dalcin   b->hd           = NULL;
2135b5df2d14SHong Zhang   b->ht_size      = 0;
2136b5df2d14SHong Zhang   b->ht_flag      = PETSC_FALSE;
2137b5df2d14SHong Zhang   b->ht_fact      = 0;
2138b5df2d14SHong Zhang   b->ht_total_ct  = 0;
2139b5df2d14SHong Zhang   b->ht_insert_ct = 0;
2140b5df2d14SHong Zhang 
21417dae84e0SHong Zhang   /* stuff for MatCreateSubMatrices_MPIBAIJ_local() */
21427a868f3eSHong Zhang   b->ijonly = PETSC_FALSE;
21437a868f3eSHong Zhang 
2144f4259b30SLisandro Dalcin   b->in_loc = NULL;
2145f4259b30SLisandro Dalcin   b->v_loc  = NULL;
214659ffdab8SBarry Smith   b->n_loc  = 0;
214794ae4db5SBarry Smith 
21489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ));
21499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ));
21509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ));
21519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ));
21526214f412SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
21539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_elemental_C",MatConvert_MPISBAIJ_Elemental));
21546214f412SHong Zhang #endif
2155d24d4204SJose E. Roman #if defined(PETSC_HAVE_SCALAPACK)
21569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_scalapack_C",MatConvert_SBAIJ_ScaLAPACK));
2157d24d4204SJose E. Roman #endif
21589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpiaij_C",MatConvert_MPISBAIJ_Basic));
21599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpibaij_C",MatConvert_MPISBAIJ_Basic));
2160aa5a9175SDahai Guo 
216123ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
216223ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
216323ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
216423ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
21659899f194SHong Zhang   B->symmetric_eternal          = PETSC_TRUE;
2166eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
216713647f61SHong Zhang   B->hermitian                  = PETSC_FALSE;
216813647f61SHong Zhang   B->hermitian_set              = PETSC_FALSE;
2169eb1ec7c1SStefano Zampini #else
2170eb1ec7c1SStefano Zampini   B->hermitian                  = PETSC_TRUE;
2171eb1ec7c1SStefano Zampini   B->hermitian_set              = PETSC_TRUE;
2172eb1ec7c1SStefano Zampini #endif
217313647f61SHong Zhang 
21749566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ));
2175d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
21769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",flg,&flg,NULL));
217794ae4db5SBarry Smith   if (flg) {
217894ae4db5SBarry Smith     PetscReal fact = 1.39;
21799566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE));
21809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL));
218194ae4db5SBarry Smith     if (fact <= 1.0) fact = 1.39;
21829566063dSJacob Faibussowitsch     PetscCall(MatMPIBAIJSetHashTableFactor(B,fact));
21839566063dSJacob Faibussowitsch     PetscCall(PetscInfo(B,"Hash table Factor used %5.2g\n",(double)fact));
218494ae4db5SBarry Smith   }
2185d0609cedSBarry Smith   PetscOptionsEnd();
2186b5df2d14SHong Zhang   PetscFunctionReturn(0);
2187b5df2d14SHong Zhang }
2188b5df2d14SHong Zhang 
2189209238afSKris Buschelman /*MC
2190002d173eSKris Buschelman    MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
2191209238afSKris Buschelman 
2192209238afSKris Buschelman    This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
2193209238afSKris Buschelman    and MATMPISBAIJ otherwise.
2194209238afSKris Buschelman 
2195209238afSKris Buschelman    Options Database Keys:
2196209238afSKris Buschelman . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
2197209238afSKris Buschelman 
2198209238afSKris Buschelman   Level: beginner
2199209238afSKris Buschelman 
2200db781477SPatrick Sanan .seealso: `MatCreateSBAIJ`, `MATSEQSBAIJ`, `MATMPISBAIJ`
2201209238afSKris Buschelman M*/
2202209238afSKris Buschelman 
2203b5df2d14SHong Zhang /*@C
2204b5df2d14SHong Zhang    MatMPISBAIJSetPreallocation - For good matrix assembly performance
2205b5df2d14SHong Zhang    the user should preallocate the matrix storage by setting the parameters
2206b5df2d14SHong Zhang    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
2207b5df2d14SHong Zhang    performance can be increased by more than a factor of 50.
2208b5df2d14SHong Zhang 
2209b5df2d14SHong Zhang    Collective on Mat
2210b5df2d14SHong Zhang 
2211b5df2d14SHong Zhang    Input Parameters:
22121c4f3114SJed Brown +  B - the matrix
2213bb7ae925SBarry 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
2214bb7ae925SBarry Smith           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2215b5df2d14SHong Zhang .  d_nz  - number of block nonzeros per block row in diagonal portion of local
2216b5df2d14SHong Zhang            submatrix  (same for all local rows)
2217b5df2d14SHong Zhang .  d_nnz - array containing the number of block nonzeros in the various block rows
22186d10fdaeSSatish Balay            in the upper triangular and diagonal part of the in diagonal portion of the local
22190298fd71SBarry Smith            (possibly different for each block row) or NULL.  If you plan to factor the matrix you must leave room
222095742e49SBarry Smith            for the diagonal entry and set a value even if it is zero.
2221b5df2d14SHong Zhang .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
2222b5df2d14SHong Zhang            submatrix (same for all local rows).
2223b5df2d14SHong Zhang -  o_nnz - array containing the number of nonzeros in the various block rows of the
2224c2fc9fa9SBarry Smith            off-diagonal portion of the local submatrix that is right of the diagonal
22250298fd71SBarry Smith            (possibly different for each block row) or NULL.
2226b5df2d14SHong Zhang 
2227b5df2d14SHong Zhang    Options Database Keys:
2228a2b725a8SWilliam Gropp +   -mat_no_unroll - uses code that does not unroll the loops in the
2229b5df2d14SHong Zhang                      block calculations (much slower)
2230a2b725a8SWilliam Gropp -   -mat_block_size - size of the blocks to use
2231b5df2d14SHong Zhang 
2232b5df2d14SHong Zhang    Notes:
2233b5df2d14SHong Zhang 
2234b5df2d14SHong Zhang    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one processor
2235b5df2d14SHong Zhang    than it must be used on all processors that share the object for that argument.
2236b5df2d14SHong Zhang 
223749a6f317SBarry Smith    If the *_nnz parameter is given then the *_nz parameter is ignored
223849a6f317SBarry Smith 
2239b5df2d14SHong Zhang    Storage Information:
2240b5df2d14SHong Zhang    For a square global matrix we define each processor's diagonal portion
2241b5df2d14SHong Zhang    to be its local rows and the corresponding columns (a square submatrix);
2242b5df2d14SHong Zhang    each processor's off-diagonal portion encompasses the remainder of the
2243b5df2d14SHong Zhang    local matrix (a rectangular submatrix).
2244b5df2d14SHong Zhang 
2245b5df2d14SHong Zhang    The user can specify preallocated storage for the diagonal part of
2246b5df2d14SHong Zhang    the local submatrix with either d_nz or d_nnz (not both).  Set
22470298fd71SBarry Smith    d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2248b5df2d14SHong Zhang    memory allocation.  Likewise, specify preallocated storage for the
2249b5df2d14SHong Zhang    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2250b5df2d14SHong Zhang 
2251aa95bbe8SBarry Smith    You can call MatGetInfo() to get information on how effective the preallocation was;
2252aa95bbe8SBarry Smith    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2253aa95bbe8SBarry Smith    You can also run with the option -info and look for messages with the string
2254aa95bbe8SBarry Smith    malloc in them to see if additional memory allocation was needed.
2255aa95bbe8SBarry Smith 
2256b5df2d14SHong Zhang    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2257b5df2d14SHong Zhang    the figure below we depict these three local rows and all columns (0-11).
2258b5df2d14SHong Zhang 
2259b5df2d14SHong Zhang .vb
2260b5df2d14SHong Zhang            0 1 2 3 4 5 6 7 8 9 10 11
2261a4b1a0f6SJed Brown           --------------------------
2262c2fc9fa9SBarry Smith    row 3  |. . . d d d o o o o  o  o
2263c2fc9fa9SBarry Smith    row 4  |. . . d d d o o o o  o  o
2264c2fc9fa9SBarry Smith    row 5  |. . . d d d o o o o  o  o
2265a4b1a0f6SJed Brown           --------------------------
2266b5df2d14SHong Zhang .ve
2267b5df2d14SHong Zhang 
2268b5df2d14SHong Zhang    Thus, any entries in the d locations are stored in the d (diagonal)
2269b5df2d14SHong Zhang    submatrix, and any entries in the o locations are stored in the
22706d10fdaeSSatish Balay    o (off-diagonal) submatrix.  Note that the d matrix is stored in
22716d10fdaeSSatish Balay    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2272b5df2d14SHong Zhang 
22736d10fdaeSSatish Balay    Now d_nz should indicate the number of block nonzeros per row in the upper triangular
22746d10fdaeSSatish Balay    plus the diagonal part of the d matrix,
2275c2fc9fa9SBarry Smith    and o_nz should indicate the number of block nonzeros per row in the o matrix
2276c2fc9fa9SBarry Smith 
2277b5df2d14SHong Zhang    In general, for PDE problems in which most nonzeros are near the diagonal,
2278b5df2d14SHong Zhang    one expects d_nz >> o_nz.   For large problems you MUST preallocate memory
2279b5df2d14SHong Zhang    or you will get TERRIBLE performance; see the users' manual chapter on
2280b5df2d14SHong Zhang    matrices.
2281b5df2d14SHong Zhang 
2282b5df2d14SHong Zhang    Level: intermediate
2283b5df2d14SHong Zhang 
2284db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateSeqSBAIJ()`, `MatSetValues()`, `MatCreateBAIJ()`, `PetscSplitOwnership()`
2285b5df2d14SHong Zhang @*/
22867087cfbeSBarry Smith PetscErrorCode  MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2287b5df2d14SHong Zhang {
2288b5df2d14SHong Zhang   PetscFunctionBegin;
22896ba663aaSJed Brown   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
22906ba663aaSJed Brown   PetscValidType(B,1);
22916ba663aaSJed Brown   PetscValidLogicalCollectiveInt(B,bs,2);
2292cac4c232SBarry Smith   PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2293b5df2d14SHong Zhang   PetscFunctionReturn(0);
2294b5df2d14SHong Zhang }
2295b5df2d14SHong Zhang 
2296a30f8f8cSSatish Balay /*@C
229769b1f4b7SBarry Smith    MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
2298a30f8f8cSSatish Balay    (block compressed row).  For good matrix assembly performance
2299a30f8f8cSSatish Balay    the user should preallocate the matrix storage by setting the parameters
2300a30f8f8cSSatish Balay    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
2301a30f8f8cSSatish Balay    performance can be increased by more than a factor of 50.
2302a30f8f8cSSatish Balay 
2303d083f849SBarry Smith    Collective
2304a30f8f8cSSatish Balay 
2305a30f8f8cSSatish Balay    Input Parameters:
2306a30f8f8cSSatish Balay +  comm - MPI communicator
2307bb7ae925SBarry 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
2308bb7ae925SBarry Smith           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2309a30f8f8cSSatish Balay .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2310a30f8f8cSSatish Balay            This value should be the same as the local size used in creating the
2311a30f8f8cSSatish Balay            y vector for the matrix-vector product y = Ax.
2312a30f8f8cSSatish Balay .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2313a30f8f8cSSatish Balay            This value should be the same as the local size used in creating the
2314a30f8f8cSSatish Balay            x vector for the matrix-vector product y = Ax.
2315a30f8f8cSSatish Balay .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2316a30f8f8cSSatish Balay .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2317a30f8f8cSSatish Balay .  d_nz  - number of block nonzeros per block row in diagonal portion of local
2318a30f8f8cSSatish Balay            submatrix  (same for all local rows)
2319a30f8f8cSSatish Balay .  d_nnz - array containing the number of block nonzeros in the various block rows
23206d10fdaeSSatish Balay            in the upper triangular portion of the in diagonal portion of the local
23210298fd71SBarry Smith            (possibly different for each block block row) or NULL.
232295742e49SBarry Smith            If you plan to factor the matrix you must leave room for the diagonal entry and
232395742e49SBarry Smith            set its value even if it is zero.
2324a30f8f8cSSatish Balay .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
2325a30f8f8cSSatish Balay            submatrix (same for all local rows).
2326a30f8f8cSSatish Balay -  o_nnz - array containing the number of nonzeros in the various block rows of the
2327a30f8f8cSSatish Balay            off-diagonal portion of the local submatrix (possibly different for
23280298fd71SBarry Smith            each block row) or NULL.
2329a30f8f8cSSatish Balay 
2330a30f8f8cSSatish Balay    Output Parameter:
2331a30f8f8cSSatish Balay .  A - the matrix
2332a30f8f8cSSatish Balay 
2333a30f8f8cSSatish Balay    Options Database Keys:
2334a2b725a8SWilliam Gropp +   -mat_no_unroll - uses code that does not unroll the loops in the
2335a30f8f8cSSatish Balay                      block calculations (much slower)
2336a30f8f8cSSatish Balay .   -mat_block_size - size of the blocks to use
2337a2b725a8SWilliam Gropp -   -mat_mpi - use the parallel matrix data structures even on one processor
2338a30f8f8cSSatish Balay                (defaults to using SeqBAIJ format on one processor)
2339a30f8f8cSSatish Balay 
2340175b88e8SBarry Smith    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2341f6f02116SRichard Tran Mills    MatXXXXSetPreallocation() paradigm instead of this routine directly.
2342175b88e8SBarry Smith    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2343175b88e8SBarry Smith 
2344a30f8f8cSSatish Balay    Notes:
2345d1be2dadSMatthew Knepley    The number of rows and columns must be divisible by blocksize.
23466d6d819aSHong Zhang    This matrix type does not support complex Hermitian operation.
2347d1be2dadSMatthew Knepley 
2348a30f8f8cSSatish Balay    The user MUST specify either the local or global matrix dimensions
2349a30f8f8cSSatish Balay    (possibly both).
2350a30f8f8cSSatish Balay 
2351a30f8f8cSSatish Balay    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one processor
2352a30f8f8cSSatish Balay    than it must be used on all processors that share the object for that argument.
2353a30f8f8cSSatish Balay 
235449a6f317SBarry Smith    If the *_nnz parameter is given then the *_nz parameter is ignored
235549a6f317SBarry Smith 
2356a30f8f8cSSatish Balay    Storage Information:
2357a30f8f8cSSatish Balay    For a square global matrix we define each processor's diagonal portion
2358a30f8f8cSSatish Balay    to be its local rows and the corresponding columns (a square submatrix);
2359a30f8f8cSSatish Balay    each processor's off-diagonal portion encompasses the remainder of the
2360a30f8f8cSSatish Balay    local matrix (a rectangular submatrix).
2361a30f8f8cSSatish Balay 
2362a30f8f8cSSatish Balay    The user can specify preallocated storage for the diagonal part of
2363a30f8f8cSSatish Balay    the local submatrix with either d_nz or d_nnz (not both).  Set
23640298fd71SBarry Smith    d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2365a30f8f8cSSatish Balay    memory allocation.  Likewise, specify preallocated storage for the
2366a30f8f8cSSatish Balay    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2367a30f8f8cSSatish Balay 
2368a30f8f8cSSatish Balay    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2369a30f8f8cSSatish Balay    the figure below we depict these three local rows and all columns (0-11).
2370a30f8f8cSSatish Balay 
2371a30f8f8cSSatish Balay .vb
2372a30f8f8cSSatish Balay            0 1 2 3 4 5 6 7 8 9 10 11
2373a4b1a0f6SJed Brown           --------------------------
2374c2fc9fa9SBarry Smith    row 3  |. . . d d d o o o o  o  o
2375c2fc9fa9SBarry Smith    row 4  |. . . d d d o o o o  o  o
2376c2fc9fa9SBarry Smith    row 5  |. . . d d d o o o o  o  o
2377a4b1a0f6SJed Brown           --------------------------
2378a30f8f8cSSatish Balay .ve
2379a30f8f8cSSatish Balay 
2380a30f8f8cSSatish Balay    Thus, any entries in the d locations are stored in the d (diagonal)
2381a30f8f8cSSatish Balay    submatrix, and any entries in the o locations are stored in the
23826d10fdaeSSatish Balay    o (off-diagonal) submatrix.  Note that the d matrix is stored in
23836d10fdaeSSatish Balay    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2384a30f8f8cSSatish Balay 
23856d10fdaeSSatish Balay    Now d_nz should indicate the number of block nonzeros per row in the upper triangular
23866d10fdaeSSatish Balay    plus the diagonal part of the d matrix,
2387a30f8f8cSSatish Balay    and o_nz should indicate the number of block nonzeros per row in the o matrix.
2388a30f8f8cSSatish Balay    In general, for PDE problems in which most nonzeros are near the diagonal,
2389a30f8f8cSSatish Balay    one expects d_nz >> o_nz.   For large problems you MUST preallocate memory
2390a30f8f8cSSatish Balay    or you will get TERRIBLE performance; see the users' manual chapter on
2391a30f8f8cSSatish Balay    matrices.
2392a30f8f8cSSatish Balay 
2393a30f8f8cSSatish Balay    Level: intermediate
2394a30f8f8cSSatish Balay 
2395db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateSeqSBAIJ()`, `MatSetValues()`, `MatCreateBAIJ()`
2396a30f8f8cSSatish Balay @*/
2397a30f8f8cSSatish Balay 
239869b1f4b7SBarry Smith PetscErrorCode  MatCreateSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
2399a30f8f8cSSatish Balay {
24001302d50aSBarry Smith   PetscMPIInt    size;
2401a30f8f8cSSatish Balay 
2402a30f8f8cSSatish Balay   PetscFunctionBegin;
24039566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,A));
24049566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A,m,n,M,N));
24059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
2406273d9f13SBarry Smith   if (size > 1) {
24079566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A,MATMPISBAIJ));
24089566063dSJacob Faibussowitsch     PetscCall(MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz));
2409273d9f13SBarry Smith   } else {
24109566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A,MATSEQSBAIJ));
24119566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz));
2412273d9f13SBarry Smith   }
2413a30f8f8cSSatish Balay   PetscFunctionReturn(0);
2414a30f8f8cSSatish Balay }
2415a30f8f8cSSatish Balay 
24166849ba73SBarry Smith static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2417a30f8f8cSSatish Balay {
2418a30f8f8cSSatish Balay   Mat            mat;
2419a30f8f8cSSatish Balay   Mat_MPISBAIJ   *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2420d0f46423SBarry Smith   PetscInt       len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2421387bc808SHong Zhang   PetscScalar    *array;
2422a30f8f8cSSatish Balay 
2423a30f8f8cSSatish Balay   PetscFunctionBegin;
2424f4259b30SLisandro Dalcin   *newmat = NULL;
242526fbe8dcSKarl Rupp 
24269566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin),&mat));
24279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N));
24289566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat,((PetscObject)matin)->type_name));
24299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->rmap,&mat->rmap));
24309566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->cmap,&mat->cmap));
2431e1b6402fSHong Zhang 
2432d5f3da31SBarry Smith   mat->factortype   = matin->factortype;
2433273d9f13SBarry Smith   mat->preallocated = PETSC_TRUE;
243482327fa8SHong Zhang   mat->assembled    = PETSC_TRUE;
24357fff6886SHong Zhang   mat->insertmode   = NOT_SET_VALUES;
24367fff6886SHong Zhang 
2437b5df2d14SHong Zhang   a      = (Mat_MPISBAIJ*)mat->data;
2438a30f8f8cSSatish Balay   a->bs2 = oldmat->bs2;
2439a30f8f8cSSatish Balay   a->mbs = oldmat->mbs;
2440a30f8f8cSSatish Balay   a->nbs = oldmat->nbs;
2441a30f8f8cSSatish Balay   a->Mbs = oldmat->Mbs;
2442a30f8f8cSSatish Balay   a->Nbs = oldmat->Nbs;
2443a30f8f8cSSatish Balay 
2444a30f8f8cSSatish Balay   a->size         = oldmat->size;
2445a30f8f8cSSatish Balay   a->rank         = oldmat->rank;
2446a30f8f8cSSatish Balay   a->donotstash   = oldmat->donotstash;
2447a30f8f8cSSatish Balay   a->roworiented  = oldmat->roworiented;
2448f4259b30SLisandro Dalcin   a->rowindices   = NULL;
2449f4259b30SLisandro Dalcin   a->rowvalues    = NULL;
2450a30f8f8cSSatish Balay   a->getrowactive = PETSC_FALSE;
2451f4259b30SLisandro Dalcin   a->barray       = NULL;
2452899cda47SBarry Smith   a->rstartbs     = oldmat->rstartbs;
2453899cda47SBarry Smith   a->rendbs       = oldmat->rendbs;
2454899cda47SBarry Smith   a->cstartbs     = oldmat->cstartbs;
2455899cda47SBarry Smith   a->cendbs       = oldmat->cendbs;
2456a30f8f8cSSatish Balay 
2457a30f8f8cSSatish Balay   /* hash table stuff */
2458f4259b30SLisandro Dalcin   a->ht           = NULL;
2459f4259b30SLisandro Dalcin   a->hd           = NULL;
2460a30f8f8cSSatish Balay   a->ht_size      = 0;
2461a30f8f8cSSatish Balay   a->ht_flag      = oldmat->ht_flag;
2462a30f8f8cSSatish Balay   a->ht_fact      = oldmat->ht_fact;
2463a30f8f8cSSatish Balay   a->ht_total_ct  = 0;
2464a30f8f8cSSatish Balay   a->ht_insert_ct = 0;
2465a30f8f8cSSatish Balay 
24669566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->rangebs,oldmat->rangebs,a->size+2));
2467a30f8f8cSSatish Balay   if (oldmat->colmap) {
2468a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
24699566063dSJacob Faibussowitsch     PetscCall(PetscTableCreateCopy(oldmat->colmap,&a->colmap));
2470a30f8f8cSSatish Balay #else
24719566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(a->Nbs,&a->colmap));
24729566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt)));
24739566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->colmap,oldmat->colmap,a->Nbs));
2474a30f8f8cSSatish Balay #endif
2475f4259b30SLisandro Dalcin   } else a->colmap = NULL;
2476387bc808SHong Zhang 
2477a30f8f8cSSatish Balay   if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
24789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len,&a->garray));
24799566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt)));
24809566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->garray,oldmat->garray,len));
2481f4259b30SLisandro Dalcin   } else a->garray = NULL;
2482a30f8f8cSSatish Balay 
24839566063dSJacob Faibussowitsch   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash));
24849566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->lvec,&a->lvec));
24859566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec));
24869566063dSJacob Faibussowitsch   PetscCall(VecScatterCopy(oldmat->Mvctx,&a->Mvctx));
24879566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx));
248882327fa8SHong Zhang 
24899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->slvec0,&a->slvec0));
24909566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0));
24919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->slvec1,&a->slvec1));
24929566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1));
2493387bc808SHong Zhang 
24949566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(a->slvec1,&nt));
24959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a->slvec1,&array));
24969566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a));
24979566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b));
24989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a->slvec1,&array));
24999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(a->slvec0,&array));
25009566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b));
25019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(a->slvec0,&array));
25029566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0));
25039566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1));
25049566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0b));
25059566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1a));
25069566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1b));
2507387bc808SHong Zhang 
2508387bc808SHong Zhang   /* ierr =  VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
25099566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)oldmat->sMvctx));
2510387bc808SHong Zhang   a->sMvctx = oldmat->sMvctx;
25119566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->sMvctx));
251282327fa8SHong Zhang 
25139566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->A,cpvalues,&a->A));
25149566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A));
25159566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->B,cpvalues,&a->B));
25169566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B));
25179566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist));
2518a30f8f8cSSatish Balay   *newmat = mat;
2519a30f8f8cSSatish Balay   PetscFunctionReturn(0);
2520a30f8f8cSSatish Balay }
2521a30f8f8cSSatish Balay 
2522618cc2edSLisandro Dalcin /* Used for both MPIBAIJ and MPISBAIJ matrices */
2523618cc2edSLisandro Dalcin #define MatLoad_MPISBAIJ_Binary MatLoad_MPIBAIJ_Binary
2524618cc2edSLisandro Dalcin 
2525618cc2edSLisandro Dalcin PetscErrorCode MatLoad_MPISBAIJ(Mat mat,PetscViewer viewer)
252695936485SShri Abhyankar {
25277f489da9SVaclav Hapla   PetscBool      isbinary;
252895936485SShri Abhyankar 
252995936485SShri Abhyankar   PetscFunctionBegin;
25309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
25315f80ce2aSJacob Faibussowitsch   PetscCheck(isbinary,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)mat)->type_name);
25329566063dSJacob Faibussowitsch   PetscCall(MatLoad_MPISBAIJ_Binary(mat,viewer));
253395936485SShri Abhyankar   PetscFunctionReturn(0);
253495936485SShri Abhyankar }
253595936485SShri Abhyankar 
2536dcf5cc72SBarry Smith /*XXXXX@
2537a30f8f8cSSatish Balay    MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2538a30f8f8cSSatish Balay 
2539a30f8f8cSSatish Balay    Input Parameters:
2540a30f8f8cSSatish Balay .  mat  - the matrix
2541a30f8f8cSSatish Balay .  fact - factor
2542a30f8f8cSSatish Balay 
2543c5eb9154SBarry Smith    Not Collective on Mat, each process can have a different hash factor
2544a30f8f8cSSatish Balay 
2545a30f8f8cSSatish Balay    Level: advanced
2546a30f8f8cSSatish Balay 
2547a30f8f8cSSatish Balay   Notes:
2548a30f8f8cSSatish Balay    This can also be set by the command line option: -mat_use_hash_table fact
2549a30f8f8cSSatish Balay 
2550db781477SPatrick Sanan .seealso: `MatSetOption()`
2551dcf5cc72SBarry Smith @XXXXX*/
2552dcf5cc72SBarry Smith 
2553985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
255424d5174aSHong Zhang {
255524d5174aSHong Zhang   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
2556f4c0e9e4SHong Zhang   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(a->B)->data;
2557ca54ac64SHong Zhang   PetscReal      atmp;
255887828ca2SBarry Smith   PetscReal      *work,*svalues,*rvalues;
25591302d50aSBarry Smith   PetscInt       i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
25601302d50aSBarry Smith   PetscMPIInt    rank,size;
25611302d50aSBarry Smith   PetscInt       *rowners_bs,dest,count,source;
256287828ca2SBarry Smith   PetscScalar    *va;
25638a1c53f2SBarry Smith   MatScalar      *ba;
2564f4c0e9e4SHong Zhang   MPI_Status     stat;
256524d5174aSHong Zhang 
256624d5174aSHong Zhang   PetscFunctionBegin;
25675f80ce2aSJacob Faibussowitsch   PetscCheck(!idx,PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
25689566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(a->A,v,NULL));
25699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v,&va));
2570f4c0e9e4SHong Zhang 
25719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size));
25729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank));
2573f4c0e9e4SHong Zhang 
2574d0f46423SBarry Smith   bs  = A->rmap->bs;
2575f4c0e9e4SHong Zhang   mbs = a->mbs;
2576f4c0e9e4SHong Zhang   Mbs = a->Mbs;
2577f4c0e9e4SHong Zhang   ba  = b->a;
2578f4c0e9e4SHong Zhang   bi  = b->i;
2579f4c0e9e4SHong Zhang   bj  = b->j;
2580f4c0e9e4SHong Zhang 
2581f4c0e9e4SHong Zhang   /* find ownerships */
2582d0f46423SBarry Smith   rowners_bs = A->rmap->range;
2583f4c0e9e4SHong Zhang 
2584f4c0e9e4SHong Zhang   /* each proc creates an array to be distributed */
25859566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(bs*Mbs,&work));
2586f4c0e9e4SHong Zhang 
2587f4c0e9e4SHong Zhang   /* row_max for B */
2588b8475685SHong Zhang   if (rank != size-1) {
2589f4c0e9e4SHong Zhang     for (i=0; i<mbs; i++) {
2590f4c0e9e4SHong Zhang       ncols = bi[1] - bi[0]; bi++;
2591f4c0e9e4SHong Zhang       brow  = bs*i;
2592f4c0e9e4SHong Zhang       for (j=0; j<ncols; j++) {
2593f4c0e9e4SHong Zhang         bcol = bs*(*bj);
2594f4c0e9e4SHong Zhang         for (kcol=0; kcol<bs; kcol++) {
2595ca54ac64SHong Zhang           col  = bcol + kcol;                /* local col index */
259604d41228SHong Zhang           col += rowners_bs[rank+1];      /* global col index */
2597f4c0e9e4SHong Zhang           for (krow=0; krow<bs; krow++) {
2598f4c0e9e4SHong Zhang             atmp = PetscAbsScalar(*ba); ba++;
2599ca54ac64SHong Zhang             row  = brow + krow;   /* local row index */
2600ca54ac64SHong Zhang             if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2601f4c0e9e4SHong Zhang             if (work[col] < atmp) work[col] = atmp;
2602f4c0e9e4SHong Zhang           }
2603f4c0e9e4SHong Zhang         }
2604f4c0e9e4SHong Zhang         bj++;
2605f4c0e9e4SHong Zhang       }
2606f4c0e9e4SHong Zhang     }
2607f4c0e9e4SHong Zhang 
2608f4c0e9e4SHong Zhang     /* send values to its owners */
2609f4c0e9e4SHong Zhang     for (dest=rank+1; dest<size; dest++) {
2610f4c0e9e4SHong Zhang       svalues = work + rowners_bs[dest];
2611ca54ac64SHong Zhang       count   = rowners_bs[dest+1]-rowners_bs[dest];
26129566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A)));
2613ca54ac64SHong Zhang     }
2614f4c0e9e4SHong Zhang   }
2615f4c0e9e4SHong Zhang 
2616f4c0e9e4SHong Zhang   /* receive values */
2617ca54ac64SHong Zhang   if (rank) {
2618f4c0e9e4SHong Zhang     rvalues = work;
2619ca54ac64SHong Zhang     count   = rowners_bs[rank+1]-rowners_bs[rank];
2620f4c0e9e4SHong Zhang     for (source=0; source<rank; source++) {
26219566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat));
2622f4c0e9e4SHong Zhang       /* process values */
2623f4c0e9e4SHong Zhang       for (i=0; i<count; i++) {
2624ca54ac64SHong Zhang         if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2625f4c0e9e4SHong Zhang       }
2626f4c0e9e4SHong Zhang     }
2627ca54ac64SHong Zhang   }
2628f4c0e9e4SHong Zhang 
26299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v,&va));
26309566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
263124d5174aSHong Zhang   PetscFunctionReturn(0);
263224d5174aSHong Zhang }
26332798e883SHong Zhang 
263441f059aeSBarry Smith PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
26352798e883SHong Zhang {
26362798e883SHong Zhang   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)matin->data;
2637d0f46423SBarry Smith   PetscInt          mbs=mat->mbs,bs=matin->rmap->bs;
26383649974fSBarry Smith   PetscScalar       *x,*ptr,*from;
2639ffe4fb16SHong Zhang   Vec               bb1;
26403649974fSBarry Smith   const PetscScalar *b;
2641ffe4fb16SHong Zhang 
2642ffe4fb16SHong Zhang   PetscFunctionBegin;
26435f80ce2aSJacob Faibussowitsch   PetscCheck(its > 0 && lits > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits);
26445f80ce2aSJacob Faibussowitsch   PetscCheck(bs <= 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2645ffe4fb16SHong Zhang 
2646a2b30743SBarry Smith   if (flag == SOR_APPLY_UPPER) {
26479566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
2648a2b30743SBarry Smith     PetscFunctionReturn(0);
2649a2b30743SBarry Smith   }
2650a2b30743SBarry Smith 
2651ffe4fb16SHong Zhang   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2652ffe4fb16SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
26539566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx));
2654ffe4fb16SHong Zhang       its--;
2655ffe4fb16SHong Zhang     }
2656ffe4fb16SHong Zhang 
26579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(bb,&bb1));
2658ffe4fb16SHong Zhang     while (its--) {
2659ffe4fb16SHong Zhang 
2660ffe4fb16SHong Zhang       /* lower triangular part: slvec0b = - B^T*xx */
26619566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b));
2662ffe4fb16SHong Zhang 
2663ffe4fb16SHong Zhang       /* copy xx into slvec0a */
26649566063dSJacob Faibussowitsch       PetscCall(VecGetArray(mat->slvec0,&ptr));
26659566063dSJacob Faibussowitsch       PetscCall(VecGetArray(xx,&x));
26669566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ptr,x,bs*mbs));
26679566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mat->slvec0,&ptr));
2668ffe4fb16SHong Zhang 
26699566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->slvec0,-1.0));
2670ffe4fb16SHong Zhang 
2671ffe4fb16SHong Zhang       /* copy bb into slvec1a */
26729566063dSJacob Faibussowitsch       PetscCall(VecGetArray(mat->slvec1,&ptr));
26739566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(bb,&b));
26749566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(ptr,b,bs*mbs));
26759566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mat->slvec1,&ptr));
2676ffe4fb16SHong Zhang 
2677ffe4fb16SHong Zhang       /* set slvec1b = 0 */
26789566063dSJacob Faibussowitsch       PetscCall(VecSet(mat->slvec1b,0.0));
2679ffe4fb16SHong Zhang 
26809566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD));
26819566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(xx,&x));
26829566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(bb,&b));
26839566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD));
2684ffe4fb16SHong Zhang 
2685ffe4fb16SHong Zhang       /* upper triangular part: bb1 = bb1 - B*x */
26869566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1));
2687ffe4fb16SHong Zhang 
2688ffe4fb16SHong Zhang       /* local diagonal sweep */
26899566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx));
2690ffe4fb16SHong Zhang     }
26919566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bb1));
2692fa22f6d0SBarry Smith   } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
26939566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
2694fa22f6d0SBarry Smith   } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
26959566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx));
2696fa22f6d0SBarry Smith   } else if (flag & SOR_EISENSTAT) {
2697fa22f6d0SBarry Smith     Vec               xx1;
2698ace3abfcSBarry Smith     PetscBool         hasop;
269920f1ed55SBarry Smith     const PetscScalar *diag;
2700887ee2caSBarry Smith     PetscScalar       *sl,scale = (omega - 2.0)/omega;
270120f1ed55SBarry Smith     PetscInt          i,n;
2702fa22f6d0SBarry Smith 
2703fa22f6d0SBarry Smith     if (!mat->xx1) {
27049566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(bb,&mat->xx1));
27059566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(bb,&mat->bb1));
2706fa22f6d0SBarry Smith     }
2707fa22f6d0SBarry Smith     xx1 = mat->xx1;
2708fa22f6d0SBarry Smith     bb1 = mat->bb1;
2709fa22f6d0SBarry Smith 
27109566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx));
2711fa22f6d0SBarry Smith 
2712fa22f6d0SBarry Smith     if (!mat->diag) {
2713effcda25SBarry Smith       /* this is wrong for same matrix with new nonzero values */
27149566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(matin,&mat->diag,NULL));
27159566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(matin,mat->diag));
2716fa22f6d0SBarry Smith     }
27179566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop));
2718fa22f6d0SBarry Smith 
2719fa22f6d0SBarry Smith     if (hasop) {
27209566063dSJacob Faibussowitsch       PetscCall(MatMultDiagonalBlock(matin,xx,bb1));
27219566063dSJacob Faibussowitsch       PetscCall(VecAYPX(mat->slvec1a,scale,bb));
272220f1ed55SBarry Smith     } else {
272320f1ed55SBarry Smith       /*
272420f1ed55SBarry Smith           These two lines are replaced by code that may be a bit faster for a good compiler
27259566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMult(mat->slvec1a,mat->diag,xx));
27269566063dSJacob Faibussowitsch       PetscCall(VecAYPX(mat->slvec1a,scale,bb));
272720f1ed55SBarry Smith       */
27289566063dSJacob Faibussowitsch       PetscCall(VecGetArray(mat->slvec1a,&sl));
27299566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(mat->diag,&diag));
27309566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(bb,&b));
27319566063dSJacob Faibussowitsch       PetscCall(VecGetArray(xx,&x));
27329566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xx,&n));
2733887ee2caSBarry Smith       if (omega == 1.0) {
273426fbe8dcSKarl Rupp         for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
27359566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0*n));
2736887ee2caSBarry Smith       } else {
273726fbe8dcSKarl Rupp         for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
27389566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(3.0*n));
2739887ee2caSBarry Smith       }
27409566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(mat->slvec1a,&sl));
27419566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(mat->diag,&diag));
27429566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(bb,&b));
27439566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(xx,&x));
274420f1ed55SBarry Smith     }
2745fa22f6d0SBarry Smith 
2746fa22f6d0SBarry Smith     /* multiply off-diagonal portion of matrix */
27479566063dSJacob Faibussowitsch     PetscCall(VecSet(mat->slvec1b,0.0));
27489566063dSJacob Faibussowitsch     PetscCall((*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b));
27499566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mat->slvec0,&from));
27509566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xx,&x));
27519566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(from,x,bs*mbs));
27529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(mat->slvec0,&from));
27539566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xx,&x));
27549566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD));
27559566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD));
27569566063dSJacob Faibussowitsch     PetscCall((*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a));
2757fa22f6d0SBarry Smith 
2758fa22f6d0SBarry Smith     /* local sweep */
27599566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1));
27609566063dSJacob Faibussowitsch     PetscCall(VecAXPY(xx,1.0,xx1));
2761f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2762ffe4fb16SHong Zhang   PetscFunctionReturn(0);
2763ffe4fb16SHong Zhang }
2764ffe4fb16SHong Zhang 
2765dfb205c3SBarry Smith /*@
2766dfb205c3SBarry Smith      MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
2767dfb205c3SBarry Smith          CSR format the local rows.
2768dfb205c3SBarry Smith 
2769d083f849SBarry Smith    Collective
2770dfb205c3SBarry Smith 
2771dfb205c3SBarry Smith    Input Parameters:
2772dfb205c3SBarry Smith +  comm - MPI communicator
2773dfb205c3SBarry Smith .  bs - the block size, only a block size of 1 is supported
2774dfb205c3SBarry Smith .  m - number of local rows (Cannot be PETSC_DECIDE)
2775dfb205c3SBarry Smith .  n - This value should be the same as the local size used in creating the
2776dfb205c3SBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2777dfb205c3SBarry Smith        calculated if N is given) For square matrices n is almost always m.
2778dfb205c3SBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2779dfb205c3SBarry Smith .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2780483a2f95SBarry 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
2781dfb205c3SBarry Smith .   j - column indices
2782dfb205c3SBarry Smith -   a - matrix values
2783dfb205c3SBarry Smith 
2784dfb205c3SBarry Smith    Output Parameter:
2785dfb205c3SBarry Smith .   mat - the matrix
2786dfb205c3SBarry Smith 
2787dfb205c3SBarry Smith    Level: intermediate
2788dfb205c3SBarry Smith 
2789dfb205c3SBarry Smith    Notes:
2790dfb205c3SBarry Smith        The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
2791dfb205c3SBarry Smith      thus you CANNOT change the matrix entries by changing the values of a[] after you have
2792dfb205c3SBarry Smith      called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
2793dfb205c3SBarry Smith 
2794dfb205c3SBarry Smith        The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
2795dfb205c3SBarry Smith 
2796db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatMPIAIJSetPreallocation()`, `MatMPIAIJSetPreallocationCSR()`,
2797db781477SPatrick Sanan           `MPIAIJ`, `MatCreateAIJ()`, `MatCreateMPIAIJWithSplitArrays()`
2798dfb205c3SBarry Smith @*/
27997087cfbeSBarry Smith PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],Mat *mat)
2800dfb205c3SBarry Smith {
2801dfb205c3SBarry Smith   PetscFunctionBegin;
28025f80ce2aSJacob Faibussowitsch   PetscCheck(!i[0],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
28035f80ce2aSJacob Faibussowitsch   PetscCheck(m >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
28049566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,mat));
28059566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*mat,m,n,M,N));
28069566063dSJacob Faibussowitsch   PetscCall(MatSetType(*mat,MATMPISBAIJ));
28079566063dSJacob Faibussowitsch   PetscCall(MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a));
2808dfb205c3SBarry Smith   PetscFunctionReturn(0);
2809dfb205c3SBarry Smith }
2810dfb205c3SBarry Smith 
2811dfb205c3SBarry Smith /*@C
2812664954b6SBarry Smith    MatMPISBAIJSetPreallocationCSR - Creates a sparse parallel matrix in SBAIJ format using the given nonzero structure and (optional) numerical values
2813dfb205c3SBarry Smith 
2814d083f849SBarry Smith    Collective
2815dfb205c3SBarry Smith 
2816dfb205c3SBarry Smith    Input Parameters:
28171c4f3114SJed Brown +  B - the matrix
2818dfb205c3SBarry Smith .  bs - the block size
2819dfb205c3SBarry Smith .  i - the indices into j for the start of each local row (starts with zero)
2820dfb205c3SBarry Smith .  j - the column indices for each local row (starts with zero) these must be sorted for each row
2821dfb205c3SBarry Smith -  v - optional values in the matrix
2822dfb205c3SBarry Smith 
2823664954b6SBarry Smith    Level: advanced
2824664954b6SBarry Smith 
2825664954b6SBarry Smith    Notes:
28260cd7f59aSBarry Smith    Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
28270cd7f59aSBarry Smith    and usually the numerical values as well
28280cd7f59aSBarry Smith 
282950c5228eSBarry Smith    Any entries below the diagonal are ignored
2830dfb205c3SBarry Smith 
2831db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateSeqAIJ()`, `MatSetValues()`, `MatMPIBAIJSetPreallocation()`, `MatCreateAIJ()`, `MPIAIJ`
2832dfb205c3SBarry Smith @*/
28337087cfbeSBarry Smith PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2834dfb205c3SBarry Smith {
2835dfb205c3SBarry Smith   PetscFunctionBegin;
2836cac4c232SBarry Smith   PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2837dfb205c3SBarry Smith   PetscFunctionReturn(0);
2838dfb205c3SBarry Smith }
2839dfb205c3SBarry Smith 
284010c56fdeSHong Zhang PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
28414dcd73b1SHong Zhang {
284210c56fdeSHong Zhang   PetscInt       m,N,i,rstart,nnz,Ii,bs,cbs;
284310c56fdeSHong Zhang   PetscInt       *indx;
284410c56fdeSHong Zhang   PetscScalar    *values;
2845dfb205c3SBarry Smith 
28464dcd73b1SHong Zhang   PetscFunctionBegin;
28479566063dSJacob Faibussowitsch   PetscCall(MatGetSize(inmat,&m,&N));
284810c56fdeSHong Zhang   if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
284910c56fdeSHong Zhang     Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inmat->data;
2850de25e9cbSPierre Jolivet     PetscInt       *dnz,*onz,mbs,Nbs,nbs;
285110c56fdeSHong Zhang     PetscInt       *bindx,rmax=a->rmax,j;
2852de25e9cbSPierre Jolivet     PetscMPIInt    rank,size;
28534dcd73b1SHong Zhang 
28549566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSizes(inmat,&bs,&cbs));
285510c56fdeSHong Zhang     mbs = m/bs; Nbs = N/cbs;
285610c56fdeSHong Zhang     if (n == PETSC_DECIDE) {
28579566063dSJacob Faibussowitsch       PetscCall(PetscSplitOwnershipBlock(comm,cbs,&n,&N));
285810c56fdeSHong Zhang     }
2859da91a574SPierre Jolivet     nbs = n/cbs;
28604dcd73b1SHong Zhang 
28619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(rmax,&bindx));
2862d0609cedSBarry Smith     MatPreallocateBegin(comm,mbs,nbs,dnz,onz); /* inline function, output __end and __rstart are used below */
2863de25e9cbSPierre Jolivet 
28649566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm,&rank));
28659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm,&size));
2866de25e9cbSPierre Jolivet     if (rank == size-1) {
2867de25e9cbSPierre Jolivet       /* Check sum(nbs) = Nbs */
28685f80ce2aSJacob Faibussowitsch       PetscCheck(__end == Nbs,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local block columns %" PetscInt_FMT " != global block columns %" PetscInt_FMT,__end,Nbs);
2869de25e9cbSPierre Jolivet     }
2870de25e9cbSPierre Jolivet 
2871d0609cedSBarry Smith     rstart = __rstart; /* block rstart of *outmat; see inline function MatPreallocateBegin */
28729566063dSJacob Faibussowitsch     PetscCall(MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE));
287310c56fdeSHong Zhang     for (i=0; i<mbs; i++) {
28749566063dSJacob Faibussowitsch       PetscCall(MatGetRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL)); /* non-blocked nnz and indx */
28754dcd73b1SHong Zhang       nnz  = nnz/bs;
28764dcd73b1SHong Zhang       for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
28779566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz));
28789566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL));
28794dcd73b1SHong Zhang     }
28809566063dSJacob Faibussowitsch     PetscCall(MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE));
28819566063dSJacob Faibussowitsch     PetscCall(PetscFree(bindx));
28824dcd73b1SHong Zhang 
28839566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm,outmat));
28849566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(*outmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE));
28859566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(*outmat,bs,cbs));
28869566063dSJacob Faibussowitsch     PetscCall(MatSetType(*outmat,MATSBAIJ));
28879566063dSJacob Faibussowitsch     PetscCall(MatSeqSBAIJSetPreallocation(*outmat,bs,0,dnz));
28889566063dSJacob Faibussowitsch     PetscCall(MatMPISBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz));
2889d0609cedSBarry Smith     MatPreallocateEnd(dnz,onz);
28904dcd73b1SHong Zhang   }
28914dcd73b1SHong Zhang 
289210c56fdeSHong Zhang   /* numeric phase */
28939566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(inmat,&bs,&cbs));
28949566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*outmat,&rstart,NULL));
28954dcd73b1SHong Zhang 
28969566063dSJacob Faibussowitsch   PetscCall(MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE));
28974dcd73b1SHong Zhang   for (i=0; i<m; i++) {
28989566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values));
28994dcd73b1SHong Zhang     Ii   = i + rstart;
29009566063dSJacob Faibussowitsch     PetscCall(MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES));
29019566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values));
29024dcd73b1SHong Zhang   }
29039566063dSJacob Faibussowitsch   PetscCall(MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE));
29049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY));
29059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY));
29064dcd73b1SHong Zhang   PetscFunctionReturn(0);
29074dcd73b1SHong Zhang }
2908